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;
39  
40  import edu.rit.pj.ParallelRegion;
41  import edu.rit.pj.ParallelSection;
42  import edu.rit.pj.ParallelTeam;
43  import ffx.crystal.Crystal;
44  import ffx.crystal.CrystalPotential;
45  import ffx.crystal.SymOp;
46  import ffx.numerics.Potential;
47  import ffx.numerics.switching.UnivariateSwitchingFunction;
48  import ffx.potential.bonded.Atom;
49  import ffx.potential.bonded.LambdaInterface;
50  import ffx.potential.openmm.OpenMMEnergy;
51  import ffx.potential.parameters.ForceField;
52  import ffx.potential.utils.EnergyException;
53  
54  import java.util.ArrayList;
55  import java.util.Arrays;
56  import java.util.List;
57  import java.util.logging.Level;
58  import java.util.logging.Logger;
59  
60  import static ffx.crystal.SymOp.applyCartesianSymOp;
61  import static ffx.crystal.SymOp.applyCartesianSymRot;
62  import static ffx.crystal.SymOp.invertSymOp;
63  import static ffx.potential.utils.Superpose.rmsd;
64  import static ffx.utilities.StringUtils.parseAtomRanges;
65  import static ffx.utilities.StringUtils.writeAtomRanges;
66  import static java.lang.Double.parseDouble;
67  import static java.lang.String.format;
68  import static java.util.Arrays.fill;
69  import static org.apache.commons.math3.util.FastMath.sqrt;
70  
71  /**
72   * Compute the potential energy and derivatives for a dual-topology system.
73   *
74   * @author Michael J. Schnieders
75   * @since 1.0
76   */
77  public class DualTopologyEnergy implements CrystalPotential, LambdaInterface {
78  
79    /**
80     * Logger for the DualTopologyEnergy class.
81     */
82    private static final Logger logger = Logger.getLogger(DualTopologyEnergy.class.getName());
83    /**
84     * Shared atoms between topologies 1 and 2.
85     */
86    private final int nShared;
87    /**
88     * Total number of variables
89     */
90    private final int nVariables;
91    /**
92     * Region for computing energy/gradient in parallel.
93     */
94    private final EnergyRegion energyRegion;
95    /**
96     * Mass array for shared and softcore atoms.
97     */
98    private final double[] mass;
99    /**
100    * VARIABLE_TYPE array for shared and softcore atoms.
101    */
102   private final VARIABLE_TYPE[] variableTypes;
103   /**
104    * Coordinates for topology 1.
105    */
106   private final double[] x1;
107   /**
108    * Coordinates for topology 2.
109    */
110   private final double[] x2;
111   /**
112    * Topology 1 coordinate gradient.
113    */
114   private final double[] g1;
115   /**
116    * Topology 2 coordinate gradient.
117    */
118   private final double[] g2;
119   /**
120    * Topology 1 restraint gradient for end state bonded terms.
121    */
122   private final double[] rg1;
123   /**
124    * Topology 2 restraint gradient for end state bonded terms.
125    */
126   private final double[] rg2;
127   /**
128    * Topology 1 derivative of the coordinate gradient with respect to lambda.
129    */
130   private final double[] gl1;
131   /**
132    * Topology 2 derivative of the coordinate gradient with respect to lambda.
133    */
134   private final double[] gl2;
135   /**
136    * Topology 1 derivative of the coordinate gradient with respect to lambda for end state bonded
137    * terms
138    */
139   private final double[] rgl1;
140   /**
141    * Topology 2 derivative of the coordinate gradient with respect to lambda for end state bonded
142    * terms
143    */
144   private final double[] rgl2;
145   /**
146    * Topology 1 Potential.
147    */
148   private final CrystalPotential potential1;
149   /**
150    * Topology 2 Potential.
151    */
152   private final CrystalPotential potential2;
153   /**
154    * Topology 1 LambdaInterface.
155    */
156   private final LambdaInterface lambdaInterface1;
157   /**
158    * Topology 2 LambdaInterface.
159    */
160   private final LambdaInterface lambdaInterface2;
161   /**
162    * Topology 1 ForceFieldEnergy.
163    */
164   private final ForceFieldEnergy forceFieldEnergy1;
165   /**
166    * Topology 2 ForceFieldEnergy.
167    */
168   private final ForceFieldEnergy forceFieldEnergy2;
169   /**
170    * Number of active atoms in Topology 1.
171    */
172   private final int nActive1;
173   /**
174    * Number of active atoms in Topology 2.
175    */
176   private final int nActive2;
177   /**
178    * Array of active atoms in Topology 1.
179    */
180   private final Atom[] activeAtoms1;
181   /**
182    * Array of active atoms in Topology 2.
183    */
184   private final Atom[] activeAtoms2;
185   /**
186    * Flag is true if the atom is shared (not alchemical) for Topology 1.
187    */
188   private final boolean[] sharedAtoms1;
189   /**
190    * Flag is true if the atom is shared (not alchemical) for Topology 1.
191    */
192   private final boolean[] sharedAtoms2;
193   /**
194    * Will default to a power-1 PowerSwitch function.
195    */
196   private final UnivariateSwitchingFunction switchFunction;
197   /**
198    * Current potential energy of topology 1 (kcal/mol).
199    */
200   private double energy1 = 0;
201   /**
202    * Current potential energy of topology 2 (kcal/mol).
203    */
204   private double energy2 = 0;
205   /**
206    * ParallelTeam to execute the EnergyRegion.
207    */
208   private ParallelTeam parallelTeam;
209   /**
210    * Include a valence restraint energy for atoms being "disappeared."
211    */
212   private final boolean doValenceRestraint1;
213   /**
214    * Include a valence restraint energy for atoms being "disappeared."
215    */
216   private final boolean doValenceRestraint2;
217   /**
218    * Use System 1 Bonded Energy.
219    */
220   private final boolean useFirstSystemBondedEnergy;
221   /**
222    * End-state restraint energy of topology 1 (kcal/mol).
223    */
224   private double restraintEnergy1 = 0;
225   /**
226    * End-state restraint energy of topology 2 (kcal/mol).
227    */
228   private double restraintEnergy2 = 0;
229   /**
230    * dEdL of topology 1 (kcal/mol).
231    */
232   private double dEdL_1 = 0;
233   /**
234    * dEdL of topology 2 (kcal/mol).
235    */
236   private double dEdL_2 = 0;
237   /**
238    * End-state restraint dEdL of topology 1 (kcal/mol).
239    */
240   private double restraintdEdL_1 = 0;
241   /**
242    * End-state restraint dEdL of topology 2 (kcal/mol).
243    */
244   private double restraintdEdL_2 = 0;
245   /**
246    * d2EdL2 of topology 1 (kcal/mol).
247    */
248   private double d2EdL2_1 = 0;
249   /**
250    * d2EdL2 of topology 2 (kcal/mol).
251    */
252   private double d2EdL2_2 = 0;
253   /**
254    * End-state restraint d2EdL2 of topology 1 (kcal/mol).
255    */
256   private double restraintd2EdL2_1 = 0;
257   /**
258    * End-state restraint d2EdL2 of topology 2 (kcal/mol).
259    */
260   private double restraintd2EdL2_2 = 0;
261   /**
262    * Total energy of the dual topology, including lambda scaling.
263    */
264   private double totalEnergy = 0;
265   /**
266    * Current lambda value.
267    */
268   private double lambda = 1.0;
269   /**
270    * Value of the switching function at lambda.
271    */
272   private double f1L = 1.0;
273   /**
274    * Value of the switching function at (1-lambda).
275    */
276   private double f2L = 0.0;
277   /**
278    * First derivative with respect to lambda of the switching function.
279    */
280   private double dF1dL = 0.0;
281   /**
282    * First derivative with respect to (1-lambda) of the switching function.
283    */
284   private double dF2dL = 0.0;
285   /**
286    * Second derivative with respect to lambda of the switching function.
287    */
288   private double d2F1dL2 = 0.0;
289   /**
290    * Second derivative with respect to (1-lambda) of the switching function.
291    */
292   private double d2F2dL2 = 0.0;
293   /**
294    * Scaling array for shared and softcore atoms.
295    */
296   private double[] scaling = null;
297   /**
298    * State for calculating energies.
299    */
300   private STATE state = STATE.BOTH;
301   /**
302    * Symmetry operator to superimpose system 2 onto system 1
303    */
304   private SymOp[] symOp;
305   /**
306    * Symmetry operator to move system 2 back to original frame
307    */
308   private SymOp[] inverse;
309   /**
310    * Atom selection for each symmetry operator.
311    */
312   private final ArrayList<List<Integer>> symOpAtoms = new ArrayList<>();
313   /**
314    * Mask to determine which atoms a sym op will be applied.
315    */
316   private boolean[][] mask = null;
317   /**
318    * Utilize a provided SymOp
319    */
320   private boolean useSymOp = false;
321 
322   /**
323    * Constructor for DualTopologyEnergy.
324    *
325    * @param topology1      a {@link ffx.potential.MolecularAssembly} object.
326    * @param topology2      a {@link ffx.potential.MolecularAssembly} object.
327    * @param switchFunction a {@link UnivariateSwitchingFunction} object.
328    */
329   public DualTopologyEnergy(
330       MolecularAssembly topology1,
331       MolecularAssembly topology2,
332       UnivariateSwitchingFunction switchFunction) {
333     forceFieldEnergy1 = topology1.getPotentialEnergy();
334     forceFieldEnergy2 = topology2.getPotentialEnergy();
335     potential1 = forceFieldEnergy1;
336     potential2 = forceFieldEnergy2;
337     lambdaInterface1 = forceFieldEnergy1;
338     lambdaInterface2 = forceFieldEnergy2;
339 
340     /* Apom array for topology 1. */
341     Atom[] atoms1 = topology1.getAtomArray();
342     /* Atom array for topology 2. */
343     Atom[] atoms2 = topology2.getAtomArray();
344 
345     ForceField forceField1 = topology1.getForceField();
346     doValenceRestraint1 = forceField1.getBoolean("LAMBDA_VALENCE_RESTRAINTS", true);
347     ForceField forceField2 = topology2.getForceField();
348     doValenceRestraint2 = forceField2.getBoolean("LAMBDA_VALENCE_RESTRAINTS", true);
349 
350     useFirstSystemBondedEnergy = forceField2.getBoolean("USE_FIRST_SYSTEM_BONDED_ENERGY", false);
351 
352     // Check that all atoms that are not undergoing alchemy are common to both topologies.
353     int shared1 = 0;
354     int shared2 = 0;
355     int activeCount1 = 0;
356     int activeCount2 = 0;
357     for (Atom a1 : atoms1) {
358       if (a1.isActive()) {
359         activeCount1++;
360         if (!a1.applyLambda()) {
361           shared1++;
362         }
363       }
364     }
365     for (Atom a2 : atoms2) {
366       if (a2.isActive()) {
367         activeCount2++;
368         if (!a2.applyLambda()) {
369           shared2++;
370         }
371       }
372     }
373     nActive1 = activeCount1;
374     nActive2 = activeCount2;
375     activeAtoms1 = new Atom[nActive1];
376     activeAtoms2 = new Atom[nActive2];
377     sharedAtoms1 = new boolean[nActive1];
378     sharedAtoms2 = new boolean[nActive2];
379     Arrays.fill(sharedAtoms1, true);
380     Arrays.fill(sharedAtoms2, true);
381     int index = 0;
382     for (Atom a1 : atoms1) {
383       if (a1.isActive()) {
384         activeAtoms1[index] = a1;
385         if (a1.applyLambda()) {
386           // This atom is softcore with independent coordinates.
387           sharedAtoms1[index] = false;
388         }
389         index++;
390       }
391     }
392     index = 0;
393     for (Atom a2 : atoms2) {
394       if (a2.isActive()) {
395         activeAtoms2[index] = a2;
396         if (a2.applyLambda()) {
397           // This atom is softcore with independent coordinates.
398           sharedAtoms2[index] = false;
399         }
400         index++;
401       }
402     }
403 
404     assert (shared1 == shared2);
405     nShared = shared1;
406     /* Topology 1 number of softcore atoms. */
407     int nSoftCore1 = nActive1 - nShared;
408     /* Topology 2 number of softcore atoms. */
409     int nSoftCore2 = nActive2 - nShared;
410     /* Total number of softcore and shared atoms: nTotal = nShared + nSoftcore1 + nSoftcore2 */
411     int nTotal = nShared + nSoftCore1 + nSoftCore2;
412     nVariables = 3 * nTotal;
413 
414     // Allocate memory for coordinates and derivatives.
415     x1 = new double[nActive1 * 3];
416     x2 = new double[nActive2 * 3];
417     g1 = new double[nActive1 * 3];
418     g2 = new double[nActive2 * 3];
419     rg1 = new double[nActive1 * 3];
420     rg2 = new double[nActive2 * 3];
421     gl1 = new double[nActive1 * 3];
422     gl2 = new double[nActive2 * 3];
423     rgl1 = new double[nActive1 * 3];
424     rgl2 = new double[nActive2 * 3];
425 
426     // All variables are coordinates.
427     index = 0;
428     variableTypes = new VARIABLE_TYPE[nVariables];
429     for (int i = 0; i < nTotal; i++) {
430       variableTypes[index++] = VARIABLE_TYPE.X;
431       variableTypes[index++] = VARIABLE_TYPE.Y;
432       variableTypes[index++] = VARIABLE_TYPE.Z;
433     }
434 
435     // Fill the mass array.
436     int commonIndex = 0;
437     int softcoreIndex = 3 * nShared;
438     mass = new double[nVariables];
439     for (int i = 0; i < nActive1; i++) {
440       Atom a = activeAtoms1[i];
441       double m = a.getMass();
442       if (sharedAtoms1[i]) {
443         mass[commonIndex++] = m;
444         mass[commonIndex++] = m;
445         mass[commonIndex++] = m;
446       } else {
447         mass[softcoreIndex++] = m;
448         mass[softcoreIndex++] = m;
449         mass[softcoreIndex++] = m;
450       }
451     }
452     commonIndex = 0;
453     for (int i = 0; i < nActive2; i++) {
454       Atom a = activeAtoms2[i];
455       double m = a.getMass();
456       if (sharedAtoms2[i]) {
457         for (int j = 0; j < 3; j++) {
458           double massCommon = mass[commonIndex];
459           massCommon = Math.max(massCommon, m);
460           mass[commonIndex++] = massCommon;
461         }
462       } else {
463         mass[softcoreIndex++] = m;
464         mass[softcoreIndex++] = m;
465         mass[softcoreIndex++] = m;
466       }
467     }
468     energyRegion = new EnergyRegion();
469     parallelTeam = new ParallelTeam(1);
470 
471     this.switchFunction = switchFunction;
472     logger.info(format("\n Dual topology using switching function:\n  %s", switchFunction));
473     logger.info(format(" Shared atoms: %d (1: %d of %d; 2: %d of %d)\n ", nShared, shared1, nActive1, shared2, nActive2));
474     String symOpString = forceField2.getString("symop", null);
475     if (symOpString != null) {
476       readSymOp(symOpString);
477     } else {
478       // Both end-states will use the same crystal object.
479       potential2.setCrystal(potential1.getCrystal());
480       useSymOp = false;
481       symOp = null;
482       inverse = null;
483     }
484 
485     // Check that all Dual-Topology atoms start with identical coordinates.
486     int i1 = 0;
487     int i2 = 0;
488     // Not true if mapping atoms via symmetry operator.
489     if (!useSymOp) {
490       for (int i = 0; i < nShared; i++) {
491         Atom a1 = atoms1[i1++];
492         while (a1.applyLambda()) {
493           a1 = atoms1[i1++];
494         }
495         Atom a2 = atoms2[i2++];
496         while (a2.applyLambda()) {
497           a2 = atoms2[i2++];
498         }
499         assert (a1.getX() == a2.getX());
500         assert (a1.getY() == a2.getY());
501         assert (a1.getZ() == a2.getZ());
502         reconcileAtoms(a1, a2, Level.INFO);
503       }
504     }
505   }
506 
507   /**
508    * Temporary method to be used as a benchmark. Move to SymOp?
509    *
510    * @param symOpString String containing sym op.
511    */
512   private void readSymOp(String symOpString) {
513     int numSymOps;
514     try {
515       String[] tokens = symOpString.split(" +");
516       int numTokens = tokens.length;
517       // Only need to move one system onto the other.
518       if (numTokens == 12) {
519         numSymOps = 1;
520         symOp = new SymOp[numSymOps];
521         inverse = new SymOp[numSymOps];
522         // Assume applies to all atoms...
523         symOpAtoms.add(parseAtomRanges("symmetryAtoms", "1-N", nActive2));
524         // Twelve values makes 1 sym op. so assign to [0][0].
525         symOp[0] = new SymOp(new double[][]{
526             {parseDouble(tokens[0]), parseDouble(tokens[1]), parseDouble(tokens[2])},
527             {parseDouble(tokens[3]), parseDouble(tokens[4]), parseDouble(tokens[5])},
528             {parseDouble(tokens[6]), parseDouble(tokens[7]), parseDouble(tokens[8])}},
529             new double[]{parseDouble(tokens[9]), parseDouble(tokens[10]), parseDouble(tokens[11])});
530       } else if (numTokens % 14 == 0) {
531         numSymOps = numTokens / 14;
532         symOp = new SymOp[numSymOps];
533         inverse = new SymOp[numSymOps];
534         if (logger.isLoggable(Level.FINER)) {
535           logger.finer(format(" Number of Tokens: %3d Number of Sym Ops: %2d", numTokens, numSymOps));
536         }
537         for (int i = 0; i < numSymOps; i++) {
538           int j = i * 14;
539           // Should be tokens[j] or tokens[j+1]?? Switches forward vs backward...
540           symOpAtoms.add(parseAtomRanges("symmetryAtoms", tokens[j], nActive2));
541           symOp[i] = new SymOp(new double[][]{
542               {parseDouble(tokens[j + 2]), parseDouble(tokens[j + 3]), parseDouble(tokens[j + 4])},
543               {parseDouble(tokens[j + 5]), parseDouble(tokens[j + 6]), parseDouble(tokens[j + 7])},
544               {parseDouble(tokens[j + 8]), parseDouble(tokens[j + 9]), parseDouble(tokens[j + 10])}},
545               new double[]{parseDouble(tokens[j + 11]), parseDouble(tokens[j + 12]), parseDouble(tokens[j + 13])});
546         }
547       } else {
548         logger.warning(format(" Provided symmetry operator is formatted incorrectly. Should have 12 or 14 entries, but found %3d.", numTokens));
549         return;
550       }
551       useSymOp = true;
552       // Define mask(s) needed to apply sym ops.
553       mask = new boolean[numSymOps][nActive2];
554       // Mask(s) used locally to check symmetry operator overlap.
555       boolean[][] sharedMask = new boolean[numSymOps][nShared];
556       for (int i = 0; i < numSymOps; i++) {
557         List<Integer> symAtoms = symOpAtoms.get(i);
558         for (int atom : symAtoms) {
559           if (sharedAtoms2[atom]) {
560             mask[i][atom] = true;
561           }
562         }
563         // TODO: Condense with loop above by replacing nActive2 with incremental index?
564         int sharedIndex = 0;
565         for (int j = 0; j < nActive2; j++) {
566           if (sharedAtoms2[j]) {
567             if (mask[i][j]) {
568               sharedMask[i][sharedIndex] = true;
569             }
570             sharedIndex++;
571           }
572         }
573       }
574 
575       // Get the coordinates for active atoms.
576       potential1.getCoordinates(x1);
577       potential2.getCoordinates(x2);
578 
579       // Collect coordinates for shared atoms (remove alchemical atoms from the superposition).
580       double[] x1Shared = new double[nShared * 3];
581       double[] x2Shared = new double[nShared * 3];
582       int sharedIndex = 0;
583       for (int i = 0; i < nActive1; i++) {
584         if (sharedAtoms1[i]) {
585           int index = i * 3;
586           x1Shared[sharedIndex++] = x1[index];
587           x1Shared[sharedIndex++] = x1[index + 1];
588           x1Shared[sharedIndex++] = x1[index + 2];
589         }
590       }
591       sharedIndex = 0;
592       for (int i = 0; i < nActive2; i++) {
593         if (sharedAtoms2[i]) {
594           int index = i * 3;
595           x2Shared[sharedIndex++] = x2[index];
596           x2Shared[sharedIndex++] = x2[index + 1];
597           x2Shared[sharedIndex++] = x2[index + 2];
598         }
599       }
600 
601       double origRMSD = rmsd(x1Shared, x2Shared, mass);
602       logger.info("\n RMSD to topology 2 coordinates from input file:");
603       logger.info(format(" RMSD: %12.3f A", origRMSD));
604 
605       // Get a fresh copy of the original coordinates.
606       double[] origX1 = Arrays.copyOf(x1Shared, nShared * 3);
607       double[] origX2 = Arrays.copyOf(x2Shared, nShared * 3);
608 
609       // Determine if first topology or second?
610       // Explicitly written (as opposed to applySymOp()) since only applied to shared atoms.
611       for (int i = 0; i < numSymOps; i++) {
612         // Transpose == inverse for orthogonal matrices (value needed later).
613         inverse[i] = invertSymOp(symOp[i]);
614         // Apply symOp to appropriate coords (validity check).
615         // Applied to topology 1 to generate approximate coordinates for topology 2.
616         applyCartesianSymOp(origX1, origX1, symOp[i], sharedMask[i]);
617         logger.info(format("\n SymOp %3d of %3d between topologies:\n Applied to atoms: %s\n %s\n Inverse SymOp:\n %s",
618             i + 1, numSymOps, writeAtomRanges(symOpAtoms.get(i).stream().mapToInt(j -> j).toArray()), symOp[i].toString(), inverse[i].toString()));
619       }
620       double symOpRMSD = rmsd(origX1, origX2, mass);
621       logger.info("\n RMSD to topology 2 coordinates via loaded symop:");
622       logger.info(format(" RMSD: %12.3f A", symOpRMSD));
623 
624       // TODO: Validate 0th order check for alchemical atoms (i.e., not bonded between sym ops).
625       logger.info("\n Checking alchemical atoms for system 1:");
626       validateAlchemicalAtoms(nActive1, sharedAtoms1, activeAtoms1);
627       logger.info(" Checking alchemical atoms for system 2:");
628       validateAlchemicalAtoms(nActive2, sharedAtoms2, activeAtoms2);
629 
630       if (logger.isLoggable(Level.FINE)) {
631         for (int i = 0; i < nShared; i++) {
632           int i3 = i * 3;
633           double rmsd = rmsd(new double[]{origX1[i3], origX1[i3 + 1], origX1[i3 + 2]}, new double[]{origX2[i3], origX2[i3 + 1], origX2[i3 + 2]}, new double[]{mass[i]});
634           if (rmsd > 0.5) {
635             logger.fine(format(" Shared atom %3d has a distance of %8.3f", i + 1, rmsd));
636           }
637         }
638       }
639     } catch (Exception ex) {
640       logger.severe(" Error parsing SymOp for Dual Topology:\n (" + symOpString + ")" + Utilities.stackTraceToString(ex));
641     }
642   }
643 
644   /**
645    * {@inheritDoc}
646    *
647    * <p>Returns true if both force field energies are zero at the ends, and the switching function's
648    * first derivative is zero at the end bounds.
649    */
650   @Override
651   public boolean dEdLZeroAtEnds() {
652     if (!forceFieldEnergy1.dEdLZeroAtEnds() || !forceFieldEnergy2.dEdLZeroAtEnds()) {
653       return false;
654     }
655     return (switchFunction.highestOrderZeroDerivativeAtZeroBound() > 0);
656   }
657 
658   /**
659    * {@inheritDoc}
660    */
661   @Override
662   public boolean destroy() {
663     boolean ffe1Destroy = forceFieldEnergy1.destroy();
664     boolean ffe2Destroy = forceFieldEnergy2.destroy();
665     try {
666       if (parallelTeam != null) {
667         parallelTeam.shutdown();
668       }
669       return ffe1Destroy && ffe2Destroy;
670     } catch (Exception ex) {
671       logger.warning(format(" Exception in shutting down DualTopologyEnergy: %s", ex));
672       logger.info(Utilities.stackTraceToString(ex));
673       return false;
674     }
675   }
676 
677   /**
678    * {@inheritDoc}
679    */
680   @Override
681   public double energy(double[] x) {
682     return energy(x, false);
683   }
684 
685   /**
686    * {@inheritDoc}
687    */
688   @Override
689   public double energy(double[] x, boolean verbose) {
690     try {
691       energyRegion.setX(x);
692       energyRegion.setVerbose(verbose);
693       parallelTeam.execute(energyRegion);
694     } catch (Exception ex) {
695       throw new EnergyException(format(" Exception in calculating dual-topology energy: %s", ex));
696     }
697     return totalEnergy;
698   }
699 
700   /**
701    * {@inheritDoc}
702    *
703    * <p>The coordinate and gradient arrays are unpacked/packed based on the dual topology.
704    */
705   @Override
706   public double energyAndGradient(double[] x, double[] g) {
707     return energyAndGradient(x, g, false);
708   }
709 
710   /**
711    * {@inheritDoc}
712    *
713    * <p>The coordinate and gradient arrays are unpacked/packed based on the dual topology.
714    */
715   @Override
716   public double energyAndGradient(double[] x, double[] g, boolean verbose) {
717     assert Arrays.stream(x).allMatch(Double::isFinite);
718     try {
719       energyRegion.setX(x);
720       energyRegion.setG(g);
721       energyRegion.setVerbose(verbose);
722       parallelTeam.execute(energyRegion);
723     } catch (Exception ex) {
724       throw new EnergyException(format(" Exception in calculating dual-topology energy: %s", ex));
725     }
726     return totalEnergy;
727   }
728 
729   /**
730    * {@inheritDoc}
731    */
732   @Override
733   public double[] getAcceleration(double[] acceleration) {
734     if (acceleration == null || acceleration.length < nVariables) {
735       acceleration = new double[nVariables];
736     }
737     int indexCommon = 0;
738     int indexUnique = nShared * 3;
739     double[] accel = new double[3];
740     for (int i = 0; i < nActive1; i++) {
741       Atom atom = activeAtoms1[i];
742       atom.getAcceleration(accel);
743       if (sharedAtoms1[i]) {
744         acceleration[indexCommon++] = accel[0];
745         acceleration[indexCommon++] = accel[1];
746         acceleration[indexCommon++] = accel[2];
747       } else {
748         acceleration[indexUnique++] = accel[0];
749         acceleration[indexUnique++] = accel[1];
750         acceleration[indexUnique++] = accel[2];
751       }
752     }
753     for (int i = 0; i < nActive2; i++) {
754       if (!sharedAtoms2[i]) {
755         Atom atom = activeAtoms2[i];
756         atom.getAcceleration(accel);
757         acceleration[indexUnique++] = accel[0];
758         acceleration[indexUnique++] = accel[1];
759         acceleration[indexUnique++] = accel[2];
760       }
761     }
762 
763     return acceleration;
764   }
765 
766   /**
767    * {@inheritDoc}
768    */
769   @Override
770   public double[] getCoordinates(double[] x) {
771     if (x == null) {
772       x = new double[nVariables];
773     }
774     int indexCommon = 0;
775     int indexUnique = nShared * 3;
776     for (int i = 0; i < nActive1; i++) {
777       Atom a = activeAtoms1[i];
778       if (sharedAtoms1[i]) {
779         x[indexCommon++] = a.getX();
780         x[indexCommon++] = a.getY();
781         x[indexCommon++] = a.getZ();
782       } else {
783         x[indexUnique++] = a.getX();
784         x[indexUnique++] = a.getY();
785         x[indexUnique++] = a.getZ();
786       }
787     }
788     for (int i = 0; i < nActive2; i++) {
789       if (!sharedAtoms2[i]) {
790         Atom a = activeAtoms2[i];
791         x[indexUnique++] = a.getX();
792         x[indexUnique++] = a.getY();
793         x[indexUnique++] = a.getZ();
794       }
795     }
796     return x;
797   }
798 
799   /**
800    * {@inheritDoc}
801    */
802   @Override
803   public Crystal getCrystal() {
804     // TODO: Handle the case where each system has a unique crystal instance (e.g., different space groups).
805     if (useSymOp && logger.isLoggable(Level.FINE)) {
806       logger.fine(" Get method only returned first crystal.");
807     }
808     return potential1.getCrystal();
809   }
810 
811   /**
812    * {@inheritDoc}
813    */
814   @Override
815   public void setCrystal(Crystal crystal) {
816     // TODO: Handle the case where each system has a unique crystal instance (e.g., different space groups).
817     if (useSymOp && logger.isLoggable(Level.FINE)) {
818       logger.fine(" Both systems set to the same crystal.");
819     }
820     potential1.setCrystal(crystal);
821     potential2.setCrystal(crystal);
822   }
823 
824   /**
825    * {@inheritDoc}
826    */
827   @Override
828   public STATE getEnergyTermState() {
829     return state;
830   }
831 
832   /**
833    * {@inheritDoc}
834    */
835   @Override
836   public void setEnergyTermState(STATE state) {
837     this.state = state;
838     potential1.setEnergyTermState(state);
839     potential2.setEnergyTermState(state);
840   }
841 
842   /**
843    * {@inheritDoc}
844    */
845   @Override
846   public double getLambda() {
847     return lambda;
848   }
849 
850   /**
851    * {@inheritDoc}
852    */
853   @Override
854   public void setLambda(double lambda) {
855     if (lambda <= 1.0 && lambda >= 0.0) {
856       this.lambda = lambda;
857       double oneMinusLambda = 1.0 - lambda;
858       lambdaInterface1.setLambda(lambda);
859       lambdaInterface2.setLambda(oneMinusLambda);
860 
861       f1L = switchFunction.valueAt(lambda);
862       dF1dL = switchFunction.firstDerivative(lambda);
863       d2F1dL2 = switchFunction.secondDerivative(lambda);
864 
865       f2L = switchFunction.valueAt(oneMinusLambda);
866       dF2dL = -1.0 * switchFunction.firstDerivative(oneMinusLambda);
867       d2F2dL2 = switchFunction.secondDerivative(oneMinusLambda);
868     } else {
869       String message = format("Lambda value %8.3f is not in the range [0..1].", lambda);
870       throw new IllegalArgumentException(message);
871     }
872   }
873 
874   /**
875    * {@inheritDoc}
876    */
877   @Override
878   public double[] getMass() {
879     return mass;
880   }
881 
882   /**
883    * Returns the number of shared variables (3 * number of shared atoms).
884    *
885    * @return An int
886    */
887   public int getNumSharedVariables() {
888     return 3 * nShared;
889   }
890 
891   /**
892    * {@inheritDoc}
893    */
894   @Override
895   public int getNumberOfVariables() {
896     return nVariables;
897   }
898 
899   /**
900    * {@inheritDoc}
901    */
902   @Override
903   public double[] getPreviousAcceleration(double[] previousAcceleration) {
904     if (previousAcceleration == null || previousAcceleration.length < nVariables) {
905       previousAcceleration = new double[nVariables];
906     }
907     int indexCommon = 0;
908     int indexUnique = nShared * 3;
909     double[] prev = new double[3];
910     for (int i = 0; i < nActive1; i++) {
911       Atom atom = activeAtoms1[i];
912       atom.getPreviousAcceleration(prev);
913       if (sharedAtoms1[i]) {
914         previousAcceleration[indexCommon++] = prev[0];
915         previousAcceleration[indexCommon++] = prev[1];
916         previousAcceleration[indexCommon++] = prev[2];
917       } else {
918         previousAcceleration[indexUnique++] = prev[0];
919         previousAcceleration[indexUnique++] = prev[1];
920         previousAcceleration[indexUnique++] = prev[2];
921       }
922     }
923     for (int i = 0; i < nActive2; i++) {
924       Atom atom = activeAtoms2[i];
925       if (!sharedAtoms2[i]) {
926         atom.getPreviousAcceleration(prev);
927         previousAcceleration[indexUnique++] = prev[0];
928         previousAcceleration[indexUnique++] = prev[1];
929         previousAcceleration[indexUnique++] = prev[2];
930       }
931     }
932 
933     return previousAcceleration;
934   }
935 
936   /**
937    * {@inheritDoc}
938    */
939   @Override
940   public double[] getScaling() {
941     return scaling;
942   }
943 
944   /**
945    * {@inheritDoc}
946    */
947   @Override
948   public void setScaling(double[] scaling) {
949     this.scaling = scaling;
950   }
951 
952   /**
953    * Returns the switching function used by this DualTopologyEnergy; presently, switching functions
954    * are immutable, and cannot be changed once a DualTopologyEnergy is constructed.
955    *
956    * @return The switching function.
957    */
958   public UnivariateSwitchingFunction getSwitch() {
959     return switchFunction;
960   }
961 
962   /**
963    * {@inheritDoc}
964    */
965   @Override
966   public double getTotalEnergy() {
967     return totalEnergy;
968   }
969 
970   @Override
971   public List<Potential> getUnderlyingPotentials() {
972     List<Potential> under = new ArrayList<>(2);
973     under.add(forceFieldEnergy1);
974     under.add(forceFieldEnergy2);
975     under.addAll(forceFieldEnergy1.getUnderlyingPotentials());
976     under.addAll(forceFieldEnergy2.getUnderlyingPotentials());
977     return under;
978   }
979 
980   /**
981    * {@inheritDoc}
982    */
983   @Override
984   public VARIABLE_TYPE[] getVariableTypes() {
985     return variableTypes;
986   }
987 
988   /**
989    * {@inheritDoc}
990    */
991   @Override
992   public double[] getVelocity(double[] velocity) {
993     if (velocity == null || velocity.length < nVariables) {
994       velocity = new double[nVariables];
995     }
996     int indexCommon = 0;
997     int indexUnique = nShared * 3;
998     double[] vel = new double[3];
999     for (int i = 0; i < nActive1; i++) {
1000       Atom atom = activeAtoms1[i];
1001       atom.getVelocity(vel);
1002       if (sharedAtoms1[i]) {
1003         velocity[indexCommon++] = vel[0];
1004         velocity[indexCommon++] = vel[1];
1005         velocity[indexCommon++] = vel[2];
1006       } else {
1007         velocity[indexUnique++] = vel[0];
1008         velocity[indexUnique++] = vel[1];
1009         velocity[indexUnique++] = vel[2];
1010       }
1011     }
1012     for (int i = 0; i < nActive2; i++) {
1013       if (!sharedAtoms2[i]) {
1014         Atom atom = activeAtoms2[i];
1015         atom.getVelocity(vel);
1016         velocity[indexUnique++] = vel[0];
1017         velocity[indexUnique++] = vel[1];
1018         velocity[indexUnique++] = vel[2];
1019       }
1020     }
1021 
1022     return velocity;
1023   }
1024 
1025   /**
1026    * {@inheritDoc}
1027    */
1028   @Override
1029   public double getd2EdL2() {
1030     double e1 =
1031         f1L * d2EdL2_1
1032             + 2.0 * dF1dL * dEdL_1
1033             + d2F1dL2 * energy1
1034             + f2L * restraintd2EdL2_1
1035             + 2.0 * dF2dL * restraintdEdL_1
1036             + d2F2dL2 * restraintEnergy1;
1037     double e2;
1038     if (!useFirstSystemBondedEnergy) {
1039       e2 = f2L * d2EdL2_2
1040           + 2.0 * dF2dL * dEdL_2
1041           + d2F2dL2 * energy2
1042           + f1L * restraintd2EdL2_2
1043           + 2.0 * dF1dL * restraintdEdL_2
1044           + d2F1dL2 * restraintEnergy2;
1045     } else {
1046       e2 = f2L * d2EdL2_2
1047           + 2.0 * dF2dL * dEdL_2
1048           + d2F2dL2 * energy2
1049           - f2L * restraintd2EdL2_2
1050           - 2.0 * dF2dL * restraintdEdL_2
1051           - d2F2dL2 * restraintEnergy2;
1052     }
1053     return e1 + e2;
1054   }
1055 
1056   /**
1057    * {@inheritDoc}
1058    */
1059   @Override
1060   public double getdEdL() {
1061     // Subtraction is implicit, as dEdL_2 and dF2dL are both multiplied by
1062     // negative 1 when set.
1063     double e1 = f1L * dEdL_1 + dF1dL * energy1 + f2L * restraintdEdL_1 + dF2dL * restraintEnergy1;
1064     double e2;
1065     if (!useFirstSystemBondedEnergy) {
1066       e2 = f2L * dEdL_2 + dF2dL * energy2 + f1L * restraintdEdL_2 + dF1dL * restraintEnergy2;
1067     } else {
1068       e2 = f2L * dEdL_2 + dF2dL * energy2 - f2L * restraintdEdL_2 - dF2dL * restraintEnergy2;
1069     }
1070     return e1 + e2;
1071   }
1072 
1073   /**
1074    * {@inheritDoc}
1075    */
1076   @Override
1077   public void getdEdXdL(double[] g) {
1078     if (g == null) {
1079       g = new double[nVariables];
1080     }
1081 
1082     int index = 0;
1083     int indexCommon = 0;
1084     int indexUnique = nShared * 3;
1085     // Coordinate Gradient from Topology 1.
1086     for (int i = 0; i < nActive1; i++) {
1087       if (sharedAtoms1[i]) {
1088         g[indexCommon++] =
1089             f1L * gl1[index] + dF1dL * g1[index] + f2L * rgl1[index] + dF2dL * rg1[index++];
1090         g[indexCommon++] =
1091             f1L * gl1[index] + dF1dL * g1[index] + f2L * rgl1[index] + dF2dL * rg1[index++];
1092         g[indexCommon++] =
1093             f1L * gl1[index] + dF1dL * g1[index] + f2L * rgl1[index] + dF2dL * rg1[index++];
1094       } else {
1095         g[indexUnique++] =
1096             f1L * gl1[index] + dF1dL * g1[index] + f2L * rgl1[index] + dF2dL * rg1[index++];
1097         g[indexUnique++] =
1098             f1L * gl1[index] + dF1dL * g1[index] + f2L * rgl1[index] + dF2dL * rg1[index++];
1099         g[indexUnique++] =
1100             f1L * gl1[index] + dF1dL * g1[index] + f2L * rgl1[index] + dF2dL * rg1[index++];
1101       }
1102     }
1103 
1104     // Coordinate Gradient from Topology 2.
1105     if (!useFirstSystemBondedEnergy) {
1106       index = 0;
1107       indexCommon = 0;
1108       for (int i = 0; i < nActive2; i++) {
1109         if (sharedAtoms2[i]) {
1110           g[indexCommon++] += (-f2L * gl2[index] + dF2dL * g2[index] - f1L * rgl2[index]
1111               + dF1dL * rg2[index++]);
1112           g[indexCommon++] += (-f2L * gl2[index] + dF2dL * g2[index] - f1L * rgl2[index]
1113               + dF1dL * rg2[index++]);
1114           g[indexCommon++] += (-f2L * gl2[index] + dF2dL * g2[index] - f1L * rgl2[index]
1115               + dF1dL * rg2[index++]);
1116         } else {
1117           g[indexUnique++] = (-f2L * gl2[index] + dF2dL * g2[index] - f1L * rgl2[index]
1118               + dF1dL * rg2[index++]);
1119           g[indexUnique++] = (-f2L * gl2[index] + dF2dL * g2[index] - f1L * rgl2[index]
1120               + dF1dL * rg2[index++]);
1121           g[indexUnique++] = (-f2L * gl2[index] + dF2dL * g2[index] - f1L * rgl2[index]
1122               + dF1dL * rg2[index++]);
1123         }
1124       }
1125     } else {
1126       index = 0;
1127       indexCommon = 0;
1128       for (int i = 0; i < nActive2; i++) {
1129         if (sharedAtoms2[i]) {
1130           g[indexCommon++] += (-f2L * gl2[index] + dF2dL * g2[index] + f2L * rgl2[index]
1131               - dF2dL * rg2[index++]);
1132           g[indexCommon++] += (-f2L * gl2[index] + dF2dL * g2[index] + f2L * rgl2[index]
1133               - dF2dL * rg2[index++]);
1134           g[indexCommon++] += (-f2L * gl2[index] + dF2dL * g2[index] + f2L * rgl2[index]
1135               - dF2dL * rg2[index++]);
1136         } else {
1137           g[indexUnique++] = (-f2L * gl2[index] + dF2dL * g2[index] + f2L * rgl2[index]
1138               - dF2dL * rg2[index++]);
1139           g[indexUnique++] = (-f2L * gl2[index] + dF2dL * g2[index] + f2L * rgl2[index]
1140               - dF2dL * rg2[index++]);
1141           g[indexUnique++] = (-f2L * gl2[index] + dF2dL * g2[index] + f2L * rgl2[index]
1142               - dF2dL * rg2[index++]);
1143         }
1144       }
1145 
1146     }
1147   }
1148 
1149   /**
1150    * {@inheritDoc}
1151    */
1152   @Override
1153   public void setAcceleration(double[] acceleration) {
1154     double[] accel = new double[3];
1155     int indexCommon = 0;
1156     int indexUnique = 3 * nShared;
1157     for (int i = 0; i < nActive1; i++) {
1158       Atom atom = activeAtoms1[i];
1159       if (sharedAtoms1[i]) {
1160         accel[0] = acceleration[indexCommon++];
1161         accel[1] = acceleration[indexCommon++];
1162         accel[2] = acceleration[indexCommon++];
1163       } else {
1164         accel[0] = acceleration[indexUnique++];
1165         accel[1] = acceleration[indexUnique++];
1166         accel[2] = acceleration[indexUnique++];
1167       }
1168       atom.setAcceleration(accel);
1169     }
1170     indexCommon = 0;
1171     for (int i = 0; i < nActive2; i++) {
1172       Atom atom = activeAtoms2[i];
1173       if (sharedAtoms2[i]) {
1174         accel[0] = acceleration[indexCommon++];
1175         accel[1] = acceleration[indexCommon++];
1176         accel[2] = acceleration[indexCommon++];
1177       } else {
1178         accel[0] = acceleration[indexUnique++];
1179         accel[1] = acceleration[indexUnique++];
1180         accel[2] = acceleration[indexUnique++];
1181       }
1182       atom.setAcceleration(accel);
1183     }
1184   }
1185 
1186   /**
1187    * setParallel.
1188    *
1189    * @param parallel a boolean.
1190    */
1191   public void setParallel(boolean parallel) {
1192     if (parallelTeam != null) {
1193       try {
1194         parallelTeam.shutdown();
1195       } catch (Exception e) {
1196         logger.severe(format(" Exception in shutting down old ParallelTeam for DualTopologyEnergy: %s", e));
1197       }
1198     }
1199     parallelTeam = parallel ? new ParallelTeam(2) : new ParallelTeam(1);
1200   }
1201 
1202   /**
1203    * {@inheritDoc}
1204    */
1205   @Override
1206   public void setPreviousAcceleration(double[] previousAcceleration) {
1207     double[] prev = new double[3];
1208     int indexCommon = 0;
1209     int indexUnique = 3 * nShared;
1210     for (int i = 0; i < nActive1; i++) {
1211       Atom atom = activeAtoms1[i];
1212       if (sharedAtoms1[i]) {
1213         prev[0] = previousAcceleration[indexCommon++];
1214         prev[1] = previousAcceleration[indexCommon++];
1215         prev[2] = previousAcceleration[indexCommon++];
1216       } else {
1217         prev[0] = previousAcceleration[indexUnique++];
1218         prev[1] = previousAcceleration[indexUnique++];
1219         prev[2] = previousAcceleration[indexUnique++];
1220       }
1221       atom.setPreviousAcceleration(prev);
1222     }
1223     indexCommon = 0;
1224     for (int i = 0; i < nActive2; i++) {
1225       Atom atom = activeAtoms2[i];
1226       if (sharedAtoms2[i]) {
1227         prev[0] = previousAcceleration[indexCommon++];
1228         prev[1] = previousAcceleration[indexCommon++];
1229         prev[2] = previousAcceleration[indexCommon++];
1230       } else {
1231         prev[0] = previousAcceleration[indexUnique++];
1232         prev[1] = previousAcceleration[indexUnique++];
1233         prev[2] = previousAcceleration[indexUnique++];
1234       }
1235       atom.setPreviousAcceleration(prev);
1236     }
1237   }
1238 
1239   /**
1240    * Sets the printOnFailure flag; if override is true, over-rides any existing property. Essentially
1241    * sets the default value of printOnFailure for an algorithm. For example, rotamer optimization
1242    * will generally run into force field issues in the normal course of execution as it tries
1243    * unphysical self and pair configurations, so the algorithm should not print out a large number of
1244    * error PDBs.
1245    *
1246    * @param onFail   To set
1247    * @param override Override properties
1248    */
1249   public void setPrintOnFailure(boolean onFail, boolean override) {
1250     forceFieldEnergy1.setPrintOnFailure(onFail, override);
1251     forceFieldEnergy2.setPrintOnFailure(onFail, override);
1252   }
1253 
1254   /**
1255    * {@inheritDoc}
1256    */
1257   @Override
1258   public void setVelocity(double[] velocity) {
1259     double[] vel = new double[3];
1260     int indexCommon = 0;
1261     int indexUnique = 3 * nShared;
1262     for (int i = 0; i < nActive1; i++) {
1263       Atom atom = activeAtoms1[i];
1264       if (sharedAtoms1[i]) {
1265         vel[0] = velocity[indexCommon++];
1266         vel[1] = velocity[indexCommon++];
1267         vel[2] = velocity[indexCommon++];
1268       } else {
1269         vel[0] = velocity[indexUnique++];
1270         vel[1] = velocity[indexUnique++];
1271         vel[2] = velocity[indexUnique++];
1272       }
1273       atom.setVelocity(vel);
1274     }
1275 
1276     indexCommon = 0;
1277     for (int i = 0; i < nActive2; i++) {
1278       Atom atom = activeAtoms2[i];
1279       if (sharedAtoms2[i]) {
1280         vel[0] = velocity[indexCommon++];
1281         vel[1] = velocity[indexCommon++];
1282         vel[2] = velocity[indexCommon++];
1283       } else {
1284         vel[0] = velocity[indexUnique++];
1285         vel[1] = velocity[indexUnique++];
1286         vel[2] = velocity[indexUnique++];
1287       }
1288       atom.setVelocity(vel);
1289     }
1290   }
1291 
1292   /**
1293    * Moves two shared atoms together if there is a small discrepancy (such as that caused by the
1294    * mutator script).
1295    *
1296    * @param a1      Atom from topology 1
1297    * @param a2      Atom from topology 2
1298    * @param warnlev Logging level to use when warning about small movements
1299    */
1300   private void reconcileAtoms(Atom a1, Atom a2, Level warnlev) {
1301     double dist = 0;
1302     double[] xyz1 = a1.getXYZ(null);
1303     double[] xyz2 = a2.getXYZ(null);
1304     double[] xyzAv = new double[3];
1305     for (int i = 0; i < 3; i++) {
1306       double dx = xyz1[i] - xyz2[i];
1307       dist += (dx * dx);
1308       xyzAv[i] = xyz1[i] + (0.5 * dx);
1309     }
1310 
1311     // Square of the maximum distance permissible between two shared atoms.
1312     double maxDist2 = 0.09;
1313 
1314     /*
1315      Square of the minimum distance between shared atoms which will cause a
1316      warning. Intended to cover anything larger than a rounding error.
1317     */
1318     double minDistWarn2 = 0.00001;
1319 
1320     if (dist > maxDist2) {
1321       logger.log(
1322           Level.SEVERE,
1323           format(
1324               " Distance between atoms %s " + "and %s is %7.4f >> maximum allowed %7.4f",
1325               a1, a2, sqrt(dist), sqrt(maxDist2)));
1326     } else if (dist > minDistWarn2) {
1327       logger.log(
1328           warnlev,
1329           format(
1330               " Distance between atoms %s " + "and %s is %7.4f; moving atoms together.",
1331               a1, a2, sqrt(dist)));
1332       a1.setXYZ(xyzAv);
1333       a2.setXYZ(xyzAv);
1334     } else if (dist > 0) {
1335       // Silently move them together; probably just a rounding error.
1336       a1.setXYZ(xyzAv);
1337       a2.setXYZ(xyzAv);
1338     }
1339   }
1340 
1341   private void packGradient(double[] x, double[] g) {
1342     if (g == null) {
1343       g = new double[nVariables];
1344     }
1345     int indexCommon = 0;
1346     int indexUnique = nShared * 3;
1347 
1348     // Coordinate Gradient from Topology 1.
1349     int index = 0;
1350     for (int i = 0; i < nActive1; i++) {
1351       if (sharedAtoms1[i]) {
1352         g[indexCommon++] = f1L * g1[index] + f2L * rg1[index++];
1353         g[indexCommon++] = f1L * g1[index] + f2L * rg1[index++];
1354         g[indexCommon++] = f1L * g1[index] + f2L * rg1[index++];
1355       } else {
1356         g[indexUnique++] = f1L * g1[index] + f2L * rg1[index++];
1357         g[indexUnique++] = f1L * g1[index] + f2L * rg1[index++];
1358         g[indexUnique++] = f1L * g1[index] + f2L * rg1[index++];
1359       }
1360     }
1361 
1362     if (!useFirstSystemBondedEnergy) {
1363       // Coordinate Gradient from Topology 2.
1364       indexCommon = 0;
1365       index = 0;
1366       for (int i = 0; i < nActive2; i++) {
1367         if (sharedAtoms2[i]) {
1368           g[indexCommon++] += f2L * g2[index] + f1L * rg2[index++];
1369           g[indexCommon++] += f2L * g2[index] + f1L * rg2[index++];
1370           g[indexCommon++] += f2L * g2[index] + f1L * rg2[index++];
1371         } else {
1372           g[indexUnique++] = f2L * g2[index] + f1L * rg2[index++];
1373           g[indexUnique++] = f2L * g2[index] + f1L * rg2[index++];
1374           g[indexUnique++] = f2L * g2[index] + f1L * rg2[index++];
1375         }
1376       }
1377     } else {
1378       // Coordinate Gradient from Topology 2.
1379       indexCommon = 0;
1380       index = 0;
1381       for (int i = 0; i < nActive2; i++) {
1382         if (sharedAtoms2[i]) {
1383           g[indexCommon++] += f2L * g2[index] - f2L * rg2[index++];
1384           g[indexCommon++] += f2L * g2[index] - f2L * rg2[index++];
1385           g[indexCommon++] += f2L * g2[index] - f2L * rg2[index++];
1386         } else {
1387           g[indexUnique++] = f2L * g2[index] - f2L * rg2[index++];
1388           g[indexUnique++] = f2L * g2[index] - f2L * rg2[index++];
1389           g[indexUnique++] = f2L * g2[index] - f2L * rg2[index++];
1390         }
1391       }
1392     }
1393 
1394     scaleCoordinatesAndGradient(x, g);
1395   }
1396 
1397   private void unpackCoordinates(double[] x) {
1398 
1399     // Unscale the coordinates.
1400     unscaleCoordinates(x);
1401 
1402     int index = 0;
1403     int indexCommon = 0;
1404     int indexUnique = 3 * nShared;
1405     for (int i = 0; i < nActive1; i++) {
1406       if (sharedAtoms1[i]) {
1407         x1[index++] = x[indexCommon++];
1408         x1[index++] = x[indexCommon++];
1409         x1[index++] = x[indexCommon++];
1410       } else {
1411         x1[index++] = x[indexUnique++];
1412         x1[index++] = x[indexUnique++];
1413         x1[index++] = x[indexUnique++];
1414       }
1415     }
1416 
1417     index = 0;
1418     indexCommon = 0;
1419     for (int i = 0; i < nActive2; i++) {
1420       if (sharedAtoms2[i]) {
1421         x2[index++] = x[indexCommon++];
1422         x2[index++] = x[indexCommon++];
1423         x2[index++] = x[indexCommon++];
1424       } else {
1425         x2[index++] = x[indexUnique++];
1426         x2[index++] = x[indexUnique++];
1427         x2[index++] = x[indexUnique++];
1428       }
1429     }
1430   }
1431 
1432   /**
1433    * Reload the common atomic masses. Intended for quad-topology dual force field corrections, where
1434    * common atoms may have slightly different masses; this was found to be the case between AMOEBA
1435    * 2013 and AMBER99SB carbons.
1436    *
1437    * @param secondTopology Load from second topology
1438    */
1439   void reloadCommonMasses(boolean secondTopology) {
1440     int commonIndex = 0;
1441     if (secondTopology) {
1442       for (int i = 0; i < nActive2; i++) {
1443         Atom a = activeAtoms2[i];
1444         double m = a.getMass();
1445         if (sharedAtoms2[i]) {
1446           mass[commonIndex++] = m;
1447           mass[commonIndex++] = m;
1448           mass[commonIndex++] = m;
1449         }
1450       }
1451     } else {
1452       for (int i = 0; i < nActive1; i++) {
1453         Atom a = activeAtoms1[i];
1454         double m = a.getMass();
1455         if (sharedAtoms1[i]) {
1456           mass[commonIndex++] = m;
1457           mass[commonIndex++] = m;
1458           mass[commonIndex++] = m;
1459         }
1460       }
1461     }
1462   }
1463 
1464   /**
1465    * Bonded terms for alchemical atoms are used to prevent them from being unconstrained.
1466    * These restraints have the potential to interact with symmetry operators applied in dual-topology simulations.
1467    *
1468    * @param nActive     Number of active atoms (shared and alchemical) in simulation for this system.
1469    * @param sharedAtoms Mask to determined which atoms are shared vs alchemical.
1470    * @param activeAtoms Active atoms for this simulation.
1471    */
1472   private void validateAlchemicalAtoms(int nActive, boolean[] sharedAtoms, Atom[] activeAtoms) {
1473     int numSymOps = symOpAtoms.size();
1474     // Loop through atoms in simulation.
1475     for (int i = 0; i < nActive; i++) {
1476       // Determine if this atom is alchemical.
1477       if (!sharedAtoms[i]) {
1478         Atom alchAtom = activeAtoms[i];
1479         int symOpForAlchemicalAtom = -1;
1480         for (int j = 0; j < numSymOps; j++) {
1481           List<Integer> symOpAtomGroup = symOpAtoms.get(j);
1482           if (symOpAtomGroup.contains(i)) {
1483             symOpForAlchemicalAtom = j;
1484             break;
1485           }
1486         }
1487 
1488         StringBuilder sb = new StringBuilder();
1489         sb.append(format("\n Alchemical atom %s\n Symmetry group %s:\n", alchAtom, symOpForAlchemicalAtom + 1));
1490         List<Integer> symOpAtomGroup = symOpAtoms.get(symOpForAlchemicalAtom);
1491         for (Integer j : symOpAtomGroup) {
1492           if (j == i) {
1493             continue;
1494           }
1495           sb.append(format("  %d %s\n", j + 1, activeAtoms[j]));
1496         }
1497         sb.append(format("\n 1-2, 1-3 or 1-4 atoms in other symmetry groups:\n"));
1498 
1499         // Collect atoms that are 1-2, 1-3 or 1-4 to the alchemical atom.
1500         ArrayList<Integer> bondedAtoms = new ArrayList<>();
1501         addUniqueBondedIndices(bondedAtoms, alchAtom.get12List());
1502         addUniqueBondedIndices(bondedAtoms, alchAtom.get13List());
1503         addUniqueBondedIndices(bondedAtoms, alchAtom.get14List());
1504         // Loop through atoms bonded to alchemical atom to determine if they're part of
1505         // a different SymOp group.
1506         boolean conflict = false;
1507         for (int j = 0; j < numSymOps; j++) {
1508           // If the bonded atom is in the same symmetry group as the alchemical atom, then there is no conflict.
1509           if (j == symOpForAlchemicalAtom) {
1510             continue;
1511           }
1512           symOpAtomGroup = symOpAtoms.get(j);
1513           for (Integer integer : bondedAtoms) {
1514             if (symOpAtomGroup.contains(integer)) {
1515               conflict = true;
1516               sb.append(format("  %d %s uses symmetry group %2d.\n", integer + 1, activeAtoms[integer], j + 1));
1517             }
1518           }
1519         }
1520         if (conflict) {
1521           logger.info(sb.toString());
1522         }
1523       }
1524     }
1525   }
1526 
1527   /**
1528    * Take list of aggregated unique indices based on input list of bonded atoms..
1529    *
1530    * @param uniqueIndices List of unique indices of bonded atoms.
1531    * @param bondedAtoms   Atoms that are bonded to the atom of interest.
1532    */
1533   private void addUniqueBondedIndices(ArrayList<Integer> uniqueIndices, List<Atom> bondedAtoms) {
1534     // Loop through bonded atoms
1535     for (Atom a : bondedAtoms) {
1536       int index = a.getIndex() - 1;
1537       // Determine if bonded atom's index is already in unique list.
1538       if (!uniqueIndices.contains(index)) {
1539         // Add to list if it is not already included.
1540         uniqueIndices.add(index);
1541       }
1542     }
1543   }
1544 
1545   private class EnergyRegion extends ParallelRegion {
1546 
1547     private final Energy1Section e1sect;
1548     private final Energy2Section e2sect;
1549     private double[] x;
1550     private double[] g;
1551     private boolean gradient = false;
1552     private boolean verbose = false;
1553 
1554     EnergyRegion() {
1555       e1sect = new Energy1Section();
1556       e2sect = new Energy2Section();
1557     }
1558 
1559     @Override
1560     public void finish() {
1561       // Apply the dual-topology scaling for the total energy.
1562       if (!useFirstSystemBondedEnergy) {
1563         totalEnergy =
1564             f1L * energy1 + f2L * restraintEnergy1 + f2L * energy2 + f1L * restraintEnergy2;
1565       } else {
1566         totalEnergy =
1567             f1L * energy1 + f2L * restraintEnergy1 + f2L * energy2 - f2L * restraintEnergy2;
1568       }
1569 
1570       if (gradient) {
1571         packGradient(x, g);
1572       } else {
1573         scaleCoordinates(x);
1574       }
1575 
1576       if (verbose) {
1577         logger.info(format(" Total dual-topology energy: %12.4f", totalEnergy));
1578       }
1579       setVerbose(false);
1580       setGradient(false);
1581     }
1582 
1583     @Override
1584     public void run() throws Exception {
1585       execute(e1sect, e2sect);
1586     }
1587 
1588     public void setG(double[] g) {
1589       this.g = g;
1590       setGradient(true);
1591     }
1592 
1593     public void setGradient(boolean grad) {
1594       this.gradient = grad;
1595       e1sect.setGradient(grad);
1596       e2sect.setGradient(grad);
1597     }
1598 
1599     public void setVerbose(boolean verb) {
1600       this.verbose = verb;
1601       e1sect.setVerbose(verb);
1602       e2sect.setVerbose(verb);
1603     }
1604 
1605     public void setX(double[] x) {
1606       this.x = x;
1607     }
1608 
1609     @Override
1610     public void start() {
1611       unpackCoordinates(x);
1612     }
1613   }
1614 
1615   private class Energy1Section extends ParallelSection {
1616 
1617     private boolean gradient = false;
1618     private boolean verbose = false;
1619 
1620     @Override
1621     public void run() {
1622       if (gradient) {
1623         fill(gl1, 0.0);
1624         fill(rgl1, 0.0);
1625         energy1 = potential1.energyAndGradient(x1, g1, verbose);
1626         dEdL_1 = lambdaInterface1.getdEdL();
1627         d2EdL2_1 = lambdaInterface1.getd2EdL2();
1628         lambdaInterface1.getdEdXdL(gl1);
1629         if (doValenceRestraint1) {
1630           if (verbose) {
1631             logger.info(" Calculating lambda bonded terms for topology 1");
1632           }
1633           forceFieldEnergy1.setLambdaBondedTerms(true, useFirstSystemBondedEnergy);
1634           if (forceFieldEnergy1 instanceof OpenMMEnergy openMMEnergy) {
1635             // Use FFX to calculate the energy and gradient for the restraints.
1636             restraintEnergy1 = openMMEnergy.energyAndGradientFFX(x1, rg1, verbose);
1637           } else {
1638             restraintEnergy1 = forceFieldEnergy1.energyAndGradient(x1, rg1, verbose);
1639           }
1640           restraintdEdL_1 = forceFieldEnergy1.getdEdL();
1641           restraintd2EdL2_1 = forceFieldEnergy1.getd2EdL2();
1642           forceFieldEnergy1.getdEdXdL(rgl1);
1643           forceFieldEnergy1.setLambdaBondedTerms(false, false);
1644         } else {
1645           restraintEnergy1 = 0.0;
1646           restraintdEdL_1 = 0.0;
1647           restraintd2EdL2_1 = 0.0;
1648         }
1649         if (logger.isLoggable(Level.FINE)) {
1650           logger.fine(format(" Topology 1 Energy & Restraints: %15.8f %15.8f\n", f1L * energy1, f2L * restraintEnergy1));
1651           logger.fine(format(" Topology 1:    %15.8f * (%.2f)", energy1, f1L));
1652           logger.fine(format(" T1 Restraints: %15.8f * (%.2f)", restraintEnergy1, f2L));
1653         }
1654       } else {
1655         energy1 = potential1.energy(x1, verbose);
1656         if (doValenceRestraint1) {
1657           if (verbose) {
1658             logger.info(" Calculating lambda bonded terms for topology 1");
1659           }
1660           forceFieldEnergy1.setLambdaBondedTerms(true, useFirstSystemBondedEnergy);
1661           if (forceFieldEnergy1 instanceof OpenMMEnergy openMMEnergy) {
1662             // Use FFX to calculate the energy and gradient for the restraints.
1663             restraintEnergy1 = openMMEnergy.energyFFX(x1, verbose);
1664           } else {
1665             restraintEnergy1 = potential1.energy(x1, verbose);
1666           }
1667           forceFieldEnergy1.setLambdaBondedTerms(false, false);
1668         } else {
1669           restraintEnergy1 = 0.0;
1670         }
1671         if (logger.isLoggable(Level.FINE)) {
1672           logger.fine(format(" Topology 1 Energy & Restraints: %15.8f %15.8f\n", f1L * energy1, f2L * restraintEnergy1));
1673         }
1674       }
1675     }
1676 
1677     public void setGradient(boolean grad) {
1678       this.gradient = grad;
1679     }
1680 
1681     public void setVerbose(boolean verb) {
1682       this.verbose = verb;
1683     }
1684   }
1685 
1686   private class Energy2Section extends ParallelSection {
1687 
1688     private boolean gradient = false;
1689     private boolean verbose = false;
1690 
1691     @Override
1692     public void run() {
1693       // Number of symmetry operators to be applied to the system.
1694       int numSymOps = 0;
1695       if (useSymOp) {
1696         numSymOps = symOpAtoms.size();
1697         for (int i = 0; i < numSymOps; i++) {
1698           // Apply sym ops to system 1 coords (duplicate of x1) to get x2.
1699           applyCartesianSymOp(x2, x2, symOp[i], mask[i]);
1700         }
1701       }
1702 
1703       if (gradient) {
1704         fill(gl2, 0.0);
1705         fill(rgl2, 0.0);
1706 
1707         // Compute the energy and gradient of topology 2.
1708         energy2 = potential2.energyAndGradient(x2, g2, verbose);
1709         dEdL_2 = -lambdaInterface2.getdEdL();
1710         d2EdL2_2 = lambdaInterface2.getd2EdL2();
1711         lambdaInterface2.getdEdXdL(gl2);
1712         if (useSymOp) {
1713           // Rotate the gradient back for appropriate atoms.
1714           for (int i = 0; i < numSymOps; i++) {
1715             applyCartesianSymRot(g2, g2, inverse[i], mask[i]);
1716             applyCartesianSymRot(gl2, gl2, inverse[i], mask[i]);
1717           }
1718         }
1719         if (doValenceRestraint2) {
1720           if (verbose) {
1721             logger.info(" Calculating lambda bonded terms for topology 2");
1722           }
1723           forceFieldEnergy2.setLambdaBondedTerms(true, useFirstSystemBondedEnergy);
1724           if (forceFieldEnergy2 instanceof OpenMMEnergy openMMEnergy) {
1725             // Use FFX to calculate the energy and gradient for the restraints.
1726             restraintEnergy2 = openMMEnergy.energyAndGradientFFX(x2, rg2, verbose);
1727           } else {
1728             restraintEnergy2 = forceFieldEnergy2.energyAndGradient(x2, rg2, verbose);
1729           }
1730           restraintdEdL_2 = -forceFieldEnergy2.getdEdL();
1731           restraintd2EdL2_2 = forceFieldEnergy2.getd2EdL2();
1732           forceFieldEnergy2.getdEdXdL(rgl2);
1733           if (useSymOp) {
1734             // Rotate the gradient back for appropriate atoms.
1735             for (int i = 0; i < numSymOps; i++) {
1736               applyCartesianSymRot(rg2, rg2, inverse[i], mask[i]);
1737               applyCartesianSymRot(rgl2, rgl2, inverse[i], mask[i]);
1738             }
1739           }
1740           forceFieldEnergy2.setLambdaBondedTerms(false, false);
1741         } else {
1742           restraintEnergy2 = 0.0;
1743           restraintdEdL_2 = 0.0;
1744           restraintd2EdL2_2 = 0.0;
1745         }
1746         if (logger.isLoggable(Level.FINE)) {
1747           logger.fine(format(" Topology 2 Energy & Restraints: %15.8f %15.8f", f2L * energy2, f1L * restraintEnergy2));
1748           logger.fine(format(" Topology 2:    %15.8f * (%.2f)", energy2, f2L));
1749           logger.fine(format(" T2 Restraints: %15.8f * (%.2f)", restraintEnergy2, f1L));
1750         }
1751       } else {
1752         energy2 = potential2.energy(x2, verbose);
1753         if (doValenceRestraint2) {
1754           if (verbose) {
1755             logger.info(" Calculating lambda bonded terms for topology 2");
1756           }
1757           forceFieldEnergy2.setLambdaBondedTerms(true, useFirstSystemBondedEnergy);
1758           if (forceFieldEnergy2 instanceof OpenMMEnergy openMMEnergy) {
1759             // Use FFX to calculate the energy and gradient for the restraints.
1760             restraintEnergy2 = openMMEnergy.energyFFX(x2, verbose);
1761           } else {
1762             restraintEnergy2 = potential2.energy(x2, verbose);
1763           }
1764           forceFieldEnergy2.setLambdaBondedTerms(false, false);
1765         } else {
1766           restraintEnergy2 = 0.0;
1767         }
1768         if (logger.isLoggable(Level.FINE)) {
1769           logger.fine(format(" Topology 2 Energy & Restraints: %15.8f %15.8f", f2L * energy2, f1L * restraintEnergy2));
1770         }
1771       }
1772     }
1773 
1774     public void setGradient(boolean grad) {
1775       this.gradient = grad;
1776     }
1777 
1778     public void setVerbose(boolean verb) {
1779       this.verbose = verb;
1780     }
1781   }
1782 }