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.ParallelTeam;
41  import ffx.crystal.Crystal;
42  import ffx.crystal.CrystalPotential;
43  import ffx.crystal.LatticeSystem;
44  import ffx.crystal.NCSCrystal;
45  import ffx.crystal.ReplicatesCrystal;
46  import ffx.crystal.SpaceGroup;
47  import ffx.crystal.SpaceGroupDefinitions;
48  import ffx.crystal.SymOp;
49  import ffx.numerics.Constraint;
50  import ffx.numerics.math.Double3;
51  import ffx.numerics.switching.ConstantSwitch;
52  import ffx.numerics.switching.UnivariateFunctionFactory;
53  import ffx.numerics.switching.UnivariateSwitchingFunction;
54  import ffx.potential.bonded.Angle;
55  import ffx.potential.bonded.AngleTorsion;
56  import ffx.potential.bonded.Atom;
57  import ffx.potential.bonded.Bond;
58  import ffx.potential.bonded.ImproperTorsion;
59  import ffx.potential.bonded.LambdaInterface;
60  import ffx.potential.bonded.MSNode;
61  import ffx.potential.bonded.OutOfPlaneBend;
62  import ffx.potential.bonded.PiOrbitalTorsion;
63  import ffx.potential.bonded.RestrainDistance;
64  import ffx.potential.bonded.RestrainPosition;
65  import ffx.potential.bonded.StretchBend;
66  import ffx.potential.bonded.StretchTorsion;
67  import ffx.potential.bonded.Torsion;
68  import ffx.potential.bonded.TorsionTorsion;
69  import ffx.potential.bonded.UreyBradley;
70  import ffx.potential.constraint.CcmaConstraint;
71  import ffx.potential.constraint.SettleConstraint;
72  import ffx.potential.extended.ExtendedSystem;
73  import ffx.potential.nonbonded.COMRestraint;
74  import ffx.potential.nonbonded.GeneralizedKirkwood;
75  import ffx.potential.nonbonded.NCSRestraint;
76  import ffx.potential.nonbonded.ParticleMeshEwald;
77  import ffx.potential.nonbonded.RestrainGroups;
78  import ffx.potential.nonbonded.VanDerWaals;
79  import ffx.potential.nonbonded.VanDerWaalsTornado;
80  import ffx.potential.openmm.OpenMMEnergy;
81  import ffx.potential.parameters.AngleType.AngleMode;
82  import ffx.potential.parameters.BondType;
83  import ffx.potential.parameters.ForceField;
84  import ffx.potential.parameters.ForceField.ELEC_FORM;
85  import ffx.potential.parameters.TorsionType;
86  import ffx.potential.terms.AnglePotentialEnergy;
87  import ffx.potential.terms.AngleTorsionPotentialEnergy;
88  import ffx.potential.terms.BondPotentialEnergy;
89  import ffx.potential.terms.EnergyTermRegion;
90  import ffx.potential.terms.ImproperTorsionPotentialEnergy;
91  import ffx.potential.terms.OutOfPlaneBendPotentialEnergy;
92  import ffx.potential.terms.PiOrbitalTorsionPotentialEnergy;
93  import ffx.potential.terms.RestrainDistancePotentialEnergy;
94  import ffx.potential.terms.RestrainPositionPotentialEnergy;
95  import ffx.potential.terms.RestrainTorsionPotentialEnergy;
96  import ffx.potential.terms.StretchBendPotentialEnergy;
97  import ffx.potential.terms.StretchTorsionPotentialEnergy;
98  import ffx.potential.terms.TorsionPotentialEnergy;
99  import ffx.potential.terms.TorsionTorsionPotentialEnergy;
100 import ffx.potential.terms.UreyBradleyPotentialEnergy;
101 import ffx.potential.utils.ConvexHullOps;
102 import ffx.potential.utils.EnergyException;
103 import ffx.potential.utils.PotentialsFunctions;
104 import ffx.potential.utils.PotentialsUtils;
105 import ffx.utilities.FFXProperty;
106 import org.apache.commons.configuration2.CompositeConfiguration;
107 
108 import javax.annotation.Nullable;
109 import java.io.File;
110 import java.time.LocalDateTime;
111 import java.time.format.DateTimeFormatter;
112 import java.util.ArrayList;
113 import java.util.Arrays;
114 import java.util.Collections;
115 import java.util.HashSet;
116 import java.util.List;
117 import java.util.Set;
118 import java.util.logging.Level;
119 import java.util.logging.Logger;
120 import java.util.stream.Collectors;
121 import java.util.stream.Stream;
122 
123 import static ffx.potential.bonded.BondedTerm.removeNeuralNetworkTerms;
124 import static ffx.potential.bonded.RestrainPosition.parseRestrainPositions;
125 import static ffx.potential.nonbonded.VanDerWaalsForm.DEFAULT_VDW_CUTOFF;
126 import static ffx.potential.nonbonded.pme.EwaldParameters.DEFAULT_EWALD_CUTOFF;
127 import static ffx.potential.parameters.ForceField.toEnumForm;
128 import static ffx.potential.parsers.XYZFileFilter.isXYZ;
129 import static ffx.utilities.PropertyGroup.NonBondedCutoff;
130 import static ffx.utilities.PropertyGroup.PotentialFunctionSelection;
131 import static java.lang.Double.isInfinite;
132 import static java.lang.Double.isNaN;
133 import static java.lang.String.format;
134 import static org.apache.commons.io.FilenameUtils.removeExtension;
135 import static org.apache.commons.math3.util.FastMath.max;
136 import static org.apache.commons.math3.util.FastMath.min;
137 import static org.apache.commons.math3.util.FastMath.sqrt;
138 
139 /**
140  * Compute the potential energy and derivatives of a molecular system described by a force field.
141  *
142  * @author Michael J. Schnieders
143  * @since 1.0
144  */
145 public class ForceFieldEnergy implements CrystalPotential, LambdaInterface {
146 
147   /**
148    * A Logger for the ForceFieldEnergy class.
149    */
150   private static final Logger logger = Logger.getLogger(ForceFieldEnergy.class.getName());
151 
152   /**
153    * Convert from nanoseconds to seconds.
154    */
155   private static final double toSeconds = 1.0e-9;
156 
157   /**
158    * RestrainMode specifies how restrain terms are applied for dual topology calculations.
159    */
160   public enum RestrainMode {
161     /**
162      * Restrain terms are applied with the energy.
163      */
164     ENERGY,
165     /**
166      * Restrain terms are applied with the alchemical restraints (Lambda Bonded Terms).
167      */
168     ALCHEMICAL
169   }
170 
171   /**
172    * The ForceFieldEnergy class implements the FFX platform.
173    */
174   private final Platform platform = Platform.FFX;
175 
176   /**
177    * The MolecularAssembly associated with this force field energy.
178    */
179   protected final MolecularAssembly molecularAssembly;
180   /**
181    * The array of Atoms being evaluated.
182    */
183   private final Atom[] atoms;
184   /**
185    * The boundary conditions used when evaluating the force field energy.
186    */
187   private Crystal crystal;
188   /**
189    * The Parallel Java ParallelTeam instance.
190    */
191   private final ParallelTeam parallelTeam;
192   /**
193    * An NCS restraint term.
194    */
195   private final NCSRestraint ncsRestraint;
196   /**
197    * A Center-of-Mass restraint term.
198    */
199   private final COMRestraint comRestraint;
200   /**
201    * Non-Bonded van der Waals energy.
202    */
203   private final VanDerWaals vanderWaals;
204   /**
205    * Particle-Mesh Ewald electrostatic energy.
206    */
207   private final ParticleMeshEwald particleMeshEwald;
208   /**
209    * ANI2x Neutral Network Potential.
210    */
211   private final ANIEnergy aniEnergy;
212 
213   /**
214    * List of constraints to be applied to the system.
215    */
216   private final List<Constraint> constraints;
217   /**
218    * Default tolerance for numerical methods of solving constraints.
219    */
220   public static final double DEFAULT_CONSTRAINT_TOLERANCE = 1E-4;
221 
222   @FFXProperty(name = "vdw-cutoff", propertyGroup = NonBondedCutoff, defaultValue = "12.0", description = """
223       Sets the cutoff distance value in Angstroms for van der Waals potential energy interactions.
224       The energy for any pair of van der Waals sites beyond the cutoff distance will be set to zero.
225       Other properties can be used to define the smoothing scheme near the cutoff distance.
226       The default cutoff distance in the absence of the vdw-cutoff keyword is infinite for nonperiodic
227       systems and 12.0 for periodic systems.
228       """)
229   private double vdwCutoff;
230 
231   @FFXProperty(name = "ewald-cutoff", propertyGroup = NonBondedCutoff, defaultValue = "7.0", description = """
232       Sets the value in Angstroms of the real-space distance cutoff for use during Ewald summation.
233       By default, in the absence of the ewald-cutoff property, a value of 7.0 is used.
234       """)
235   private double ewaldCutoff;
236 
237   @FFXProperty(name = "gk-cutoff", propertyGroup = NonBondedCutoff, defaultValue = "12.0", description = """
238       Sets the value in Angstroms of the generalized Kirkwood distance cutoff for use
239       during implicit solvent simulations. By default, in the absence of the gk-cutoff property,
240       no cutoff is used under aperiodic boundary conditions and the vdw-cutoff is used under PBC.
241       """)
242   private double gkCutoff;
243 
244   private static final double DEFAULT_LIST_BUFFER = 2.0;
245 
246   @FFXProperty(name = "list-buffer", propertyGroup = NonBondedCutoff, defaultValue = "2.0", description = """
247       Sets the size of the neighbor list buffer in Angstroms for potential energy functions.
248       This value is added to the actual cutoff distance to determine which pairs will be kept on the neighbor list.
249       This buffer value is used for all potential function neighbor lists.
250       The default value in the absence of the list-buffer keyword is 2.0 Angstroms.
251       """)
252   private double listBuffer;
253 
254   /**
255    * 2.0 times the neighbor list cutoff.
256    */
257   private final double cutOff2;
258   /**
259    * The non-bonded cut-off plus buffer distance (Angstroms).
260    */
261   private final double cutoffPlusBuffer;
262   /**
263    * The Electrostatic functional form.
264    */
265   private final ELEC_FORM elecForm;
266 
267   /**
268    * Indicates use of the Lambda state variable.
269    */
270   @FFXProperty(name = "lambdaterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
271       defaultValue = "false", description = "Specifies use of the Lambda state variable.")
272   protected boolean lambdaTerm;
273   /**
274    * Current value of the Lambda state variable.
275    */
276   private double lambda = 1.0;
277   /**
278    * Indicates only bonded energy terms effected by Lambda should be evaluated.
279    */
280   public boolean lambdaBondedTerms = false;
281   /**
282    * Indicates all bonded energy terms should be evaluated if lambdaBondedTerms is true.
283    */
284   boolean lambdaAllBondedTerms = false;
285   /**
286    * Indicates application of lambda scaling to all Torsion based energy terms.
287    */
288   private final boolean lambdaTorsions;
289 
290   /**
291    * An instance of the STATE enumeration to specify calculation of slowly varying energy terms, fast
292    * varying or both.
293    */
294   private STATE state = STATE.BOTH;
295   /**
296    * Optimization scaling value to use for each degree of freedom.
297    */
298   protected double[] optimizationScaling = null;
299 
300   /**
301    * Evaluate bonded energy terms in parallel.
302    */
303   private final EnergyTermRegion forceFieldBondedEnergyRegion;
304   /**
305    * Implements the Bond potential.
306    */
307   private final BondPotentialEnergy bondPotentialEnergy;
308   /**
309    * Implements the Angle potential.
310    */
311   private final AnglePotentialEnergy anglePotentialEnergy;
312   /**
313    * Implements the Stretch-Bend potential.
314    */
315   private final StretchBendPotentialEnergy stretchBendPotentialEnergy;
316   /**
317    * Implements the Urey-Bradley potential.
318    */
319   private final UreyBradleyPotentialEnergy ureyBradleyPotentialEnergy;
320   /**
321    * Implements the Out-of-Plane Bend potential.
322    */
323   private final OutOfPlaneBendPotentialEnergy outOfPlaneBendPotentialEnergy;
324   /**
325    * Implements the Torsion potential.
326    */
327   private final TorsionPotentialEnergy torsionPotentialEnergy;
328   /**
329    * Implements the Stretch-Torsion potential.
330    */
331   private final StretchTorsionPotentialEnergy stretchTorsionPotentialEnergy;
332   /**
333    * Implements the Angle-Torsion potential.
334    */
335   private final AngleTorsionPotentialEnergy angleTorsionPotentialEnergy;
336   /**
337    * Implements the Improper Torsion potential.
338    */
339   private final ImproperTorsionPotentialEnergy improperTorsionPotentialEnergy;
340   /**
341    * Implements the Pi-Orbital Torsion potential.
342    */
343   private final PiOrbitalTorsionPotentialEnergy piOrbitalTorsionPotentialEnergy;
344   /**
345    * Implements the Torsion-Torsion potential.
346    */
347   private final TorsionTorsionPotentialEnergy torsionTorsionPotentialEnergy;
348 
349   /**
350    * Evaluate bonded energy terms in parallel.
351    */
352   private final EnergyTermRegion restrainEnergyRegion;
353   /**
354    * RestrainMode specifies how restrain terms are applied for dual topology calculations.
355    */
356   private final RestrainMode restrainMode;
357   /**
358    * Implements the Restrain-Position potential.
359    */
360   private final RestrainPositionPotentialEnergy restrainPositionPotentialEnergy;
361   /**
362    * Implements the Restrain-Distance potential.
363    */
364   private final RestrainDistancePotentialEnergy restrainDistancePotentialEnergy;
365   /**
366    * Implements the Restrain-Torsion potential.
367    */
368   private final RestrainTorsionPotentialEnergy restrainTorsionPotentialEnergy;
369   /**
370    * The Restrain Groups energy term.
371    */
372   private final RestrainGroups restrainGroups;
373 
374   /**
375    * Number of atoms in the system.
376    */
377   private final int nAtoms;
378   /**
379    * Number of Restrain Groups in the system.
380    */
381   private int nRestrainGroups;
382   /**
383    * Number of van der Waals interactions evaluated.
384    */
385   private int nVanDerWaalInteractions;
386   /**
387    * Number of electrostatic interactions evaluated.
388    */
389   private int nPermanentInteractions;
390   /**
391    * Number of implicit solvent interactions evaluated.
392    */
393   private int nGKInteractions;
394 
395   /**
396    * Specifies use of an intramolecular neural network.
397    */
398   private boolean nnTerm;
399 
400   /**
401    * Specifies use of the bond stretch potential.
402    */
403   @FFXProperty(name = "bondterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
404       defaultValue = "true", description = "Specifies use of the bond stretch potential.")
405   private final boolean bondTerm;
406 
407   /**
408    * Specifies use of the angle bend potential.
409    */
410   @FFXProperty(name = "angleterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
411       defaultValue = "true", description = "Specifies use of the angle bend potential.")
412   private final boolean angleTerm;
413 
414   /**
415    * Evaluate Stretch-Bend energy terms.
416    */
417   @FFXProperty(name = "strbndterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
418       defaultValue = "true", description = "Specifies use of the stretch-bend potential.")
419   private final boolean stretchBendTerm;
420 
421   /**
422    * Evaluate Urey-Bradley energy terms.
423    */
424   @FFXProperty(name = "ureyterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
425       defaultValue = "true", description = "Specifies use of the Urey-Bradley potential.")
426   private final boolean ureyBradleyTerm;
427 
428   /**
429    * Evaluate Out of Plane Bend energy terms.
430    */
431   @FFXProperty(name = "opbendterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
432       defaultValue = "true", description = "Specifies use of the out-of-plane potential.")
433   private final boolean outOfPlaneBendTerm;
434 
435   /**
436    * Evaluate Torsion energy terms.
437    */
438   @FFXProperty(name = "torsionterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
439       defaultValue = "true", description = "Specifies use of the torsional potential.")
440   private final boolean torsionTerm;
441 
442   /**
443    * Evaluate Stretch-Torsion energy terms.
444    */
445   @FFXProperty(name = "strtorterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
446       defaultValue = "true", description = "Specifies use of the stretch-torsion potential.")
447   private final boolean stretchTorsionTerm;
448 
449   /**
450    * Evaluate Angle-Torsion energy terms.
451    */
452   @FFXProperty(name = "angtorterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
453       defaultValue = "true", description = "Specifies use of the angle-torsion potential.")
454   private final boolean angleTorsionTerm;
455 
456   /**
457    * Evaluate Improper Torsion energy terms.
458    */
459   @FFXProperty(name = "imptorterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
460       defaultValue = "true", description = "Specifies use of the improper torsion potential.")
461   private final boolean improperTorsionTerm;
462 
463   /**
464    * Evaluate Pi-Orbital Torsion energy terms.
465    */
466   @FFXProperty(name = "pitorsterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
467       defaultValue = "true", description = "Specifies use of the pi-system torsion potential.")
468   private final boolean piOrbitalTorsionTerm;
469 
470   /**
471    * Evaluate Torsion-Torsion energy terms.
472    */
473   @FFXProperty(name = "tortorterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
474       defaultValue = "true", description = "Specifies use of the pi-system torsion potential.")
475   private final boolean torsionTorsionTerm;
476 
477   /**
478    * Evaluate RestrainPosition term.
479    */
480   private final boolean restrainPositionTerm;
481   /**
482    * Evaluate RestrainDistance terms.
483    */
484   private final boolean restrainDistanceTerm;
485   /**
486    * Evaluate RestrainTorsion terms.
487    */
488   private final boolean restrainTorsionTerm;
489   /**
490    * Evaluate RestrainGroup terms.
491    */
492   private boolean restrainGroupTerm;
493 
494   /**
495    * Evaluate van der Waals energy term.
496    */
497   @FFXProperty(name = "vdwterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
498       defaultValue = "true", description = """
499       Specifies use of the vdw der Waals potential.
500       If set to false, all non-bonded terms are turned off.
501       """)
502   private boolean vanderWaalsTerm;
503 
504   /**
505    * Evaluate permanent multipole electrostatics energy term.
506    */
507   @FFXProperty(name = "mpoleterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
508       defaultValue = "true", description = """
509       Specifies use of the fixed charge electrostatic potential.
510       Setting mpoleterm to false also turns off polarization and generalized Kirkwood,
511       overriding the polarizeterm and gkterm properties.
512       """)
513   private boolean multipoleTerm;
514 
515   /**
516    * Evaluate polarization energy term.
517    */
518   @FFXProperty(name = "polarizeterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
519       defaultValue = "true", description = """
520       Specifies use of the polarizable electrostatic potential.
521       Setting polarizeterm to false overrides the polarization property.
522       """)
523   private boolean polarizationTerm;
524 
525   /**
526    * Evaluate generalized Kirkwood energy term.
527    */
528   @FFXProperty(name = "gkterm", clazz = Boolean.class, propertyGroup = PotentialFunctionSelection,
529       defaultValue = "false", description = "Specifies use of generalized Kirkwood electrostatics.")
530   private boolean generalizedKirkwoodTerm;
531 
532   /**
533    * Evaluate COM energy term.
534    */
535   private boolean comTerm;
536   /**
537    * Evaluate NCS energy term.
538    */
539   private boolean ncsTerm;
540 
541   /**
542    * Original state of the neural network energy term flag.
543    */
544   private final boolean nnTermOrig;
545   /**
546    * Original state of the RestrainGroup term flag.
547    */
548   private final boolean restrainGroupTermOrig;
549   /**
550    * Original state of the van der Waals energy term flag.
551    */
552   private final boolean vanderWaalsTermOrig;
553   /**
554    * Original state of the multipole energy term flag.
555    */
556   private final boolean multipoleTermOrig;
557   /**
558    * Original state of the polarization energy term flag.
559    */
560   private final boolean polarizationTermOrig;
561   /**
562    * Original state of the GK energy term flag.
563    */
564   private final boolean generalizedKirkwoodTermOrig;
565   /**
566    * Original state of the ESV energy term flag.
567    */
568   private boolean esvTermOrig;
569   /**
570    * Original state of the NCS energy term flag.
571    */
572   private boolean ncsTermOrig;
573   /**
574    * Original state of the COM energy term flag.
575    */
576   private boolean comTermOrig;
577 
578   /**
579    * The total neutral network energy.
580    */
581   private double nnEnergy;
582   /**
583    * The RestrainGroup Energy.
584    */
585   private double restrainGroupEnergy;
586   /**
587    * The total van der Waals energy.
588    */
589   private double vanDerWaalsEnergy;
590   /**
591    * The permanent Multipole energy.
592    */
593   private double permanentMultipoleEnergy;
594   /**
595    * The total polarization energy.
596    */
597   private double polarizationEnergy;
598   /**
599    * The solvation energy (GK).
600    */
601   private double solvationEnergy;
602   /**
603    * The total NCS Energy.
604    */
605   private double ncsEnergy;
606   /**
607    * The total COM Restraint Energy.
608    */
609   private double comRestraintEnergy;
610   /**
611    * The total system energy.
612    */
613   private double totalEnergy;
614 
615   /**
616    * Time to evaluate NCS term.
617    */
618   private long ncsTime;
619   /**
620    * Time to evaluate the neutral network.
621    */
622   private long nnTime;
623   /**
624    * Time to evaluate restrain group term.
625    */
626   private long restrainGroupTime;
627   /**
628    * Time to evaluate Center of Mass restraint term.
629    */
630   private long comRestraintTime;
631   /**
632    * Time to evaluate van der Waals term.
633    */
634   private long vanDerWaalsTime;
635   /**
636    * Time to evaluate electrostatics term.
637    */
638   private long electrostaticTime;
639   /**
640    * Time to evaluate all energy terms.
641    */
642   private long totalTime;
643 
644   /**
645    * Constant pH extended system (TODO: needs further testing).
646    */
647   private ExtendedSystem esvSystem = null;
648   private boolean esvTerm;
649   private double esvBias;
650 
651   /**
652    * Flag to indicate proper shutdown of the ForceFieldEnergy.
653    */
654   boolean destroyed = false;
655   /**
656    * If the absolute value of a gradient component is greater than "maxDebugGradient", verbose
657    * logging results.
658    */
659   public final double maxDebugGradient;
660   /**
661    * Enable verbose printing if large energy gradient components are observed.
662    */
663   private boolean printOnFailure;
664 
665   /**
666    * Constructor for ForceFieldEnergy.
667    *
668    * @param molecularAssembly a {@link ffx.potential.MolecularAssembly} object.
669    */
670   protected ForceFieldEnergy(MolecularAssembly molecularAssembly) {
671     this(molecularAssembly, ParallelTeam.getDefaultThreadCount());
672   }
673 
674   /**
675    * Constructor for ForceFieldEnergy.
676    *
677    * @param molecularAssembly a {@link ffx.potential.MolecularAssembly} object.
678    * @param numThreads        a int.
679    */
680   protected ForceFieldEnergy(MolecularAssembly molecularAssembly, int numThreads) {
681     // Get a reference to the sorted atom array.
682     this.molecularAssembly = molecularAssembly;
683     atoms = molecularAssembly.getAtomArray();
684     nAtoms = atoms.length;
685 
686     // Check that atom ordering is correct and count the number of active atoms.
687     for (int i = 0; i < nAtoms; i++) {
688       int index = atoms[i].getXyzIndex() - 1;
689       assert (i == index);
690     }
691 
692     // Enforce that the number of threads be less than or equal to the number of atoms.
693     int nThreads = min(nAtoms, numThreads);
694     parallelTeam = new ParallelTeam(nThreads);
695 
696     ForceField forceField = molecularAssembly.getForceField();
697     CompositeConfiguration properties = forceField.getProperties();
698     String name = forceField.toString().toUpperCase();
699 
700     logger.info(format(" Constructing Force Field %s", name));
701     logger.info(format("\n SMP threads:                        %10d", nThreads));
702 
703     // Rename water protons if requested.
704     boolean standardizeAtomNames = properties.getBoolean("standardizeAtomNames", true);
705     if (standardizeAtomNames) {
706       molecularAssembly.renameWaterProtons();
707     }
708 
709     nnTerm = forceField.getBoolean("NNTERM", false);
710     // For now, all atoms are neural network atoms, or none are.
711     if (nnTerm) {
712       for (Atom atom : atoms) {
713         atom.setNeuralNetwork(true);
714       }
715     }
716     bondTerm = forceField.getBoolean("BONDTERM", true);
717     angleTerm = forceField.getBoolean("ANGLETERM", true);
718     stretchBendTerm = forceField.getBoolean("STRBNDTERM", true);
719     ureyBradleyTerm = forceField.getBoolean("UREYTERM", true);
720     outOfPlaneBendTerm = forceField.getBoolean("OPBENDTERM", true);
721     torsionTerm = forceField.getBoolean("TORSIONTERM", true);
722     stretchTorsionTerm = forceField.getBoolean("STRTORTERM", true);
723     angleTorsionTerm = forceField.getBoolean("ANGTORTERM", true);
724     piOrbitalTorsionTerm = forceField.getBoolean("PITORSTERM", true);
725     torsionTorsionTerm = forceField.getBoolean("TORTORTERM", true);
726     improperTorsionTerm = forceField.getBoolean("IMPROPERTERM", true);
727     vanderWaalsTerm = forceField.getBoolean("VDWTERM", true);
728 
729     boolean tornadoVM = forceField.getBoolean("tornado", false);
730     if (vanderWaalsTerm && !tornadoVM) {
731       multipoleTerm = forceField.getBoolean("MPOLETERM", true);
732       if (multipoleTerm) {
733         String polarizeString = forceField.getString("POLARIZATION", "NONE");
734         boolean defaultPolarizeTerm = !polarizeString.equalsIgnoreCase("NONE");
735         polarizationTerm = forceField.getBoolean("POLARIZETERM", defaultPolarizeTerm);
736         generalizedKirkwoodTerm = forceField.getBoolean("GKTERM", false);
737       } else {
738         // If multipole electrostatics is turned off, turn off all electrostatics.
739         polarizationTerm = false;
740         generalizedKirkwoodTerm = false;
741       }
742     } else {
743       // If van der Waals is turned off, turn off all non-bonded terms.
744       multipoleTerm = false;
745       polarizationTerm = false;
746       generalizedKirkwoodTerm = false;
747     }
748 
749     lambdaTerm = forceField.getBoolean("LAMBDATERM", false);
750     comTerm = forceField.getBoolean("COMRESTRAINTERM", false);
751     lambdaTorsions = forceField.getBoolean("TORSION_LAMBDATERM", false);
752     printOnFailure = forceField.getBoolean("PRINT_ON_FAILURE", false);
753 
754     String mode = forceField.getString("RESTRAIN_MODE", "ENERGY");
755     if (mode.equalsIgnoreCase("ALCHEMICAL")) {
756       restrainMode = RestrainMode.ALCHEMICAL;
757     } else {
758       restrainMode = RestrainMode.ENERGY;
759     }
760 
761     // Detect "restrain-groups" property.
762     if (properties.containsKey("restrain-groups")) {
763       restrainGroupTerm = true;
764     } else {
765       restrainGroupTerm = false;
766     }
767 
768     // For RESPA
769     nnTermOrig = nnTerm;
770     restrainGroupTermOrig = restrainGroupTerm;
771     vanderWaalsTermOrig = vanderWaalsTerm;
772     multipoleTermOrig = multipoleTerm;
773     polarizationTermOrig = polarizationTerm;
774     generalizedKirkwoodTermOrig = generalizedKirkwoodTerm;
775     ncsTermOrig = ncsTerm;
776     comTermOrig = comTerm;
777 
778     // Determine the unit cell dimensions and space group.
779     String spacegroup;
780     double a, b, c, alpha, beta, gamma;
781     boolean aperiodic;
782     try {
783       spacegroup = forceField.getString("SPACEGROUP", "P 1");
784       SpaceGroup sg = SpaceGroupDefinitions.spaceGroupFactory(spacegroup);
785       LatticeSystem latticeSystem = sg.latticeSystem;
786       a = forceField.getDouble("A_AXIS");
787       aperiodic = false;
788       b = forceField.getDouble("B_AXIS", latticeSystem.getDefaultBAxis(a));
789       c = forceField.getDouble("C_AXIS", latticeSystem.getDefaultCAxis(a, b));
790       alpha = forceField.getDouble("ALPHA", latticeSystem.getDefaultAlpha());
791       beta = forceField.getDouble("BETA", latticeSystem.getDefaultBeta());
792       gamma = forceField.getDouble("GAMMA", latticeSystem.getDefaultGamma());
793       if (!sg.latticeSystem.validParameters(a, b, c, alpha, beta, gamma)) {
794         logger.severe(" Check lattice parameters.");
795       }
796       if (a == 1.0 && b == 1.0 && c == 1.0) {
797         String message = " A-, B-, and C-axis values equal to 1.0.";
798         logger.info(message);
799         throw new Exception(message);
800       }
801       if (a <= 0.0 || b <= 0.0 || c <= 0.0 || alpha <= 0.0 || beta <= 0.0 || gamma <= 0.0) {
802         // Parameters are not valid. Switch to aperiodic.
803         String message = " Crystal parameters are not valid due to negative or zero value.";
804         logger.warning(message);
805         throw new Exception(message);
806       }
807     } catch (Exception e) {
808       aperiodic = true;
809 
810       // Determine the maximum separation between atoms.
811       double maxR = 0.0;
812       if (nAtoms < 10) {
813         for (int i = 0; i < nAtoms - 1; i++) {
814           Double3 xi = atoms[i].getXYZ();
815           for (int j = 1; j < nAtoms; j++) {
816             double r = atoms[j].getXYZ().dist(xi);
817             maxR = max(r, maxR);
818           }
819         }
820       } else {
821         try {
822           maxR = ConvexHullOps.maxDist(ConvexHullOps.constructHull(atoms));
823         } catch (Exception ex) {
824           // If Convex Hull approach fails (e.g., coplanar input, brute force...
825           logger.info(" Convex Hull operation failed with message " + ex + "\n Trying brute force approach...");
826           if (logger.isLoggable(Level.FINE)) {
827             logger.fine(Utilities.stackTraceToString(ex));
828           }
829           for (int i = 0; i < nAtoms - 1; i++) {
830             Double3 xi = atoms[i].getXYZ();
831             for (int j = 1; j < nAtoms; j++) {
832               double r = atoms[j].getXYZ().dist(xi);
833               maxR = max(r, maxR);
834             }
835           }
836         }
837       }
838       maxR = max(10.0, maxR);
839 
840       logger.info(format(" The system will be treated as aperiodic (max separation = %6.1f A).", maxR));
841 
842       // Turn off reciprocal space calculations.
843       forceField.addProperty("EWALD_ALPHA", "0.0");
844 
845       // Specify some dummy values for the crystal.
846       spacegroup = "P1";
847       a = 2.0 * maxR;
848       b = 2.0 * maxR;
849       c = 2.0 * maxR;
850       alpha = 90.0;
851       beta = 90.0;
852       gamma = 90.0;
853     }
854     Crystal unitCell = new Crystal(a, b, c, alpha, beta, gamma, spacegroup);
855     unitCell.setAperiodic(aperiodic);
856 
857     // First set cutoffs assuming PBC.
858     vdwCutoff = DEFAULT_VDW_CUTOFF;
859     ewaldCutoff = DEFAULT_EWALD_CUTOFF;
860     gkCutoff = vdwCutoff;
861     if (unitCell.aperiodic()) {
862       // Default to no cutoffs for aperiodic systems.
863       vdwCutoff = Double.POSITIVE_INFINITY;
864       ewaldCutoff = Double.POSITIVE_INFINITY;
865       gkCutoff = Double.POSITIVE_INFINITY;
866     }
867     // Load user specified cutoffs.
868     vdwCutoff = forceField.getDouble("VDW_CUTOFF", vdwCutoff);
869     ewaldCutoff = forceField.getDouble("EWALD_CUTOFF", ewaldCutoff);
870     gkCutoff = forceField.getDouble("GK_CUTOFF", gkCutoff);
871 
872     // Neighbor list cutoff is the maximum of the vdW, Ewald and GK cutoffs (i.e. to use a single list).
873     double nlistCutoff = vdwCutoff;
874     if (multipoleTerm) {
875       nlistCutoff = max(vdwCutoff, ewaldCutoff);
876     }
877     if (generalizedKirkwoodTerm) {
878       nlistCutoff = max(nlistCutoff, gkCutoff);
879     }
880 
881     // Check for a frozen neighbor list.
882     boolean disabledNeighborUpdates = forceField.getBoolean("DISABLE_NEIGHBOR_UPDATES", false);
883     if (disabledNeighborUpdates) {
884       logger.info(format(" Neighbor list updates disabled; interactions will "
885           + "only be calculated between atoms that started the simulation "
886           + "within a radius of %9.3g Angstroms of each other.", nlistCutoff));
887       // The individual cutoffs are infinity (i.e. they are limited by interactions in the list).
888       vdwCutoff = Double.POSITIVE_INFINITY;
889       ewaldCutoff = Double.POSITIVE_INFINITY;
890       gkCutoff = Double.POSITIVE_INFINITY;
891     }
892 
893     listBuffer = forceField.getDouble("LIST_BUFFER", DEFAULT_LIST_BUFFER);
894     cutoffPlusBuffer = nlistCutoff + listBuffer;
895     cutOff2 = 2.0 * cutoffPlusBuffer;
896     unitCell = configureNCS(forceField, unitCell);
897 
898     // If necessary, create a ReplicatesCrystal.
899     if (!aperiodic) {
900       this.crystal = ReplicatesCrystal.replicatesCrystalFactory(unitCell, cutOff2);
901       logger.info(format("\n Density:                                %6.3f (g/cc)",
902           crystal.getDensity(molecularAssembly.getMass())));
903       logger.info(crystal.toString());
904     } else {
905       this.crystal = unitCell;
906     }
907 
908     if (!unitCell.aperiodic() && unitCell.spaceGroup.number == 1) {
909       ncsTerm = forceField.getBoolean("NCSTERM", false);
910       ncsTermOrig = ncsTerm;
911     } else {
912       ncsTerm = false;
913       ncsTermOrig = false;
914     }
915 
916     // Now that the crystal is set, check for special positions.
917     int nSpecial = checkForSpecialPositions(forceField);
918 
919     // Flag to indicate bonded terms with hydrogen should be scaled up.
920     boolean rigidHydrogens = forceField.getBoolean("RIGID_HYDROGENS", false);
921     // Scale factor for increasing the strength of bonded terms involving hydrogen atoms.
922     double rigidScale = forceField.getDouble("RIGID_SCALE", 10.0);
923 
924     boolean checkAllNodeCharges = forceField.getBoolean("CHECK_ALL_NODE_CHARGES", false);
925 
926     if (rigidScale <= 1.0) {
927       rigidScale = 1.0;
928     }
929 
930     logger.info("\n Bonded Terms");
931     if (rigidHydrogens && rigidScale > 1.0) {
932       logger.info(format("  Rigid hydrogens:                   %10.2f", rigidScale));
933     }
934 
935     forceFieldBondedEnergyRegion = new EnergyTermRegion(parallelTeam, molecularAssembly, lambdaTerm);
936     restrainEnergyRegion = new EnergyTermRegion(parallelTeam, molecularAssembly, lambdaTerm);
937 
938     // Collect bonds.
939     if (bondTerm) {
940       List<Bond> bondList = molecularAssembly.getBondList();
941       if (nnTerm) {
942         removeNeuralNetworkTerms(bondList);
943       }
944       if (!bondList.isEmpty()) {
945         int forceGroup = forceField.getInteger("BOND_FORCE_GROUP", 0);
946         bondPotentialEnergy = new BondPotentialEnergy("Bonds", forceGroup, bondList);
947         forceFieldBondedEnergyRegion.addEnergyTerm(bondPotentialEnergy);
948       } else {
949         bondPotentialEnergy = null;
950       }
951     } else {
952       bondPotentialEnergy = null;
953     }
954 
955     // Collect angles.
956     if (angleTerm) {
957       List<Angle> angleList = molecularAssembly.getAngleList();
958       if (nnTerm) {
959         removeNeuralNetworkTerms(angleList);
960       }
961       if (!angleList.isEmpty()) {
962         int forceGroup = forceField.getInteger("ANGLE_FORCE_GROUP", 0);
963         anglePotentialEnergy = new AnglePotentialEnergy("Angles", forceGroup, angleList);
964         forceFieldBondedEnergyRegion.addEnergyTerm(anglePotentialEnergy);
965       } else {
966         anglePotentialEnergy = null;
967       }
968     } else {
969       anglePotentialEnergy = null;
970     }
971 
972     // Collect stretch-bends.
973     if (stretchBendTerm) {
974       List<StretchBend> stretchBendList = molecularAssembly.getStretchBendList();
975       if (nnTerm) {
976         removeNeuralNetworkTerms(stretchBendList);
977       }
978       if (!stretchBendList.isEmpty()) {
979         int forceGroup = forceField.getInteger("STRETCH_BEND_FORCE_GROUP", 0);
980         stretchBendPotentialEnergy = new StretchBendPotentialEnergy("Stretch-Bends", forceGroup, stretchBendList);
981         forceFieldBondedEnergyRegion.addEnergyTerm(stretchBendPotentialEnergy);
982       } else {
983         stretchBendPotentialEnergy = null;
984       }
985     } else {
986       stretchBendPotentialEnergy = null;
987     }
988 
989     // Collect Urey-Bradleys.
990     if (ureyBradleyTerm) {
991       List<UreyBradley> ureyBradleyList = molecularAssembly.getUreyBradleyList();
992       if (nnTerm) {
993         removeNeuralNetworkTerms(ureyBradleyList);
994       }
995       if (!ureyBradleyList.isEmpty()) {
996         int forceGroup = forceField.getInteger("UREY_BRADLEY_FORCE_GROUP", 0);
997         ureyBradleyPotentialEnergy = new UreyBradleyPotentialEnergy("Urey-Bradleys", forceGroup, ureyBradleyList);
998         forceFieldBondedEnergyRegion.addEnergyTerm(ureyBradleyPotentialEnergy);
999       } else {
1000         ureyBradleyPotentialEnergy = null;
1001       }
1002     } else {
1003       ureyBradleyPotentialEnergy = null;
1004     }
1005 
1006     // Set a multiplier on the force constants of bonded terms containing hydrogen atoms.
1007     if (rigidHydrogens) {
1008       if (bondPotentialEnergy != null) {
1009         for (Bond bond : bondPotentialEnergy.getBonds()) {
1010           if (bond.containsHydrogen()) {
1011             bond.setRigidScale(rigidScale);
1012           }
1013         }
1014       }
1015       if (anglePotentialEnergy != null) {
1016         for (Angle angle : anglePotentialEnergy.getAngles()) {
1017           if (angle.containsHydrogen()) {
1018             angle.setRigidScale(rigidScale);
1019           }
1020         }
1021       }
1022       if (stretchBendPotentialEnergy != null) {
1023         for (StretchBend stretchBend : stretchBendPotentialEnergy.getStretchBends()) {
1024           if (stretchBend.containsHydrogen()) {
1025             stretchBend.setRigidScale(rigidScale);
1026           }
1027         }
1028       }
1029       if (ureyBradleyPotentialEnergy != null) {
1030         for (UreyBradley ureyBradley : ureyBradleyPotentialEnergy.getUreyBradleys()) {
1031           if (ureyBradley.containsHydrogen()) {
1032             ureyBradley.setRigidScale(rigidScale);
1033           }
1034         }
1035       }
1036     }
1037 
1038     // Collect out-of-plane bends.
1039     if (outOfPlaneBendTerm) {
1040       List<OutOfPlaneBend> outOfPlaneBendList = molecularAssembly.getOutOfPlaneBendList();
1041       if (nnTerm) {
1042         removeNeuralNetworkTerms(outOfPlaneBendList);
1043       }
1044       if (!outOfPlaneBendList.isEmpty()) {
1045         int forceGroup = forceField.getInteger("OUT_OF_PLANE_BEND_FORCE_GROUP", 0);
1046         outOfPlaneBendPotentialEnergy = new OutOfPlaneBendPotentialEnergy("Out-0f-Plane Bends", forceGroup, outOfPlaneBendList);
1047         forceFieldBondedEnergyRegion.addEnergyTerm(outOfPlaneBendPotentialEnergy);
1048       } else {
1049         outOfPlaneBendPotentialEnergy = null;
1050       }
1051     } else {
1052       outOfPlaneBendPotentialEnergy = null;
1053     }
1054 
1055     // Collect torsions.
1056     if (torsionTerm) {
1057       List<Torsion> torsionList = molecularAssembly.getTorsionList();
1058       if (nnTerm) {
1059         removeNeuralNetworkTerms(torsionList);
1060       }
1061       if (!torsionList.isEmpty()) {
1062         int forceGroup = forceField.getInteger("TORSION_FORCE_GROUP", 0);
1063         torsionPotentialEnergy = new TorsionPotentialEnergy("Torsions", forceGroup, torsionList);
1064         forceFieldBondedEnergyRegion.addEnergyTerm(torsionPotentialEnergy);
1065       } else
1066         torsionPotentialEnergy = null;
1067     } else {
1068       torsionPotentialEnergy = null;
1069     }
1070 
1071     // Collect stretch-torsions.
1072     if (stretchTorsionTerm) {
1073       List<StretchTorsion> stretchTorsionList = molecularAssembly.getStretchTorsionList();
1074       if (nnTerm) {
1075         removeNeuralNetworkTerms(stretchTorsionList);
1076       }
1077       if (!stretchTorsionList.isEmpty()) {
1078         int forceGroup = forceField.getInteger("STRETCH_TORSION_FORCE_GROUP", 0);
1079         stretchTorsionPotentialEnergy = new StretchTorsionPotentialEnergy("Stretch-Torsions", forceGroup, stretchTorsionList);
1080         forceFieldBondedEnergyRegion.addEnergyTerm(stretchTorsionPotentialEnergy);
1081       } else {
1082         stretchTorsionPotentialEnergy = null;
1083       }
1084     } else {
1085       stretchTorsionPotentialEnergy = null;
1086     }
1087 
1088     // Collect angle-torsions.
1089     if (angleTorsionTerm) {
1090       List<AngleTorsion> angleTorsionList = molecularAssembly.getAngleTorsionList();
1091       if (nnTerm) {
1092         removeNeuralNetworkTerms(angleTorsionList);
1093       }
1094       if (!angleTorsionList.isEmpty()) {
1095         int forceGroup = forceField.getInteger("ANGLE_TORSION_FORCE_GROUP", 0);
1096         angleTorsionPotentialEnergy = new AngleTorsionPotentialEnergy("Angle-Torsions", forceGroup, angleTorsionList);
1097         forceFieldBondedEnergyRegion.addEnergyTerm(angleTorsionPotentialEnergy);
1098       } else {
1099         angleTorsionPotentialEnergy = null;
1100       }
1101     } else {
1102       angleTorsionPotentialEnergy = null;
1103     }
1104 
1105     // Collect improper torsions.
1106     if (improperTorsionTerm) {
1107       List<ImproperTorsion> improperTorsionList = molecularAssembly.getImproperTorsionList();
1108       if (nnTerm) {
1109         removeNeuralNetworkTerms(improperTorsionList);
1110       }
1111       if (!improperTorsionList.isEmpty()) {
1112         int forceGroup = forceField.getInteger("IMPROPER_TORSION_FORCE_GROUP", 0);
1113         improperTorsionPotentialEnergy = new ImproperTorsionPotentialEnergy("Improper Torsions", forceGroup, improperTorsionList);
1114         forceFieldBondedEnergyRegion.addEnergyTerm(improperTorsionPotentialEnergy);
1115       } else {
1116         improperTorsionPotentialEnergy = null;
1117       }
1118     } else {
1119       improperTorsionPotentialEnergy = null;
1120     }
1121 
1122     // Collect pi-orbital torsions.
1123     if (piOrbitalTorsionTerm) {
1124       List<PiOrbitalTorsion> piOrbitalTorsionList = molecularAssembly.getPiOrbitalTorsionList();
1125       if (nnTerm) {
1126         removeNeuralNetworkTerms(piOrbitalTorsionList);
1127       }
1128       if (!piOrbitalTorsionList.isEmpty()) {
1129         int forceGroup = forceField.getInteger("PI_ORBITAL_TORSION_FORCE_GROUP", 0);
1130         piOrbitalTorsionPotentialEnergy = new PiOrbitalTorsionPotentialEnergy("Pi-Orbital Torsions", forceGroup, piOrbitalTorsionList);
1131         forceFieldBondedEnergyRegion.addEnergyTerm(piOrbitalTorsionPotentialEnergy);
1132       } else {
1133         piOrbitalTorsionPotentialEnergy = null;
1134       }
1135     } else {
1136       piOrbitalTorsionPotentialEnergy = null;
1137     }
1138 
1139     // Collect torsion-torsions.
1140     if (torsionTorsionTerm) {
1141       List<TorsionTorsion> torsionTorsionList = molecularAssembly.getTorsionTorsionList();
1142       if (nnTerm) {
1143         removeNeuralNetworkTerms(torsionTorsionList);
1144       }
1145       if (!torsionTorsionList.isEmpty()) {
1146         int forceGroup = forceField.getInteger("TORSION_TORSION_FORCE_GROUP", 0);
1147         torsionTorsionPotentialEnergy = new TorsionTorsionPotentialEnergy("Torsion-Torsions", forceGroup, torsionTorsionList);
1148         forceFieldBondedEnergyRegion.addEnergyTerm(torsionTorsionPotentialEnergy);
1149       } else {
1150         torsionTorsionPotentialEnergy = null;
1151       }
1152     } else {
1153       torsionTorsionPotentialEnergy = null;
1154     }
1155 
1156     // Collect restrain positions.
1157     RestrainPosition[] restrainPositions = parseRestrainPositions(molecularAssembly);
1158     if (restrainPositions != null && restrainPositions.length > 0) {
1159       restrainPositionTerm = true;
1160     } else {
1161       restrainPositionTerm = false;
1162     }
1163     if (restrainPositionTerm) {
1164       int forceGroup = forceField.getInteger("RESTRAIN_POSITION_FORCE_GROUP", 0);
1165       restrainPositionPotentialEnergy = new RestrainPositionPotentialEnergy("Restrain Positions", forceGroup, java.util.Arrays.asList(restrainPositions));
1166       restrainEnergyRegion.addEnergyTerm(restrainPositionPotentialEnergy);
1167     } else {
1168       restrainPositionPotentialEnergy = null;
1169     }
1170 
1171     // Collect restrain distances.
1172     RestrainDistance[] restrainDistances = configureRestrainDistances(properties);
1173     if (restrainDistances != null && restrainDistances.length > 0) {
1174       restrainDistanceTerm = true;
1175     } else {
1176       restrainDistanceTerm = false;
1177     }
1178     if (restrainDistanceTerm) {
1179       int forceGroup = forceField.getInteger("RESTRAIN_DISTANCE_FORCE_GROUP", 0);
1180       restrainDistancePotentialEnergy = new RestrainDistancePotentialEnergy("Restrain Distances", forceGroup, java.util.Arrays.asList(restrainDistances));
1181       restrainEnergyRegion.addEnergyTerm(restrainDistancePotentialEnergy);
1182     } else {
1183       restrainDistancePotentialEnergy = null;
1184     }
1185 
1186     // Collect restrain torsions.
1187     Torsion[] restrainTorsions = configureRestrainTorsions(properties, forceField);
1188     if (restrainTorsions != null && restrainTorsions.length > 0) {
1189       restrainTorsionTerm = true;
1190     } else {
1191       restrainTorsionTerm = false;
1192     }
1193     if (restrainTorsionTerm) {
1194       logger.info(" Restrain Torsion Mode: " + restrainMode);
1195       int forceGroup = forceField.getInteger("RESTRAIN_TORSION_FORCE_GROUP", 0);
1196       restrainTorsionPotentialEnergy = new RestrainTorsionPotentialEnergy("Restrain Torsions", forceGroup, java.util.Arrays.asList(restrainTorsions));
1197       restrainEnergyRegion.addEnergyTerm(restrainTorsionPotentialEnergy);
1198     } else {
1199       restrainTorsionPotentialEnergy = null;
1200     }
1201 
1202     int[] molecule = molecularAssembly.getMoleculeNumbers();
1203     if (vanderWaalsTerm) {
1204       logger.info("\n Non-Bonded Terms");
1205       boolean[] nn = null;
1206       if (nnTerm) {
1207         nn = molecularAssembly.getNeuralNetworkIdentity();
1208       } else {
1209         nn = new boolean[nAtoms];
1210         Arrays.fill(nn, false);
1211       }
1212 
1213       if (!tornadoVM) {
1214         vanderWaals = new VanDerWaals(atoms, molecule, nn, crystal, forceField, parallelTeam,
1215             vdwCutoff, nlistCutoff);
1216       } else {
1217         vanderWaals = new VanDerWaalsTornado(atoms, crystal, forceField, vdwCutoff);
1218       }
1219     } else {
1220       vanderWaals = null;
1221     }
1222 
1223     if (nnTerm) {
1224       aniEnergy = new ANIEnergy(molecularAssembly);
1225     } else {
1226       aniEnergy = null;
1227     }
1228 
1229     if (multipoleTerm) {
1230       if (name.contains("OPLS") || name.contains("AMBER") || name.contains("CHARMM")) {
1231         elecForm = ELEC_FORM.FIXED_CHARGE;
1232       } else {
1233         elecForm = ELEC_FORM.PAM;
1234       }
1235       particleMeshEwald = new ParticleMeshEwald(atoms, molecule, forceField, crystal,
1236           vanderWaals.getNeighborList(), elecForm, ewaldCutoff, gkCutoff, parallelTeam);
1237       double charge = molecularAssembly.getCharge(checkAllNodeCharges);
1238       logger.info(format("\n  Overall system charge:             %10.3f", charge));
1239     } else {
1240       elecForm = null;
1241       particleMeshEwald = null;
1242     }
1243 
1244     if (ncsTerm) {
1245       String sg = forceField.getString("NCSGROUP", "P 1");
1246       Crystal ncsCrystal = new Crystal(a, b, c, alpha, beta, gamma, sg);
1247       ncsRestraint = new NCSRestraint(atoms, forceField, ncsCrystal);
1248     } else {
1249       ncsRestraint = null;
1250     }
1251 
1252     if (restrainGroupTerm) {
1253       restrainGroups = new RestrainGroups(molecularAssembly);
1254       nRestrainGroups = restrainGroups.getNumberOfRestraints();
1255     } else {
1256       restrainGroups = null;
1257     }
1258 
1259     if (comTerm) {
1260       comRestraint = new COMRestraint(molecularAssembly);
1261     } else {
1262       comRestraint = null;
1263     }
1264 
1265     // bondedRegion = new BondedRegion();
1266     maxDebugGradient = forceField.getDouble("MAX_DEBUG_GRADIENT", Double.POSITIVE_INFINITY);
1267 
1268     molecularAssembly.setPotential(this);
1269 
1270     // Configure constraints.
1271     constraints = configureConstraints(forceField);
1272 
1273     if (lambdaTerm) {
1274       this.setLambda(1.0);
1275       if (nSpecial == 0) {
1276         // For lambda calculations (e.g., polymorph searches), turn-off special position checks.
1277         // In this case, as the crystal lattice is alchemically turned on, special position overlaps are expected.
1278         crystal.setSpecialPositionCutoff(0.0);
1279       } else {
1280         // If special positions have already been identified, leave checks on.
1281         logger.info(" Special positions checking will be performed during a lambda simulation.");
1282       }
1283     }
1284   }
1285 
1286   /**
1287    * Get the MolecularAssembly associated with this ForceFieldEnergy.
1288    *
1289    * @return a {@link ffx.potential.MolecularAssembly} object.
1290    */
1291   public MolecularAssembly getMolecularAssembly() {
1292     return molecularAssembly;
1293   }
1294 
1295   /**
1296    * Set the lambdaTerm flag.
1297    *
1298    * @param lambdaTerm The value to set.
1299    */
1300   public void setLambdaTerm(boolean lambdaTerm) {
1301     this.lambdaTerm = lambdaTerm;
1302   }
1303 
1304   /**
1305    * Get the lambdaTerm flag.
1306    *
1307    * @return lambdaTerm.
1308    */
1309   public boolean getLambdaTerm() {
1310     return lambdaTerm;
1311   }
1312 
1313   private int checkForSpecialPositions(ForceField forceField) {
1314     // Check for atoms at special positions. These should normally be set to inactive.
1315     boolean specialPositionsInactive = forceField.getBoolean("SPECIAL_POSITIONS_INACTIVE", true);
1316     int nSpecial = 0;
1317     if (specialPositionsInactive) {
1318       int nSymm = crystal.getNumSymOps();
1319       if (nSymm > 1) {
1320         SpaceGroup spaceGroup = crystal.spaceGroup;
1321         double sp2 = crystal.getSpecialPositionCutoff2();
1322         double[] mate = new double[3];
1323         StringBuilder sb = new StringBuilder("\n Atoms at Special Positions set to Inactive:\n");
1324         for (int i = 0; i < nAtoms; i++) {
1325           Atom atom = atoms[i];
1326           double[] xyz = atom.getXYZ(null);
1327           for (int iSymm = 1; iSymm < nSymm; iSymm++) {
1328             SymOp symOp = spaceGroup.getSymOp(iSymm);
1329             crystal.applySymOp(xyz, mate, symOp);
1330             double dr2 = crystal.image(xyz, mate);
1331             if (dr2 < sp2) {
1332               sb.append(
1333                   format("  %s separation with SymOp %d at %8.6f A.\n", atoms[i], iSymm, sqrt(dr2)));
1334               atom.setActive(false);
1335               nSpecial++;
1336               break;
1337             }
1338           }
1339         }
1340         if (nSpecial > 0) {
1341           logger.info(sb.toString());
1342         }
1343       }
1344     }
1345     return nSpecial;
1346   }
1347 
1348   /**
1349    * Method to parse the restrain-torsion-cos records.
1350    *
1351    * @param properties Configuration properties.
1352    * @param forceField Force field properties.
1353    * @return An array of restrain torsions, or null if none were found.
1354    */
1355   private Torsion[] configureRestrainTorsions(CompositeConfiguration properties, ForceField forceField) {
1356     StringBuilder restrainLog = new StringBuilder("\n  Restrain-Torsions");
1357 
1358     String[] restrainTorsions = properties.getStringArray("restrain-torsion");
1359     double torsionUnits = forceField.getDouble("TORSIONUNIT", TorsionType.DEFAULT_TORSION_UNIT);
1360     List<Torsion> restrainTorsionList = new ArrayList<>(restrainTorsions.length);
1361     for (String restrainString : restrainTorsions) {
1362       // Add the key back to the input line.
1363       restrainString = "restrain-torsion " + restrainString;
1364       // Split the line on the pound symbol to remove comments.
1365       String input = restrainString.split("#+")[0];
1366       // Split the line on whitespace.
1367       String[] tokens = input.trim().split(" +");
1368       // Restrain torsion records have a similar form as torsion records.
1369       // The first four tokens are atom indices instead of atom classes.
1370       TorsionType torsionType = TorsionType.parse(input, tokens);
1371       torsionType.torsionUnit = torsionUnits;
1372 
1373       // Collect the atom indices.
1374       int[] atomIndices = torsionType.atomClasses;
1375       int ai1 = atomIndices[0] - 1;
1376       int ai2 = atomIndices[1] - 1;
1377       int ai3 = atomIndices[2] - 1;
1378       int ai4 = atomIndices[3] - 1;
1379       Atom a1 = atoms[ai1];
1380       Atom a2 = atoms[ai2];
1381       Atom a3 = atoms[ai3];
1382       Atom a4 = atoms[ai4];
1383 
1384       // Collect the bonds between the atoms making up the restrain torsion.
1385       Bond firstBond = a1.getBond(a2);
1386       Bond middleBond = a2.getBond(a3);
1387       Bond lastBond = a3.getBond(a4);
1388       Torsion torsion = new Torsion(firstBond, middleBond, lastBond);
1389       torsion.setTorsionType(torsionType);
1390       restrainTorsionList.add(torsion);
1391       restrainLog.append("\n   ").append(torsion);
1392     }
1393 
1394     if (!restrainTorsionList.isEmpty()) {
1395       logger.info(restrainLog.toString());
1396       return restrainTorsionList.toArray(new Torsion[0]);
1397     } else {
1398       return null;
1399     }
1400   }
1401 
1402   /**
1403    * Method to parse the restrain-distance records.
1404    *
1405    * @param properties Configuration properties.
1406    */
1407   private RestrainDistance[] configureRestrainDistances(CompositeConfiguration properties) {
1408     List<RestrainDistance> restrainDistanceList = new ArrayList<>();
1409     String[] bondRestraints = properties.getStringArray("restrain-distance");
1410     for (String bondRest : bondRestraints) {
1411       try {
1412         String[] toks = bondRest.split("\\s+");
1413         if (toks.length < 2) {
1414           throw new IllegalArgumentException(
1415               format(" restrain-distance value %s could not be parsed!", bondRest));
1416         }
1417         // Internally, everything starts with 0, but restrain distance starts at 1, so that 1 has to
1418         // be subtracted
1419         int at1 = Integer.parseInt(toks[0]) - 1;
1420         int at2 = Integer.parseInt(toks[1]) - 1;
1421 
1422         double forceConst = 100.0;
1423         double flatBottomRadius = 0;
1424         Atom a1 = atoms[at1];
1425         Atom a2 = atoms[at2];
1426 
1427         if (toks.length > 2) {
1428           forceConst = Double.parseDouble(toks[2]);
1429         }
1430         double dist;
1431         switch (toks.length) {
1432           case 2:
1433           case 3:
1434             double[] xyz1 = new double[3];
1435             xyz1 = a1.getXYZ(xyz1);
1436             double[] xyz2 = new double[3];
1437             xyz2 = a2.getXYZ(xyz2);
1438             // Current distance between restrained atoms
1439             dist = crystal.minDistOverSymOps(xyz1, xyz2);
1440             break;
1441           case 4:
1442             dist = Double.parseDouble(toks[3]);
1443             break;
1444           case 5:
1445           default:
1446             double minDist = Double.parseDouble(toks[3]);
1447             double maxDist = Double.parseDouble(toks[4]);
1448             dist = 0.5 * (minDist + maxDist);
1449             flatBottomRadius = 0.5 * Math.abs(maxDist - minDist);
1450             break;
1451         }
1452 
1453         UnivariateSwitchingFunction switchF;
1454         double lamStart = RestrainDistance.DEFAULT_RB_LAM_START;
1455         double lamEnd = RestrainDistance.DEFAULT_RB_LAM_END;
1456         if (toks.length > 5) {
1457           int offset = 5;
1458           if (toks[5].matches("^[01](?:\\.[0-9]*)?")) {
1459             offset = 6;
1460             lamStart = Double.parseDouble(toks[5]);
1461             if (toks[6].matches("^[01](?:\\.[0-9]*)?")) {
1462               offset = 7;
1463               lamEnd = Double.parseDouble(toks[6]);
1464             }
1465           }
1466           switchF = UnivariateFunctionFactory.parseUSF(toks, offset);
1467         } else {
1468           switchF = new ConstantSwitch();
1469         }
1470 
1471         RestrainDistance restrainDistance = createRestrainDistance(a1, a2, dist,
1472             forceConst, flatBottomRadius, lamStart, lamEnd, switchF);
1473         restrainDistanceList.add(restrainDistance);
1474       } catch (Exception ex) {
1475         logger.info(format(" Exception in parsing restrain-distance: %s", ex));
1476       }
1477     }
1478     if (!restrainDistanceList.isEmpty()) {
1479       return restrainDistanceList.toArray(new RestrainDistance[0]);
1480     } else {
1481       return null;
1482     }
1483   }
1484 
1485   /**
1486    * Configure bonded constraints.
1487    *
1488    * @param forceField Force field properties.
1489    * @return List of constraints.
1490    */
1491   private List<Constraint> configureConstraints(ForceField forceField) {
1492     String constraintStrings = forceField.getString("CONSTRAIN", forceField.getString("RATTLE", null));
1493     if (constraintStrings == null) {
1494       return Collections.emptyList();
1495     }
1496 
1497     ArrayList<Constraint> constraints = new ArrayList<>();
1498     logger.info(format(" Experimental: parsing constraints option %s", constraintStrings));
1499 
1500     Set<Bond> numericBonds = new HashSet<>(1);
1501     Set<Angle> numericAngles = new HashSet<>(1);
1502 
1503     // Totally empty constrain option: constrain only X-H bonds. No other options applied.
1504     if (constraintStrings.isEmpty() || constraintStrings.matches("^\\s*$")) {
1505       // Assume constraining only X-H bonds (i.e. RIGID-HYDROGEN).
1506       logger.info(" Constraining X-H bonds.");
1507       Bond[] bonds = bondPotentialEnergy.getBondArray();
1508       numericBonds = Arrays.stream(bonds).filter(
1509           (Bond bond) -> bond.getAtom(0).getAtomicNumber() == 1
1510               || bond.getAtom(1).getAtomicNumber() == 1).collect(Collectors.toSet());
1511     } else {
1512       String[] constraintToks = constraintStrings.split("\\s+");
1513 
1514       // First, accumulate SETTLE constraints.
1515       for (String tok : constraintToks) {
1516         if (tok.equalsIgnoreCase("WATER")) {
1517           logger.info(" Constraining waters to be rigid based on angle & bonds.");
1518           // XYZ files, in particular, have waters mislabeled as generic Molecules.
1519           // First, find any such mislabeled water.
1520           Stream<MSNode> settleStream = molecularAssembly.getMolecules().stream()
1521               .filter((MSNode m) -> m.getAtomList().size() == 3).filter((MSNode m) -> {
1522                 List<Atom> atoms = m.getAtomList();
1523                 Atom O = null;
1524                 List<Atom> H = new ArrayList<>(2);
1525                 for (Atom at : atoms) {
1526                   int atN = at.getAtomicNumber();
1527                   if (atN == 8) {
1528                     O = at;
1529                   } else if (atN == 1) {
1530                     H.add(at);
1531                   }
1532                 }
1533                 return O != null && H.size() == 2;
1534               });
1535           // Now concatenate the stream with the properly labeled waters.
1536           settleStream = Stream.concat(settleStream, molecularAssembly.getWater().stream());
1537           // Map them into new Settle constraints and collect.
1538           List<SettleConstraint> settleConstraints = settleStream.map(
1539               (MSNode m) -> m.getAngleList().getFirst()).map(SettleConstraint::settleFactory).toList();
1540           constraints.addAll(settleConstraints);
1541 
1542         } else if (tok.equalsIgnoreCase("DIATOMIC")) {
1543           logger.severe(" Diatomic distance constraints not yet implemented properly.");
1544         } else if (tok.equalsIgnoreCase("TRIATOMIC")) {
1545           logger.severe(
1546               " Triatomic SETTLE constraints for non-water molecules not yet implemented properly.");
1547         }
1548       }
1549 
1550       // Second, accumulate bond/angle constraints.
1551       for (String tok : constraintToks) {
1552         if (tok.equalsIgnoreCase("BONDS") && bondPotentialEnergy != null) {
1553           Bond[] bonds = bondPotentialEnergy.getBondArray();
1554           numericBonds = new HashSet<>(Arrays.asList(bonds));
1555         } else if (tok.equalsIgnoreCase("ANGLES") && anglePotentialEnergy != null) {
1556           Angle[] angles = anglePotentialEnergy.getAngleArray();
1557           numericAngles = new HashSet<>(Arrays.asList(angles));
1558         }
1559       }
1560     }
1561 
1562     // Remove bonds that are already dealt with via angles.
1563     for (Angle angle : numericAngles) {
1564       angle.getBondList().forEach(numericBonds::remove);
1565     }
1566 
1567     // Remove already-constrained angles and bonds (e.g. SETTLE-constrained ones).
1568     List<Angle> ccmaAngles = numericAngles.stream().filter((Angle ang) -> !ang.isConstrained()).collect(Collectors.toList());
1569     List<Bond> ccmaBonds = numericBonds.stream().filter((Bond bond) -> !bond.isConstrained()).collect(Collectors.toList());
1570 
1571     CcmaConstraint ccmaConstraint = CcmaConstraint.ccmaFactory(ccmaBonds, ccmaAngles, atoms,
1572         getMass(), CcmaConstraint.DEFAULT_CCMA_NONZERO_CUTOFF);
1573     constraints.add(ccmaConstraint);
1574     logger.info(format(" Added %d constraints.", constraints.size()));
1575 
1576     return constraints;
1577   }
1578 
1579   /**
1580    * Static factory method to create a ForceFieldEnergy, possibly via FFX or OpenMM implementations.
1581    *
1582    * @param assembly To create FFE over
1583    * @return a {@link ffx.potential.ForceFieldEnergy} object.
1584    */
1585   public static ForceFieldEnergy energyFactory(MolecularAssembly assembly) {
1586     return energyFactory(assembly, ParallelTeam.getDefaultThreadCount());
1587   }
1588 
1589   /**
1590    * Static factory method to create a ForceFieldEnergy, possibly via FFX or OpenMM implementations.
1591    *
1592    * @param assembly   To create FFE over
1593    * @param numThreads Number of threads to use for FFX energy
1594    * @return A ForceFieldEnergy on some Platform
1595    */
1596   @SuppressWarnings("fallthrough")
1597   public static ForceFieldEnergy energyFactory(MolecularAssembly assembly, int numThreads) {
1598     ForceField forceField = assembly.getForceField();
1599     String platformString = toEnumForm(forceField.getString("PLATFORM", "FFX"));
1600     Platform platform;
1601     try {
1602       platform = Platform.valueOf(platformString);
1603     } catch (IllegalArgumentException e) {
1604       logger.warning(format(" String %s did not match a known energy implementation", platformString));
1605       platform = Platform.FFX;
1606     }
1607 
1608     // Check if the dual-topology platform flag is set.
1609     String dtPlatformString = toEnumForm(forceField.getString("PLATFORM_DT", "FFX"));
1610     try {
1611       Platform dtPlatform = Platform.valueOf(dtPlatformString);
1612       // If the dtPlatform uses OpenMM, then single topology energy will use FFX.
1613       if (dtPlatform != Platform.FFX) {
1614         // If the dtPlatform platform uses OpenMM, then use FFX for single topology.
1615         platform = Platform.FFX;
1616       }
1617     } catch (IllegalArgumentException e) {
1618       // If the dtPlatform is not recognized, ignore it.
1619     }
1620 
1621     switch (platform) {
1622       case OMM, OMM_REF, OMM_CUDA, OMM_OPENCL:
1623         try {
1624           return new OpenMMEnergy(assembly, platform, numThreads);
1625         } catch (Exception ex) {
1626           logger.warning(format(" Exception creating OpenMMEnergy: %s", ex));
1627           ForceFieldEnergy ffxEnergy = assembly.getPotentialEnergy();
1628           if (ffxEnergy == null) {
1629             ffxEnergy = new ForceFieldEnergy(assembly, numThreads);
1630             assembly.setPotential(ffxEnergy);
1631           }
1632           return ffxEnergy;
1633         }
1634       case OMM_CPU:
1635         logger.warning(format(" Platform %s not supported; defaulting to FFX", platform));
1636       default:
1637         ForceFieldEnergy ffxEnergy = assembly.getPotentialEnergy();
1638         if (ffxEnergy == null) {
1639           ffxEnergy = new ForceFieldEnergy(assembly, numThreads);
1640           assembly.setPotential(ffxEnergy);
1641         }
1642         return ffxEnergy;
1643     }
1644   }
1645 
1646   /**
1647    * Get all atoms that make up this ForceFieldEnergy.
1648    *
1649    * @return An array of Atoms.
1650    */
1651   public Atom[] getAtomArray() {
1652     return atoms;
1653   }
1654 
1655   /**
1656    * Applies constraints to positions
1657    *
1658    * @param xPrior Prior coodinates.
1659    * @param xNew   New coordinates.
1660    */
1661   public void applyAllConstraintPositions(double[] xPrior, double[] xNew) {
1662     applyAllConstraintPositions(xPrior, xNew, DEFAULT_CONSTRAINT_TOLERANCE);
1663   }
1664 
1665   /**
1666    * Applies constraints to positions
1667    *
1668    * @param xPrior Prior coodinates.
1669    * @param xNew   New coordinates.
1670    * @param tol    Constraint Tolerance.
1671    */
1672   public void applyAllConstraintPositions(double[] xPrior, double[] xNew, double tol) {
1673     if (xPrior == null) {
1674       xPrior = Arrays.copyOf(xNew, xNew.length);
1675     }
1676     for (Constraint constraint : constraints) {
1677       constraint.applyConstraintToStep(xPrior, xNew, getMass(), tol);
1678     }
1679   }
1680 
1681   /**
1682    * Overwrites current esvSystem if present. Multiple ExtendedSystems is possible but unnecessary;
1683    * add all ESVs to one system (per FFE, at least).
1684    *
1685    * @param system a {@link ffx.potential.extended.ExtendedSystem} object.
1686    */
1687   public void attachExtendedSystem(ExtendedSystem system) {
1688     if (system == null) {
1689       throw new IllegalArgumentException();
1690     }
1691     esvTerm = true;
1692     esvTermOrig = esvTerm;
1693     esvSystem = system;
1694     if (vanderWaalsTerm) {
1695       if (vanderWaals == null) {
1696         logger.warning("Null VdW during ESV setup.");
1697       } else {
1698         vanderWaals.attachExtendedSystem(system);
1699       }
1700     }
1701     if (multipoleTerm) {
1702       if (particleMeshEwald == null) {
1703         logger.warning("Null PME during ESV setup.");
1704       } else {
1705         particleMeshEwald.attachExtendedSystem(system);
1706       }
1707     }
1708   }
1709 
1710   /**
1711    * {@inheritDoc}
1712    *
1713    * <p>Returns true if lambda term is not enabled for this ForceFieldEnergy.
1714    */
1715   @Override
1716   public boolean dEdLZeroAtEnds() {
1717     // This may actually be true even with softcored atoms.
1718     // For now, serves the purpose of reporting true when nothing is softcored.
1719     return !lambdaTerm;
1720   }
1721 
1722   /**
1723    * Frees up assets associated with this ForceFieldEnergy, such as worker Threads.
1724    *
1725    * @return If successful in freeing up assets.
1726    */
1727   public boolean destroy() {
1728     if (destroyed) {
1729       // This regularly occurs with Repex OST, as multiple OrthogonalSpaceTempering objects wrap a
1730       // single FFE.
1731       logger.fine(format(" This ForceFieldEnergy is already destroyed: %s", this));
1732       return true;
1733     } else {
1734       try {
1735         if (parallelTeam != null) {
1736           parallelTeam.shutdown();
1737         }
1738         if (vanderWaals != null) {
1739           vanderWaals.destroy();
1740         }
1741         if (particleMeshEwald != null) {
1742           particleMeshEwald.destroy();
1743         }
1744         molecularAssembly.finishDestruction();
1745         destroyed = true;
1746         return true;
1747       } catch (Exception ex) {
1748         logger.warning(format(" Exception in shutting down a ForceFieldEnergy: %s", ex));
1749         logger.info(Utilities.stackTraceToString(ex));
1750         return false;
1751       }
1752     }
1753   }
1754 
1755   /**
1756    * energy.
1757    *
1758    * @return a double.
1759    */
1760   public double energy() {
1761     return energy(false, false);
1762   }
1763 
1764   /**
1765    * Compute the potential energy of the system.
1766    *
1767    * @param gradient If true, compute the Cartesian coordinate gradient.
1768    * @param print    If true, print the energy terms.
1769    * @return the energy in kcal/mol.
1770    */
1771   public double energy(boolean gradient, boolean print) {
1772     try {
1773       totalTime = System.nanoTime();
1774       nnTime = 0;
1775       restrainGroupTime = 0;
1776       vanDerWaalsTime = 0;
1777       electrostaticTime = 0;
1778       ncsTime = 0;
1779 
1780       // Zero out the potential energy of each bonded term.
1781       double forceFieldBondedEnergy = 0.0;
1782       nnEnergy = 0.0;
1783 
1784       // Zero out potential energy of restraint terms
1785       double restrainEnergy = 0.0;
1786       restrainGroupEnergy = 0.0;
1787       ncsEnergy = 0.0;
1788 
1789       // Zero out the potential energy of each non-bonded term.
1790       double totalNonBondedEnergy = 0.0;
1791       double totalMultipoleEnergy = 0.0;
1792       vanDerWaalsEnergy = 0.0;
1793       permanentMultipoleEnergy = 0.0;
1794       polarizationEnergy = 0.0;
1795 
1796       // Zero out the solvation energy.
1797       solvationEnergy = 0.0;
1798       esvBias = 0.0;
1799 
1800       // Zero out the total potential energy.
1801       totalEnergy = 0.0;
1802 
1803       // BondedEnergyRegion to compute bonded terms in parallel.
1804       try {
1805         forceFieldBondedEnergyRegion.setGradient(gradient);
1806         forceFieldBondedEnergyRegion.setLambdaBondedTerms(lambdaBondedTerms);
1807         forceFieldBondedEnergyRegion.setLambdaAllBondedTerms(lambdaAllBondedTerms);
1808         forceFieldBondedEnergyRegion.setInitAtomGradients(true);
1809         parallelTeam.execute(forceFieldBondedEnergyRegion);
1810         forceFieldBondedEnergy = forceFieldBondedEnergyRegion.getEnergy();
1811       } catch (RuntimeException ex) {
1812         logger.warning("Runtime exception during bonded term calculation.");
1813         throw ex;
1814       } catch (Exception ex) {
1815         logger.info(Utilities.stackTraceToString(ex));
1816         logger.severe(ex.toString());
1817       }
1818 
1819       // Restrain Energy Region to compute restrain terms in parallel.
1820       try {
1821         if ((restrainMode == RestrainMode.ENERGY) ||
1822             (restrainMode == RestrainMode.ALCHEMICAL && lambdaBondedTerms)) {
1823           boolean checkAlchemicalAtoms = (restrainMode == RestrainMode.ENERGY);
1824           restrainEnergyRegion.setGradient(gradient);
1825           restrainEnergyRegion.setLambdaBondedTerms(lambdaBondedTerms);
1826           restrainEnergyRegion.setLambdaAllBondedTerms(lambdaAllBondedTerms);
1827           restrainEnergyRegion.setCheckAlchemicalAtoms(checkAlchemicalAtoms);
1828           restrainEnergyRegion.setInitAtomGradients(false);
1829           parallelTeam.execute(restrainEnergyRegion);
1830           restrainEnergy = restrainEnergyRegion.getEnergy();
1831         }
1832       } catch (RuntimeException ex) {
1833         logger.warning("Runtime exception during restrain term calculation.");
1834         throw ex;
1835       } catch (Exception ex) {
1836         logger.info(Utilities.stackTraceToString(ex));
1837         logger.severe(ex.toString());
1838       }
1839 
1840       if (!lambdaBondedTerms) {
1841         // Compute restraint terms.
1842         if (ncsTerm) {
1843           ncsTime = -System.nanoTime();
1844           ncsEnergy = ncsRestraint.residual(gradient, print);
1845           ncsTime += System.nanoTime();
1846         }
1847 
1848         /*
1849         if (restrainPositionTerm) {
1850           restrainPositionTime = -System.nanoTime();
1851           for (RestrainPosition restrainPosition : restrainPositions) {
1852             restrainPositionEnergy += restrainPosition.residual(gradient);
1853           }
1854           restrainPositionTime += System.nanoTime();
1855         }
1856         */
1857 
1858         if (restrainGroupTerm) {
1859           restrainGroupTime = -System.nanoTime();
1860           restrainGroupEnergy = restrainGroups.energy(gradient);
1861           restrainGroupTime += System.nanoTime();
1862         }
1863 
1864         if (comTerm) {
1865           comRestraintTime = -System.nanoTime();
1866           comRestraintEnergy = comRestraint.residual(gradient, print);
1867           comRestraintTime += System.nanoTime();
1868         }
1869 
1870 
1871         // Compute the neural network term.
1872         if (nnTerm) {
1873           nnTime = -System.nanoTime();
1874           nnEnergy = aniEnergy.energy(gradient, print);
1875           nnTime += System.nanoTime();
1876         }
1877 
1878         // Compute non-bonded terms.
1879         if (vanderWaalsTerm) {
1880           vanDerWaalsTime = -System.nanoTime();
1881           vanDerWaalsEnergy = vanderWaals.energy(gradient, print);
1882           nVanDerWaalInteractions = this.vanderWaals.getInteractions();
1883           vanDerWaalsTime += System.nanoTime();
1884         }
1885 
1886         if (multipoleTerm) {
1887           electrostaticTime = -System.nanoTime();
1888           totalMultipoleEnergy = particleMeshEwald.energy(gradient, print);
1889           permanentMultipoleEnergy = particleMeshEwald.getPermanentEnergy();
1890           polarizationEnergy = particleMeshEwald.getPolarizationEnergy();
1891           nPermanentInteractions = particleMeshEwald.getInteractions();
1892           solvationEnergy = particleMeshEwald.getSolvationEnergy();
1893           nGKInteractions = particleMeshEwald.getGKInteractions();
1894           electrostaticTime += System.nanoTime();
1895         }
1896       }
1897 
1898       totalTime = System.nanoTime() - totalTime;
1899 
1900       double totalBondedEnergy = forceFieldBondedEnergy + restrainEnergy
1901           + nnEnergy + ncsEnergy + comRestraintEnergy + restrainGroupEnergy;
1902 
1903       totalNonBondedEnergy = vanDerWaalsEnergy + totalMultipoleEnergy;
1904       totalEnergy = totalBondedEnergy + totalNonBondedEnergy + solvationEnergy;
1905       if (esvTerm) {
1906         esvBias = esvSystem.getBiasEnergy();
1907         totalEnergy += esvBias;
1908       }
1909 
1910       if (print) {
1911         StringBuilder sb = new StringBuilder();
1912         if (gradient) {
1913           sb.append("\n Computed Potential Energy and Atomic Coordinate Gradients\n");
1914         } else {
1915           sb.append("\n Computed Potential Energy\n");
1916         }
1917         sb.append(this);
1918         logger.info(sb.toString());
1919       }
1920       return totalEnergy;
1921     } catch (EnergyException ex) {
1922       if (printOnFailure) {
1923         printFailure();
1924       }
1925       if (ex.doCauseSevere()) {
1926         logger.info(Utilities.stackTraceToString(ex));
1927         logger.log(Level.SEVERE, " Error in calculating energies or gradients", ex);
1928       } else {
1929         logger.log(Level.INFO, format(" Exception in energy calculation:\n %s", ex));
1930       }
1931       throw ex;
1932     }
1933   }
1934 
1935   /**
1936    * {@inheritDoc}
1937    */
1938   @Override
1939   public double energy(double[] x) {
1940     return energy(x, false);
1941   }
1942 
1943   /**
1944    * {@inheritDoc}
1945    */
1946   @Override
1947   public double energy(double[] x, boolean verbose) {
1948     assert Arrays.stream(x).allMatch(Double::isFinite);
1949 
1950     // Unscale the coordinates.
1951     unscaleCoordinates(x);
1952 
1953     // Set coordinates.
1954     setCoordinates(x);
1955 
1956     double e = this.energy(false, verbose);
1957 
1958     // Rescale the coordinates.
1959     scaleCoordinates(x);
1960 
1961     return e;
1962   }
1963 
1964   /**
1965    * {@inheritDoc}
1966    */
1967   @Override
1968   public double energyAndGradient(double[] x, double[] g) {
1969     return energyAndGradient(x, g, false);
1970   }
1971 
1972   /**
1973    * {@inheritDoc}
1974    */
1975   @Override
1976   public double energyAndGradient(double[] x, double[] g, boolean verbose) {
1977     // Un-scale the coordinates.
1978     unscaleCoordinates(x);
1979     // Set coordinates.
1980     setCoordinates(x);
1981     double e = energy(true, verbose);
1982 
1983     // Try block already exists inside energy(boolean, boolean), so only
1984     // need to try-catch fillGradient.
1985     try {
1986       fillGradient(g);
1987 
1988       // Scale the coordinates and gradients.
1989       scaleCoordinatesAndGradient(x, g);
1990 
1991       if (maxDebugGradient < Double.MAX_VALUE) {
1992         boolean extremeGrad = Arrays.stream(g)
1993             .anyMatch((double gi) -> (gi > maxDebugGradient || gi < -maxDebugGradient));
1994         if (extremeGrad) {
1995           File origFile = molecularAssembly.getFile();
1996           String timeString = LocalDateTime.now()
1997               .format(DateTimeFormatter.ofPattern("yyyy_MM_dd-HH_mm_ss"));
1998 
1999           String filename = format("%s-LARGEGRAD-%s.pdb",
2000               removeExtension(molecularAssembly.getFile().getName()), timeString);
2001           PotentialsFunctions ef = new PotentialsUtils();
2002           filename = ef.versionFile(filename);
2003 
2004           logger.warning(format(" Excessively large gradient detected; printing snapshot to file %s",
2005               filename));
2006           ef.saveAsPDB(molecularAssembly, new File(filename));
2007           molecularAssembly.setFile(origFile);
2008         }
2009       }
2010       return e;
2011     } catch (EnergyException ex) {
2012       if (printOnFailure) {
2013         printFailure();
2014       }
2015       if (ex.doCauseSevere()) {
2016         logger.info(Utilities.stackTraceToString(ex));
2017         logger.log(Level.SEVERE, " Error in calculating energies or gradients", ex);
2018       } else {
2019         logger.log(Level.INFO, format(" Exception in energy calculation:\n %s", ex));
2020       }
2021       throw ex;
2022     }
2023   }
2024 
2025   /**
2026    * Save coordinates when an EnergyException is caught.
2027    */
2028   private void printFailure() {
2029     String timeString = LocalDateTime.now()
2030         .format(DateTimeFormatter.ofPattern("yyyy_MM_dd-HH_mm_ss"));
2031     File file = molecularAssembly.getFile();
2032     String ext = "pdb";
2033     if (isXYZ(file)) {
2034       ext = "xyz";
2035     }
2036     String filename = format("%s-ERROR-%s.%s", removeExtension(file.getName()), timeString, ext);
2037     PotentialsFunctions ef = new PotentialsUtils();
2038     filename = ef.versionFile(filename);
2039     logger.info(format(" Writing on-error snapshot to file %s", filename));
2040     ef.save(molecularAssembly, new File(filename));
2041   }
2042 
2043   /**
2044    * {@inheritDoc}
2045    *
2046    * <p>Returns an array of acceleration values for active atoms.
2047    */
2048   @Override
2049   public double[] getAcceleration(double[] acceleration) {
2050     int n = getNumberOfVariables();
2051     if (acceleration == null || acceleration.length < n) {
2052       acceleration = new double[n];
2053     }
2054     int index = 0;
2055     double[] a = new double[3];
2056     for (int i = 0; i < nAtoms; i++) {
2057       if (atoms[i].isActive()) {
2058         atoms[i].getAcceleration(a);
2059         acceleration[index++] = a[0];
2060         acceleration[index++] = a[1];
2061         acceleration[index++] = a[2];
2062       }
2063     }
2064     return acceleration;
2065   }
2066 
2067   /**
2068    * Getter for the field <code>angles</code> with only <code>AngleMode</code> angles.
2069    *
2070    * @param angleMode Only angles of this mode will be returned.
2071    * @return an array of {@link ffx.potential.bonded.Angle} objects.
2072    */
2073   public Angle[] getAngles(AngleMode angleMode) {
2074     if (anglePotentialEnergy == null) {
2075       return null;
2076     }
2077     Angle[] angles = anglePotentialEnergy.getAngleArray();
2078     int nAngles = angles.length;
2079     List<Angle> angleList = new ArrayList<>();
2080     // Sort all normal angles from in-plane angles
2081     for (int i = 0; i < nAngles; i++) {
2082       if (angles[i].getAngleMode() == angleMode) {
2083         angleList.add(angles[i]);
2084       }
2085     }
2086     nAngles = angleList.size();
2087     if (nAngles < 1) {
2088       return null;
2089     }
2090     return angleList.toArray(new Angle[0]);
2091   }
2092 
2093   /**
2094    * Returns a copy of the list of constraints this ForceFieldEnergy has.
2095    *
2096    * @return Copied list of constraints.
2097    */
2098   @Override
2099   public List<Constraint> getConstraints() {
2100     return constraints.isEmpty() ? Collections.emptyList() : new ArrayList<>(constraints);
2101   }
2102 
2103   /**
2104    * {@inheritDoc}
2105    */
2106   @Override
2107   public double[] getCoordinates(double[] x) {
2108     int n = getNumberOfVariables();
2109     if (x == null || x.length < n) {
2110       x = new double[n];
2111     }
2112     int index = 0;
2113     for (int i = 0; i < nAtoms; i++) {
2114       Atom a = atoms[i];
2115       if (a.isActive()) {
2116         x[index++] = a.getX();
2117         x[index++] = a.getY();
2118         x[index++] = a.getZ();
2119       }
2120     }
2121     return x;
2122   }
2123 
2124   /**
2125    * {@inheritDoc}
2126    *
2127    * <p>Getter for the field <code>crystal</code>.
2128    */
2129   @Override
2130   public Crystal getCrystal() {
2131     return crystal;
2132   }
2133 
2134   /**
2135    * {@inheritDoc}
2136    *
2137    * <p>Set the boundary conditions for this calculation.
2138    */
2139   @Override
2140   public void setCrystal(Crystal crystal) {
2141     setCrystal(crystal, false);
2142   }
2143 
2144   /**
2145    * Getter for the field <code>cutoffPlusBuffer</code>.
2146    *
2147    * @return a double.
2148    */
2149   public double getCutoffPlusBuffer() {
2150     return cutoffPlusBuffer;
2151   }
2152 
2153   /**
2154    * {@inheritDoc}
2155    */
2156   @Override
2157   public STATE getEnergyTermState() {
2158     return state;
2159   }
2160 
2161   /**
2162    * {@inheritDoc}
2163    *
2164    * <p>This method is for the RESPA integrator only.
2165    */
2166   @Override
2167   public void setEnergyTermState(STATE state) {
2168     this.state = state;
2169     if (forceFieldBondedEnergyRegion != null) {
2170       forceFieldBondedEnergyRegion.setState(state);
2171     }
2172     if (restrainEnergyRegion != null) {
2173       restrainEnergyRegion.setState(state);
2174     }
2175     switch (state) {
2176       case FAST:
2177         nnTerm = nnTermOrig;
2178         restrainGroupTerm = restrainGroupTermOrig;
2179         ncsTerm = ncsTermOrig;
2180         comTerm = comTermOrig;
2181         vanderWaalsTerm = false;
2182         multipoleTerm = false;
2183         polarizationTerm = false;
2184         generalizedKirkwoodTerm = false;
2185         esvTerm = false;
2186         break;
2187       case SLOW:
2188         vanderWaalsTerm = vanderWaalsTermOrig;
2189         multipoleTerm = multipoleTermOrig;
2190         polarizationTerm = polarizationTermOrig;
2191         generalizedKirkwoodTerm = generalizedKirkwoodTermOrig;
2192         esvTerm = esvTermOrig;
2193         nnTerm = false;
2194         restrainGroupTerm = false;
2195         ncsTerm = false;
2196         comTerm = false;
2197         break;
2198       default:
2199         nnTerm = nnTermOrig;
2200         restrainGroupTerm = restrainGroupTermOrig;
2201         ncsTerm = ncsTermOrig;
2202         comTermOrig = comTerm;
2203         vanderWaalsTerm = vanderWaalsTermOrig;
2204         multipoleTerm = multipoleTermOrig;
2205         polarizationTerm = polarizationTermOrig;
2206         generalizedKirkwoodTerm = generalizedKirkwoodTermOrig;
2207     }
2208   }
2209 
2210   /**
2211    * getEsvBiasEnergy.
2212    *
2213    * @return a double.
2214    */
2215   public double getEsvBiasEnergy() {
2216     return esvBias;
2217   }
2218 
2219   /**
2220    * getExtendedSystem.
2221    *
2222    * @return a {@link ffx.potential.extended.ExtendedSystem} object.
2223    */
2224   public ExtendedSystem getExtendedSystem() {
2225     return esvSystem;
2226   }
2227 
2228   /**
2229    * getGK.
2230    *
2231    * @return a {@link ffx.potential.nonbonded.GeneralizedKirkwood} object.
2232    */
2233   public GeneralizedKirkwood getGK() {
2234     if (particleMeshEwald != null) {
2235       return particleMeshEwald.getGK();
2236     } else {
2237       return null;
2238     }
2239   }
2240 
2241   /**
2242    * Returns the gradient array for this ForceFieldEnergy.
2243    *
2244    * @param g an array of double.
2245    * @return the gradient array.
2246    */
2247   public double[] getGradient(double[] g) {
2248     return fillGradient(g);
2249   }
2250 
2251   /**
2252    * {@inheritDoc}
2253    */
2254   @Override
2255   public double getLambda() {
2256     return lambda;
2257   }
2258 
2259   /**
2260    * {@inheritDoc}
2261    */
2262   @Override
2263   public void setLambda(double lambda) {
2264     if (lambdaTerm) {
2265       if (lambda <= 1.0 && lambda >= 0.0) {
2266         this.lambda = lambda;
2267         if (vanderWaalsTerm) {
2268           vanderWaals.setLambda(lambda);
2269         }
2270         if (multipoleTerm) {
2271           particleMeshEwald.setLambda(lambda);
2272         }
2273 
2274         /*
2275         if (restrainPositionTerm && nRestrainPositions > 0) {
2276           for (RestrainPosition restrainPosition : restrainPositions) {
2277             restrainPosition.setLambda(lambda);
2278           }
2279         } */
2280 
2281         if (restrainDistanceTerm && restrainDistancePotentialEnergy != null) {
2282           for (RestrainDistance restrainDistance : restrainDistancePotentialEnergy.getRestrainDistances()) {
2283             restrainDistance.setLambda(lambda);
2284           }
2285         }
2286         if (restrainTorsionTerm && restrainTorsionPotentialEnergy != null) {
2287           for (Torsion restrainTorsion : restrainTorsionPotentialEnergy.getRestrainTorsions()) {
2288             restrainTorsion.setLambda(lambda);
2289           }
2290         }
2291         if (ncsTerm && ncsRestraint != null) {
2292           ncsRestraint.setLambda(lambda);
2293         }
2294         if (comTerm && comRestraint != null) {
2295           comRestraint.setLambda(lambda);
2296         }
2297         if (lambdaTorsions) {
2298           if (torsionPotentialEnergy != null) {
2299             torsionPotentialEnergy.setLambda(lambda);
2300           }
2301           if (piOrbitalTorsionPotentialEnergy != null) {
2302             piOrbitalTorsionPotentialEnergy.setLambda(lambda);
2303           }
2304           if (torsionTorsionPotentialEnergy != null) {
2305             torsionTorsionPotentialEnergy.setLambda(lambda);
2306           }
2307         }
2308       } else {
2309         String message = format("Lambda value %8.3f is not in the range [0..1].", lambda);
2310         logger.warning(message);
2311       }
2312     } else {
2313       logger.fine(" Attempting to set a lambda value on a ForceFieldEnergy with lambdaterm false.");
2314     }
2315   }
2316 
2317   /**
2318    * {@inheritDoc}
2319    */
2320   @Override
2321   public double[] getMass() {
2322     int n = getNumberOfVariables();
2323     double[] mass = new double[n];
2324     int index = 0;
2325     for (int i = 0; i < nAtoms; i++) {
2326       Atom a = atoms[i];
2327       if (a.isActive()) {
2328         double m = a.getMass();
2329         mass[index++] = m;
2330         mass[index++] = m;
2331         mass[index++] = m;
2332       }
2333     }
2334     return mass;
2335   }
2336 
2337   /**
2338    * {@inheritDoc}
2339    */
2340   @Override
2341   public int getNumberOfVariables() {
2342     int nActive = 0;
2343     for (int i = 0; i < nAtoms; i++) {
2344       if (atoms[i].isActive()) {
2345         nActive++;
2346       }
2347     }
2348     return nActive * 3;
2349   }
2350 
2351   /**
2352    * Create a PDB REMARK 3 string containing the potential energy terms.
2353    *
2354    * @return a String containing the PDB REMARK 3 formatted energy terms.
2355    */
2356   public String getPDBHeaderString() {
2357     energy(false, false);
2358     StringBuilder sb = new StringBuilder();
2359     sb.append("REMARK   3  CALCULATED POTENTIAL ENERGY\n");
2360     sb.append(forceFieldBondedEnergyRegion.toPDBString());
2361     sb.append(restrainEnergyRegion.toPDBString());
2362     if (nnTerm) {
2363       sb.append(
2364           format("REMARK   3   %s %g (%d)\n", "NEUTRAL NETWORK            : ", nnEnergy, nAtoms));
2365     }
2366     if (ncsTerm) {
2367       sb.append(
2368           format("REMARK   3   %s %g (%d)\n", "NCS RESTRAINT              : ", ncsEnergy, nAtoms));
2369     }
2370     if (comTerm) {
2371       sb.append(
2372           format("REMARK   3   %s %g (%d)\n", "COM RESTRAINT              : ", comRestraintEnergy,
2373               nAtoms));
2374     }
2375     if (vanderWaalsTerm) {
2376       sb.append(
2377           format("REMARK   3   %s %g (%d)\n", "VAN DER WAALS              : ", vanDerWaalsEnergy,
2378               nVanDerWaalInteractions));
2379     }
2380     if (multipoleTerm) {
2381       sb.append(format("REMARK   3   %s %g (%d)\n", "ATOMIC MULTIPOLES          : ",
2382           permanentMultipoleEnergy, nPermanentInteractions));
2383     }
2384     if (polarizationTerm) {
2385       sb.append(format("REMARK   3   %s %g (%d)\n", "POLARIZATION               : ",
2386           polarizationEnergy, nPermanentInteractions));
2387     }
2388     sb.append(format("REMARK   3   %s %g\n", "TOTAL POTENTIAL (KCAL/MOL) : ", totalEnergy));
2389     int nsymm = crystal.getUnitCell().spaceGroup.getNumberOfSymOps();
2390     if (nsymm > 1) {
2391       sb.append(format("REMARK   3   %s %g\n", "UNIT CELL POTENTIAL        : ", totalEnergy * nsymm));
2392     }
2393     if (crystal.getUnitCell() != crystal) {
2394       nsymm = crystal.spaceGroup.getNumberOfSymOps();
2395       if (nsymm > 1) {
2396         sb.append(format("REMARK   3   %s %g\n", "REPLICATES CELL POTENTIAL  : ", totalEnergy * nsymm));
2397       }
2398     }
2399     sb.append("REMARK   3\n");
2400 
2401     return sb.toString();
2402   }
2403 
2404   /**
2405    * Getter for the field <code>parallelTeam</code>.
2406    *
2407    * @return a {@link edu.rit.pj.ParallelTeam} object.
2408    */
2409   public ParallelTeam getParallelTeam() {
2410     return parallelTeam;
2411   }
2412 
2413   /**
2414    * getPermanentInteractions.
2415    *
2416    * @return a int.
2417    */
2418   public int getPermanentInteractions() {
2419     return nPermanentInteractions;
2420   }
2421 
2422   /**
2423    * Getter for the field <code>permanentMultipoleEnergy</code>.
2424    *
2425    * @return a double.
2426    */
2427   public double getPermanentMultipoleEnergy() {
2428     return permanentMultipoleEnergy;
2429   }
2430 
2431   /**
2432    * Gets the Platform associated with this force field energy. For the reference platform, always
2433    * returns FFX.
2434    *
2435    * @return A Platform.
2436    */
2437   public Platform getPlatform() {
2438     return platform;
2439   }
2440 
2441   /**
2442    * getPmeNode.
2443    *
2444    * @return a {@link ParticleMeshEwald} object.
2445    */
2446   public ParticleMeshEwald getPmeNode() {
2447     return particleMeshEwald;
2448   }
2449 
2450   /**
2451    * Getter for the field <code>polarizationEnergy</code>.
2452    *
2453    * @return a double.
2454    */
2455   public double getPolarizationEnergy() {
2456     return polarizationEnergy;
2457   }
2458 
2459   /**
2460    * {@inheritDoc}
2461    *
2462    * <p>Returns an array of previous acceleration values for active atoms.
2463    */
2464   @Override
2465   public double[] getPreviousAcceleration(double[] previousAcceleration) {
2466     int n = getNumberOfVariables();
2467     if (previousAcceleration == null || previousAcceleration.length < n) {
2468       previousAcceleration = new double[n];
2469     }
2470     int index = 0;
2471     double[] a = new double[3];
2472     for (int i = 0; i < nAtoms; i++) {
2473       if (atoms[i].isActive()) {
2474         atoms[i].getPreviousAcceleration(a);
2475         previousAcceleration[index++] = a[0];
2476         previousAcceleration[index++] = a[1];
2477         previousAcceleration[index++] = a[2];
2478       }
2479     }
2480     return previousAcceleration;
2481   }
2482 
2483   /**
2484    * {@inheritDoc}
2485    */
2486   @Override
2487   public double[] getScaling() {
2488     return optimizationScaling;
2489   }
2490 
2491   /**
2492    * {@inheritDoc}
2493    */
2494   @Override
2495   public void setScaling(double[] scaling) {
2496     optimizationScaling = scaling;
2497   }
2498 
2499   /**
2500    * Getter for the field <code>solvationEnergy</code>.
2501    *
2502    * @return a double.
2503    */
2504   public double getSolvationEnergy() {
2505     return solvationEnergy;
2506   }
2507 
2508   /**
2509    * getSolvationInteractions.
2510    *
2511    * @return a int.
2512    */
2513   public int getSolvationInteractions() {
2514     return nGKInteractions;
2515   }
2516 
2517   /**
2518    * {@inheritDoc}
2519    */
2520   @Override
2521   public double getTotalEnergy() {
2522     return totalEnergy;
2523   }
2524 
2525   /**
2526    * Getter for the field <code>vanDerWaalsEnergy</code>.
2527    *
2528    * @return a double.
2529    */
2530   public double getVanDerWaalsEnergy() {
2531     return vanDerWaalsEnergy;
2532   }
2533 
2534   /**
2535    * getVanDerWaalsInteractions.
2536    *
2537    * @return a int.
2538    */
2539   public int getVanDerWaalsInteractions() {
2540     return nVanDerWaalInteractions;
2541   }
2542 
2543   /**
2544    * {@inheritDoc}
2545    *
2546    * <p>Return a reference to each variables type.
2547    */
2548   @Override
2549   public VARIABLE_TYPE[] getVariableTypes() {
2550     int n = getNumberOfVariables();
2551     VARIABLE_TYPE[] type = new VARIABLE_TYPE[n];
2552     int i = 0;
2553     for (int j = 0; j < nAtoms; j++) {
2554       if (atoms[j].isActive()) {
2555         type[i++] = VARIABLE_TYPE.X;
2556         type[i++] = VARIABLE_TYPE.Y;
2557         type[i++] = VARIABLE_TYPE.Z;
2558       }
2559     }
2560     return type;
2561   }
2562 
2563   /**
2564    * getVdwNode.
2565    *
2566    * @return a {@link ffx.potential.nonbonded.VanDerWaals} object.
2567    */
2568   public VanDerWaals getVdwNode() {
2569     return vanderWaals;
2570   }
2571 
2572   /**
2573    * {@inheritDoc}
2574    *
2575    * <p>Returns an array of velocity values for active atoms.
2576    */
2577   @Override
2578   public double[] getVelocity(double[] velocity) {
2579     int n = getNumberOfVariables();
2580     if (velocity == null || velocity.length < n) {
2581       velocity = new double[n];
2582     }
2583     int index = 0;
2584     double[] v = new double[3];
2585     for (int i = 0; i < nAtoms; i++) {
2586       Atom a = atoms[i];
2587       if (a.isActive()) {
2588         a.getVelocity(v);
2589         velocity[index++] = v[0];
2590         velocity[index++] = v[1];
2591         velocity[index++] = v[2];
2592       }
2593     }
2594     return velocity;
2595   }
2596 
2597   /**
2598    * {@inheritDoc}
2599    */
2600   @Override
2601   public double getd2EdL2() {
2602     double d2EdLambda2 = 0.0;
2603     if (!lambdaBondedTerms) {
2604       if (vanderWaalsTerm) {
2605         d2EdLambda2 = vanderWaals.getd2EdL2();
2606       }
2607       if (multipoleTerm) {
2608         d2EdLambda2 += particleMeshEwald.getd2EdL2();
2609       }
2610       /*
2611       if (restrainPositionTerm && nRestrainPositions > 0) {
2612         for (RestrainPosition restrainPosition : restrainPositions) {
2613           d2EdLambda2 += restrainPosition.getd2EdL2();
2614         }
2615       }
2616       */
2617       if (restrainDistanceTerm && restrainDistancePotentialEnergy != null) {
2618         for (RestrainDistance restrainDistance : restrainDistancePotentialEnergy.getRestrainDistanceArray()) {
2619           d2EdLambda2 += restrainDistance.getd2EdL2();
2620         }
2621       }
2622       if (ncsTerm && ncsRestraint != null) {
2623         d2EdLambda2 += ncsRestraint.getd2EdL2();
2624       }
2625       if (comTerm && comRestraint != null) {
2626         d2EdLambda2 += comRestraint.getd2EdL2();
2627       }
2628       if (lambdaTorsions) {
2629         if (torsionPotentialEnergy != null) {
2630           d2EdLambda2 += torsionPotentialEnergy.getd2EdL2();
2631         }
2632         if (piOrbitalTorsionPotentialEnergy != null) {
2633           d2EdLambda2 += piOrbitalTorsionPotentialEnergy.getd2EdL2();
2634         }
2635         if (torsionTorsionPotentialEnergy != null) {
2636           d2EdLambda2 += torsionTorsionPotentialEnergy.getd2EdL2();
2637         }
2638       }
2639     }
2640     return d2EdLambda2;
2641   }
2642 
2643   /**
2644    * {@inheritDoc}
2645    */
2646   @Override
2647   public double getdEdL() {
2648     double dEdLambda = 0.0;
2649     if (!lambdaBondedTerms) {
2650       if (vanderWaalsTerm) {
2651         dEdLambda = vanderWaals.getdEdL();
2652       }
2653       if (multipoleTerm) {
2654         dEdLambda += particleMeshEwald.getdEdL();
2655       }
2656       if (restrainDistanceTerm && restrainDistancePotentialEnergy != null) {
2657         for (RestrainDistance restrainDistance : restrainDistancePotentialEnergy.getRestrainDistanceArray()) {
2658           dEdLambda += restrainDistance.getdEdL();
2659         }
2660       }
2661       /*
2662       if (restrainPositionTerm && nRestrainPositions > 0) {
2663         for (RestrainPosition restrainPosition : restrainPositions) {
2664           dEdLambda += restrainPosition.getdEdL();
2665         }
2666       }
2667       */
2668       if (ncsTerm && ncsRestraint != null) {
2669         dEdLambda += ncsRestraint.getdEdL();
2670       }
2671       if (comTerm && comRestraint != null) {
2672         dEdLambda += comRestraint.getdEdL();
2673       }
2674       if (lambdaTorsions) {
2675         if (torsionPotentialEnergy != null) {
2676           dEdLambda += torsionPotentialEnergy.getdEdL();
2677         }
2678         if (piOrbitalTorsionPotentialEnergy != null) {
2679           dEdLambda += piOrbitalTorsionPotentialEnergy.getdEdL();
2680         }
2681         if (torsionTorsionPotentialEnergy != null) {
2682           dEdLambda += torsionTorsionPotentialEnergy.getdEdL();
2683         }
2684       }
2685     }
2686     return dEdLambda;
2687   }
2688 
2689   /**
2690    * {@inheritDoc}
2691    */
2692   @Override
2693   public void getdEdXdL(double[] gradients) {
2694     if (!lambdaBondedTerms) {
2695       if (vanderWaalsTerm) {
2696         vanderWaals.getdEdXdL(gradients);
2697       }
2698       if (multipoleTerm) {
2699         particleMeshEwald.getdEdXdL(gradients);
2700       }
2701       if (restrainDistanceTerm && restrainDistancePotentialEnergy != null) {
2702         for (RestrainDistance restrainDistance : restrainDistancePotentialEnergy.getRestrainDistanceArray()) {
2703           restrainDistance.getdEdXdL(gradients);
2704         }
2705       }
2706       /*
2707       if (restrainPositionTerm && nRestrainPositions > 0) {
2708         for (RestrainPosition restrainPosition : restrainPositions) {
2709           restrainPosition.getdEdXdL(gradients);
2710         }
2711       }
2712       */
2713       if (ncsTerm && ncsRestraint != null) {
2714         ncsRestraint.getdEdXdL(gradients);
2715       }
2716       if (comTerm && comRestraint != null) {
2717         comRestraint.getdEdXdL(gradients);
2718       }
2719       if (lambdaTorsions) {
2720         double[] grad = new double[3];
2721         int index = 0;
2722         for (int i = 0; i < nAtoms; i++) {
2723           Atom a = atoms[i];
2724           if (a.isActive()) {
2725             a.getLambdaXYZGradient(grad);
2726             gradients[index++] += grad[0];
2727             gradients[index++] += grad[1];
2728             gradients[index++] += grad[2];
2729           }
2730         }
2731       }
2732     }
2733   }
2734 
2735   /**
2736    * {@inheritDoc}
2737    *
2738    * <p>The acceleration array should only contain acceleration data for active atoms.
2739    */
2740   @Override
2741   public void setAcceleration(double[] acceleration) {
2742     if (acceleration == null) {
2743       return;
2744     }
2745     int index = 0;
2746     double[] accel = new double[3];
2747     for (int i = 0; i < nAtoms; i++) {
2748       if (atoms[i].isActive()) {
2749         accel[0] = acceleration[index++];
2750         accel[1] = acceleration[index++];
2751         accel[2] = acceleration[index++];
2752         atoms[i].setAcceleration(accel);
2753       }
2754     }
2755   }
2756 
2757   /**
2758    * The coordinate array should only contain active atoms.
2759    *
2760    * @param coords the coordinates to set.
2761    */
2762   public void setCoordinates(@Nullable double[] coords) {
2763     if (coords == null) {
2764       return;
2765     }
2766     int index = 0;
2767     for (Atom a : atoms) {
2768       if (a.isActive()) {
2769         double x = coords[index++];
2770         double y = coords[index++];
2771         double z = coords[index++];
2772         a.moveTo(x, y, z);
2773       }
2774     }
2775   }
2776 
2777   /**
2778    * Set the boundary conditions for this calculation.
2779    *
2780    * @param crystal             Crystal to set.
2781    * @param checkReplicatesCell Check if a replicates cell must be created.
2782    */
2783   public void setCrystal(Crystal crystal, boolean checkReplicatesCell) {
2784     if (checkReplicatesCell) {
2785       this.crystal = ReplicatesCrystal.replicatesCrystalFactory(crystal.getUnitCell(), cutOff2);
2786     } else {
2787       this.crystal = crystal;
2788     }
2789     /*
2790      Update VanDerWaals first, in case the NeighborList needs to be
2791      re-allocated to include a larger number of replicated cells.
2792     */
2793     if (vanderWaalsTerm) {
2794       vanderWaals.setCrystal(this.crystal);
2795     }
2796     if (multipoleTerm) {
2797       particleMeshEwald.setCrystal(this.crystal);
2798     }
2799   }
2800 
2801   /**
2802    * {@inheritDoc}
2803    *
2804    * <p>The previousAcceleration array should only contain previous acceleration data for active
2805    * atoms.
2806    */
2807   @Override
2808   public void setPreviousAcceleration(double[] previousAcceleration) {
2809     if (previousAcceleration == null) {
2810       return;
2811     }
2812     int index = 0;
2813     double[] prev = new double[3];
2814     for (int i = 0; i < nAtoms; i++) {
2815       if (atoms[i].isActive()) {
2816         prev[0] = previousAcceleration[index++];
2817         prev[1] = previousAcceleration[index++];
2818         prev[2] = previousAcceleration[index++];
2819         atoms[i].setPreviousAcceleration(prev);
2820       }
2821     }
2822   }
2823 
2824   /**
2825    * Sets the printOnFailure flag; if override is true, over-rides any existing property. Essentially
2826    * sets the default value of printOnFailure for an algorithm. For example, rotamer optimization
2827    * will generally run into force field issues in the normal course of execution as it tries
2828    * unphysical self and pair configurations, so the algorithm should not print out a large number of
2829    * error PDBs.
2830    *
2831    * @param onFail   To set
2832    * @param override Override properties
2833    */
2834   public void setPrintOnFailure(boolean onFail, boolean override) {
2835     if (override) {
2836       // Ignore any pre-existing value
2837       printOnFailure = onFail;
2838     } else {
2839       try {
2840         molecularAssembly.getForceField().getBoolean("PRINT_ON_FAILURE");
2841         /*
2842          If the call was successful, the property was explicitly set
2843          somewhere and should be kept. If an exception was thrown, the
2844          property was never set explicitly, so over-write.
2845         */
2846       } catch (Exception ex) {
2847         printOnFailure = onFail;
2848       }
2849     }
2850   }
2851 
2852   /**
2853    * {@inheritDoc}
2854    *
2855    * <p>The velocity array should only contain velocity data for active atoms.
2856    */
2857   @Override
2858   public void setVelocity(double[] velocity) {
2859     if (velocity == null) {
2860       return;
2861     }
2862     int index = 0;
2863     double[] vel = new double[3];
2864     for (Atom a : atoms) {
2865       if (a.isActive()) {
2866         vel[0] = velocity[index++];
2867         vel[1] = velocity[index++];
2868         vel[2] = velocity[index++];
2869       } else {
2870         // If the atom is not active, set the velocity to zero.
2871         vel[0] = 0.0;
2872         vel[1] = 0.0;
2873         vel[2] = 0.0;
2874       }
2875       a.setVelocity(vel);
2876     }
2877   }
2878 
2879   /**
2880    * {@inheritDoc}
2881    */
2882   @Override
2883   public String toString() {
2884     StringBuilder sb = new StringBuilder();
2885     sb.append(forceFieldBondedEnergyRegion.toString());
2886     sb.append(restrainEnergyRegion.toString());
2887     if (restrainGroupTerm && nRestrainGroups > 0) {
2888       sb.append(format("  %s %20.8f %12d %12.3f\n", "Restrain Groups   ", restrainGroupEnergy,
2889           nRestrainGroups, restrainGroupTime * toSeconds));
2890     }
2891     if (ncsTerm) {
2892       sb.append(format("  %s %20.8f %12d %12.3f\n", "NCS Restraint     ", ncsEnergy, nAtoms,
2893           ncsTime * toSeconds));
2894     }
2895     if (comTerm) {
2896       sb.append(format("  %s %20.8f %12d %12.3f\n", "COM Restraint     ", comRestraintEnergy, nAtoms,
2897           comRestraintTime * toSeconds));
2898     }
2899     if (vanderWaalsTerm && nVanDerWaalInteractions > 0) {
2900       sb.append(format("  %s %20.8f %12d %12.3f\n", "Van der Waals     ", vanDerWaalsEnergy,
2901           nVanDerWaalInteractions, vanDerWaalsTime * toSeconds));
2902     }
2903     if (multipoleTerm && nPermanentInteractions > 0) {
2904       if (polarizationTerm) {
2905         sb.append(format("  %s %20.8f %12d\n", "Atomic Multipoles ", permanentMultipoleEnergy,
2906             nPermanentInteractions));
2907       } else {
2908         if (elecForm == ELEC_FORM.FIXED_CHARGE) {
2909           sb.append(format("  %s %20.8f %12d %12.3f\n", "Atomic Charges    ", permanentMultipoleEnergy,
2910               nPermanentInteractions, electrostaticTime * toSeconds));
2911         } else
2912           sb.append(format("  %s %20.8f %12d %12.3f\n", "Atomic Multipoles ", permanentMultipoleEnergy,
2913               nPermanentInteractions, electrostaticTime * toSeconds));
2914       }
2915     }
2916     if (polarizationTerm && nPermanentInteractions > 0) {
2917       sb.append(format("  %s %20.8f %12d %12.3f\n", "Polarization      ", polarizationEnergy,
2918           nPermanentInteractions, electrostaticTime * toSeconds));
2919     }
2920     if (generalizedKirkwoodTerm && nGKInteractions > 0) {
2921       sb.append(
2922           format("  %s %20.8f %12d\n", "Solvation         ", solvationEnergy, nGKInteractions));
2923     }
2924     if (esvTerm) {
2925       sb.append(
2926           format("  %s %20.8f  %s\n", "ExtendedSystemBias", esvBias, esvSystem.getLambdaList()));
2927       sb.append(esvSystem.getBiasDecomposition());
2928     }
2929     if (nnTerm) {
2930       sb.append(format("  %s %20.8f %12d %12.3f\n", "Neural Network    ", nnEnergy, nAtoms,
2931           nnTime * toSeconds));
2932     }
2933     sb.append(format("  %s %20.8f  %s %12.3f (sec)",
2934         "Total Potential   ", totalEnergy, "(Kcal/mole)", totalTime * toSeconds));
2935     int nsymm = crystal.getUnitCell().spaceGroup.getNumberOfSymOps();
2936     if (nsymm > 1) {
2937       sb.append(format("\n  %s %20.8f", "Unit Cell         ", totalEnergy * nsymm));
2938     }
2939     if (crystal.getUnitCell() != crystal) {
2940       nsymm = crystal.spaceGroup.getNumberOfSymOps();
2941       sb.append(format("\n  %s %20.8f", "Replicates Cell   ", totalEnergy * nsymm));
2942     }
2943 
2944     return sb.toString();
2945   }
2946 
2947   /**
2948    * Log bonded energy terms and restraints.
2949    */
2950   public void logBondedTermsAndRestraints() {
2951     if (forceFieldBondedEnergyRegion != null) {
2952       forceFieldBondedEnergyRegion.log();
2953     }
2954     if (restrainEnergyRegion != null) {
2955       restrainEnergyRegion.log();
2956     }
2957     if (restrainGroupTerm && nRestrainGroups > 0) {
2958       logger.info("\n Restrain Group Interactions:");
2959       logger.info(restrainGroups.toString());
2960     }
2961   }
2962 
2963   /**
2964    * Setter for the field <code>lambdaBondedTerms</code>.
2965    *
2966    * @param lambdaBondedTerms    a boolean.
2967    * @param lambdaAllBondedTerms a boolean.
2968    */
2969   void setLambdaBondedTerms(boolean lambdaBondedTerms, boolean lambdaAllBondedTerms) {
2970     this.lambdaBondedTerms = lambdaBondedTerms;
2971     this.lambdaAllBondedTerms = lambdaAllBondedTerms;
2972   }
2973 
2974   /**
2975    * Private method for internal use, so we don't have subclasses calling super.energy, and this
2976    * class delegating to the subclass's getGradient method.
2977    *
2978    * @param g Gradient array to fill.
2979    * @return Gradient array.
2980    */
2981   private double[] fillGradient(double[] g) {
2982     assert (g != null);
2983     double[] grad = new double[3];
2984     int n = getNumberOfVariables();
2985     if (g == null || g.length < n) {
2986       g = new double[n];
2987     }
2988     int index = 0;
2989     for (int i = 0; i < nAtoms; i++) {
2990       Atom a = atoms[i];
2991       if (a.isActive()) {
2992         a.getXYZGradient(grad);
2993         double gx = grad[0];
2994         double gy = grad[1];
2995         double gz = grad[2];
2996         if (isNaN(gx) || isInfinite(gx) || isNaN(gy) || isInfinite(gy) || isNaN(gz) || isInfinite(
2997             gz)) {
2998           StringBuilder sb = new StringBuilder(
2999               format("The gradient of atom %s is (%8.3f,%8.3f,%8.3f).", a, gx, gy, gz));
3000           double[] vals = new double[3];
3001           a.getVelocity(vals);
3002           sb.append(format("\n Velocities: %8.3g %8.3g %8.3g", vals[0], vals[1], vals[2]));
3003           a.getAcceleration(vals);
3004           sb.append(format("\n Accelerations: %8.3g %8.3g %8.3g", vals[0], vals[1], vals[2]));
3005           a.getPreviousAcceleration(vals);
3006           sb.append(
3007               format("\n Previous accelerations: %8.3g %8.3g %8.3g", vals[0], vals[1], vals[2]));
3008 
3009           throw new EnergyException(sb.toString());
3010         }
3011         g[index++] = gx;
3012         g[index++] = gy;
3013         g[index++] = gz;
3014       }
3015     }
3016     return g;
3017   }
3018 
3019   /**
3020    * setRestrainDistance
3021    *
3022    * @param a1                a {@link ffx.potential.bonded.Atom} object.
3023    * @param a2                a {@link ffx.potential.bonded.Atom} object.
3024    * @param distance          a double.
3025    * @param forceConstant     the force constant in kcal/mole.
3026    * @param flatBottom        Radius of a flat-bottom potential in Angstroms.
3027    * @param lamStart          At what lambda does the restraint begin to take effect?
3028    * @param lamEnd            At what lambda does the restraint hit full strength?
3029    * @param switchingFunction Switching function to use as a lambda dependence.
3030    */
3031   private RestrainDistance createRestrainDistance(Atom a1, Atom a2, double distance, double forceConstant,
3032                                                   double flatBottom, double lamStart, double lamEnd,
3033                                                   UnivariateSwitchingFunction switchingFunction) {
3034     boolean rbLambda = !(switchingFunction instanceof ConstantSwitch) && lambdaTerm;
3035     RestrainDistance restrainDistance = new RestrainDistance(a1, a2, crystal, rbLambda, lamStart, lamEnd, switchingFunction);
3036     int[] classes = {a1.getAtomType().atomClass, a2.getAtomType().atomClass};
3037     if (flatBottom != 0) {
3038       BondType bondType = new BondType(classes, forceConstant, distance,
3039           BondType.BondFunction.FLAT_BOTTOM_HARMONIC, flatBottom);
3040       restrainDistance.setBondType(bondType);
3041     } else {
3042       BondType bondType = new BondType(classes, forceConstant, distance,
3043           BondType.BondFunction.HARMONIC);
3044       restrainDistance.setBondType(bondType);
3045     }
3046     return restrainDistance;
3047   }
3048 
3049   private Crystal configureNCS(ForceField forceField, Crystal unitCell) {
3050     // MTRIXn to be permuted with standard space group in NCSCrystal.java for experimental
3051     // refinement.
3052     if (forceField.getProperties().containsKey("MTRIX1")
3053         && forceField.getProperties().containsKey("MTRIX2") && forceField.getProperties()
3054         .containsKey("MTRIX3")) {
3055       Crystal unitCell2 = new Crystal(unitCell.a, unitCell.b, unitCell.c, unitCell.alpha,
3056           unitCell.beta, unitCell.gamma, unitCell.spaceGroup.pdbName);
3057       SpaceGroup spaceGroup = unitCell2.spaceGroup;
3058       // Separate string list MTRIXn into Double matricies then pass into symops
3059       CompositeConfiguration properties = forceField.getProperties();
3060       String[] MTRX1List = properties.getStringArray("MTRIX1");
3061       String[] MTRX2List = properties.getStringArray("MTRIX2");
3062       String[] MTRX3List = properties.getStringArray("MTRIX3");
3063       spaceGroup.symOps.clear();
3064       double number1;
3065       double number2;
3066       double number3;
3067       for (int i = 0; i < MTRX1List.length; i++) {
3068         double[][] Rot_MTRX = {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}};
3069         double[] Tr_MTRX = {0, 0, 0};
3070         String[] tokens1 = MTRX1List[i].trim().split(" +"); // 4 items: rot [0][1-3] * trans[0]
3071         String[] tokens2 = MTRX2List[i].trim().split(" +"); // 4 items: rot [1][1-3] * trans[1]
3072         String[] tokens3 = MTRX3List[i].trim().split(" +"); // 4 items: rot [2][1-3] * trans[2]
3073         for (int k = 0; k < 4; k++) {
3074           number1 = Double.parseDouble(tokens1[k]);
3075           number2 = Double.parseDouble(tokens2[k]);
3076           number3 = Double.parseDouble(tokens3[k]);
3077           if (k != 3) {
3078             Rot_MTRX[0][k] = number1;
3079             Rot_MTRX[1][k] = number2;
3080             Rot_MTRX[2][k] = number3;
3081           } else {
3082             Tr_MTRX[0] = number1;
3083             Tr_MTRX[1] = number2;
3084             Tr_MTRX[2] = number3;
3085           }
3086         }
3087         SymOp symOp = new SymOp(Rot_MTRX, Tr_MTRX);
3088         if (logger.isLoggable(Level.FINEST)) {
3089           logger.info(
3090               format(" MTRIXn SymOp: %d of %d\n" + symOp, i + 1, MTRX1List.length));
3091         }
3092         spaceGroup.symOps.add(symOp);
3093       }
3094       unitCell = NCSCrystal.NCSCrystalFactory(unitCell, spaceGroup.symOps);
3095       unitCell.updateCrystal();
3096     }
3097     return unitCell;
3098   }
3099 
3100   public RestrainMode getRestrainMode() {
3101     return restrainMode;
3102   }
3103 
3104   /**
3105    * Getter for the RestrainGroup field.
3106    *
3107    * @return Returns RestrainGroup.
3108    */
3109   public RestrainGroups getRestrainGroups() {
3110     return restrainGroups;
3111   }
3112 
3113   /**
3114    * Get the TorsionTorsionPotentialEnergy.
3115    *
3116    * @return The TorsionTorsionPotentialEnergy.
3117    */
3118   public TorsionTorsionPotentialEnergy getTorsionTorsionPotentialEnergy() {
3119     return torsionTorsionPotentialEnergy;
3120   }
3121 
3122   /**
3123    * Get the AnglePotentialEnergy.
3124    *
3125    * @return The AnglePotentialEnergy.
3126    */
3127   public AnglePotentialEnergy getAnglePotentialEnergy() {
3128     return anglePotentialEnergy;
3129   }
3130 
3131   /**
3132    * Get the BondPotentialEnergy.
3133    *
3134    * @return The BondPotentialEnergy.
3135    */
3136   public BondPotentialEnergy getBondPotentialEnergy() {
3137     return bondPotentialEnergy;
3138   }
3139 
3140   /**
3141    * Get the StretchBendPotentialEnergy.
3142    *
3143    * @return The StretchBendPotentialEnergy.
3144    */
3145   public StretchBendPotentialEnergy getStretchBendPotentialEnergy() {
3146     return stretchBendPotentialEnergy;
3147   }
3148 
3149   /**
3150    * Get the UreyBradleyPotentialEnergy.
3151    *
3152    * @return The UreyBradleyPotentialEnergy.
3153    */
3154   public UreyBradleyPotentialEnergy getUreyBradleyPotentialEnergy() {
3155     return ureyBradleyPotentialEnergy;
3156   }
3157 
3158   /**
3159    * Get the OutOfPlaneBendPotentialEnergy.
3160    *
3161    * @return The OutOfPlaneBendPotentialEnergy.
3162    */
3163   public OutOfPlaneBendPotentialEnergy getOutOfPlaneBendPotentialEnergy() {
3164     return outOfPlaneBendPotentialEnergy;
3165   }
3166 
3167   /**
3168    * Get the TorsionPotentialEnergy.
3169    *
3170    * @return The TorsionPotentialEnergy.
3171    */
3172   public TorsionPotentialEnergy getTorsionPotentialEnergy() {
3173     return torsionPotentialEnergy;
3174   }
3175 
3176   /**
3177    * Get the StretchTorsionPotentialEnergy.
3178    *
3179    * @return The StretchTorsionPotentialEnergy.
3180    */
3181   public StretchTorsionPotentialEnergy getStretchTorsionPotentialEnergy() {
3182     return stretchTorsionPotentialEnergy;
3183   }
3184 
3185   /**
3186    * Get the AngleTorsionPotentialEnergy.
3187    *
3188    * @return The AngleTorsionPotentialEnergy.
3189    */
3190   public AngleTorsionPotentialEnergy getAngleTorsionPotentialEnergy() {
3191     return angleTorsionPotentialEnergy;
3192   }
3193 
3194   /**
3195    * Get the ImproperTorsionPotentialEnergy.
3196    *
3197    * @return The ImproperTorsionPotentialEnergy.
3198    */
3199   public ImproperTorsionPotentialEnergy getImproperTorsionPotentialEnergy() {
3200     return improperTorsionPotentialEnergy;
3201   }
3202 
3203   /**
3204    * Get the PiOrbitalTorsionPotentialEnergy.
3205    *
3206    * @return The PiOrbitalTorsionPotentialEnergy.
3207    */
3208   public PiOrbitalTorsionPotentialEnergy getPiOrbitalTorsionPotentialEnergy() {
3209     return piOrbitalTorsionPotentialEnergy;
3210   }
3211 
3212   /**
3213    * Get the RestrainPositionPotentialEnergy.
3214    *
3215    * @return The RestrainPositionPotentialEnergy.
3216    */
3217   public RestrainPositionPotentialEnergy getRestrainPositionPotentialEnergy() {
3218     return restrainPositionPotentialEnergy;
3219   }
3220 
3221   /**
3222    * Get the RestrainDistancePotentialEnergy.
3223    *
3224    * @return The RestrainDistancePotentialEnergy.
3225    */
3226   public RestrainDistancePotentialEnergy getRestrainDistancePotentialEnergy() {
3227     return restrainDistancePotentialEnergy;
3228   }
3229 
3230   /**
3231    * Get the RestrainTorsionPotentialEnergy.
3232    *
3233    * @return The RestrainTorsionPotentialEnergy.
3234    */
3235   public RestrainTorsionPotentialEnergy getRestrainTorsionPotentialEnergy() {
3236     return restrainTorsionPotentialEnergy;
3237   }
3238 
3239 }