View Javadoc
1   // ******************************************************************************
2   //
3   // Title:       Force Field X.
4   // Description: Force Field X - Software for Molecular Biophysics.
5   // Copyright:   Copyright (c) Michael J. Schnieders 2001-2025.
6   //
7   // This file is part of Force Field X.
8   //
9   // Force Field X is free software; you can redistribute it and/or modify it
10  // under the terms of the GNU General Public License version 3 as published by
11  // the Free Software Foundation.
12  //
13  // Force Field X is distributed in the hope that it will be useful, but WITHOUT
14  // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
15  // FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
16  // details.
17  //
18  // You should have received a copy of the GNU General Public License along with
19  // Force Field X; if not, write to the Free Software Foundation, Inc., 59 Temple
20  // Place, Suite 330, Boston, MA 02111-1307 USA
21  //
22  // Linking this library statically or dynamically with other modules is making a
23  // combined work based on this library. Thus, the terms and conditions of the
24  // GNU General Public License cover the whole combination.
25  //
26  // As a special exception, the copyright holders of this library give you
27  // permission to link this library with independent modules to produce an
28  // executable, regardless of the license terms of these independent modules, and
29  // to copy and distribute the resulting executable under terms of your choice,
30  // provided that you also meet, for each linked independent module, the terms
31  // and conditions of the license of that module. An independent module is a
32  // module which is not derived from or based on this library. If you modify this
33  // library, you may extend this exception to your version of the library, but
34  // you are not obligated to do so. If you do not wish to do so, delete this
35  // exception statement from your version.
36  //
37  // ******************************************************************************
38  package ffx.numerics.multipole;
39  
40  import jdk.incubator.vector.DoubleVector;
41  import jdk.incubator.vector.VectorOperators;
42  
43  import static ffx.numerics.multipole.EwaldTensorGlobal.initEwaldSource;
44  import static ffx.numerics.special.Erf.erfc;
45  import static java.lang.Math.PI;
46  import static org.apache.commons.math3.util.FastMath.sqrt;
47  
48  /**
49   * The EwaldMultipoleTensorGlobal class computes derivatives of erfc(<b>r</b>)/|<b>r</b>| via
50   * recursion to arbitrary order for Cartesian multipoles in the global frame.
51   *
52   * @author Michael J. Schnieders
53   * @see <a href="http://doi.org/10.1142/9789812830364_0002" target="_blank"> Matt Challacombe, Eric
54   * Schwegler and Jan Almlof, Modern developments in Hartree-Fock theory: Fast methods for
55   * computing the Coulomb matrix. Computational Chemistry: Review of Current Trends. pp. 53-107,
56   * Ed. J. Leczszynski, World Scientifc, 1996. </a>
57   * @since 1.0
58   */
59  public class EwaldTensorGlobalSIMD extends CoulombTensorGlobalSIMD {
60  
61    /**
62     * Constant <code>sqrtPI = sqrt(PI)</code>
63     */
64    private static final double sqrtPI = sqrt(PI);
65  
66    /**
67     * These are the "source" terms for the recursion for the screened Coulomb operator erfc(R)/R.
68     */
69    private final double[] ewaldSource;
70  
71    /**
72     * A work array for generation of source terms that cannot be vectorized (exp and erfc).
73     */
74    private final double[] work;
75  
76    /**
77     * The Ewald convergence parameter.
78     */
79    private final double beta;
80  
81    /**
82     * Constructor for EwaldMultipoleTensorGlobal.
83     *
84     * @param order Tensor order.
85     * @param beta  The Ewald convergence parameter.
86     */
87    public EwaldTensorGlobalSIMD(int order, double beta) {
88      super(order);
89      this.beta = beta;
90      operator = Operator.SCREENED_COULOMB;
91  
92      // Auxiliary terms for screened Coulomb (Sagui et al. Eq. 2.28)
93      ewaldSource = new double[o1];
94      work = new double[o1];
95      initEwaldSource(order, beta, ewaldSource);
96    }
97  
98    /**
99     * Generate source terms for the Ewald Challacombe et al. recursion.
100    *
101    * @param T000 Location to store the source terms.
102    */
103   @Override
104   protected void source(DoubleVector[] T000) {
105     // Generate source terms for real space Ewald summation.
106     if (beta > 0.0) {
107       fillEwaldSource(order, beta, ewaldSource, R, T000, work);
108     } else {
109       // For beta = 0, generate tensors for the Coulomb operator.
110       super.source(T000);
111     }
112   }
113 
114   /**
115    * Fill the Ewald source terms.
116    *
117    * @param order       The order plus one.
118    * @param beta        The Ewald convergence parameter.
119    * @param ewaldSource The source terms.
120    * @param R           The separation distance.
121    * @param T000        The location to store the source terms.
122    * @param work        A work array for generation of source terms that cannot be vectorized.
123    */
124   protected static void fillEwaldSource(int order, double beta, double[] ewaldSource,
125                                         DoubleVector R, DoubleVector[] T000, double[] work) {
126     // Sagui et al. Eq. 2.22
127     DoubleVector betaR = R.mul(beta);
128     DoubleVector betaR2 = betaR.mul(betaR);
129     DoubleVector iBetaR2 = DoubleVector.broadcast(R.species(), 1.0);
130     iBetaR2 = iBetaR2.div(betaR2.mul(2.0));
131     DoubleVector expBR2 = betaR2.neg().lanewise(VectorOperators.EXP);
132     // Fnc(x^2) = Sqrt(PI) * erfc(x) / (2*x)
133     // where x = Beta*R
134     // Serial portion to handle the erfc.
135     betaR.intoArray(work, 0);
136     for (int i = 0; i < R.length(); i++) {
137       work[i] = erfc(work[i]);
138     }
139     DoubleVector Fnc = DoubleVector.fromArray(R.species(), work, 0);
140 
141     Fnc = Fnc.mul(sqrtPI).div(betaR.mul(2.0));
142     for (int n = 0; n <= order; n++) {
143       T000[n] = Fnc.mul(ewaldSource[n]);
144       // Generate F(n+1)c from Fnc (Eq. 2.24 in Sagui et al.)
145       // F(n+1)c = [(2*n+1) Fnc(x) + exp(-x)] / 2x
146       // where x = (Beta*R)^2
147       Fnc = ((Fnc.mul(2.0 * n + 1.0)).add(expBR2)).mul(iBetaR2);
148     }
149   }
150 
151 }