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 static ffx.numerics.special.Erf.erfc;
41  import static java.lang.Math.PI;
42  import static org.apache.commons.math3.util.FastMath.exp;
43  import static org.apache.commons.math3.util.FastMath.pow;
44  import static org.apache.commons.math3.util.FastMath.sqrt;
45  
46  /**
47   * The EwaldMultipoleTensorGlobal class computes derivatives of erfc(<b>r</b>)/|<b>r</b>| via
48   * recursion to arbitrary order for Cartesian multipoles in the global frame.
49   *
50   * @author Michael J. Schnieders
51   * @see <a href="http://doi.org/10.1142/9789812830364_0002" target="_blank"> Matt Challacombe, Eric
52   * Schwegler and Jan Almlof, Modern developments in Hartree-Fock theory: Fast methods for
53   * computing the Coulomb matrix. Computational Chemistry: Review of Current Trends. pp. 53-107,
54   * Ed. J. Leczszynski, World Scientifc, 1996. </a>
55   * @since 1.0
56   */
57  public class EwaldTensorGlobal extends CoulombTensorGlobal {
58  
59    /**
60     * Constant <code>sqrtPI = sqrt(PI)</code>
61     */
62    private static final double sqrtPI = sqrt(PI);
63  
64    /**
65     * These are the "source" terms for the recursion for the screened Coulomb operator erfc(R)/R.
66     */
67    private final double[] ewaldSource;
68  
69    /**
70     * The Ewald convergence parameter.
71     */
72    private final double beta;
73  
74    /**
75     * Constructor for EwaldMultipoleTensorGlobal.
76     *
77     * @param order Tensor order.
78     * @param beta  The Ewald convergence parameter.
79     */
80    public EwaldTensorGlobal(int order, double beta) {
81      super(order);
82      this.beta = beta;
83      operator = Operator.SCREENED_COULOMB;
84  
85      // Auxiliary terms for screened Coulomb (Sagui et al. Eq. 2.28)
86      ewaldSource = new double[o1];
87      initEwaldSource(order, beta, ewaldSource);
88    }
89  
90    /**
91     * Initialize the Ewald source terms.
92     *
93     * @param order       Tensor order.
94     * @param beta        The Ewald convergence parameter.
95     * @param ewaldSource Location to store the source terms.
96     * @return The source terms.
97     */
98    protected static double[] initEwaldSource(int order, double beta, double[] ewaldSource) {
99      double prefactor = 2.0 * beta / sqrtPI;
100     double twoBeta2 = -2.0 * beta * beta;
101     for (int n = 0; n <= order; n++) {
102       ewaldSource[n] = prefactor * pow(twoBeta2, n);
103     }
104     return ewaldSource;
105   }
106 
107   /**
108    * Generate source terms for the Ewald Challacombe et al. recursion.
109    *
110    * @param T000 Location to store the source terms.
111    */
112   @Override
113   protected void source(double[] T000) {
114     // Generate source terms for real space Ewald summation.
115     if (beta > 0.0) {
116       fillEwaldSource(order, beta, ewaldSource, R, T000);
117     } else {
118       // For beta = 0, generate tensors for the Coulomb operator.
119       super.source(T000);
120     }
121   }
122 
123   /**
124    * Fill the Ewald source terms.
125    *
126    * @param order       The order plus one.
127    * @param beta        The Ewald convergence parameter.
128    * @param ewaldSource The source terms.
129    * @param R           The separation distance.
130    * @param T000        The location to store the source terms.
131    */
132   protected static void fillEwaldSource(int order, double beta, double[] ewaldSource, double R, double[] T000) {
133     // Sagui et al. Eq. 2.22
134     double betaR = beta * R;
135     double betaR2 = betaR * betaR;
136     double iBetaR2 = 1.0 / (2.0 * betaR2);
137     double expBR2 = exp(-betaR2);
138     // Fnc(x^2) = Sqrt(PI) * erfc(x) / (2*x)
139     // where x = Beta*R
140     double Fnc = sqrtPI * erfc(betaR) / (2.0 * betaR);
141     for (int n = 0; n <= order; n++) {
142       T000[n] = ewaldSource[n] * Fnc;
143       // Generate F(n+1)c from Fnc (Eq. 2.24 in Sagui et al.)
144       // F(n+1)c = [(2*n+1) Fnc(x) + exp(-x)] / 2x
145       // where x = (Beta*R)^2
146       Fnc = ((2.0 * n + 1.0) * Fnc + expBR2) * iBetaR2;
147     }
148   }
149 
150 }