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.potential.nonbonded;
39  
40  import ffx.potential.nonbonded.pme.LambdaMode;
41  import ffx.potential.parameters.ForceField;
42  import org.apache.commons.math3.fitting.leastsquares.EvaluationRmsChecker;
43  import org.apache.commons.math3.fitting.leastsquares.LeastSquaresFactory;
44  import org.apache.commons.math3.fitting.leastsquares.LeastSquaresOptimizer;
45  import org.apache.commons.math3.fitting.leastsquares.LeastSquaresProblem;
46  import org.apache.commons.math3.fitting.leastsquares.LevenbergMarquardtOptimizer;
47  import org.apache.commons.math3.fitting.leastsquares.MultivariateJacobianFunction;
48  import org.apache.commons.math3.linear.Array2DRowRealMatrix;
49  import org.apache.commons.math3.linear.ArrayRealVector;
50  import org.apache.commons.math3.linear.RealMatrix;
51  import org.apache.commons.math3.linear.RealVector;
52  import org.apache.commons.math3.optim.ConvergenceChecker;
53  import org.apache.commons.math3.util.Pair;
54  
55  import java.util.logging.Level;
56  import java.util.logging.Logger;
57  
58  import static ffx.numerics.math.ScalarMath.binomial;
59  import static java.util.Arrays.fill;
60  
61  /**
62   * Predict Mutual Induced Dipoles based on previous steps.
63   *
64   * @author Stephen LuCore
65   * @since 1.0
66   */
67  public class ScfPredictor {
68  
69    private static final Logger logger = Logger.getLogger(ScfPredictor.class.getName());
70    /** Convergence tolerance of the LeastSquares optimizer. */
71    private static final double eps = 1.0e-4;
72  
73    private final PredictorMode predictorMode;
74    /** Induced dipole predictor order. */
75    private final int predictorOrder;
76    /** Dimensions of [nsymm][nAtoms][3] */
77    protected double[][][] inducedDipole;
78  
79    protected double[][][] inducedDipoleCR;
80    /** Number of atoms. */
81    private int nAtoms;
82    /** Maps LambdaMode to array indices: OFF/CONDENSED=0, CONDENSED_NOLIGAND=1, VAPOR=2 */
83    private int mode;
84    /** Dimensions of [mode][predictorOrder][nAtoms][3] */
85    private double[][][][] predictorInducedDipole;
86  
87    private double[][][][] predictorInducedDipoleCR;
88    /** Induced dipole predictor index. */
89    private int predictorStartIndex;
90    /** Induced dipole predictor count. */
91    private int predictorCount;
92    /**
93     * Predicts induced dipoles by locally minimizing and testing eps against squared change in
94     * parameters.
95     */
96    private LeastSquaresPredictor leastSquaresPredictor;
97    /**
98     * Constructor for ScfPredictor.
99     *
100    * @param mode a {@link ffx.potential.nonbonded.ScfPredictor.PredictorMode} object.
101    * @param order a int.
102    * @param ff a {@link ffx.potential.parameters.ForceField} object.
103    */
104   public ScfPredictor(PredictorMode mode, int order, ForceField ff) {
105     predictorMode = mode;
106     predictorOrder = order;
107     predictorCount = 0;
108     predictorStartIndex = 0;
109     if (predictorMode != PredictorMode.NONE) {
110       if (predictorMode == PredictorMode.LS) {
111         leastSquaresPredictor = new LeastSquaresPredictor(eps);
112       }
113       if (ff.getBoolean("LAMBDATERM", false) ) {
114         predictorInducedDipole = new double[3][predictorOrder][nAtoms][3];
115         predictorInducedDipoleCR = new double[3][predictorOrder][nAtoms][3];
116       } else {
117         predictorInducedDipole = new double[1][predictorOrder][nAtoms][3];
118         predictorInducedDipoleCR = new double[1][predictorOrder][nAtoms][3];
119       }
120     }
121   }
122 
123   /**
124    * run.
125    *
126    * @param lambdaMode a {@link LambdaMode} object.
127    */
128   public void run(LambdaMode lambdaMode) {
129     if (predictorMode == PredictorMode.NONE) {
130       return;
131     }
132     switch (lambdaMode) {
133       case CONDENSED_NO_LIGAND:
134         mode = 1;
135         break;
136       case VAPOR:
137         mode = 2;
138         break;
139       default:
140         mode = 0;
141     }
142     if (predictorMode != PredictorMode.NONE) {
143       switch (predictorMode) {
144         case ASPC:
145           aspcPredictor();
146           break;
147         case LS:
148           leastSquaresPredictor();
149           break;
150         case POLY:
151           polynomialPredictor();
152           break;
153         case NONE:
154         default:
155           break;
156       }
157     }
158   }
159 
160   /**
161    * Save the current converged mutual induced dipoles.
162    *
163    * @param inducedDipole an array of induced dipoles.
164    * @param inducedDipoleCR an array of induced dipoles chain rule terms.
165    * @param directDipole an array of direct dipoles.
166    * @param directDipoleCR an array of direct dipoles chain rule terms.
167    */
168   public void saveMutualInducedDipoles(
169       double[][][] inducedDipole,
170       double[][][] inducedDipoleCR,
171       double[][] directDipole,
172       double[][] directDipoleCR) {
173 
174     // Current induced dipoles are saved before those from the previous step.
175     predictorStartIndex--;
176     if (predictorStartIndex < 0) {
177       predictorStartIndex = predictorOrder - 1;
178     }
179 
180     if (predictorCount < predictorOrder) {
181       predictorCount++;
182     }
183 
184     for (int i = 0; i < nAtoms; i++) {
185       for (int j = 0; j < 3; j++) {
186         predictorInducedDipole[mode][predictorStartIndex][i][j] =
187             inducedDipole[0][i][j] - directDipole[i][j];
188         predictorInducedDipoleCR[mode][predictorStartIndex][i][j] =
189             inducedDipoleCR[0][i][j] - directDipoleCR[i][j];
190       }
191     }
192   }
193 
194   /**
195    * To be called upon initialization and update of inducedDipole arrays in parent.
196    *
197    * @param inducedDipole an array of induced dipoles.
198    * @param inducedDipoleCR an array of induced dipoles chain rule terms.
199    * @param lambdaTerm a boolean.
200    */
201   public void setInducedDipoleReferences(
202       double[][][] inducedDipole, double[][][] inducedDipoleCR, boolean lambdaTerm) {
203     this.inducedDipole = inducedDipole;
204     this.inducedDipoleCR = inducedDipoleCR;
205     if (lambdaTerm) {
206       predictorInducedDipole = new double[3][predictorOrder][nAtoms][3];
207       predictorInducedDipoleCR = new double[3][predictorOrder][nAtoms][3];
208     } else {
209       predictorInducedDipole = new double[1][predictorOrder][nAtoms][3];
210       predictorInducedDipoleCR = new double[1][predictorOrder][nAtoms][3];
211     }
212   }
213 
214   /** {@inheritDoc} */
215   @Override
216   public String toString() {
217     return predictorMode.toString();
218   }
219 
220   /**
221    * The least-squares predictor with induced dipole information from 8-10 previous steps reduces
222    * the number SCF iterations by ~50%.
223    */
224   private void leastSquaresPredictor() {
225     if (predictorCount < 2) {
226       return;
227     }
228     try {
229       /**
230        * The Jacobian and target do not change during the LS optimization, so it's most efficient to
231        * update them once before the Least-Squares optimizer starts.
232        */
233       leastSquaresPredictor.updateJacobianAndTarget();
234       int maxEvals, maxIter;
235       maxEvals = maxIter = 1000;
236       LeastSquaresOptimizer.Optimum optimum = leastSquaresPredictor.predict(maxEvals, maxIter);
237       double[] optimalValues = optimum.getPoint().toArray();
238       if (logger.isLoggable(Level.FINEST)) {
239         logger.finest(String.format("\n LS RMS:            %10.6f", optimum.getRMS()));
240         logger.finest(String.format(" LS Iterations:     %10d", optimum.getIterations()));
241         logger.finest(String.format(" Jacobian Evals:    %10d", optimum.getEvaluations()));
242         logger.finest(String.format(" Root Mean Square:  %10.6f", optimum.getRMS()));
243         logger.finest(" LS Coefficients");
244         for (int i = 0; i < predictorOrder - 1; i++) {
245           logger.finest(String.format(" %2d  %10.6f", i + 1, optimalValues[i]));
246         }
247       }
248 
249       /** Initialize a pointer into predictor induced dipole array. */
250       int index = predictorStartIndex;
251       /**
252        * Apply the LS coefficients in order to provide an initial guess at the converged induced
253        * dipoles.
254        */
255       for (int k = 0; k < predictorOrder - 1; k++) {
256         /** Set the current coefficient. */
257         double c = optimalValues[k];
258         for (int i = 0; i < nAtoms; i++) {
259           for (int j = 0; j < 3; j++) {
260             inducedDipole[0][i][j] += c * predictorInducedDipole[mode][index][i][j];
261             inducedDipoleCR[0][i][j] += c * predictorInducedDipoleCR[mode][index][i][j];
262           }
263         }
264         index++;
265         if (index >= predictorOrder) {
266           index = 0;
267         }
268       }
269     } catch (Exception e) {
270       logger.log(Level.WARNING, " Exception computing predictor coefficients", e);
271     }
272   }
273 
274   /** Always-stable predictor-corrector for the mutual induced dipoles. */
275   private void aspcPredictor() {
276 
277     if (predictorCount < 6) {
278       return;
279     }
280 
281     final double[] aspc = {
282       22.0 / 7.0, -55.0 / 14.0, 55.0 / 21.0, -22.0 / 21.0, 5.0 / 21.0, -1.0 / 42.0
283     };
284     /* Initialize a pointer into predictor induced dipole array. */
285     int index = predictorStartIndex;
286     /* Expansion loop. */
287     for (int k = 0; k < 6; k++) {
288       /* Set the current predictor coefficient. */
289       double c = aspc[k];
290       for (int i = 0; i < nAtoms; i++) {
291         for (int j = 0; j < 3; j++) {
292           inducedDipole[0][i][j] += c * predictorInducedDipole[mode][index][i][j];
293           inducedDipoleCR[0][i][j] += c * predictorInducedDipoleCR[mode][index][i][j];
294         }
295       }
296       index++;
297       if (index >= predictorOrder) {
298         index = 0;
299       }
300     }
301   }
302 
303   /** Polynomial predictor for the mutual induced dipoles. */
304   private void polynomialPredictor() {
305 
306     if (predictorCount == 0) {
307       return;
308     }
309 
310     /** Check the number of previous induced dipole vectors available. */
311     int n = predictorOrder;
312     if (predictorCount < predictorOrder) {
313       n = predictorCount;
314     }
315     /** Initialize a pointer into predictor induced dipole array. */
316     int index = predictorStartIndex;
317     /** Initialize the sign of the polynomial expansion. */
318     double sign = -1.0;
319     /** Expansion loop. */
320     for (int k = 0; k < n; k++) {
321       /** Set the current predictor sign and coefficient. */
322       sign *= -1.0;
323       double c = sign * binomial(n, k);
324       for (int i = 0; i < nAtoms; i++) {
325         for (int j = 0; j < 3; j++) {
326           inducedDipole[0][i][j] += c * predictorInducedDipole[mode][index][i][j];
327           inducedDipoleCR[0][i][j] += c * predictorInducedDipoleCR[mode][index][i][j];
328         }
329       }
330       index++;
331       if (index >= predictorOrder) {
332         index = 0;
333       }
334     }
335   }
336 
337   public enum PredictorMode {
338     NONE,
339     LS,
340     POLY,
341     ASPC
342   }
343 
344   private class LeastSquaresPredictor {
345 
346     double[] weights;
347     double[] target;
348     double[][] jacobian;
349     double[] initialSolution;
350     double tolerance;
351     RealVector valuesVector;
352     RealVector targetVector;
353     RealMatrix jacobianMatrix;
354     LeastSquaresOptimizer optimizer;
355     ConvergenceChecker<LeastSquaresProblem.Evaluation> checker;
356     Pair<RealVector, RealMatrix> test = new Pair<>(targetVector, jacobianMatrix);
357     MultivariateJacobianFunction function =
358         new MultivariateJacobianFunction() {
359           @Override
360           public Pair<RealVector, RealMatrix> value(RealVector point) {
361             return new Pair<>(targetVector, jacobianMatrix);
362           }
363         };
364 
365     public LeastSquaresPredictor(double eps) {
366       tolerance = eps;
367       weights = new double[2 * nAtoms * 3];
368       target = new double[2 * nAtoms * 3];
369       jacobian = new double[2 * nAtoms * 3][predictorOrder - 1];
370       initialSolution = new double[predictorOrder - 1];
371       fill(weights, 1.0);
372       fill(initialSolution, 0.0);
373       initialSolution[0] = 1.0;
374       optimizer = new LevenbergMarquardtOptimizer().withParameterRelativeTolerance(eps);
375       checker = new EvaluationRmsChecker(eps);
376     }
377 
378     public LeastSquaresOptimizer.Optimum predict(int maxEval, int maxIter) {
379       RealVector start = new ArrayRealVector(initialSolution);
380       LeastSquaresProblem lsp =
381           LeastSquaresFactory.create(function, targetVector, start, checker, maxEval, maxIter);
382 
383       LeastSquaresOptimizer.Optimum optimum = optimizer.optimize(lsp);
384       if (true)
385         logger.info(
386             String.format(
387                 " LS Optimization parameters:" + "  %s %s\n" + "  %s %s\n" + "  %d %d",
388                 function,
389                 targetVector.toString(),
390                 start,
391                 checker.toString(),
392                 maxIter,
393                 maxEval));
394       return optimum;
395     }
396 
397     public void updateJacobianAndTarget() {
398 
399       // Update the target.
400       int index = 0;
401       for (int i = 0; i < nAtoms; i++) {
402         target[index++] = predictorInducedDipole[mode][predictorStartIndex][i][0];
403         target[index++] = predictorInducedDipole[mode][predictorStartIndex][i][1];
404         target[index++] = predictorInducedDipole[mode][predictorStartIndex][i][2];
405         target[index++] = predictorInducedDipoleCR[mode][predictorStartIndex][i][0];
406         target[index++] = predictorInducedDipoleCR[mode][predictorStartIndex][i][1];
407         target[index++] = predictorInducedDipoleCR[mode][predictorStartIndex][i][2];
408       }
409       targetVector = new ArrayRealVector(target);
410 
411       // Update the Jacobian.
412       index = predictorStartIndex + 1;
413       if (index >= predictorOrder) {
414         index = 0;
415       }
416       for (int j = 0; j < predictorOrder - 1; j++) {
417         int ji = 0;
418         for (int i = 0; i < nAtoms; i++) {
419           jacobian[ji++][j] = predictorInducedDipole[mode][index][i][0];
420           jacobian[ji++][j] = predictorInducedDipole[mode][index][i][1];
421           jacobian[ji++][j] = predictorInducedDipole[mode][index][i][2];
422           jacobian[ji++][j] = predictorInducedDipoleCR[mode][index][i][0];
423           jacobian[ji++][j] = predictorInducedDipoleCR[mode][index][i][1];
424           jacobian[ji++][j] = predictorInducedDipoleCR[mode][index][i][2];
425         }
426         index++;
427         if (index >= predictorOrder) {
428           index = 0;
429         }
430       }
431       jacobianMatrix = new Array2DRowRealMatrix(jacobian);
432     }
433 
434     private RealVector value(double[] variables) {
435 
436       double[] values = new double[2 * nAtoms * 3];
437       for (int i = 0; i < nAtoms; i++) {
438         int index = 6 * i;
439         values[index] = 0;
440         values[index + 1] = 0;
441         values[index + 2] = 0;
442         values[index + 3] = 0;
443         values[index + 4] = 0;
444         values[index + 5] = 0;
445         int pi = predictorStartIndex + 1;
446         if (pi >= predictorOrder) {
447           pi = 0;
448         }
449         for (int j = 0; j < predictorOrder - 1; j++) {
450           values[index] += variables[j] * predictorInducedDipole[mode][pi][i][0];
451           values[index + 1] += variables[j] * predictorInducedDipole[mode][pi][i][1];
452           values[index + 2] += variables[j] * predictorInducedDipole[mode][pi][i][2];
453           values[index + 3] += variables[j] * predictorInducedDipoleCR[mode][pi][i][0];
454           values[index + 4] += variables[j] * predictorInducedDipoleCR[mode][pi][i][1];
455           values[index + 5] += variables[j] * predictorInducedDipoleCR[mode][pi][i][2];
456           pi++;
457           if (pi >= predictorOrder) {
458             pi = 0;
459           }
460         }
461       }
462       return new ArrayRealVector(values);
463     }
464   }
465 }