1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
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.UnivariateSwitchingFunction;
53 import ffx.potential.bonded.Angle;
54 import ffx.potential.bonded.AngleTorsion;
55 import ffx.potential.bonded.Atom;
56 import ffx.potential.bonded.Bond;
57 import ffx.potential.bonded.ImproperTorsion;
58 import ffx.potential.bonded.LambdaInterface;
59 import ffx.potential.bonded.MSNode;
60 import ffx.potential.bonded.OutOfPlaneBend;
61 import ffx.potential.bonded.PiOrbitalTorsion;
62 import ffx.potential.bonded.RestrainDistance;
63 import ffx.potential.bonded.RestrainPosition;
64 import ffx.potential.bonded.StretchBend;
65 import ffx.potential.bonded.StretchTorsion;
66 import ffx.potential.bonded.Torsion;
67 import ffx.potential.bonded.TorsionTorsion;
68 import ffx.potential.bonded.UreyBradley;
69 import ffx.potential.constraint.CcmaConstraint;
70 import ffx.potential.constraint.SettleConstraint;
71 import ffx.potential.extended.ExtendedSystem;
72 import ffx.potential.nonbonded.COMRestraint;
73 import ffx.potential.nonbonded.GeneralizedKirkwood;
74 import ffx.potential.nonbonded.NCSRestraint;
75 import ffx.potential.nonbonded.ParticleMeshEwald;
76 import ffx.potential.nonbonded.RestrainGroups;
77 import ffx.potential.nonbonded.VanDerWaals;
78 import ffx.potential.nonbonded.VanDerWaalsTornado;
79 import ffx.potential.openmm.OpenMMEnergy;
80 import ffx.potential.parameters.AngleType.AngleMode;
81 import ffx.potential.parameters.BondType;
82 import ffx.potential.parameters.ForceField;
83 import ffx.potential.parameters.ForceField.ELEC_FORM;
84 import ffx.potential.terms.AnglePotentialEnergy;
85 import ffx.potential.terms.AngleTorsionPotentialEnergy;
86 import ffx.potential.terms.BondPotentialEnergy;
87 import ffx.potential.terms.EnergyTermRegion;
88 import ffx.potential.terms.ImproperTorsionPotentialEnergy;
89 import ffx.potential.terms.OutOfPlaneBendPotentialEnergy;
90 import ffx.potential.terms.PiOrbitalTorsionPotentialEnergy;
91 import ffx.potential.terms.RestrainDistancePotentialEnergy;
92 import ffx.potential.terms.RestrainPositionPotentialEnergy;
93 import ffx.potential.terms.RestrainTorsionPotentialEnergy;
94 import ffx.potential.terms.StretchBendPotentialEnergy;
95 import ffx.potential.terms.StretchTorsionPotentialEnergy;
96 import ffx.potential.terms.TorsionPotentialEnergy;
97 import ffx.potential.terms.TorsionTorsionPotentialEnergy;
98 import ffx.potential.terms.UreyBradleyPotentialEnergy;
99 import ffx.potential.utils.ConvexHullOps;
100 import ffx.potential.utils.EnergyException;
101 import ffx.potential.utils.PotentialsFunctions;
102 import ffx.potential.utils.PotentialsUtils;
103 import ffx.utilities.FFXProperty;
104 import org.apache.commons.configuration2.CompositeConfiguration;
105
106 import javax.annotation.Nullable;
107 import java.io.File;
108 import java.time.LocalDateTime;
109 import java.time.format.DateTimeFormatter;
110 import java.util.ArrayList;
111 import java.util.Arrays;
112 import java.util.Collections;
113 import java.util.HashSet;
114 import java.util.List;
115 import java.util.Set;
116 import java.util.logging.Level;
117 import java.util.logging.Logger;
118 import java.util.stream.Collectors;
119 import java.util.stream.Stream;
120
121 import static ffx.potential.bonded.BondedTerm.removeNeuralNetworkTerms;
122 import static ffx.potential.bonded.RestrainPosition.parseRestrainPositions;
123 import static ffx.potential.nonbonded.VanDerWaalsForm.DEFAULT_VDW_CUTOFF;
124 import static ffx.potential.nonbonded.pme.EwaldParameters.DEFAULT_EWALD_CUTOFF;
125 import static ffx.potential.parameters.ForceField.toEnumForm;
126 import static ffx.potential.parsers.XYZFileFilter.isXYZ;
127 import static ffx.potential.terms.RestrainDistancePotentialEnergy.configureRestrainDistances;
128 import static ffx.potential.terms.RestrainTorsionPotentialEnergy.configureRestrainTorsions;
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
141
142
143
144
145 public class ForceFieldEnergy implements CrystalPotential, LambdaInterface {
146
147
148
149
150 private static final Logger logger = Logger.getLogger(ForceFieldEnergy.class.getName());
151
152
153
154
155 private static final double toSeconds = 1.0e-9;
156
157
158
159
160 public enum RestrainMode {
161
162
163
164 ENERGY,
165
166
167
168 ALCHEMICAL
169 }
170
171
172
173
174 private final Platform platform = Platform.FFX;
175
176
177
178
179 protected final MolecularAssembly molecularAssembly;
180
181
182
183 private final Atom[] atoms;
184
185
186
187 private Crystal crystal;
188
189
190
191 private final ParallelTeam parallelTeam;
192
193
194
195 private final NCSRestraint ncsRestraint;
196
197
198
199 private final COMRestraint comRestraint;
200
201
202
203 private final VanDerWaals vanderWaals;
204
205
206
207 private final ParticleMeshEwald particleMeshEwald;
208
209
210
211 private final ANIEnergy aniEnergy;
212
213
214
215
216 private final List<Constraint> constraints;
217
218
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
256
257 private final double cutOff2;
258
259
260
261 private final double cutoffPlusBuffer;
262
263
264
265 private final ELEC_FORM elecForm;
266
267
268
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
275
276 private double lambda = 1.0;
277
278
279
280 public boolean lambdaBondedTerms = false;
281
282
283
284 boolean lambdaAllBondedTerms = false;
285
286
287
288 private final boolean lambdaTorsions;
289
290
291
292
293
294 private STATE state = STATE.BOTH;
295
296
297
298 protected double[] optimizationScaling = null;
299
300
301
302
303 private final EnergyTermRegion forceFieldBondedEnergyRegion;
304
305
306
307 private final BondPotentialEnergy bondPotentialEnergy;
308
309
310
311 private final AnglePotentialEnergy anglePotentialEnergy;
312
313
314
315 private final StretchBendPotentialEnergy stretchBendPotentialEnergy;
316
317
318
319 private final UreyBradleyPotentialEnergy ureyBradleyPotentialEnergy;
320
321
322
323 private final OutOfPlaneBendPotentialEnergy outOfPlaneBendPotentialEnergy;
324
325
326
327 private final TorsionPotentialEnergy torsionPotentialEnergy;
328
329
330
331 private final StretchTorsionPotentialEnergy stretchTorsionPotentialEnergy;
332
333
334
335 private final AngleTorsionPotentialEnergy angleTorsionPotentialEnergy;
336
337
338
339 private final ImproperTorsionPotentialEnergy improperTorsionPotentialEnergy;
340
341
342
343 private final PiOrbitalTorsionPotentialEnergy piOrbitalTorsionPotentialEnergy;
344
345
346
347 private final TorsionTorsionPotentialEnergy torsionTorsionPotentialEnergy;
348
349
350
351
352 private final EnergyTermRegion restrainEnergyRegion;
353
354
355
356 private final RestrainMode restrainMode;
357
358
359
360 private final RestrainPositionPotentialEnergy restrainPositionPotentialEnergy;
361
362
363
364 private final RestrainDistancePotentialEnergy restrainDistancePotentialEnergy;
365
366
367
368 private final RestrainTorsionPotentialEnergy restrainTorsionPotentialEnergy;
369
370
371
372 private final RestrainGroups restrainGroups;
373
374
375
376
377 private final int nAtoms;
378
379
380
381 private int nRestrainGroups;
382
383
384
385 private int nVanDerWaalInteractions;
386
387
388
389 private int nPermanentInteractions;
390
391
392
393 private int nGKInteractions;
394
395
396
397
398 private boolean nnTerm;
399
400
401
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
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
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
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
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
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
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
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
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
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
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
479
480 private final boolean restrainPositionTerm;
481
482
483
484 private final boolean restrainDistanceTerm;
485
486
487
488 private final boolean restrainTorsionTerm;
489
490
491
492 private boolean restrainGroupTerm;
493
494
495
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
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
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
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
534
535 private boolean comTerm;
536
537
538
539 private boolean ncsTerm;
540
541
542
543
544 private final boolean nnTermOrig;
545
546
547
548 private final boolean restrainGroupTermOrig;
549
550
551
552 private final boolean vanderWaalsTermOrig;
553
554
555
556 private final boolean multipoleTermOrig;
557
558
559
560 private final boolean polarizationTermOrig;
561
562
563
564 private final boolean generalizedKirkwoodTermOrig;
565
566
567
568 private boolean esvTermOrig;
569
570
571
572 private boolean ncsTermOrig;
573
574
575
576 private boolean comTermOrig;
577
578
579
580
581 private double nnEnergy;
582
583
584
585 private double restrainGroupEnergy;
586
587
588
589 private double vanDerWaalsEnergy;
590
591
592
593 private double permanentMultipoleEnergy;
594
595
596
597 private double polarizationEnergy;
598
599
600
601 private double solvationEnergy;
602
603
604
605 private double ncsEnergy;
606
607
608
609 private double comRestraintEnergy;
610
611
612
613 private double totalEnergy;
614
615
616
617
618 private long ncsTime;
619
620
621
622 private long nnTime;
623
624
625
626 private long restrainGroupTime;
627
628
629
630 private long comRestraintTime;
631
632
633
634 private long vanDerWaalsTime;
635
636
637
638 private long electrostaticTime;
639
640
641
642 private long totalTime;
643
644
645
646
647 private ExtendedSystem esvSystem = null;
648 private boolean esvTerm;
649 private double esvBias;
650
651
652
653
654 boolean destroyed = false;
655
656
657
658
659 public final double maxDebugGradient;
660
661
662
663 private boolean printOnFailure;
664
665
666
667
668
669
670 protected ForceFieldEnergy(MolecularAssembly molecularAssembly) {
671 this(molecularAssembly, ParallelTeam.getDefaultThreadCount());
672 }
673
674
675
676
677
678
679
680 protected ForceFieldEnergy(MolecularAssembly molecularAssembly, int numThreads) {
681
682 this.molecularAssembly = molecularAssembly;
683 atoms = molecularAssembly.getAtomArray();
684 nAtoms = atoms.length;
685
686
687 for (int i = 0; i < nAtoms; i++) {
688 int index = atoms[i].getXyzIndex() - 1;
689 assert (i == index);
690 }
691
692
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
704 boolean standardizeAtomNames = properties.getBoolean("standardizeAtomNames", true);
705 if (standardizeAtomNames) {
706 molecularAssembly.renameWaterProtons();
707 }
708
709 nnTerm = forceField.getBoolean("NNTERM", false);
710
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
739 polarizationTerm = false;
740 generalizedKirkwoodTerm = false;
741 }
742 } else {
743
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
762 if (properties.containsKey("restrain-groups")) {
763 restrainGroupTerm = true;
764 } else {
765 restrainGroupTerm = false;
766 }
767
768
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
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
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
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
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
843 forceField.addProperty("EWALD_ALPHA", "0.0");
844
845
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
858 vdwCutoff = DEFAULT_VDW_CUTOFF;
859 ewaldCutoff = DEFAULT_EWALD_CUTOFF;
860 gkCutoff = vdwCutoff;
861 if (unitCell.aperiodic()) {
862
863 vdwCutoff = Double.POSITIVE_INFINITY;
864 ewaldCutoff = Double.POSITIVE_INFINITY;
865 gkCutoff = Double.POSITIVE_INFINITY;
866 }
867
868 vdwCutoff = forceField.getDouble("VDW_CUTOFF", vdwCutoff);
869 ewaldCutoff = forceField.getDouble("EWALD_CUTOFF", ewaldCutoff);
870 gkCutoff = forceField.getDouble("GK_CUTOFF", gkCutoff);
871
872
873 double nlistCutoff = vdwCutoff;
874 if (multipoleTerm) {
875 nlistCutoff = max(vdwCutoff, ewaldCutoff);
876 }
877 if (generalizedKirkwoodTerm) {
878 nlistCutoff = max(nlistCutoff, gkCutoff);
879 }
880
881
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
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
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
917 int nSpecial = checkForSpecialPositions(forceField);
918
919
920 boolean rigidHydrogens = forceField.getBoolean("RIGID_HYDROGENS", false);
921
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
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
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
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
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
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
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
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
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
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
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
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
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
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
1172 RestrainDistance[] restrainDistances = configureRestrainDistances(properties, atoms, crystal, lambdaTerm);
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
1187 Torsion[] torsions = null;
1188 if (torsionPotentialEnergy != null) {
1189 torsions = torsionPotentialEnergy.getTorsionArray();
1190 }
1191 Torsion[] restrainTorsions = configureRestrainTorsions(properties, forceField, atoms, torsions);
1192 if (restrainTorsions != null && restrainTorsions.length > 0) {
1193 restrainTorsionTerm = true;
1194 } else {
1195 restrainTorsionTerm = false;
1196 }
1197 if (restrainTorsionTerm) {
1198 logger.info(" Restrain Torsion Mode: " + restrainMode);
1199 int forceGroup = forceField.getInteger("RESTRAIN_TORSION_FORCE_GROUP", 0);
1200 restrainTorsionPotentialEnergy = new RestrainTorsionPotentialEnergy("Restrain Torsions", forceGroup, java.util.Arrays.asList(restrainTorsions));
1201 restrainEnergyRegion.addEnergyTerm(restrainTorsionPotentialEnergy);
1202 } else {
1203 restrainTorsionPotentialEnergy = null;
1204 }
1205
1206 int[] molecule = molecularAssembly.getMoleculeNumbers();
1207 if (vanderWaalsTerm) {
1208 logger.info("\n Non-Bonded Terms");
1209 boolean[] nn = null;
1210 if (nnTerm) {
1211 nn = molecularAssembly.getNeuralNetworkIdentity();
1212 } else {
1213 nn = new boolean[nAtoms];
1214 Arrays.fill(nn, false);
1215 }
1216
1217 if (!tornadoVM) {
1218 vanderWaals = new VanDerWaals(atoms, molecule, nn, crystal, forceField, parallelTeam,
1219 vdwCutoff, nlistCutoff);
1220 } else {
1221 vanderWaals = new VanDerWaalsTornado(atoms, crystal, forceField, vdwCutoff);
1222 }
1223 } else {
1224 vanderWaals = null;
1225 }
1226
1227 if (nnTerm) {
1228 aniEnergy = new ANIEnergy(molecularAssembly);
1229 } else {
1230 aniEnergy = null;
1231 }
1232
1233 if (multipoleTerm) {
1234 if (name.contains("OPLS") || name.contains("AMBER") || name.contains("CHARMM")) {
1235 elecForm = ELEC_FORM.FIXED_CHARGE;
1236 } else {
1237 elecForm = ELEC_FORM.PAM;
1238 }
1239 particleMeshEwald = new ParticleMeshEwald(atoms, molecule, forceField, crystal,
1240 vanderWaals.getNeighborList(), elecForm, ewaldCutoff, gkCutoff, parallelTeam);
1241 double charge = molecularAssembly.getCharge(checkAllNodeCharges);
1242 logger.info(format("\n Overall system charge: %10.3f", charge));
1243 } else {
1244 elecForm = null;
1245 particleMeshEwald = null;
1246 }
1247
1248 if (ncsTerm) {
1249 String sg = forceField.getString("NCSGROUP", "P 1");
1250 Crystal ncsCrystal = new Crystal(a, b, c, alpha, beta, gamma, sg);
1251 ncsRestraint = new NCSRestraint(atoms, forceField, ncsCrystal);
1252 } else {
1253 ncsRestraint = null;
1254 }
1255
1256 if (restrainGroupTerm) {
1257 restrainGroups = new RestrainGroups(molecularAssembly);
1258 nRestrainGroups = restrainGroups.getNumberOfRestraints();
1259 } else {
1260 restrainGroups = null;
1261 }
1262
1263 if (comTerm) {
1264 comRestraint = new COMRestraint(molecularAssembly);
1265 } else {
1266 comRestraint = null;
1267 }
1268
1269
1270 maxDebugGradient = forceField.getDouble("MAX_DEBUG_GRADIENT", Double.POSITIVE_INFINITY);
1271
1272 molecularAssembly.setPotential(this);
1273
1274
1275 constraints = configureConstraints(forceField);
1276
1277 if (lambdaTerm) {
1278 this.setLambda(1.0);
1279 if (nSpecial == 0) {
1280
1281
1282 crystal.setSpecialPositionCutoff(0.0);
1283 } else {
1284
1285 logger.info(" Special positions checking will be performed during a lambda simulation.");
1286 }
1287 }
1288 }
1289
1290
1291
1292
1293
1294
1295 public MolecularAssembly getMolecularAssembly() {
1296 return molecularAssembly;
1297 }
1298
1299
1300
1301
1302
1303
1304 public void setLambdaTerm(boolean lambdaTerm) {
1305 this.lambdaTerm = lambdaTerm;
1306 }
1307
1308
1309
1310
1311
1312
1313 public boolean getLambdaTerm() {
1314 return lambdaTerm;
1315 }
1316
1317 private int checkForSpecialPositions(ForceField forceField) {
1318
1319 boolean specialPositionsInactive = forceField.getBoolean("SPECIAL_POSITIONS_INACTIVE", true);
1320 int nSpecial = 0;
1321 if (specialPositionsInactive) {
1322 int nSymm = crystal.getNumSymOps();
1323 if (nSymm > 1) {
1324 SpaceGroup spaceGroup = crystal.spaceGroup;
1325 double sp2 = crystal.getSpecialPositionCutoff2();
1326 double[] mate = new double[3];
1327 StringBuilder sb = new StringBuilder("\n Atoms at Special Positions set to Inactive:\n");
1328 for (int i = 0; i < nAtoms; i++) {
1329 Atom atom = atoms[i];
1330 double[] xyz = atom.getXYZ(null);
1331 for (int iSymm = 1; iSymm < nSymm; iSymm++) {
1332 SymOp symOp = spaceGroup.getSymOp(iSymm);
1333 crystal.applySymOp(xyz, mate, symOp);
1334 double dr2 = crystal.image(xyz, mate);
1335 if (dr2 < sp2) {
1336 sb.append(
1337 format(" %s separation with SymOp %d at %8.6f A.\n", atoms[i], iSymm, sqrt(dr2)));
1338 atom.setActive(false);
1339 nSpecial++;
1340 break;
1341 }
1342 }
1343 }
1344 if (nSpecial > 0) {
1345 logger.info(sb.toString());
1346 }
1347 }
1348 }
1349 return nSpecial;
1350 }
1351
1352
1353
1354
1355
1356
1357
1358 private List<Constraint> configureConstraints(ForceField forceField) {
1359 String constraintStrings = forceField.getString("CONSTRAIN", forceField.getString("RATTLE", null));
1360 if (constraintStrings == null) {
1361 return Collections.emptyList();
1362 }
1363
1364 ArrayList<Constraint> constraints = new ArrayList<>();
1365 logger.info(format(" Experimental: parsing constraints option %s", constraintStrings));
1366
1367 Set<Bond> numericBonds = new HashSet<>(1);
1368 Set<Angle> numericAngles = new HashSet<>(1);
1369
1370
1371 if (constraintStrings.isEmpty() || constraintStrings.matches("^\\s*$")) {
1372
1373 logger.info(" Constraining X-H bonds.");
1374 Bond[] bonds = bondPotentialEnergy.getBondArray();
1375 numericBonds = Arrays.stream(bonds).filter(
1376 (Bond bond) -> bond.getAtom(0).getAtomicNumber() == 1
1377 || bond.getAtom(1).getAtomicNumber() == 1).collect(Collectors.toSet());
1378 } else {
1379 String[] constraintToks = constraintStrings.split("\\s+");
1380
1381
1382 for (String tok : constraintToks) {
1383 if (tok.equalsIgnoreCase("WATER")) {
1384 logger.info(" Constraining waters to be rigid based on angle & bonds.");
1385
1386
1387 Stream<MSNode> settleStream = molecularAssembly.getMolecules().stream()
1388 .filter((MSNode m) -> m.getAtomList().size() == 3).filter((MSNode m) -> {
1389 List<Atom> atoms = m.getAtomList();
1390 Atom O = null;
1391 List<Atom> H = new ArrayList<>(2);
1392 for (Atom at : atoms) {
1393 int atN = at.getAtomicNumber();
1394 if (atN == 8) {
1395 O = at;
1396 } else if (atN == 1) {
1397 H.add(at);
1398 }
1399 }
1400 return O != null && H.size() == 2;
1401 });
1402
1403 settleStream = Stream.concat(settleStream, molecularAssembly.getWater().stream());
1404
1405 List<SettleConstraint> settleConstraints = settleStream.map(
1406 (MSNode m) -> m.getAngleList().getFirst()).map(SettleConstraint::settleFactory).toList();
1407 constraints.addAll(settleConstraints);
1408
1409 } else if (tok.equalsIgnoreCase("DIATOMIC")) {
1410 logger.severe(" Diatomic distance constraints not yet implemented properly.");
1411 } else if (tok.equalsIgnoreCase("TRIATOMIC")) {
1412 logger.severe(
1413 " Triatomic SETTLE constraints for non-water molecules not yet implemented properly.");
1414 }
1415 }
1416
1417
1418 for (String tok : constraintToks) {
1419 if (tok.equalsIgnoreCase("BONDS") && bondPotentialEnergy != null) {
1420 Bond[] bonds = bondPotentialEnergy.getBondArray();
1421 numericBonds = new HashSet<>(Arrays.asList(bonds));
1422 } else if (tok.equalsIgnoreCase("ANGLES") && anglePotentialEnergy != null) {
1423 Angle[] angles = anglePotentialEnergy.getAngleArray();
1424 numericAngles = new HashSet<>(Arrays.asList(angles));
1425 }
1426 }
1427 }
1428
1429
1430 for (Angle angle : numericAngles) {
1431 angle.getBondList().forEach(numericBonds::remove);
1432 }
1433
1434
1435 List<Angle> ccmaAngles = numericAngles.stream().filter((Angle ang) -> !ang.isConstrained()).collect(Collectors.toList());
1436 List<Bond> ccmaBonds = numericBonds.stream().filter((Bond bond) -> !bond.isConstrained()).collect(Collectors.toList());
1437
1438 CcmaConstraint ccmaConstraint = CcmaConstraint.ccmaFactory(ccmaBonds, ccmaAngles, atoms,
1439 getMass(), CcmaConstraint.DEFAULT_CCMA_NONZERO_CUTOFF);
1440 constraints.add(ccmaConstraint);
1441 logger.info(format(" Added %d constraints.", constraints.size()));
1442
1443 return constraints;
1444 }
1445
1446
1447
1448
1449
1450
1451
1452 public static ForceFieldEnergy energyFactory(MolecularAssembly assembly) {
1453 return energyFactory(assembly, ParallelTeam.getDefaultThreadCount());
1454 }
1455
1456
1457
1458
1459
1460
1461
1462
1463 @SuppressWarnings("fallthrough")
1464 public static ForceFieldEnergy energyFactory(MolecularAssembly assembly, int numThreads) {
1465 ForceField forceField = assembly.getForceField();
1466 String platformString = toEnumForm(forceField.getString("PLATFORM", "FFX"));
1467 Platform platform;
1468 try {
1469 platform = Platform.valueOf(platformString);
1470 } catch (IllegalArgumentException e) {
1471 logger.warning(format(" String %s did not match a known energy implementation", platformString));
1472 platform = Platform.FFX;
1473 }
1474
1475
1476 String dtPlatformString = toEnumForm(forceField.getString("PLATFORM_DT", "FFX"));
1477 try {
1478 Platform dtPlatform = Platform.valueOf(dtPlatformString);
1479
1480 if (dtPlatform != Platform.FFX) {
1481
1482 platform = Platform.FFX;
1483 }
1484 } catch (IllegalArgumentException e) {
1485
1486 }
1487
1488 switch (platform) {
1489 case OMM, OMM_REF, OMM_CUDA, OMM_OPENCL:
1490 try {
1491 return new OpenMMEnergy(assembly, platform, numThreads);
1492 } catch (Exception ex) {
1493 logger.warning(format(" Exception creating OpenMMEnergy: %s", ex));
1494 ForceFieldEnergy ffxEnergy = assembly.getPotentialEnergy();
1495 if (ffxEnergy == null) {
1496 ffxEnergy = new ForceFieldEnergy(assembly, numThreads);
1497 assembly.setPotential(ffxEnergy);
1498 }
1499 return ffxEnergy;
1500 }
1501 case OMM_CPU:
1502 logger.warning(format(" Platform %s not supported; defaulting to FFX", platform));
1503 default:
1504 ForceFieldEnergy ffxEnergy = assembly.getPotentialEnergy();
1505 if (ffxEnergy == null) {
1506 ffxEnergy = new ForceFieldEnergy(assembly, numThreads);
1507 assembly.setPotential(ffxEnergy);
1508 }
1509 return ffxEnergy;
1510 }
1511 }
1512
1513
1514
1515
1516
1517
1518 public Atom[] getAtomArray() {
1519 return atoms;
1520 }
1521
1522
1523
1524
1525
1526
1527
1528 public void applyAllConstraintPositions(double[] xPrior, double[] xNew) {
1529 applyAllConstraintPositions(xPrior, xNew, DEFAULT_CONSTRAINT_TOLERANCE);
1530 }
1531
1532
1533
1534
1535
1536
1537
1538
1539 public void applyAllConstraintPositions(double[] xPrior, double[] xNew, double tol) {
1540 if (xPrior == null) {
1541 xPrior = Arrays.copyOf(xNew, xNew.length);
1542 }
1543 for (Constraint constraint : constraints) {
1544 constraint.applyConstraintToStep(xPrior, xNew, getMass(), tol);
1545 }
1546 }
1547
1548
1549
1550
1551
1552
1553
1554 public void attachExtendedSystem(ExtendedSystem system) {
1555 if (system == null) {
1556 throw new IllegalArgumentException();
1557 }
1558 esvTerm = true;
1559 esvTermOrig = esvTerm;
1560 esvSystem = system;
1561 if (vanderWaalsTerm) {
1562 if (vanderWaals == null) {
1563 logger.warning("Null VdW during ESV setup.");
1564 } else {
1565 vanderWaals.attachExtendedSystem(system);
1566 }
1567 }
1568 if (multipoleTerm) {
1569 if (particleMeshEwald == null) {
1570 logger.warning("Null PME during ESV setup.");
1571 } else {
1572 particleMeshEwald.attachExtendedSystem(system);
1573 }
1574 }
1575 }
1576
1577
1578
1579
1580
1581
1582 @Override
1583 public boolean dEdLZeroAtEnds() {
1584
1585
1586 return !lambdaTerm;
1587 }
1588
1589
1590
1591
1592
1593
1594 public boolean destroy() {
1595 if (destroyed) {
1596
1597
1598 logger.fine(format(" This ForceFieldEnergy is already destroyed: %s", this));
1599 return true;
1600 } else {
1601 try {
1602 if (parallelTeam != null) {
1603 parallelTeam.shutdown();
1604 }
1605 if (vanderWaals != null) {
1606 vanderWaals.destroy();
1607 }
1608 if (particleMeshEwald != null) {
1609 particleMeshEwald.destroy();
1610 }
1611 molecularAssembly.finishDestruction();
1612 destroyed = true;
1613 return true;
1614 } catch (Exception ex) {
1615 logger.warning(format(" Exception in shutting down a ForceFieldEnergy: %s", ex));
1616 logger.info(Utilities.stackTraceToString(ex));
1617 return false;
1618 }
1619 }
1620 }
1621
1622
1623
1624
1625
1626
1627 public double energy() {
1628 return energy(false, false);
1629 }
1630
1631
1632
1633
1634
1635
1636
1637
1638 public double energy(boolean gradient, boolean print) {
1639 try {
1640 totalTime = System.nanoTime();
1641 nnTime = 0;
1642 restrainGroupTime = 0;
1643 vanDerWaalsTime = 0;
1644 electrostaticTime = 0;
1645 ncsTime = 0;
1646
1647
1648 double forceFieldBondedEnergy = 0.0;
1649 nnEnergy = 0.0;
1650
1651
1652 double restrainEnergy = 0.0;
1653 restrainGroupEnergy = 0.0;
1654 ncsEnergy = 0.0;
1655
1656
1657 double totalNonBondedEnergy = 0.0;
1658 double totalMultipoleEnergy = 0.0;
1659 vanDerWaalsEnergy = 0.0;
1660 permanentMultipoleEnergy = 0.0;
1661 polarizationEnergy = 0.0;
1662
1663
1664 solvationEnergy = 0.0;
1665 esvBias = 0.0;
1666
1667
1668 totalEnergy = 0.0;
1669
1670
1671 try {
1672 forceFieldBondedEnergyRegion.setGradient(gradient);
1673 forceFieldBondedEnergyRegion.setLambdaBondedTerms(lambdaBondedTerms);
1674 forceFieldBondedEnergyRegion.setLambdaAllBondedTerms(lambdaAllBondedTerms);
1675 forceFieldBondedEnergyRegion.setInitAtomGradients(true);
1676 parallelTeam.execute(forceFieldBondedEnergyRegion);
1677 forceFieldBondedEnergy = forceFieldBondedEnergyRegion.getEnergy();
1678 } catch (RuntimeException ex) {
1679 logger.warning("Runtime exception during bonded term calculation.");
1680 throw ex;
1681 } catch (Exception ex) {
1682 logger.info(Utilities.stackTraceToString(ex));
1683 logger.severe(ex.toString());
1684 }
1685
1686
1687 try {
1688 if ((restrainMode == RestrainMode.ENERGY) ||
1689 (restrainMode == RestrainMode.ALCHEMICAL && lambdaBondedTerms)) {
1690 boolean checkAlchemicalAtoms = (restrainMode == RestrainMode.ENERGY);
1691 restrainEnergyRegion.setGradient(gradient);
1692 restrainEnergyRegion.setLambdaBondedTerms(lambdaBondedTerms);
1693 restrainEnergyRegion.setLambdaAllBondedTerms(lambdaAllBondedTerms);
1694 restrainEnergyRegion.setCheckAlchemicalAtoms(checkAlchemicalAtoms);
1695 restrainEnergyRegion.setInitAtomGradients(false);
1696 parallelTeam.execute(restrainEnergyRegion);
1697 restrainEnergy = restrainEnergyRegion.getEnergy();
1698 }
1699 } catch (RuntimeException ex) {
1700 logger.warning("Runtime exception during restrain term calculation.");
1701 throw ex;
1702 } catch (Exception ex) {
1703 logger.info(Utilities.stackTraceToString(ex));
1704 logger.severe(ex.toString());
1705 }
1706
1707 if (!lambdaBondedTerms) {
1708
1709 if (ncsTerm) {
1710 ncsTime = -System.nanoTime();
1711 ncsEnergy = ncsRestraint.residual(gradient, print);
1712 ncsTime += System.nanoTime();
1713 }
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725 if (restrainGroupTerm) {
1726 restrainGroupTime = -System.nanoTime();
1727 restrainGroupEnergy = restrainGroups.energy(gradient);
1728 restrainGroupTime += System.nanoTime();
1729 }
1730
1731 if (comTerm) {
1732 comRestraintTime = -System.nanoTime();
1733 comRestraintEnergy = comRestraint.residual(gradient, print);
1734 comRestraintTime += System.nanoTime();
1735 }
1736
1737
1738
1739 if (nnTerm) {
1740 nnTime = -System.nanoTime();
1741 nnEnergy = aniEnergy.energy(gradient, print);
1742 nnTime += System.nanoTime();
1743 }
1744
1745
1746 if (vanderWaalsTerm) {
1747 vanDerWaalsTime = -System.nanoTime();
1748 vanDerWaalsEnergy = vanderWaals.energy(gradient, print);
1749 nVanDerWaalInteractions = this.vanderWaals.getInteractions();
1750 vanDerWaalsTime += System.nanoTime();
1751 }
1752
1753 if (multipoleTerm) {
1754 electrostaticTime = -System.nanoTime();
1755 totalMultipoleEnergy = particleMeshEwald.energy(gradient, print);
1756 permanentMultipoleEnergy = particleMeshEwald.getPermanentEnergy();
1757 polarizationEnergy = particleMeshEwald.getPolarizationEnergy();
1758 nPermanentInteractions = particleMeshEwald.getInteractions();
1759 solvationEnergy = particleMeshEwald.getSolvationEnergy();
1760 nGKInteractions = particleMeshEwald.getGKInteractions();
1761 electrostaticTime += System.nanoTime();
1762 }
1763 }
1764
1765 totalTime = System.nanoTime() - totalTime;
1766
1767 double totalBondedEnergy = forceFieldBondedEnergy + restrainEnergy
1768 + nnEnergy + ncsEnergy + comRestraintEnergy + restrainGroupEnergy;
1769
1770 totalNonBondedEnergy = vanDerWaalsEnergy + totalMultipoleEnergy;
1771 totalEnergy = totalBondedEnergy + totalNonBondedEnergy + solvationEnergy;
1772 if (esvTerm) {
1773 esvBias = esvSystem.getBiasEnergy();
1774 totalEnergy += esvBias;
1775 }
1776
1777 if (print) {
1778 StringBuilder sb = new StringBuilder();
1779 if (gradient) {
1780 sb.append("\n Computed Potential Energy and Atomic Coordinate Gradients\n");
1781 } else {
1782 sb.append("\n Computed Potential Energy\n");
1783 }
1784 sb.append(this);
1785 logger.info(sb.toString());
1786 }
1787 return totalEnergy;
1788 } catch (EnergyException ex) {
1789 if (printOnFailure) {
1790 printFailure();
1791 }
1792 if (ex.doCauseSevere()) {
1793 logger.info(Utilities.stackTraceToString(ex));
1794 logger.log(Level.SEVERE, " Error in calculating energies or gradients", ex);
1795 } else {
1796 logger.log(Level.INFO, format(" Exception in energy calculation:\n %s", ex));
1797 }
1798 throw ex;
1799 }
1800 }
1801
1802
1803
1804
1805 @Override
1806 public double energy(double[] x) {
1807 return energy(x, false);
1808 }
1809
1810
1811
1812
1813 @Override
1814 public double energy(double[] x, boolean verbose) {
1815 assert Arrays.stream(x).allMatch(Double::isFinite);
1816
1817
1818 unscaleCoordinates(x);
1819
1820
1821 setCoordinates(x);
1822
1823 double e = this.energy(false, verbose);
1824
1825
1826 scaleCoordinates(x);
1827
1828 return e;
1829 }
1830
1831
1832
1833
1834 @Override
1835 public double energyAndGradient(double[] x, double[] g) {
1836 return energyAndGradient(x, g, false);
1837 }
1838
1839
1840
1841
1842 @Override
1843 public double energyAndGradient(double[] x, double[] g, boolean verbose) {
1844
1845 unscaleCoordinates(x);
1846
1847 setCoordinates(x);
1848 double e = energy(true, verbose);
1849
1850
1851
1852 try {
1853 fillGradient(g);
1854
1855
1856 scaleCoordinatesAndGradient(x, g);
1857
1858 if (maxDebugGradient < Double.MAX_VALUE) {
1859 boolean extremeGrad = Arrays.stream(g)
1860 .anyMatch((double gi) -> (gi > maxDebugGradient || gi < -maxDebugGradient));
1861 if (extremeGrad) {
1862 File origFile = molecularAssembly.getFile();
1863 String timeString = LocalDateTime.now()
1864 .format(DateTimeFormatter.ofPattern("yyyy_MM_dd-HH_mm_ss"));
1865
1866 String filename = format("%s-LARGEGRAD-%s.pdb",
1867 removeExtension(molecularAssembly.getFile().getName()), timeString);
1868 PotentialsFunctions ef = new PotentialsUtils();
1869 filename = ef.versionFile(filename);
1870
1871 logger.warning(format(" Excessively large gradient detected; printing snapshot to file %s",
1872 filename));
1873 ef.saveAsPDB(molecularAssembly, new File(filename));
1874 molecularAssembly.setFile(origFile);
1875 }
1876 }
1877 return e;
1878 } catch (EnergyException ex) {
1879 if (printOnFailure) {
1880 printFailure();
1881 }
1882 if (ex.doCauseSevere()) {
1883 logger.info(Utilities.stackTraceToString(ex));
1884 logger.log(Level.SEVERE, " Error in calculating energies or gradients", ex);
1885 } else {
1886 logger.log(Level.INFO, format(" Exception in energy calculation:\n %s", ex));
1887 }
1888 throw ex;
1889 }
1890 }
1891
1892
1893
1894
1895 private void printFailure() {
1896 String timeString = LocalDateTime.now()
1897 .format(DateTimeFormatter.ofPattern("yyyy_MM_dd-HH_mm_ss"));
1898 File file = molecularAssembly.getFile();
1899 String ext = "pdb";
1900 if (isXYZ(file)) {
1901 ext = "xyz";
1902 }
1903 String filename = format("%s-ERROR-%s.%s", removeExtension(file.getName()), timeString, ext);
1904 PotentialsFunctions ef = new PotentialsUtils();
1905 filename = ef.versionFile(filename);
1906 logger.info(format(" Writing on-error snapshot to file %s", filename));
1907 ef.save(molecularAssembly, new File(filename));
1908 }
1909
1910
1911
1912
1913
1914
1915 @Override
1916 public double[] getAcceleration(double[] acceleration) {
1917 int n = getNumberOfVariables();
1918 if (acceleration == null || acceleration.length < n) {
1919 acceleration = new double[n];
1920 }
1921 int index = 0;
1922 double[] a = new double[3];
1923 for (int i = 0; i < nAtoms; i++) {
1924 if (atoms[i].isActive()) {
1925 atoms[i].getAcceleration(a);
1926 acceleration[index++] = a[0];
1927 acceleration[index++] = a[1];
1928 acceleration[index++] = a[2];
1929 }
1930 }
1931 return acceleration;
1932 }
1933
1934
1935
1936
1937
1938
1939
1940 public Angle[] getAngles(AngleMode angleMode) {
1941 if (anglePotentialEnergy == null) {
1942 return null;
1943 }
1944 Angle[] angles = anglePotentialEnergy.getAngleArray();
1945 int nAngles = angles.length;
1946 List<Angle> angleList = new ArrayList<>();
1947
1948 for (int i = 0; i < nAngles; i++) {
1949 if (angles[i].getAngleMode() == angleMode) {
1950 angleList.add(angles[i]);
1951 }
1952 }
1953 nAngles = angleList.size();
1954 if (nAngles < 1) {
1955 return null;
1956 }
1957 return angleList.toArray(new Angle[0]);
1958 }
1959
1960
1961
1962
1963
1964
1965 @Override
1966 public List<Constraint> getConstraints() {
1967 return constraints.isEmpty() ? Collections.emptyList() : new ArrayList<>(constraints);
1968 }
1969
1970
1971
1972
1973 @Override
1974 public double[] getCoordinates(double[] x) {
1975 int n = getNumberOfVariables();
1976 if (x == null || x.length < n) {
1977 x = new double[n];
1978 }
1979 int index = 0;
1980 for (int i = 0; i < nAtoms; i++) {
1981 Atom a = atoms[i];
1982 if (a.isActive()) {
1983 x[index++] = a.getX();
1984 x[index++] = a.getY();
1985 x[index++] = a.getZ();
1986 }
1987 }
1988 return x;
1989 }
1990
1991
1992
1993
1994
1995
1996 @Override
1997 public Crystal getCrystal() {
1998 return crystal;
1999 }
2000
2001
2002
2003
2004
2005
2006 @Override
2007 public void setCrystal(Crystal crystal) {
2008 setCrystal(crystal, false);
2009 }
2010
2011
2012
2013
2014
2015
2016 public double getCutoffPlusBuffer() {
2017 return cutoffPlusBuffer;
2018 }
2019
2020
2021
2022
2023 @Override
2024 public STATE getEnergyTermState() {
2025 return state;
2026 }
2027
2028
2029
2030
2031
2032
2033 @Override
2034 public void setEnergyTermState(STATE state) {
2035 this.state = state;
2036 if (forceFieldBondedEnergyRegion != null) {
2037 forceFieldBondedEnergyRegion.setState(state);
2038 }
2039 if (restrainEnergyRegion != null) {
2040 restrainEnergyRegion.setState(state);
2041 }
2042 switch (state) {
2043 case FAST:
2044 nnTerm = nnTermOrig;
2045 restrainGroupTerm = restrainGroupTermOrig;
2046 ncsTerm = ncsTermOrig;
2047 comTerm = comTermOrig;
2048 vanderWaalsTerm = false;
2049 multipoleTerm = false;
2050 polarizationTerm = false;
2051 generalizedKirkwoodTerm = false;
2052 esvTerm = false;
2053 break;
2054 case SLOW:
2055 vanderWaalsTerm = vanderWaalsTermOrig;
2056 multipoleTerm = multipoleTermOrig;
2057 polarizationTerm = polarizationTermOrig;
2058 generalizedKirkwoodTerm = generalizedKirkwoodTermOrig;
2059 esvTerm = esvTermOrig;
2060 nnTerm = false;
2061 restrainGroupTerm = false;
2062 ncsTerm = false;
2063 comTerm = false;
2064 break;
2065 default:
2066 nnTerm = nnTermOrig;
2067 restrainGroupTerm = restrainGroupTermOrig;
2068 ncsTerm = ncsTermOrig;
2069 comTermOrig = comTerm;
2070 vanderWaalsTerm = vanderWaalsTermOrig;
2071 multipoleTerm = multipoleTermOrig;
2072 polarizationTerm = polarizationTermOrig;
2073 generalizedKirkwoodTerm = generalizedKirkwoodTermOrig;
2074 }
2075 }
2076
2077
2078
2079
2080
2081
2082 public double getEsvBiasEnergy() {
2083 return esvBias;
2084 }
2085
2086
2087
2088
2089
2090
2091 public ExtendedSystem getExtendedSystem() {
2092 return esvSystem;
2093 }
2094
2095
2096
2097
2098
2099
2100 public GeneralizedKirkwood getGK() {
2101 if (particleMeshEwald != null) {
2102 return particleMeshEwald.getGK();
2103 } else {
2104 return null;
2105 }
2106 }
2107
2108
2109
2110
2111
2112
2113
2114 public double[] getGradient(double[] g) {
2115 return fillGradient(g);
2116 }
2117
2118
2119
2120
2121 @Override
2122 public double getLambda() {
2123 return lambda;
2124 }
2125
2126
2127
2128
2129 @Override
2130 public void setLambda(double lambda) {
2131 if (lambdaTerm) {
2132 if (lambda <= 1.0 && lambda >= 0.0) {
2133 this.lambda = lambda;
2134 if (vanderWaalsTerm) {
2135 vanderWaals.setLambda(lambda);
2136 }
2137 if (multipoleTerm) {
2138 particleMeshEwald.setLambda(lambda);
2139 }
2140
2141
2142
2143
2144
2145
2146
2147
2148 if (restrainDistanceTerm && restrainDistancePotentialEnergy != null) {
2149 for (RestrainDistance restrainDistance : restrainDistancePotentialEnergy.getRestrainDistances()) {
2150 restrainDistance.setLambda(lambda);
2151 }
2152 }
2153 if (restrainTorsionTerm && restrainTorsionPotentialEnergy != null) {
2154 for (Torsion restrainTorsion : restrainTorsionPotentialEnergy.getRestrainTorsions()) {
2155 restrainTorsion.setLambda(lambda);
2156 }
2157 }
2158 if (ncsTerm && ncsRestraint != null) {
2159 ncsRestraint.setLambda(lambda);
2160 }
2161 if (comTerm && comRestraint != null) {
2162 comRestraint.setLambda(lambda);
2163 }
2164 if (lambdaTorsions) {
2165 if (torsionPotentialEnergy != null) {
2166 torsionPotentialEnergy.setLambda(lambda);
2167 }
2168 if (piOrbitalTorsionPotentialEnergy != null) {
2169 piOrbitalTorsionPotentialEnergy.setLambda(lambda);
2170 }
2171 if (torsionTorsionPotentialEnergy != null) {
2172 torsionTorsionPotentialEnergy.setLambda(lambda);
2173 }
2174 }
2175 } else {
2176 String message = format("Lambda value %8.3f is not in the range [0..1].", lambda);
2177 logger.warning(message);
2178 }
2179 } else {
2180 logger.fine(" Attempting to set a lambda value on a ForceFieldEnergy with lambdaterm false.");
2181 }
2182 }
2183
2184
2185
2186
2187 @Override
2188 public double[] getMass() {
2189 int n = getNumberOfVariables();
2190 double[] mass = new double[n];
2191 int index = 0;
2192 for (int i = 0; i < nAtoms; i++) {
2193 Atom a = atoms[i];
2194 if (a.isActive()) {
2195 double m = a.getMass();
2196 mass[index++] = m;
2197 mass[index++] = m;
2198 mass[index++] = m;
2199 }
2200 }
2201 return mass;
2202 }
2203
2204
2205
2206
2207 @Override
2208 public int getNumberOfVariables() {
2209 int nActive = 0;
2210 for (int i = 0; i < nAtoms; i++) {
2211 if (atoms[i].isActive()) {
2212 nActive++;
2213 }
2214 }
2215 return nActive * 3;
2216 }
2217
2218
2219
2220
2221
2222
2223 public String getPDBHeaderString() {
2224 energy(false, false);
2225 StringBuilder sb = new StringBuilder();
2226 sb.append("REMARK 3 CALCULATED POTENTIAL ENERGY\n");
2227 sb.append(forceFieldBondedEnergyRegion.toPDBString());
2228 sb.append(restrainEnergyRegion.toPDBString());
2229 if (nnTerm) {
2230 sb.append(
2231 format("REMARK 3 %s %g (%d)\n", "NEUTRAL NETWORK : ", nnEnergy, nAtoms));
2232 }
2233 if (ncsTerm) {
2234 sb.append(
2235 format("REMARK 3 %s %g (%d)\n", "NCS RESTRAINT : ", ncsEnergy, nAtoms));
2236 }
2237 if (comTerm) {
2238 sb.append(
2239 format("REMARK 3 %s %g (%d)\n", "COM RESTRAINT : ", comRestraintEnergy,
2240 nAtoms));
2241 }
2242 if (vanderWaalsTerm) {
2243 sb.append(
2244 format("REMARK 3 %s %g (%d)\n", "VAN DER WAALS : ", vanDerWaalsEnergy,
2245 nVanDerWaalInteractions));
2246 }
2247 if (multipoleTerm) {
2248 sb.append(format("REMARK 3 %s %g (%d)\n", "ATOMIC MULTIPOLES : ",
2249 permanentMultipoleEnergy, nPermanentInteractions));
2250 }
2251 if (polarizationTerm) {
2252 sb.append(format("REMARK 3 %s %g (%d)\n", "POLARIZATION : ",
2253 polarizationEnergy, nPermanentInteractions));
2254 }
2255 sb.append(format("REMARK 3 %s %g\n", "TOTAL POTENTIAL (KCAL/MOL) : ", totalEnergy));
2256 int nsymm = crystal.getUnitCell().spaceGroup.getNumberOfSymOps();
2257 if (nsymm > 1) {
2258 sb.append(format("REMARK 3 %s %g\n", "UNIT CELL POTENTIAL : ", totalEnergy * nsymm));
2259 }
2260 if (crystal.getUnitCell() != crystal) {
2261 nsymm = crystal.spaceGroup.getNumberOfSymOps();
2262 if (nsymm > 1) {
2263 sb.append(format("REMARK 3 %s %g\n", "REPLICATES CELL POTENTIAL : ", totalEnergy * nsymm));
2264 }
2265 }
2266 sb.append("REMARK 3\n");
2267
2268 return sb.toString();
2269 }
2270
2271
2272
2273
2274
2275
2276 public ParallelTeam getParallelTeam() {
2277 return parallelTeam;
2278 }
2279
2280
2281
2282
2283
2284
2285 public int getPermanentInteractions() {
2286 return nPermanentInteractions;
2287 }
2288
2289
2290
2291
2292
2293
2294 public double getPermanentMultipoleEnergy() {
2295 return permanentMultipoleEnergy;
2296 }
2297
2298
2299
2300
2301
2302
2303
2304 public Platform getPlatform() {
2305 return platform;
2306 }
2307
2308
2309
2310
2311
2312
2313 public ParticleMeshEwald getPmeNode() {
2314 return particleMeshEwald;
2315 }
2316
2317
2318
2319
2320
2321
2322 public double getPolarizationEnergy() {
2323 return polarizationEnergy;
2324 }
2325
2326
2327
2328
2329
2330
2331 @Override
2332 public double[] getPreviousAcceleration(double[] previousAcceleration) {
2333 int n = getNumberOfVariables();
2334 if (previousAcceleration == null || previousAcceleration.length < n) {
2335 previousAcceleration = new double[n];
2336 }
2337 int index = 0;
2338 double[] a = new double[3];
2339 for (int i = 0; i < nAtoms; i++) {
2340 if (atoms[i].isActive()) {
2341 atoms[i].getPreviousAcceleration(a);
2342 previousAcceleration[index++] = a[0];
2343 previousAcceleration[index++] = a[1];
2344 previousAcceleration[index++] = a[2];
2345 }
2346 }
2347 return previousAcceleration;
2348 }
2349
2350
2351
2352
2353 @Override
2354 public double[] getScaling() {
2355 return optimizationScaling;
2356 }
2357
2358
2359
2360
2361 @Override
2362 public void setScaling(double[] scaling) {
2363 optimizationScaling = scaling;
2364 }
2365
2366
2367
2368
2369
2370
2371 public double getSolvationEnergy() {
2372 return solvationEnergy;
2373 }
2374
2375
2376
2377
2378
2379
2380 public int getSolvationInteractions() {
2381 return nGKInteractions;
2382 }
2383
2384
2385
2386
2387 @Override
2388 public double getTotalEnergy() {
2389 return totalEnergy;
2390 }
2391
2392
2393
2394
2395
2396
2397 public double getVanDerWaalsEnergy() {
2398 return vanDerWaalsEnergy;
2399 }
2400
2401
2402
2403
2404
2405
2406 public int getVanDerWaalsInteractions() {
2407 return nVanDerWaalInteractions;
2408 }
2409
2410
2411
2412
2413
2414
2415 @Override
2416 public VARIABLE_TYPE[] getVariableTypes() {
2417 int n = getNumberOfVariables();
2418 VARIABLE_TYPE[] type = new VARIABLE_TYPE[n];
2419 int i = 0;
2420 for (int j = 0; j < nAtoms; j++) {
2421 if (atoms[j].isActive()) {
2422 type[i++] = VARIABLE_TYPE.X;
2423 type[i++] = VARIABLE_TYPE.Y;
2424 type[i++] = VARIABLE_TYPE.Z;
2425 }
2426 }
2427 return type;
2428 }
2429
2430
2431
2432
2433
2434
2435 public VanDerWaals getVdwNode() {
2436 return vanderWaals;
2437 }
2438
2439
2440
2441
2442
2443
2444 @Override
2445 public double[] getVelocity(double[] velocity) {
2446 int n = getNumberOfVariables();
2447 if (velocity == null || velocity.length < n) {
2448 velocity = new double[n];
2449 }
2450 int index = 0;
2451 double[] v = new double[3];
2452 for (int i = 0; i < nAtoms; i++) {
2453 Atom a = atoms[i];
2454 if (a.isActive()) {
2455 a.getVelocity(v);
2456 velocity[index++] = v[0];
2457 velocity[index++] = v[1];
2458 velocity[index++] = v[2];
2459 }
2460 }
2461 return velocity;
2462 }
2463
2464
2465
2466
2467 @Override
2468 public double getd2EdL2() {
2469 double d2EdLambda2 = 0.0;
2470 if (!lambdaBondedTerms) {
2471 if (vanderWaalsTerm) {
2472 d2EdLambda2 = vanderWaals.getd2EdL2();
2473 }
2474 if (multipoleTerm) {
2475 d2EdLambda2 += particleMeshEwald.getd2EdL2();
2476 }
2477
2478
2479
2480
2481
2482
2483
2484 if (restrainDistanceTerm && restrainDistancePotentialEnergy != null) {
2485 for (RestrainDistance restrainDistance : restrainDistancePotentialEnergy.getRestrainDistanceArray()) {
2486 d2EdLambda2 += restrainDistance.getd2EdL2();
2487 }
2488 }
2489 if (ncsTerm && ncsRestraint != null) {
2490 d2EdLambda2 += ncsRestraint.getd2EdL2();
2491 }
2492 if (comTerm && comRestraint != null) {
2493 d2EdLambda2 += comRestraint.getd2EdL2();
2494 }
2495 if (lambdaTorsions) {
2496 if (torsionPotentialEnergy != null) {
2497 d2EdLambda2 += torsionPotentialEnergy.getd2EdL2();
2498 }
2499 if (piOrbitalTorsionPotentialEnergy != null) {
2500 d2EdLambda2 += piOrbitalTorsionPotentialEnergy.getd2EdL2();
2501 }
2502 if (torsionTorsionPotentialEnergy != null) {
2503 d2EdLambda2 += torsionTorsionPotentialEnergy.getd2EdL2();
2504 }
2505 }
2506 }
2507 return d2EdLambda2;
2508 }
2509
2510
2511
2512
2513 @Override
2514 public double getdEdL() {
2515 double dEdLambda = 0.0;
2516 if (!lambdaBondedTerms) {
2517 if (vanderWaalsTerm) {
2518 dEdLambda = vanderWaals.getdEdL();
2519 }
2520 if (multipoleTerm) {
2521 dEdLambda += particleMeshEwald.getdEdL();
2522 }
2523 if (restrainDistanceTerm && restrainDistancePotentialEnergy != null) {
2524 for (RestrainDistance restrainDistance : restrainDistancePotentialEnergy.getRestrainDistanceArray()) {
2525 dEdLambda += restrainDistance.getdEdL();
2526 }
2527 }
2528
2529
2530
2531
2532
2533
2534
2535 if (ncsTerm && ncsRestraint != null) {
2536 dEdLambda += ncsRestraint.getdEdL();
2537 }
2538 if (comTerm && comRestraint != null) {
2539 dEdLambda += comRestraint.getdEdL();
2540 }
2541 if (lambdaTorsions) {
2542 if (torsionPotentialEnergy != null) {
2543 dEdLambda += torsionPotentialEnergy.getdEdL();
2544 }
2545 if (piOrbitalTorsionPotentialEnergy != null) {
2546 dEdLambda += piOrbitalTorsionPotentialEnergy.getdEdL();
2547 }
2548 if (torsionTorsionPotentialEnergy != null) {
2549 dEdLambda += torsionTorsionPotentialEnergy.getdEdL();
2550 }
2551 }
2552 }
2553 return dEdLambda;
2554 }
2555
2556
2557
2558
2559 @Override
2560 public void getdEdXdL(double[] gradients) {
2561 if (!lambdaBondedTerms) {
2562 if (vanderWaalsTerm) {
2563 vanderWaals.getdEdXdL(gradients);
2564 }
2565 if (multipoleTerm) {
2566 particleMeshEwald.getdEdXdL(gradients);
2567 }
2568 if (restrainDistanceTerm && restrainDistancePotentialEnergy != null) {
2569 for (RestrainDistance restrainDistance : restrainDistancePotentialEnergy.getRestrainDistanceArray()) {
2570 restrainDistance.getdEdXdL(gradients);
2571 }
2572 }
2573
2574
2575
2576
2577
2578
2579
2580 if (ncsTerm && ncsRestraint != null) {
2581 ncsRestraint.getdEdXdL(gradients);
2582 }
2583 if (comTerm && comRestraint != null) {
2584 comRestraint.getdEdXdL(gradients);
2585 }
2586 if (lambdaTorsions) {
2587 double[] grad = new double[3];
2588 int index = 0;
2589 for (int i = 0; i < nAtoms; i++) {
2590 Atom a = atoms[i];
2591 if (a.isActive()) {
2592 a.getLambdaXYZGradient(grad);
2593 gradients[index++] += grad[0];
2594 gradients[index++] += grad[1];
2595 gradients[index++] += grad[2];
2596 }
2597 }
2598 }
2599 }
2600 }
2601
2602
2603
2604
2605
2606
2607 @Override
2608 public void setAcceleration(double[] acceleration) {
2609 if (acceleration == null) {
2610 return;
2611 }
2612 int index = 0;
2613 double[] accel = new double[3];
2614 for (int i = 0; i < nAtoms; i++) {
2615 if (atoms[i].isActive()) {
2616 accel[0] = acceleration[index++];
2617 accel[1] = acceleration[index++];
2618 accel[2] = acceleration[index++];
2619 atoms[i].setAcceleration(accel);
2620 }
2621 }
2622 }
2623
2624
2625
2626
2627
2628
2629 public void setCoordinates(@Nullable double[] coords) {
2630 if (coords == null) {
2631 return;
2632 }
2633 int index = 0;
2634 for (Atom a : atoms) {
2635 if (a.isActive()) {
2636 double x = coords[index++];
2637 double y = coords[index++];
2638 double z = coords[index++];
2639 a.moveTo(x, y, z);
2640 }
2641 }
2642 }
2643
2644
2645
2646
2647
2648
2649
2650 public void setCrystal(Crystal crystal, boolean checkReplicatesCell) {
2651 if (checkReplicatesCell) {
2652 this.crystal = ReplicatesCrystal.replicatesCrystalFactory(crystal.getUnitCell(), cutOff2);
2653 } else {
2654 this.crystal = crystal;
2655 }
2656
2657
2658
2659
2660 if (vanderWaalsTerm) {
2661 vanderWaals.setCrystal(this.crystal);
2662 }
2663 if (multipoleTerm) {
2664 particleMeshEwald.setCrystal(this.crystal);
2665 }
2666 }
2667
2668
2669
2670
2671
2672
2673
2674 @Override
2675 public void setPreviousAcceleration(double[] previousAcceleration) {
2676 if (previousAcceleration == null) {
2677 return;
2678 }
2679 int index = 0;
2680 double[] prev = new double[3];
2681 for (int i = 0; i < nAtoms; i++) {
2682 if (atoms[i].isActive()) {
2683 prev[0] = previousAcceleration[index++];
2684 prev[1] = previousAcceleration[index++];
2685 prev[2] = previousAcceleration[index++];
2686 atoms[i].setPreviousAcceleration(prev);
2687 }
2688 }
2689 }
2690
2691
2692
2693
2694
2695
2696
2697
2698
2699
2700
2701 public void setPrintOnFailure(boolean onFail, boolean override) {
2702 if (override) {
2703
2704 printOnFailure = onFail;
2705 } else {
2706 try {
2707 molecularAssembly.getForceField().getBoolean("PRINT_ON_FAILURE");
2708
2709
2710
2711
2712
2713 } catch (Exception ex) {
2714 printOnFailure = onFail;
2715 }
2716 }
2717 }
2718
2719
2720
2721
2722
2723
2724 @Override
2725 public void setVelocity(double[] velocity) {
2726 if (velocity == null) {
2727 return;
2728 }
2729 int index = 0;
2730 double[] vel = new double[3];
2731 for (Atom a : atoms) {
2732 if (a.isActive()) {
2733 vel[0] = velocity[index++];
2734 vel[1] = velocity[index++];
2735 vel[2] = velocity[index++];
2736 } else {
2737
2738 vel[0] = 0.0;
2739 vel[1] = 0.0;
2740 vel[2] = 0.0;
2741 }
2742 a.setVelocity(vel);
2743 }
2744 }
2745
2746
2747
2748
2749 @Override
2750 public String toString() {
2751 StringBuilder sb = new StringBuilder();
2752 sb.append(forceFieldBondedEnergyRegion.toString());
2753 sb.append(restrainEnergyRegion.toString());
2754 if (restrainGroupTerm && nRestrainGroups > 0) {
2755 sb.append(format(" %s %20.8f %12d %12.3f\n", "Restrain Groups ", restrainGroupEnergy,
2756 nRestrainGroups, restrainGroupTime * toSeconds));
2757 }
2758 if (ncsTerm) {
2759 sb.append(format(" %s %20.8f %12d %12.3f\n", "NCS Restraint ", ncsEnergy, nAtoms,
2760 ncsTime * toSeconds));
2761 }
2762 if (comTerm) {
2763 sb.append(format(" %s %20.8f %12d %12.3f\n", "COM Restraint ", comRestraintEnergy, nAtoms,
2764 comRestraintTime * toSeconds));
2765 }
2766 if (vanderWaalsTerm && nVanDerWaalInteractions > 0) {
2767 sb.append(format(" %s %20.8f %12d %12.3f\n", "Van der Waals ", vanDerWaalsEnergy,
2768 nVanDerWaalInteractions, vanDerWaalsTime * toSeconds));
2769 }
2770 if (multipoleTerm && nPermanentInteractions > 0) {
2771 if (polarizationTerm) {
2772 sb.append(format(" %s %20.8f %12d\n", "Atomic Multipoles ", permanentMultipoleEnergy,
2773 nPermanentInteractions));
2774 } else {
2775 if (elecForm == ELEC_FORM.FIXED_CHARGE) {
2776 sb.append(format(" %s %20.8f %12d %12.3f\n", "Atomic Charges ", permanentMultipoleEnergy,
2777 nPermanentInteractions, electrostaticTime * toSeconds));
2778 } else
2779 sb.append(format(" %s %20.8f %12d %12.3f\n", "Atomic Multipoles ", permanentMultipoleEnergy,
2780 nPermanentInteractions, electrostaticTime * toSeconds));
2781 }
2782 }
2783 if (polarizationTerm && nPermanentInteractions > 0) {
2784 sb.append(format(" %s %20.8f %12d %12.3f\n", "Polarization ", polarizationEnergy,
2785 nPermanentInteractions, electrostaticTime * toSeconds));
2786 }
2787 if (generalizedKirkwoodTerm && nGKInteractions > 0) {
2788 sb.append(
2789 format(" %s %20.8f %12d\n", "Solvation ", solvationEnergy, nGKInteractions));
2790 }
2791 if (esvTerm) {
2792 sb.append(
2793 format(" %s %20.8f %s\n", "ExtendedSystemBias", esvBias, esvSystem.getLambdaList()));
2794 sb.append(esvSystem.getBiasDecomposition());
2795 }
2796 if (nnTerm) {
2797 sb.append(format(" %s %20.8f %12d %12.3f\n", "Neural Network ", nnEnergy, nAtoms,
2798 nnTime * toSeconds));
2799 }
2800 sb.append(format(" %s %20.8f %s %12.3f (sec)",
2801 "Total Potential ", totalEnergy, "(Kcal/mole)", totalTime * toSeconds));
2802 int nsymm = crystal.getUnitCell().spaceGroup.getNumberOfSymOps();
2803 if (nsymm > 1) {
2804 sb.append(format("\n %s %20.8f", "Unit Cell ", totalEnergy * nsymm));
2805 }
2806 if (crystal.getUnitCell() != crystal) {
2807 nsymm = crystal.spaceGroup.getNumberOfSymOps();
2808 sb.append(format("\n %s %20.8f", "Replicates Cell ", totalEnergy * nsymm));
2809 }
2810
2811 return sb.toString();
2812 }
2813
2814
2815
2816
2817 public void logBondedTermsAndRestraints() {
2818 if (forceFieldBondedEnergyRegion != null) {
2819 forceFieldBondedEnergyRegion.log();
2820 }
2821 if (restrainEnergyRegion != null) {
2822 restrainEnergyRegion.log();
2823 }
2824 if (restrainGroupTerm && nRestrainGroups > 0) {
2825 logger.info("\n Restrain Group Interactions:");
2826 logger.info(restrainGroups.toString());
2827 }
2828 }
2829
2830
2831
2832
2833
2834
2835
2836 void setLambdaBondedTerms(boolean lambdaBondedTerms, boolean lambdaAllBondedTerms) {
2837 this.lambdaBondedTerms = lambdaBondedTerms;
2838 this.lambdaAllBondedTerms = lambdaAllBondedTerms;
2839 }
2840
2841
2842
2843
2844
2845
2846
2847
2848 private double[] fillGradient(double[] g) {
2849 assert (g != null);
2850 double[] grad = new double[3];
2851 int n = getNumberOfVariables();
2852 if (g == null || g.length < n) {
2853 g = new double[n];
2854 }
2855 int index = 0;
2856 for (int i = 0; i < nAtoms; i++) {
2857 Atom a = atoms[i];
2858 if (a.isActive()) {
2859 a.getXYZGradient(grad);
2860 double gx = grad[0];
2861 double gy = grad[1];
2862 double gz = grad[2];
2863 if (isNaN(gx) || isInfinite(gx) || isNaN(gy) || isInfinite(gy) || isNaN(gz) || isInfinite(
2864 gz)) {
2865 StringBuilder sb = new StringBuilder(
2866 format("The gradient of atom %s is (%8.3f,%8.3f,%8.3f).", a, gx, gy, gz));
2867 double[] vals = new double[3];
2868 a.getVelocity(vals);
2869 sb.append(format("\n Velocities: %8.3g %8.3g %8.3g", vals[0], vals[1], vals[2]));
2870 a.getAcceleration(vals);
2871 sb.append(format("\n Accelerations: %8.3g %8.3g %8.3g", vals[0], vals[1], vals[2]));
2872 a.getPreviousAcceleration(vals);
2873 sb.append(
2874 format("\n Previous accelerations: %8.3g %8.3g %8.3g", vals[0], vals[1], vals[2]));
2875
2876 throw new EnergyException(sb.toString());
2877 }
2878 g[index++] = gx;
2879 g[index++] = gy;
2880 g[index++] = gz;
2881 }
2882 }
2883 return g;
2884 }
2885
2886 private Crystal configureNCS(ForceField forceField, Crystal unitCell) {
2887
2888
2889 if (forceField.getProperties().containsKey("MTRIX1")
2890 && forceField.getProperties().containsKey("MTRIX2") && forceField.getProperties()
2891 .containsKey("MTRIX3")) {
2892 Crystal unitCell2 = new Crystal(unitCell.a, unitCell.b, unitCell.c, unitCell.alpha,
2893 unitCell.beta, unitCell.gamma, unitCell.spaceGroup.pdbName);
2894 SpaceGroup spaceGroup = unitCell2.spaceGroup;
2895
2896 CompositeConfiguration properties = forceField.getProperties();
2897 String[] MTRX1List = properties.getStringArray("MTRIX1");
2898 String[] MTRX2List = properties.getStringArray("MTRIX2");
2899 String[] MTRX3List = properties.getStringArray("MTRIX3");
2900 spaceGroup.symOps.clear();
2901 double number1;
2902 double number2;
2903 double number3;
2904 for (int i = 0; i < MTRX1List.length; i++) {
2905 double[][] Rot_MTRX = {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}};
2906 double[] Tr_MTRX = {0, 0, 0};
2907 String[] tokens1 = MTRX1List[i].trim().split(" +");
2908 String[] tokens2 = MTRX2List[i].trim().split(" +");
2909 String[] tokens3 = MTRX3List[i].trim().split(" +");
2910 for (int k = 0; k < 4; k++) {
2911 number1 = Double.parseDouble(tokens1[k]);
2912 number2 = Double.parseDouble(tokens2[k]);
2913 number3 = Double.parseDouble(tokens3[k]);
2914 if (k != 3) {
2915 Rot_MTRX[0][k] = number1;
2916 Rot_MTRX[1][k] = number2;
2917 Rot_MTRX[2][k] = number3;
2918 } else {
2919 Tr_MTRX[0] = number1;
2920 Tr_MTRX[1] = number2;
2921 Tr_MTRX[2] = number3;
2922 }
2923 }
2924 SymOp symOp = new SymOp(Rot_MTRX, Tr_MTRX);
2925 if (logger.isLoggable(Level.FINEST)) {
2926 logger.info(
2927 format(" MTRIXn SymOp: %d of %d\n" + symOp, i + 1, MTRX1List.length));
2928 }
2929 spaceGroup.symOps.add(symOp);
2930 }
2931 unitCell = NCSCrystal.NCSCrystalFactory(unitCell, spaceGroup.symOps);
2932 unitCell.updateCrystal();
2933 }
2934 return unitCell;
2935 }
2936
2937 public RestrainMode getRestrainMode() {
2938 return restrainMode;
2939 }
2940
2941
2942
2943
2944
2945
2946 public RestrainGroups getRestrainGroups() {
2947 return restrainGroups;
2948 }
2949
2950
2951
2952
2953
2954
2955 public TorsionTorsionPotentialEnergy getTorsionTorsionPotentialEnergy() {
2956 return torsionTorsionPotentialEnergy;
2957 }
2958
2959
2960
2961
2962
2963
2964 public AnglePotentialEnergy getAnglePotentialEnergy() {
2965 return anglePotentialEnergy;
2966 }
2967
2968
2969
2970
2971
2972
2973 public BondPotentialEnergy getBondPotentialEnergy() {
2974 return bondPotentialEnergy;
2975 }
2976
2977
2978
2979
2980
2981
2982 public StretchBendPotentialEnergy getStretchBendPotentialEnergy() {
2983 return stretchBendPotentialEnergy;
2984 }
2985
2986
2987
2988
2989
2990
2991 public UreyBradleyPotentialEnergy getUreyBradleyPotentialEnergy() {
2992 return ureyBradleyPotentialEnergy;
2993 }
2994
2995
2996
2997
2998
2999
3000 public OutOfPlaneBendPotentialEnergy getOutOfPlaneBendPotentialEnergy() {
3001 return outOfPlaneBendPotentialEnergy;
3002 }
3003
3004
3005
3006
3007
3008
3009 public TorsionPotentialEnergy getTorsionPotentialEnergy() {
3010 return torsionPotentialEnergy;
3011 }
3012
3013
3014
3015
3016
3017
3018 public StretchTorsionPotentialEnergy getStretchTorsionPotentialEnergy() {
3019 return stretchTorsionPotentialEnergy;
3020 }
3021
3022
3023
3024
3025
3026
3027 public AngleTorsionPotentialEnergy getAngleTorsionPotentialEnergy() {
3028 return angleTorsionPotentialEnergy;
3029 }
3030
3031
3032
3033
3034
3035
3036 public ImproperTorsionPotentialEnergy getImproperTorsionPotentialEnergy() {
3037 return improperTorsionPotentialEnergy;
3038 }
3039
3040
3041
3042
3043
3044
3045 public PiOrbitalTorsionPotentialEnergy getPiOrbitalTorsionPotentialEnergy() {
3046 return piOrbitalTorsionPotentialEnergy;
3047 }
3048
3049
3050
3051
3052
3053
3054 public RestrainPositionPotentialEnergy getRestrainPositionPotentialEnergy() {
3055 return restrainPositionPotentialEnergy;
3056 }
3057
3058
3059
3060
3061
3062
3063 public RestrainDistancePotentialEnergy getRestrainDistancePotentialEnergy() {
3064 return restrainDistancePotentialEnergy;
3065 }
3066
3067
3068
3069
3070
3071
3072 public RestrainTorsionPotentialEnergy getRestrainTorsionPotentialEnergy() {
3073 return restrainTorsionPotentialEnergy;
3074 }
3075
3076 }