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.nonbonded;
39
40 import edu.rit.pj.IntegerSchedule;
41 import edu.rit.pj.ParallelTeam;
42 import edu.rit.pj.reduction.SharedDouble;
43 import edu.rit.util.Range;
44 import ffx.crystal.Crystal;
45 import ffx.numerics.atomic.AtomicDoubleArray.AtomicDoubleArrayImpl;
46 import ffx.numerics.atomic.AtomicDoubleArray3D;
47 import ffx.numerics.multipole.MultipoleUtilities;
48 import ffx.potential.Platform;
49 import ffx.potential.bonded.Atom;
50 import ffx.potential.bonded.Bond;
51 import ffx.potential.bonded.LambdaInterface;
52 import ffx.potential.extended.ExtendedSystem;
53 import ffx.potential.nonbonded.pme.AlchemicalParameters;
54 import ffx.potential.nonbonded.pme.DirectRegion;
55 import ffx.potential.nonbonded.pme.EwaldParameters;
56 import ffx.potential.nonbonded.pme.ExpandInducedDipolesRegion;
57 import ffx.potential.nonbonded.pme.InducedDipoleFieldReduceRegion;
58 import ffx.potential.nonbonded.pme.InducedDipoleFieldRegion;
59 import ffx.potential.nonbonded.pme.InitializationRegion;
60 import ffx.potential.nonbonded.pme.LambdaMode;
61 import ffx.potential.nonbonded.pme.OPTRegion;
62 import ffx.potential.nonbonded.pme.PCGSolver;
63 import ffx.potential.nonbonded.pme.PMETimings;
64 import ffx.potential.nonbonded.pme.PermanentFieldRegion;
65 import ffx.potential.nonbonded.pme.Polarization;
66 import ffx.potential.nonbonded.pme.RealSpaceEnergyRegion;
67 import ffx.potential.nonbonded.pme.RealSpaceNeighborParameters;
68 import ffx.potential.nonbonded.pme.ReciprocalEnergyRegion;
69 import ffx.potential.nonbonded.pme.ReduceRegion;
70 import ffx.potential.nonbonded.pme.SCFAlgorithm;
71 import ffx.potential.nonbonded.pme.SCFPredictor;
72 import ffx.potential.nonbonded.pme.SCFPredictorParameters;
73 import ffx.potential.nonbonded.pme.SORRegion;
74 import ffx.potential.nonbonded.pme.ScaleParameters;
75 import ffx.potential.parameters.AtomType;
76 import ffx.potential.parameters.ForceField;
77 import ffx.potential.parameters.ForceField.ELEC_FORM;
78 import ffx.potential.parameters.ForceField.ForceFieldType;
79 import ffx.potential.parameters.MultipoleType;
80 import ffx.potential.parameters.PolarizeType;
81 import ffx.potential.utils.EnergyException;
82 import ffx.utilities.Constants;
83 import ffx.utilities.FFXProperty;
84 import org.apache.commons.math3.linear.Array2DRowRealMatrix;
85 import org.apache.commons.math3.linear.EigenDecomposition;
86
87 import java.util.ArrayList;
88 import java.util.Arrays;
89 import java.util.List;
90 import java.util.logging.Level;
91 import java.util.logging.Logger;
92
93 import static ffx.potential.nonbonded.pme.EwaldParameters.DEFAULT_EWALD_COEFFICIENT;
94 import static ffx.potential.nonbonded.pme.SCFAlgorithm.SOR;
95 import static ffx.potential.parameters.ForceField.ELEC_FORM.PAM;
96 import static ffx.potential.parameters.ForceField.toEnumForm;
97 import static ffx.potential.parameters.MultipoleType.MultipoleFrameDefinition;
98 import static ffx.potential.parameters.MultipoleType.assignMultipole;
99 import static ffx.potential.parameters.MultipoleType.t000;
100 import static ffx.potential.parameters.MultipoleType.t001;
101 import static ffx.potential.parameters.MultipoleType.t002;
102 import static ffx.potential.parameters.MultipoleType.t010;
103 import static ffx.potential.parameters.MultipoleType.t011;
104 import static ffx.potential.parameters.MultipoleType.t020;
105 import static ffx.potential.parameters.MultipoleType.t100;
106 import static ffx.potential.parameters.MultipoleType.t101;
107 import static ffx.potential.parameters.MultipoleType.t110;
108 import static ffx.potential.parameters.MultipoleType.t200;
109 import static ffx.utilities.Constants.ELEC_ANG_TO_DEBYE;
110 import static ffx.utilities.Constants.NS2SEC;
111 import static ffx.utilities.PropertyGroup.ElectrostaticsFunctionalForm;
112 import static ffx.utilities.PropertyGroup.LocalGeometryFunctionalForm;
113 import static java.lang.String.format;
114 import static java.lang.System.arraycopy;
115 import static java.util.Arrays.fill;
116 import static java.util.Collections.sort;
117 import static org.apache.commons.math3.util.FastMath.max;
118 import static org.apache.commons.math3.util.FastMath.pow;
119 import static org.apache.commons.math3.util.FastMath.sqrt;
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149 @SuppressWarnings("deprecation")
150 public class ParticleMeshEwald implements LambdaInterface {
151
152 private static final Logger logger = Logger.getLogger(ParticleMeshEwald.class.getName());
153
154
155
156 private static final int tensorCount = MultipoleUtilities.tensorCount(3);
157
158
159
160
161 private final boolean lambdaTerm;
162
163
164
165 private final boolean nnTerm;
166 private final boolean reciprocalSpaceTerm;
167
168
169
170 private final ForceField forceField;
171 private static final double DEFAULT_POLAR_EPS = 1.0e-6;
172 @FFXProperty(name = "polar-eps", propertyGroup = ElectrostaticsFunctionalForm, defaultValue = "1.0e-6",
173 description = """
174 Sets the convergence criterion applied during computation of self-consistent induced dipoles.
175 The calculation is deemed to have converged when the rms change in Debyes in
176 the induced dipole at all polarizable sites is less than the value specified with this property.
177 The default value in the absence of the keyword is 1.0e-6 Debyes.
178 """)
179 private final double poleps;
180
181
182
183
184 private final boolean sorFallback;
185
186
187
188 private final double sorFallBackPolarSOR = 0.65;
189
190
191
192 private final boolean directFallback;
193
194
195
196
197 private final SCFPredictorParameters scfPredictorParameters;
198 private final EwaldParameters ewaldParameters;
199 private final ScaleParameters scaleParameters;
200 private final AlchemicalParameters alchemicalParameters;
201 private final RealSpaceNeighborParameters realSpaceNeighborParameters;
202 private final PMETimings pmeTimings;
203
204
205
206 private final int maxThreads;
207
208
209
210
211 private final ParallelTeam parallelTeam;
212
213
214
215
216
217
218
219
220 private final ParallelTeam sectionTeam;
221
222
223
224
225
226
227
228 private final ParallelTeam realSpaceTeam;
229
230
231
232
233
234
235
236 private final ParallelTeam fftTeam;
237 private final NeighborList neighborList;
238 private final InitializationRegion initializationRegion;
239 private final PermanentFieldRegion permanentFieldRegion;
240 private final InducedDipoleFieldRegion inducedDipoleFieldRegion;
241 private final InducedDipoleFieldReduceRegion inducedDipoleFieldReduceRegion;
242 private final ExpandInducedDipolesRegion expandInducedDipolesRegion;
243 private final DirectRegion directRegion;
244 private final SORRegion sorRegion;
245 private final OPTRegion optRegion;
246 private final PCGSolver pcgSolver;
247 private final ReciprocalSpace reciprocalSpace;
248 private final ReciprocalEnergyRegion reciprocalEnergyRegion;
249 private final RealSpaceEnergyRegion realSpaceEnergyRegion;
250 private final ReduceRegion reduceRegion;
251 private final GeneralizedKirkwood generalizedKirkwood;
252
253
254
255
256 public Polarization polarization;
257
258
259
260 public double[][][] coordinates;
261
262
263
264 public int[][][] neighborLists;
265
266
267
268 public double[][][] globalMultipole;
269
270
271
272 public double[][][] fractionalMultipole;
273
274
275
276 public double[][][] inducedDipole;
277 public double[][][] inducedDipoleCR;
278
279
280
281 public double[][] directDipole;
282 public double[][] directDipoleCR;
283 public double[][] directField;
284 public double[][] directFieldCR;
285
286
287
288 public double[][][] vacuumInducedDipole;
289 public double[][][] vacuumInducedDipoleCR;
290
291
292
293 public double[][] vacuumDirectDipole;
294 public double[][] vacuumDirectDipoleCR;
295
296
297
298 @FFXProperty(name = "electric", propertyGroup = LocalGeometryFunctionalForm, defaultValue = "332.063713",
299 description = """
300 Specifies a value for the so-called "electric constant" allowing conversion unit of electrostatic
301 potential energy values from electrons^2/Angstrom to kcal/mol. Internally, FFX stores a default value
302 for this constant as 332.063713 based on CODATA reference values. Since different force fields are
303 intended for use with slightly different values, this keyword allows overriding the default value.
304 """)
305 public double electric;
306
307
308
309 public double soluteDielectric;
310
311
312
313 protected Atom[] atoms;
314
315
316
317 protected int nAtoms;
318
319
320
321 protected int[][] ip11;
322 protected int[][] ip12;
323 protected int[][] ip13;
324
325
326
327
328 protected double totalMultipoleEnergy;
329
330
331
332
333 protected double permanentMultipoleEnergy;
334 protected double permanentRealSpaceEnergy;
335 protected double permanentSelfEnergy;
336 protected double permanentReciprocalEnergy;
337 protected double permanentChargeCorrectionEnergy;
338
339
340
341 protected double polarizationEnergy;
342 protected double inducedRealSpaceEnergy;
343 protected double inducedSelfEnergy;
344 protected double inducedReciprocalEnergy;
345 protected SCFAlgorithm scfAlgorithm;
346
347
348
349 private final SharedDouble shareddEdLambda;
350
351
352
353 private final SharedDouble sharedd2EdLambda2;
354
355
356
357 private final ELEC_FORM elecForm;
358
359
360
361 private Crystal crystal;
362
363
364
365 private int nSymm;
366
367
368
369 private boolean generalizedKirkwoodTerm;
370
371
372
373 private boolean gradient = false;
374
375
376
377 private int interactions;
378
379
380
381 private int gkInteractions;
382
383
384
385 private double solvationEnergy;
386
387
388
389 private LambdaMode lambdaMode = LambdaMode.OFF;
390
391
392
393 private double lambda = 1.0;
394
395
396
397 private double[][] localMultipole;
398 private MultipoleFrameDefinition[] frame;
399 private int[][] axisAtom;
400 private double[] ipdamp;
401 private double[] thole;
402 private double[] polarizability;
403
404
405
406 private int[][] mask12;
407 private int[][] mask13;
408 private int[][] mask14;
409 private int[][] mask15;
410
411
412
413 private boolean[] isSoft;
414
415
416
417 private int[] molecule;
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436 private boolean[] use;
437 private IntegerSchedule permanentSchedule;
438 private double[][] cartesianMultipolePhi;
439 private double[][] fracMultipolePhi;
440 private double[][] cartesianInducedDipolePhi;
441 private double[][] cartesianInducedDipolePhiCR;
442 private double[][] fractionalInducedDipolePhi;
443 private double[][] fractionalInducedDipolePhiCR;
444 private double[][] cartesianVacuumDipolePhi;
445 private double[][] cartesianVacuumDipolePhiCR;
446 private double[][] fractionalVacuumDipolePhi;
447 private double[][] fractionalVacuumDipolePhiCR;
448
449
450
451 private AtomicDoubleArrayImpl atomicDoubleArrayImpl;
452
453
454
455 private AtomicDoubleArray3D field;
456
457
458
459 private AtomicDoubleArray3D fieldCR;
460
461
462
463 private AtomicDoubleArray3D grad;
464
465
466
467 private AtomicDoubleArray3D torque;
468
469
470
471 private AtomicDoubleArray3D lambdaGrad;
472
473
474
475 private AtomicDoubleArray3D lambdaTorque;
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490 public ParticleMeshEwald(
491 Atom[] atoms,
492 int[] molecule,
493 ForceField forceField,
494 Crystal crystal,
495 NeighborList neighborList,
496 ELEC_FORM elecForm,
497 double ewaldCutoff,
498 double gkCutoff,
499 ParallelTeam parallelTeam) {
500 this.atoms = atoms;
501 this.molecule = molecule;
502 this.forceField = forceField;
503 this.crystal = crystal;
504 this.parallelTeam = parallelTeam;
505 this.neighborList = neighborList;
506 this.elecForm = elecForm;
507 neighborLists = neighborList.getNeighborList();
508 permanentSchedule = neighborList.getPairwiseSchedule();
509 nAtoms = atoms.length;
510 nSymm = crystal.spaceGroup.getNumberOfSymOps();
511 maxThreads = parallelTeam.getThreadCount();
512
513 electric = forceField.getDouble("ELECTRIC", Constants.DEFAULT_ELECTRIC);
514
515
516 soluteDielectric = forceField.getDouble("SOLUTE_DIELECTRIC", 1.0);
517 if (soluteDielectric > 1.0) {
518 electric = electric / soluteDielectric;
519 } else {
520 soluteDielectric = 1.0;
521 }
522
523 poleps = forceField.getDouble("POLAR_EPS", DEFAULT_POLAR_EPS);
524
525
526 lambdaTerm = forceField.getBoolean("ELEC_LAMBDATERM",
527 forceField.getBoolean("LAMBDATERM", false));
528
529 nnTerm = forceField.getBoolean("NNTERM", false);
530
531 if (nnTerm && lambdaTerm) {
532 logger.severe(" Use a neural network potential with alchemical simulations is not yet supported.");
533 }
534
535 double aewald = forceField.getDouble("EWALD_ALPHA", DEFAULT_EWALD_COEFFICIENT);
536 ewaldParameters = new EwaldParameters(ewaldCutoff, aewald);
537 scaleParameters = new ScaleParameters(elecForm, forceField);
538 reciprocalSpaceTerm = forceField.getBoolean("RECIPTERM", true);
539
540 SCFPredictor scfPredictor;
541 try {
542 String predictor = forceField.getString("SCF_PREDICTOR", "NONE");
543 predictor = predictor.replaceAll("-", "_").toUpperCase();
544 scfPredictor = SCFPredictor.valueOf(predictor);
545 } catch (Exception e) {
546 scfPredictor = SCFPredictor.NONE;
547 }
548
549 scfPredictorParameters = new SCFPredictorParameters(scfPredictor, nAtoms);
550 if (scfPredictor != SCFPredictor.NONE) {
551 scfPredictorParameters.init(forceField);
552 }
553
554 String algorithm = forceField.getString("SCF_ALGORITHM", "CG");
555 try {
556 algorithm = algorithm.replaceAll("-", "_").toUpperCase();
557 scfAlgorithm = SCFAlgorithm.valueOf(algorithm);
558 } catch (Exception e) {
559 scfAlgorithm = SCFAlgorithm.CG;
560 }
561
562
563 sorFallback = forceField.getBoolean("SOR_SCF_FALLBACK", true);
564
565 directFallback = forceField.getBoolean("DIRECT_SCF_FALLBACK", !sorFallback);
566
567
568 String value = forceField.getString("ARRAY_REDUCTION", "MULTI");
569 try {
570 atomicDoubleArrayImpl = AtomicDoubleArrayImpl.valueOf(toEnumForm(value));
571 } catch (Exception e) {
572 atomicDoubleArrayImpl = AtomicDoubleArrayImpl.MULTI;
573 logger.info(
574 format(
575 " Unrecognized ARRAY-REDUCTION %s; defaulting to %s", value, atomicDoubleArrayImpl));
576 }
577 logger.fine(format(" PME using %s arrays.", atomicDoubleArrayImpl.toString()));
578
579 if (!scfAlgorithm.isSupported(Platform.FFX)) {
580
581
582 logger.fine(format(
583 " SCF algorithm %s is not supported by FFX reference implementation; falling back to CG!",
584 scfAlgorithm));
585 scfAlgorithm = SCFAlgorithm.CG;
586 }
587
588 pcgSolver = new PCGSolver(maxThreads, poleps, forceField, nAtoms);
589
590 alchemicalParameters = new AlchemicalParameters(forceField, lambdaTerm, nnTerm, polarization);
591 if (nnTerm) {
592 initNN();
593 }
594
595 String polar = forceField.getString("POLARIZATION", "MUTUAL");
596 if (elecForm == ELEC_FORM.FIXED_CHARGE) {
597 polar = "NONE";
598 }
599 boolean polarizationTerm = forceField.getBoolean("POLARIZETERM", true);
600 if (!polarizationTerm || polar.equalsIgnoreCase("NONE")) {
601 polarization = Polarization.NONE;
602 } else if (polar.equalsIgnoreCase("DIRECT")) {
603 polarization = Polarization.DIRECT;
604 } else {
605 polarization = Polarization.MUTUAL;
606 }
607
608 if (lambdaTerm) {
609 shareddEdLambda = new SharedDouble();
610 sharedd2EdLambda2 = new SharedDouble();
611 } else {
612 shareddEdLambda = null;
613 sharedd2EdLambda2 = null;
614 lambdaGrad = null;
615 lambdaTorque = null;
616 }
617
618 directRegion = new DirectRegion(maxThreads);
619 sorRegion = new SORRegion(maxThreads, forceField);
620 optRegion = new OPTRegion(maxThreads);
621
622 if (logger.isLoggable(Level.INFO)) {
623 StringBuilder sb = new StringBuilder();
624 sb.append("\n Electrostatics\n");
625 sb.append(format(" Polarization: %8s\n", polarization.toString()));
626 if (polarization == Polarization.MUTUAL) {
627 sb.append(format(" SCF Convergence Criteria: %8.3e\n", poleps));
628 sb.append(format(" SCF Predictor: %8s\n",
629 scfPredictorParameters.scfPredictor));
630 sb.append(format(" SCF Algorithm: %8s\n", scfAlgorithm));
631 if (scfAlgorithm == SOR) {
632 sb.append(format(" SOR Parameter: %8.3f\n", sorRegion.getSOR()));
633 } else {
634 sb.append(format(" CG Preconditioner Cut-Off: %8.3f\n",
635 pcgSolver.getPreconditionerCutoff()));
636 sb.append(format(" CG Preconditioner Ewald Coeff.: %8.3f\n",
637 pcgSolver.getPreconditionerEwald()));
638 sb.append(format(" CG Preconditioner Scale: %8.3f\n",
639 pcgSolver.getPreconditionerScale()));
640 sb.append(
641 format(" CG Preconditioner Mode: %15s\n", pcgSolver.getPreconditionerMode()));
642 }
643 }
644 if (ewaldParameters.aewald > 0.0) {
645 sb.append(" Particle-mesh Ewald\n");
646 sb.append(format(" Ewald Coefficient: %8.3f\n", ewaldParameters.aewald));
647 sb.append(format(" Particle Cut-Off: %8.3f (A)", ewaldParameters.off));
648 } else if (ewaldParameters.off != Double.POSITIVE_INFINITY) {
649 sb.append(format(" Electrostatics Cut-Off: %8.3f (A)\n",
650 ewaldParameters.off));
651 } else {
652 sb.append(" Electrostatics Cut-Off: NONE\n");
653 }
654
655 logger.info(sb.toString());
656 }
657
658
659 int sectionThreads;
660
661
662
663
664
665 int realSpaceThreads;
666
667
668
669
670
671 int reciprocalThreads;
672
673 boolean concurrent;
674 int realThreads = 1;
675 try {
676 realThreads = forceField.getInteger("PME_REAL_THREADS");
677 if (realThreads >= maxThreads || realThreads < 1) {
678 throw new Exception("pme-real-threads must be < ffx.nt and greater than 0");
679 }
680 concurrent = true;
681 } catch (Exception e) {
682 concurrent = false;
683 }
684 if (concurrent) {
685 sectionThreads = 2;
686 realSpaceThreads = realThreads;
687 reciprocalThreads = maxThreads - realThreads;
688 sectionTeam = new ParallelTeam(sectionThreads);
689 realSpaceTeam = new ParallelTeam(realSpaceThreads);
690 fftTeam = new ParallelTeam(reciprocalThreads);
691 } else {
692
693 sectionThreads = 1;
694 sectionTeam = new ParallelTeam(sectionThreads);
695 realSpaceTeam = parallelTeam;
696 fftTeam = parallelTeam;
697 }
698
699 realSpaceNeighborParameters = new RealSpaceNeighborParameters(maxThreads);
700 initializationRegion = new InitializationRegion(this, maxThreads, forceField);
701 expandInducedDipolesRegion = new ExpandInducedDipolesRegion(maxThreads);
702 initAtomArrays();
703
704
705
706
707
708
709 if (ewaldParameters.aewald > 0.0 && reciprocalSpaceTerm) {
710 reciprocalSpace = new ReciprocalSpace(this,
711 crystal.getUnitCell(), forceField, atoms,
712 ewaldParameters.aewald, fftTeam, parallelTeam);
713 reciprocalEnergyRegion = new ReciprocalEnergyRegion(maxThreads, ewaldParameters.aewald, electric);
714 } else {
715 reciprocalSpace = null;
716 reciprocalEnergyRegion = null;
717 }
718 permanentFieldRegion = new PermanentFieldRegion(realSpaceTeam, forceField, lambdaTerm);
719 inducedDipoleFieldRegion = new InducedDipoleFieldRegion(realSpaceTeam, forceField, lambdaTerm);
720 inducedDipoleFieldReduceRegion = new InducedDipoleFieldReduceRegion(maxThreads);
721 realSpaceEnergyRegion = new RealSpaceEnergyRegion(elecForm, maxThreads, forceField, lambdaTerm, electric);
722 reduceRegion = new ReduceRegion(maxThreads, forceField);
723
724 pmeTimings = new PMETimings(maxThreads);
725
726 if (lambdaTerm) {
727 logger.info(alchemicalParameters.toString());
728 }
729
730
731
732 generalizedKirkwoodTerm = forceField.getBoolean("GKTERM", false);
733 if (generalizedKirkwoodTerm || alchemicalParameters.doLigandGKElec) {
734 generalizedKirkwood = new GeneralizedKirkwood(forceField, atoms, this, crystal, parallelTeam, gkCutoff);
735 } else {
736 generalizedKirkwood = null;
737 }
738 }
739
740
741
742
743 public ELEC_FORM getElecForm() {
744 return elecForm;
745 }
746
747
748
749
750
751
752 public boolean getLambdaTerm() {
753 return lambdaTerm;
754 }
755
756
757
758
759 public AlchemicalParameters getAlchemicalParameters() {
760 return alchemicalParameters;
761 }
762
763 public void computeInduceDipoleField() {
764 expandInducedDipoles();
765
766 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
767 reciprocalSpace.splineInducedDipoles(inducedDipole, inducedDipoleCR, use);
768 }
769 field.reset(parallelTeam);
770 fieldCR.reset(parallelTeam);
771 inducedDipoleFieldRegion.init(
772 atoms, crystal, use, molecule, ipdamp, thole, coordinates, realSpaceNeighborParameters,
773 inducedDipole, inducedDipoleCR, reciprocalSpaceTerm, reciprocalSpace, lambdaMode,
774 ewaldParameters, field, fieldCR, pmeTimings);
775 inducedDipoleFieldRegion.executeWith(sectionTeam);
776
777 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
778 reciprocalSpace.computeInducedPhi(
779 cartesianInducedDipolePhi, cartesianInducedDipolePhiCR,
780 fractionalInducedDipolePhi, fractionalInducedDipolePhiCR);
781 }
782
783 if (generalizedKirkwoodTerm) {
784 pmeTimings.gkEnergyTotal = -System.nanoTime();
785 generalizedKirkwood.computeInducedGKField();
786 pmeTimings.gkEnergyTotal += System.nanoTime();
787 logger.fine(
788 format(" Computed GK induced field %8.3f (sec)", pmeTimings.gkEnergyTotal * 1.0e-9));
789 }
790
791 inducedDipoleFieldReduceRegion.init(
792 atoms,
793 inducedDipole,
794 inducedDipoleCR,
795 generalizedKirkwoodTerm,
796 generalizedKirkwood,
797 ewaldParameters,
798 soluteDielectric,
799 cartesianInducedDipolePhi,
800 cartesianInducedDipolePhiCR,
801 field,
802 fieldCR);
803 inducedDipoleFieldReduceRegion.executeWith(parallelTeam);
804 }
805
806 public void destroy() {
807 if (fftTeam != null) {
808 try {
809 fftTeam.shutdown();
810 } catch (Exception ex) {
811 logger.warning(" Exception in shutting down fftTeam");
812 }
813 }
814 if (sectionTeam != null) {
815 try {
816 sectionTeam.shutdown();
817 } catch (Exception ex) {
818 logger.warning(" Exception in shutting down sectionTeam");
819 }
820 }
821 if (realSpaceTeam != null) {
822 try {
823 realSpaceTeam.shutdown();
824 } catch (Exception ex) {
825 logger.warning(" Exception in shutting down realSpaceTeam");
826 }
827 }
828 }
829
830
831
832
833
834
835
836
837 public double energy(boolean gradient, boolean print) {
838
839 this.gradient = gradient;
840
841
842 totalMultipoleEnergy = 0.0;
843 permanentMultipoleEnergy = 0.0;
844 permanentRealSpaceEnergy = 0.0;
845 permanentSelfEnergy = 0.0;
846 permanentReciprocalEnergy = 0.0;
847 permanentChargeCorrectionEnergy = 0.0;
848 polarizationEnergy = 0.0;
849 inducedRealSpaceEnergy = 0.0;
850 inducedSelfEnergy = 0.0;
851 inducedReciprocalEnergy = 0.0;
852 solvationEnergy = 0.0;
853
854
855 interactions = 0;
856 gkInteractions = 0;
857
858
859 pmeTimings.init();
860 permanentFieldRegion.initTimings();
861
862 if (reciprocalSpace != null) {
863 reciprocalSpace.initTimings();
864 }
865
866
867 if (lambdaTerm) {
868 shareddEdLambda.set(0.0);
869 sharedd2EdLambda2.set(0.0);
870 }
871
872 if (esvTerm) {
873 extendedSystem.initEsvPermElec();
874 extendedSystem.initEsvIndElec();
875 }
876
877 alchemicalParameters.doPermanentRealSpace = true;
878 alchemicalParameters.permanentScale = 1.0;
879 alchemicalParameters.doPolarization = true;
880 alchemicalParameters.polarizationScale = 1.0;
881
882
883 initializationRegion.init(lambdaTerm, alchemicalParameters, extendedSystem,
884 atoms, coordinates, crystal, frame, axisAtom, globalMultipole,
885 dMultipoledTirationESV, dMultipoledTautomerESV, polarizability, dPolardTitrationESV, dPolardTautomerESV,
886 thole, ipdamp, use, neighborLists, realSpaceNeighborParameters.realSpaceLists,
887 grad, torque, lambdaGrad, lambdaTorque);
888 initializationRegion.executeWith(parallelTeam);
889
890
891 if (generalizedKirkwoodTerm || alchemicalParameters.doLigandGKElec) {
892 generalizedKirkwood.init();
893 }
894
895
896 double energy;
897 if (!lambdaTerm) {
898 lambdaMode = LambdaMode.OFF;
899 energy = computeEnergy(print);
900 if (nnTerm) {
901 lambdaMode = LambdaMode.VAPOR;
902 double temp = energy;
903 energy = nnElec();
904 logger.fine(format(" Vacuum energy: %20.8f", energy - temp));
905 }
906 } else {
907
908 lambdaMode = LambdaMode.CONDENSED;
909 energy = condensedEnergy();
910 if (logger.isLoggable(Level.FINE)) {
911 logger.fine(format(" Condensed energy: %20.8f", energy));
912 }
913
914
915 if (alchemicalParameters.mode == AlchemicalParameters.AlchemicalMode.OST) {
916
917 lambdaMode = LambdaMode.CONDENSED_NO_LIGAND;
918 double temp = energy;
919 energy = condensedNoLigandSCF();
920 if (logger.isLoggable(Level.FINE)) {
921 logger.fine(format(" Condensed no ligand energy: %20.8f", energy - temp));
922 }
923
924
925 if (alchemicalParameters.doLigandVaporElec) {
926 lambdaMode = LambdaMode.VAPOR;
927 temp = energy;
928 energy = ligandElec();
929 if (logger.isLoggable(Level.FINE)) {
930 logger.fine(format(" Vacuum energy: %20.8f", energy - temp));
931 }
932 }
933 }
934 }
935
936
937
938
939
940 if (gradient || lambdaTerm) {
941 reduceRegion.init(lambdaTerm, gradient, atoms, coordinates, frame, axisAtom, grad, torque, lambdaGrad, lambdaTorque);
942 reduceRegion.executeWith(parallelTeam);
943 }
944
945
946 if (logger.isLoggable(Level.FINE)) {
947 pmeTimings.printRealSpaceTimings(maxThreads, permanentFieldRegion, realSpaceEnergyRegion);
948 if (ewaldParameters.aewald > 0.0 && reciprocalSpaceTerm) {
949 reciprocalSpace.printTimings();
950 }
951 }
952
953 return permanentMultipoleEnergy + polarizationEnergy;
954 }
955
956 public void expandInducedDipoles() {
957 if (nSymm > 1) {
958 expandInducedDipolesRegion.init(atoms, crystal, inducedDipole, inducedDipoleCR);
959 expandInducedDipolesRegion.executeWith(parallelTeam);
960 }
961 }
962
963 public int[][] getAxisAtoms() {
964 return axisAtom;
965 }
966
967 public double getCavitationEnergy() {
968 return generalizedKirkwood.getCavitationEnergy();
969 }
970
971 public double[][][] getCoordinates() {
972 return coordinates;
973 }
974
975 public double getDispersionEnergy() {
976 return generalizedKirkwood.getDispersionEnergy();
977 }
978
979 public double getEwaldCoefficient() {
980 return ewaldParameters.aewald;
981 }
982
983 public double getEwaldCutoff() {
984 return ewaldParameters.off;
985 }
986
987 public GeneralizedKirkwood getGK() {
988 return generalizedKirkwood;
989 }
990
991
992
993
994
995
996 public double getGKEnergy() {
997 return generalizedKirkwood.getGeneralizedKirkwoordEnergy();
998 }
999
1000
1001
1002
1003
1004
1005 public int getGKInteractions() {
1006 return gkInteractions;
1007 }
1008
1009
1010
1011
1012
1013
1014 public int getInteractions() {
1015 return interactions;
1016 }
1017
1018
1019
1020
1021
1022
1023 @Override
1024 public double getLambda() {
1025 return lambda;
1026 }
1027
1028
1029
1030
1031
1032
1033 @Override
1034 public void setLambda(double lambda) {
1035 assert (lambda >= 0.0 && lambda <= 1.0);
1036 if (!lambdaTerm) {
1037 return;
1038 }
1039 this.lambda = lambda;
1040
1041 initSoftCore();
1042 alchemicalParameters.update(lambda);
1043 if (generalizedKirkwoodTerm) {
1044 generalizedKirkwood.setLambda(alchemicalParameters.polLambda);
1045 generalizedKirkwood.setLambdaFunction(
1046 alchemicalParameters.lPowPol,
1047 alchemicalParameters.dlPowPol,
1048 alchemicalParameters.d2lPowPol);
1049 }
1050 }
1051
1052 public double getPermanentEnergy() {
1053 return permanentMultipoleEnergy;
1054 }
1055
1056
1057
1058
1059
1060
1061 public double getIndRealEnergy() {
1062 return inducedRealSpaceEnergy;
1063 }
1064
1065
1066
1067
1068
1069
1070 public double getIndRecipEnergy() {
1071 return inducedReciprocalEnergy;
1072 }
1073
1074
1075
1076
1077
1078
1079 public double getIndSelfEnergy() {
1080 return inducedSelfEnergy;
1081 }
1082
1083
1084
1085
1086
1087
1088 public double getPermSelfEnergy() {
1089 return permanentSelfEnergy;
1090 }
1091
1092 public double getPermRealEnergy() {
1093 return permanentRealSpaceEnergy;
1094 }
1095
1096 public double getPermRecipEnergy() {
1097 return permanentReciprocalEnergy;
1098 }
1099
1100 public double getPermanentChargeCorrectionEnergy() {
1101 return permanentChargeCorrectionEnergy;
1102 }
1103
1104 public double getPolarEps() {
1105 return poleps;
1106 }
1107
1108 public int[][] getPolarization11() {
1109 return ip11;
1110 }
1111
1112 public int[][] getPolarization12() {
1113 return ip12;
1114 }
1115
1116 public int[][] getPolarization13() {
1117 return ip13;
1118 }
1119
1120
1121
1122
1123
1124
1125 public double getPolarizationEnergy() {
1126 return polarizationEnergy;
1127 }
1128
1129 public Polarization getPolarizationType() {
1130 return polarization;
1131 }
1132
1133 public ReciprocalSpace getReciprocalSpace() {
1134 return reciprocalSpace;
1135 }
1136
1137 public double getScale14() {
1138 return scaleParameters.m14scale;
1139 }
1140
1141
1142
1143
1144
1145
1146 public double getSolvationEnergy() {
1147 return solvationEnergy;
1148 }
1149
1150
1151
1152
1153
1154
1155
1156 public MultipoleType getMultipoleType(int i) {
1157 Atom atom = atoms[i];
1158 MultipoleType multipoleType = atom.getMultipoleType();
1159 if (esvTerm && extendedSystem.isTitrating(i)) {
1160 double[] multipole = multipoleType.getMultipole();
1161 double[] esvMultipole = new double[10];
1162 arraycopy(multipole, 0, esvMultipole, 0, multipole.length);
1163 double titrationLambda = extendedSystem.getTitrationLambda(i);
1164 double tautomerLambda = extendedSystem.getTautomerLambda(i);
1165 esvMultipole = extendedSystem.getTitrationUtils()
1166 .getMultipole(atom, titrationLambda, tautomerLambda, esvMultipole);
1167
1168 multipoleType = new MultipoleType(multipoleType, esvMultipole);
1169 }
1170 return multipoleType;
1171 }
1172
1173
1174
1175
1176
1177
1178
1179 public PolarizeType getPolarizeType(int i) {
1180 Atom atom = atoms[i];
1181 PolarizeType polarizeType = atom.getPolarizeType();
1182 if (polarizeType != null) {
1183 if (esvTerm && extendedSystem.isTitrating(i) && (extendedSystem.isTitratingHydrogen(i)
1184 || extendedSystem.isTitratingHeavy(i))) {
1185 double titrationLambda = extendedSystem.getTitrationLambda(i);
1186 double tautomerLambda = extendedSystem.getTautomerLambda(i);
1187 double esvPolarizability = extendedSystem.getTitrationUtils()
1188 .getPolarizability(atom, titrationLambda, tautomerLambda, polarizeType.polarizability);
1189 polarizeType = new PolarizeType(polarizeType, esvPolarizability);
1190 }
1191 }
1192 return polarizeType;
1193 }
1194
1195
1196
1197
1198 @Override
1199 public double getd2EdL2() {
1200 if (sharedd2EdLambda2 == null || !lambdaTerm) {
1201 return 0.0;
1202 }
1203 double d2EdL2 = sharedd2EdLambda2.get();
1204 if (generalizedKirkwoodTerm || alchemicalParameters.doLigandGKElec) {
1205 d2EdL2 += generalizedKirkwood.getd2EdL2();
1206 }
1207 return d2EdL2;
1208 }
1209
1210
1211
1212
1213 @Override
1214 public double getdEdL() {
1215 if (shareddEdLambda == null || !lambdaTerm) {
1216 return 0.0;
1217 }
1218 double dEdL = shareddEdLambda.get();
1219 if (generalizedKirkwoodTerm || alchemicalParameters.doLigandGKElec) {
1220 dEdL += generalizedKirkwood.getdEdL();
1221 }
1222 return dEdL;
1223 }
1224
1225
1226
1227
1228 @Override
1229 public void getdEdXdL(double[] gradient) {
1230 if (lambdaGrad == null || !lambdaTerm) {
1231 return;
1232 }
1233
1234 int index = 0;
1235 for (int i = 0; i < nAtoms; i++) {
1236 if (atoms[i].isActive()) {
1237 gradient[index++] += lambdaGrad.getX(i);
1238 gradient[index++] += lambdaGrad.getY(i);
1239 gradient[index++] += lambdaGrad.getZ(i);
1240 }
1241 }
1242 }
1243
1244 public void setAtoms(Atom[] atoms, int[] molecule) {
1245 if (lambdaTerm && atoms.length != nAtoms) {
1246 logger.severe(" Changing the number of atoms is not compatible with use of Lambda.");
1247 }
1248 this.atoms = atoms;
1249 this.molecule = molecule;
1250 nAtoms = atoms.length;
1251 initAtomArrays();
1252
1253 if (reciprocalSpace != null) {
1254 reciprocalSpace.setAtoms(atoms);
1255 }
1256
1257 if (generalizedKirkwood != null) {
1258 generalizedKirkwood.setAtoms(atoms);
1259 }
1260 }
1261
1262 public void setCrystal(Crystal crystal) {
1263
1264 int nSymmNew = crystal.spaceGroup.getNumberOfSymOps();
1265 if (nSymm < nSymmNew) {
1266 coordinates = new double[nSymmNew][3][nAtoms];
1267 globalMultipole = new double[nSymmNew][nAtoms][10];
1268 fractionalMultipole = new double[nSymmNew][nAtoms][10];
1269 inducedDipole = new double[nSymmNew][nAtoms][3];
1270 inducedDipoleCR = new double[nSymmNew][nAtoms][3];
1271 if (generalizedKirkwood != null) {
1272 generalizedKirkwood.setAtoms(atoms);
1273 }
1274 realSpaceNeighborParameters.allocate(nAtoms, nSymmNew);
1275 pcgSolver.allocateLists(nSymmNew, nAtoms);
1276 }
1277 nSymm = nSymmNew;
1278 neighborLists = neighborList.getNeighborList();
1279 this.crystal = crystal;
1280 if (reciprocalSpace != null) {
1281 reciprocalSpace.setCrystal(crystal.getUnitCell());
1282 }
1283 }
1284
1285 private void initAtomArrays() {
1286 if (localMultipole == null || localMultipole.length < nAtoms) {
1287 localMultipole = new double[nAtoms][10];
1288 frame = new MultipoleFrameDefinition[nAtoms];
1289 axisAtom = new int[nAtoms][];
1290 cartesianMultipolePhi = new double[nAtoms][tensorCount];
1291 fracMultipolePhi = new double[nAtoms][tensorCount];
1292 directDipole = new double[nAtoms][3];
1293 directDipoleCR = new double[nAtoms][3];
1294 directField = new double[nAtoms][3];
1295 directFieldCR = new double[nAtoms][3];
1296 vacuumDirectDipole = new double[nAtoms][3];
1297 vacuumDirectDipoleCR = new double[nAtoms][3];
1298 if (optRegion != null) {
1299 int optOrder = optRegion.optOrder;
1300 optRegion.optDipole = new double[optOrder + 1][nAtoms][3];
1301 optRegion.optDipoleCR = new double[optOrder + 1][nAtoms][3];
1302 }
1303
1304 cartesianInducedDipolePhi = new double[nAtoms][tensorCount];
1305 cartesianInducedDipolePhiCR = new double[nAtoms][tensorCount];
1306 fractionalInducedDipolePhi = new double[nAtoms][tensorCount];
1307 fractionalInducedDipolePhiCR = new double[nAtoms][tensorCount];
1308
1309 cartesianVacuumDipolePhi = new double[nAtoms][tensorCount];
1310 cartesianVacuumDipolePhiCR = new double[nAtoms][tensorCount];
1311 fractionalVacuumDipolePhi = new double[nAtoms][tensorCount];
1312 fractionalVacuumDipolePhiCR = new double[nAtoms][tensorCount];
1313
1314 mask12 = new int[nAtoms][];
1315 mask13 = new int[nAtoms][];
1316 mask14 = new int[nAtoms][];
1317 mask15 = new int[nAtoms][];
1318 ip11 = new int[nAtoms][];
1319 ip12 = new int[nAtoms][];
1320 ip13 = new int[nAtoms][];
1321 thole = new double[nAtoms];
1322 ipdamp = new double[nAtoms];
1323 polarizability = new double[nAtoms];
1324
1325 if (scfAlgorithm == SCFAlgorithm.CG) {
1326 pcgSolver.allocateVectors(nAtoms);
1327 }
1328 pcgSolver.allocateLists(nSymm, nAtoms);
1329
1330 if (scfPredictorParameters.scfPredictor != SCFPredictor.NONE) {
1331 int predictorOrder = scfPredictorParameters.predictorOrder;
1332 if (lambdaTerm) {
1333 scfPredictorParameters.predictorInducedDipole = new double[3][predictorOrder][nAtoms][3];
1334 scfPredictorParameters.predictorInducedDipoleCR =
1335 new double[3][predictorOrder][nAtoms][3];
1336 } else {
1337 scfPredictorParameters.predictorInducedDipole = new double[1][predictorOrder][nAtoms][3];
1338 scfPredictorParameters.predictorInducedDipoleCR =
1339 new double[1][predictorOrder][nAtoms][3];
1340 }
1341 }
1342
1343
1344
1345 grad = new AtomicDoubleArray3D(atomicDoubleArrayImpl, nAtoms, maxThreads);
1346 torque = new AtomicDoubleArray3D(atomicDoubleArrayImpl, nAtoms, maxThreads);
1347 field = new AtomicDoubleArray3D(atomicDoubleArrayImpl, nAtoms, maxThreads);
1348 fieldCR = new AtomicDoubleArray3D(atomicDoubleArrayImpl, nAtoms, maxThreads);
1349 lambdaGrad = new AtomicDoubleArray3D(atomicDoubleArrayImpl, nAtoms, maxThreads);
1350 lambdaTorque = new AtomicDoubleArray3D(atomicDoubleArrayImpl, nAtoms, maxThreads);
1351 isSoft = new boolean[nAtoms];
1352 use = new boolean[nAtoms];
1353
1354 coordinates = new double[nSymm][3][nAtoms];
1355 globalMultipole = new double[nSymm][nAtoms][10];
1356 fractionalMultipole = new double[nSymm][nAtoms][10];
1357 inducedDipole = new double[nSymm][nAtoms][3];
1358 inducedDipoleCR = new double[nSymm][nAtoms][3];
1359 vacuumInducedDipole = new double[nSymm][nAtoms][3];
1360 vacuumInducedDipoleCR = new double[nSymm][nAtoms][3];
1361
1362
1363 realSpaceNeighborParameters.allocate(nAtoms, nSymm);
1364
1365
1366 lambdaFactors = new LambdaFactors[maxThreads];
1367 for (int i = 0; i < maxThreads; i++) {
1368
1369
1370
1371
1372 if (lambdaTerm) {
1373
1374 lambdaFactors[i] = new LambdaFactorsOST();
1375 } else {
1376
1377 lambdaFactors[i] = LambdaDefaults;
1378 }
1379 }
1380 }
1381
1382
1383 fill(isSoft, false);
1384
1385
1386 fill(use, true);
1387
1388
1389 assignMultipoles();
1390
1391
1392 if (elecForm == PAM) {
1393 assignPolarizationGroups();
1394 }
1395
1396
1397 for (int i = 0; i < nAtoms; i++) {
1398 Atom ai = atoms[i];
1399
1400 if (elecForm == PAM) {
1401 PolarizeType polarizeType = ai.getPolarizeType();
1402 int index = ai.getIndex() - 1;
1403 thole[index] = polarizeType.thole;
1404 ipdamp[index] = polarizeType.pdamp;
1405 if (!(ipdamp[index] > 0.0)) {
1406 ipdamp[index] = Double.POSITIVE_INFINITY;
1407 } else {
1408 ipdamp[index] = 1.0 / ipdamp[index];
1409 }
1410 polarizability[index] = polarizeType.polarizability;
1411 }
1412
1413
1414 List<Atom> n12 = ai.get12List();
1415 mask12[i] = new int[n12.size()];
1416 int j = 0;
1417 for (Atom a12 : n12) {
1418 mask12[i][j++] = a12.getIndex() - 1;
1419 }
1420
1421
1422 List<Atom> n13 = ai.get13List();
1423 mask13[i] = new int[n13.size()];
1424 j = 0;
1425 for (Atom a13 : n13) {
1426 mask13[i][j++] = a13.getIndex() - 1;
1427 }
1428
1429
1430 List<Atom> n14 = ai.get14List();
1431 mask14[i] = new int[n14.size()];
1432 j = 0;
1433 for (Atom a14 : n14) {
1434 mask14[i][j++] = a14.getIndex() - 1;
1435 }
1436
1437
1438 List<Atom> n15 = ai.get15List();
1439 mask15[i] = new int[n15.size()];
1440 j = 0;
1441 for (Atom a15 : n15) {
1442 mask15[i][j++] = a15.getIndex() - 1;
1443 }
1444 }
1445 }
1446
1447
1448
1449
1450 private void initSoftCore() {
1451
1452 boolean rebuild = false;
1453
1454
1455 StringBuilder sb = new StringBuilder("\n Softcore Atoms:\n");
1456 int count = 0;
1457 for (int i = 0; i < nAtoms; i++) {
1458 Atom ai = atoms[i];
1459 boolean soft = ai.applyLambda();
1460 if (soft != isSoft[i]) {
1461 rebuild = true;
1462 }
1463 isSoft[i] = soft;
1464 if (soft) {
1465 count++;
1466 sb.append(ai).append("\n");
1467 }
1468 }
1469
1470
1471 if (alchemicalParameters.doLigandVaporElec
1472 && alchemicalParameters.vaporCrystal == null) {
1473 rebuild = true;
1474 }
1475
1476 if (!rebuild) {
1477 return;
1478 }
1479
1480 if (count > 0 && logger.isLoggable(Level.FINE)) {
1481 logger.fine(format(" Softcore atom count: %d", count));
1482 logger.fine(sb.toString());
1483 }
1484
1485
1486
1487 if (alchemicalParameters.doLigandVaporElec) {
1488 double vacuumOff = getVacuumOff();
1489 alchemicalParameters.vaporCrystal =
1490 new Crystal(3 * vacuumOff, 3 * vacuumOff, 3 * vacuumOff, 90.0, 90.0, 90.0, "P1");
1491 alchemicalParameters.vaporCrystal.setAperiodic(true);
1492 NeighborList vacuumNeighborList = new NeighborList(alchemicalParameters.vaporCrystal,
1493 atoms, vacuumOff, 2.0, parallelTeam);
1494 vacuumNeighborList.setIntermolecular(false);
1495
1496 alchemicalParameters.vaporLists = new int[1][nAtoms][];
1497 double[][] coords = new double[1][nAtoms * 3];
1498 for (int i = 0; i < nAtoms; i++) {
1499 coords[0][i * 3] = atoms[i].getX();
1500 coords[0][i * 3 + 1] = atoms[i].getY();
1501 coords[0][i * 3 + 2] = atoms[i].getZ();
1502 }
1503 boolean print = logger.isLoggable(Level.FINE);
1504 vacuumNeighborList.buildList(coords, alchemicalParameters.vaporLists, isSoft,
1505 true, print);
1506 alchemicalParameters.vaporPermanentSchedule = vacuumNeighborList.getPairwiseSchedule();
1507 alchemicalParameters.vaporEwaldSchedule = alchemicalParameters.vaporPermanentSchedule;
1508 alchemicalParameters.vacuumRanges = new Range[maxThreads];
1509 try {
1510 vacuumNeighborList.destroy();
1511 } catch (Exception ex) {
1512 logger.warning(" Exception in destroying vacuumNeighborList");
1513 }
1514 } else {
1515 alchemicalParameters.vaporCrystal = null;
1516 alchemicalParameters.vaporLists = null;
1517 alchemicalParameters.vaporPermanentSchedule = null;
1518 alchemicalParameters.vaporEwaldSchedule = null;
1519 alchemicalParameters.vacuumRanges = null;
1520 }
1521 }
1522
1523 private double getVacuumOff() {
1524 double maxr = 10.0;
1525 for (int i = 0; i < nAtoms; i++) {
1526 Atom ai = atoms[i];
1527 if (ai.applyLambda()) {
1528
1529
1530 for (int j = i + 1; j < nAtoms; j++) {
1531 Atom aj = atoms[j];
1532 if (aj.applyLambda()) {
1533 double dx = ai.getX() - aj.getX();
1534 double dy = ai.getY() - aj.getY();
1535 double dz = ai.getZ() - aj.getZ();
1536 double r = sqrt(dx * dx + dy * dy + dz * dz);
1537 maxr = max(r, maxr);
1538 }
1539 }
1540 }
1541 }
1542
1543 double vacuumOff = 2 * maxr;
1544 return vacuumOff;
1545 }
1546
1547
1548
1549
1550 private void initNN() {
1551
1552
1553 StringBuilder sb = new StringBuilder("\n NN Atoms:\n");
1554 boolean[] nnAtoms = new boolean[nAtoms];
1555
1556 int count = 0;
1557 for (int i = 0; i < nAtoms; i++) {
1558 Atom ai = atoms[i];
1559 nnAtoms[i] = ai.isNeuralNetwork();
1560 if (nnAtoms[i]) {
1561 count++;
1562 sb.append(ai).append("\n");
1563 }
1564 }
1565
1566 if (count > 0 && logger.isLoggable(Level.FINE)) {
1567 logger.fine(format(" Neural network atom count: %d", count));
1568 logger.fine(sb.toString());
1569 }
1570
1571
1572
1573 if (alchemicalParameters.doLigandVaporElec) {
1574 alchemicalParameters.vaporCrystal = new Crystal(10.0, 10.0, 10.0,
1575 90.0, 90.0, 90.0, "P1");
1576 alchemicalParameters.vaporCrystal.setAperiodic(true);
1577 NeighborList vacuumNeighborList = new NeighborList(alchemicalParameters.vaporCrystal,
1578 atoms, Double.POSITIVE_INFINITY, 2.0, parallelTeam);
1579 vacuumNeighborList.setIntermolecular(false);
1580
1581 alchemicalParameters.vaporLists = new int[1][nAtoms][];
1582 double[][] coords = new double[1][nAtoms * 3];
1583 for (int i = 0; i < nAtoms; i++) {
1584 coords[0][i * 3] = atoms[i].getX();
1585 coords[0][i * 3 + 1] = atoms[i].getY();
1586 coords[0][i * 3 + 2] = atoms[i].getZ();
1587 }
1588 boolean print = logger.isLoggable(Level.FINE);
1589 vacuumNeighborList.buildList(coords, alchemicalParameters.vaporLists, nnAtoms, true, print);
1590 alchemicalParameters.vaporPermanentSchedule = vacuumNeighborList.getPairwiseSchedule();
1591 alchemicalParameters.vaporEwaldSchedule = alchemicalParameters.vaporPermanentSchedule;
1592 alchemicalParameters.vacuumRanges = new Range[maxThreads];
1593 vacuumNeighborList.setDisableUpdates(forceField.getBoolean("DISABLE_NEIGHBOR_UPDATES", false));
1594 } else {
1595 alchemicalParameters.vaporCrystal = null;
1596 alchemicalParameters.vaporLists = null;
1597 alchemicalParameters.vaporPermanentSchedule = null;
1598 alchemicalParameters.vaporEwaldSchedule = null;
1599 alchemicalParameters.vacuumRanges = null;
1600 }
1601 }
1602
1603
1604
1605
1606
1607
1608
1609 private double condensedEnergy() {
1610
1611 if (alchemicalParameters.mode == AlchemicalParameters.AlchemicalMode.SCALE) {
1612 alchemicalParameters.polarizationScale = 1.0;
1613 alchemicalParameters.doPolarization = true;
1614 } else if (lambda < alchemicalParameters.polLambdaStart) {
1615
1616
1617
1618
1619
1620 alchemicalParameters.polarizationScale = 0.0;
1621 alchemicalParameters.doPolarization = false;
1622 } else if (lambda <= alchemicalParameters.polLambdaEnd) {
1623 alchemicalParameters.polarizationScale = alchemicalParameters.lPowPol;
1624 alchemicalParameters.doPolarization = true;
1625 } else {
1626 alchemicalParameters.polarizationScale = 1.0;
1627 alchemicalParameters.doPolarization = true;
1628 }
1629
1630
1631 if (alchemicalParameters.mode == AlchemicalParameters.AlchemicalMode.SCALE) {
1632 alchemicalParameters.doPermanentRealSpace = true;
1633 alchemicalParameters.permanentScale = 1.0;
1634 alchemicalParameters.dEdLSign = 1.0;
1635 } else {
1636 alchemicalParameters.doPermanentRealSpace = true;
1637 alchemicalParameters.permanentScale = alchemicalParameters.lPowPerm;
1638 alchemicalParameters.dEdLSign = 1.0;
1639 }
1640
1641 return computeEnergy(false);
1642 }
1643
1644
1645
1646
1647
1648
1649
1650 private double condensedNoLigandSCF() {
1651
1652 boolean skip = true;
1653 for (int i = 0; i < nAtoms; i++) {
1654 if (atoms[i].applyLambda()) {
1655 use[i] = false;
1656 } else {
1657 use[i] = true;
1658 skip = false;
1659 }
1660 }
1661
1662
1663 alchemicalParameters.doPermanentRealSpace = false;
1664 alchemicalParameters.permanentScale =
1665 1.0 - alchemicalParameters.lPowPerm;
1666 alchemicalParameters.dEdLSign = -1.0;
1667
1668
1669
1670 if (lambda <= alchemicalParameters.polLambdaEnd
1671 && alchemicalParameters.doNoLigandCondensedSCF) {
1672 alchemicalParameters.doPolarization = true;
1673 alchemicalParameters.polarizationScale =
1674 1.0 - alchemicalParameters.lPowPol;
1675 } else {
1676 alchemicalParameters.doPolarization = false;
1677 alchemicalParameters.polarizationScale = 0.0;
1678 }
1679
1680
1681 boolean gkBack = generalizedKirkwoodTerm;
1682 generalizedKirkwoodTerm = false;
1683
1684
1685
1686 double energy;
1687 if (skip) {
1688 energy = permanentMultipoleEnergy + polarizationEnergy + solvationEnergy;
1689 } else {
1690 energy = computeEnergy(false);
1691 Arrays.fill(use, true);
1692 }
1693
1694 generalizedKirkwoodTerm = gkBack;
1695
1696 return energy;
1697 }
1698
1699
1700
1701
1702
1703
1704 private double ligandElec() {
1705 for (int i = 0; i < nAtoms; i++) {
1706 use[i] = atoms[i].applyLambda();
1707 }
1708
1709
1710
1711 alchemicalParameters.doPermanentRealSpace = true;
1712 alchemicalParameters.permanentScale = 1.0 - alchemicalParameters.lPowPerm;
1713 alchemicalParameters.dEdLSign = -1.0;
1714 double lAlphaBack = alchemicalParameters.lAlpha;
1715 double dlAlphaBack = alchemicalParameters.dlAlpha;
1716 double d2lAlphaBack = alchemicalParameters.d2lAlpha;
1717 alchemicalParameters.lAlpha = 0.0;
1718 alchemicalParameters.dlAlpha = 0.0;
1719 alchemicalParameters.d2lAlpha = 0.0;
1720
1721
1722
1723 if (lambda <= alchemicalParameters.polLambdaEnd) {
1724 alchemicalParameters.doPolarization = true;
1725 alchemicalParameters.polarizationScale = 1.0 - alchemicalParameters.lPowPol;
1726 } else {
1727 alchemicalParameters.doPolarization = false;
1728 alchemicalParameters.polarizationScale = 0.0;
1729 }
1730
1731
1732 double offBack = ewaldParameters.off;
1733 double aewaldBack = ewaldParameters.aewald;
1734 ewaldParameters.setEwaldParameters(Double.MAX_VALUE, 0.0);
1735
1736 IntegerSchedule permanentScheduleBack = permanentSchedule;
1737 IntegerSchedule ewaldScheduleBack = realSpaceNeighborParameters.realSpaceSchedule;
1738 Range[] rangesBack = realSpaceNeighborParameters.realSpaceRanges;
1739 permanentSchedule = alchemicalParameters.vaporPermanentSchedule;
1740 realSpaceNeighborParameters.realSpaceSchedule = alchemicalParameters.vaporEwaldSchedule;
1741 realSpaceNeighborParameters.realSpaceRanges = alchemicalParameters.vacuumRanges;
1742
1743
1744 Crystal crystalBack = crystal;
1745 int nSymmBack = nSymm;
1746 int[][][] listsBack = neighborLists;
1747 neighborLists = alchemicalParameters.vaporLists;
1748 crystal = alchemicalParameters.vaporCrystal;
1749 nSymm = 1;
1750
1751
1752 boolean gkBack = generalizedKirkwoodTerm;
1753
1754
1755 SCFAlgorithm scfBack = scfAlgorithm;
1756 scfAlgorithm = SOR;
1757
1758 if (alchemicalParameters.doLigandGKElec) {
1759 generalizedKirkwoodTerm = true;
1760 generalizedKirkwood.setNeighborList(alchemicalParameters.vaporLists);
1761 generalizedKirkwood.setLambda(lambda);
1762 generalizedKirkwood.setCutoff(ewaldParameters.off);
1763 generalizedKirkwood.setCrystal(alchemicalParameters.vaporCrystal);
1764 generalizedKirkwood.setLambdaFunction(alchemicalParameters.polarizationScale,
1765 alchemicalParameters.dEdLSign * alchemicalParameters.dlPowPol,
1766 alchemicalParameters.dEdLSign * alchemicalParameters.d2lPowPol);
1767 } else {
1768 generalizedKirkwoodTerm = false;
1769 }
1770
1771 double energy = computeEnergy(false);
1772
1773
1774 ewaldParameters.setEwaldParameters(offBack, aewaldBack);
1775 neighborLists = listsBack;
1776 crystal = crystalBack;
1777 nSymm = nSymmBack;
1778 permanentSchedule = permanentScheduleBack;
1779 realSpaceNeighborParameters.realSpaceSchedule = ewaldScheduleBack;
1780 realSpaceNeighborParameters.realSpaceRanges = rangesBack;
1781 alchemicalParameters.lAlpha = lAlphaBack;
1782 alchemicalParameters.dlAlpha = dlAlphaBack;
1783 alchemicalParameters.d2lAlpha = d2lAlphaBack;
1784 generalizedKirkwoodTerm = gkBack;
1785 scfAlgorithm = scfBack;
1786
1787 fill(use, true);
1788
1789 return energy;
1790 }
1791
1792
1793
1794
1795
1796
1797 private double nnElec() {
1798 for (int i = 0; i < nAtoms; i++) {
1799 use[i] = atoms[i].isNeuralNetwork();
1800 }
1801
1802
1803
1804 alchemicalParameters.doPermanentRealSpace = true;
1805 alchemicalParameters.permanentScale = -1.0;
1806 alchemicalParameters.dEdLSign = 1.0;
1807 double lAlphaBack = alchemicalParameters.lAlpha;
1808 double dlAlphaBack = alchemicalParameters.dlAlpha;
1809 double d2lAlphaBack = alchemicalParameters.d2lAlpha;
1810 alchemicalParameters.lAlpha = 0.0;
1811 alchemicalParameters.dlAlpha = 0.0;
1812 alchemicalParameters.d2lAlpha = 0.0;
1813 alchemicalParameters.doPolarization = true;
1814 alchemicalParameters.polarizationScale = -1.0;
1815
1816
1817 double offBack = ewaldParameters.off;
1818 double aewaldBack = ewaldParameters.aewald;
1819 ewaldParameters.setEwaldParameters(Double.POSITIVE_INFINITY, 0.0);
1820
1821 IntegerSchedule permanentScheduleBack = permanentSchedule;
1822 IntegerSchedule ewaldScheduleBack = realSpaceNeighborParameters.realSpaceSchedule;
1823 Range[] rangesBack = realSpaceNeighborParameters.realSpaceRanges;
1824 permanentSchedule = alchemicalParameters.vaporPermanentSchedule;
1825 realSpaceNeighborParameters.realSpaceSchedule = alchemicalParameters.vaporEwaldSchedule;
1826 realSpaceNeighborParameters.realSpaceRanges = alchemicalParameters.vacuumRanges;
1827
1828
1829 Crystal crystalBack = crystal;
1830 int nSymmBack = nSymm;
1831 int[][][] listsBack = neighborLists;
1832 neighborLists = alchemicalParameters.vaporLists;
1833 crystal = alchemicalParameters.vaporCrystal;
1834 nSymm = 1;
1835
1836
1837 boolean gkBack = generalizedKirkwoodTerm;
1838
1839
1840 SCFAlgorithm scfBack = scfAlgorithm;
1841 scfAlgorithm = SOR;
1842 generalizedKirkwoodTerm = false;
1843
1844 double energy = computeEnergy(false);
1845
1846
1847 ewaldParameters.setEwaldParameters(offBack, aewaldBack);
1848 neighborLists = listsBack;
1849 crystal = crystalBack;
1850 nSymm = nSymmBack;
1851 permanentSchedule = permanentScheduleBack;
1852 realSpaceNeighborParameters.realSpaceSchedule = ewaldScheduleBack;
1853 realSpaceNeighborParameters.realSpaceRanges = rangesBack;
1854 alchemicalParameters.lAlpha = lAlphaBack;
1855 alchemicalParameters.dlAlpha = dlAlphaBack;
1856 alchemicalParameters.d2lAlpha = d2lAlphaBack;
1857 generalizedKirkwoodTerm = gkBack;
1858 scfAlgorithm = scfBack;
1859
1860 fill(use, true);
1861
1862 return energy;
1863 }
1864
1865
1866
1867
1868
1869
1870
1871 private double computeEnergy(boolean print) {
1872
1873 permanentMultipoleField();
1874
1875
1876 if (generalizedKirkwoodTerm) {
1877 pmeTimings.bornRadiiTotal -= System.nanoTime();
1878 generalizedKirkwood.setUse(use);
1879 generalizedKirkwood.computeBornRadii();
1880 pmeTimings.bornRadiiTotal += System.nanoTime();
1881 }
1882
1883
1884 if (polarization != Polarization.NONE && alchemicalParameters.doPolarization) {
1885
1886 if (generalizedKirkwoodTerm) {
1887 generalizedKirkwoodTerm = false;
1888
1889 selfConsistentField(logger.isLoggable(Level.FINE));
1890
1891 saveInducedDipolesToVacuumDipoles();
1892
1893 generalizedKirkwoodTerm = true;
1894
1895 permanentMultipoleField();
1896 }
1897
1898
1899 selfConsistentField(logger.isLoggable(Level.FINE));
1900
1901 if (esvTerm && polarization != Polarization.NONE) {
1902 for (int i = 0; i < nAtoms; i++) {
1903 if (extendedSystem.isTitrating(i)
1904 && (extendedSystem.isTitratingHydrogen(i) || extendedSystem.isTitratingHeavy(i))) {
1905 double dx = field.getX(i);
1906 double dy = field.getY(i);
1907 double dz = field.getZ(i);
1908 double dxCR = fieldCR.getX(i);
1909 double dyCR = fieldCR.getY(i);
1910 double dzCR = fieldCR.getZ(i);
1911
1912
1913 if (polarization == Polarization.MUTUAL) {
1914 dx += directField[i][0];
1915 dy += directField[i][1];
1916 dz += directField[i][2];
1917 dxCR += directFieldCR[i][0];
1918 dyCR += directFieldCR[i][1];
1919 dzCR += directFieldCR[i][2];
1920 }
1921 double fix = dx * dPolardTitrationESV[i] * dxCR;
1922 double fiy = dy * dPolardTitrationESV[i] * dyCR;
1923 double fiz = dz * dPolardTitrationESV[i] * dzCR;
1924 double titrdUdL = fix + fiy + fiz;
1925 double tautdUdL = 0.0;
1926 if (extendedSystem.isTautomerizing(i)) {
1927 fix = dx * dPolardTautomerESV[i] * dxCR;
1928 fiy = dy * dPolardTautomerESV[i] * dyCR;
1929 fiz = dz * dPolardTautomerESV[i] * dzCR;
1930 tautdUdL = fix + fiy + fiz;
1931 }
1932 extendedSystem.addIndElecDeriv(i, titrdUdL * electric * -0.5, tautdUdL * electric * -0.5);
1933 }
1934 }
1935 }
1936
1937 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
1938 if (gradient && polarization == Polarization.DIRECT) {
1939 reciprocalSpace.splineInducedDipoles(inducedDipole, inducedDipoleCR, use);
1940 field.reset(parallelTeam);
1941 fieldCR.reset(parallelTeam);
1942 inducedDipoleFieldRegion.init(atoms, crystal, use, molecule, ipdamp, thole, coordinates,
1943 realSpaceNeighborParameters, inducedDipole, inducedDipoleCR, reciprocalSpaceTerm,
1944 reciprocalSpace, lambdaMode, ewaldParameters, field, fieldCR, pmeTimings);
1945 inducedDipoleFieldRegion.executeWith(sectionTeam);
1946 reciprocalSpace.computeInducedPhi(
1947 cartesianInducedDipolePhi, cartesianInducedDipolePhiCR,
1948 fractionalInducedDipolePhi, fractionalInducedDipolePhiCR);
1949 }
1950 }
1951
1952 if (scfPredictorParameters.scfPredictor != SCFPredictor.NONE) {
1953 scfPredictorParameters.saveMutualInducedDipoles(lambdaMode,
1954 inducedDipole, inducedDipoleCR, directDipole, directDipoleCR);
1955 }
1956 }
1957
1958
1959
1960
1961
1962
1963
1964
1965
1966
1967 double[][][] inputDipole = inducedDipole;
1968 double[][][] inputDipoleCR = inducedDipoleCR;
1969 double[][] inputPhi = cartesianInducedDipolePhi;
1970 double[][] inputPhiCR = cartesianInducedDipolePhiCR;
1971 double[][] fracInputPhi = fractionalInducedDipolePhi;
1972 double[][] fracInputPhiCR = fractionalInducedDipolePhiCR;
1973 if (generalizedKirkwoodTerm) {
1974 inputDipole = vacuumInducedDipole;
1975 inputDipoleCR = vacuumInducedDipoleCR;
1976 inputPhi = cartesianVacuumDipolePhi;
1977 inputPhiCR = cartesianVacuumDipolePhiCR;
1978 fracInputPhi = fractionalVacuumDipolePhi;
1979 fracInputPhiCR = fractionalVacuumDipolePhiCR;
1980 }
1981
1982 double eself = 0.0;
1983 double erecip = 0.0;
1984 double ecorrect = 0.0;
1985 double eselfi = 0.0;
1986 double erecipi = 0.0;
1987 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
1988 reciprocalEnergyRegion.init(atoms, crystal, gradient, lambdaTerm, esvTerm, use,
1989 globalMultipole, fractionalMultipole, dMultipoledTirationESV, dMultipoledTautomerESV,
1990 cartesianMultipolePhi, fracMultipolePhi,
1991 polarization, inputDipole, inputDipoleCR, inputPhi, inputPhiCR, fracInputPhi,
1992 fracInputPhiCR, reciprocalSpace, alchemicalParameters, extendedSystem,
1993
1994 grad, torque, lambdaGrad, lambdaTorque, shareddEdLambda, sharedd2EdLambda2);
1995 reciprocalEnergyRegion.executeWith(parallelTeam);
1996 eself = reciprocalEnergyRegion.getPermanentSelfEnergy();
1997 erecip = reciprocalEnergyRegion.getPermanentReciprocalEnergy();
1998 ecorrect = reciprocalEnergyRegion.getPermanentChargeCorrectionEnergy();
1999 eselfi = reciprocalEnergyRegion.getInducedDipoleSelfEnergy();
2000 erecipi = reciprocalEnergyRegion.getInducedDipoleReciprocalEnergy();
2001 interactions += nAtoms;
2002 }
2003
2004 pmeTimings.realSpaceEnergyTotal -= System.nanoTime();
2005 realSpaceEnergyRegion.init(atoms, crystal, extendedSystem, esvTerm, coordinates, frame, axisAtom,
2006 globalMultipole, dMultipoledTirationESV, dMultipoledTautomerESV,
2007 inputDipole, inputDipoleCR, use, molecule,
2008 ip11, mask12, mask13, mask14, mask15, isSoft, ipdamp, thole, realSpaceNeighborParameters,
2009 gradient, lambdaTerm, nnTerm, lambdaMode, polarization,
2010 ewaldParameters, scaleParameters, alchemicalParameters,
2011 pmeTimings.realSpaceEnergyTime,
2012
2013 grad, torque, lambdaGrad, lambdaTorque, shareddEdLambda, sharedd2EdLambda2);
2014 realSpaceEnergyRegion.executeWith(parallelTeam);
2015 double ereal = realSpaceEnergyRegion.getPermanentEnergy();
2016 double ereali = realSpaceEnergyRegion.getPolarizationEnergy();
2017 interactions += realSpaceEnergyRegion.getInteractions();
2018 pmeTimings.realSpaceEnergyTotal += System.nanoTime();
2019
2020 if (generalizedKirkwoodTerm) {
2021
2022
2023 double eGK = 0.0;
2024
2025
2026 AlchemicalParameters alchemicalParametersGK = new AlchemicalParameters(
2027 forceField, false, false, polarization);
2028
2029 alchemicalParametersGK.permanentScale = 0.0;
2030 alchemicalParametersGK.doPermanentRealSpace = false;
2031
2032 alchemicalParametersGK.polarizationScale = -1.0;
2033
2034
2035 AtomicDoubleArray3D gradGK = generalizedKirkwood.getGrad();
2036 AtomicDoubleArray3D torqueGK = generalizedKirkwood.getTorque();
2037
2038 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
2039 reciprocalEnergyRegion.init(
2040 atoms, crystal, gradient, false, esvTerm,
2041 use, globalMultipole, fractionalMultipole,
2042 dMultipoledTirationESV, dMultipoledTautomerESV,
2043 cartesianMultipolePhi, fracMultipolePhi,
2044 polarization,
2045 vacuumInducedDipole, vacuumInducedDipoleCR,
2046 cartesianVacuumDipolePhi, cartesianVacuumDipolePhiCR,
2047 fractionalVacuumDipolePhi, fractionalVacuumDipolePhiCR,
2048 reciprocalSpace, alchemicalParametersGK,
2049 extendedSystem,
2050 gradGK, torqueGK, null, null,
2051 shareddEdLambda, sharedd2EdLambda2);
2052 reciprocalEnergyRegion.executeWith(parallelTeam);
2053 eGK = reciprocalEnergyRegion.getInducedDipoleSelfEnergy()
2054 + reciprocalEnergyRegion.getInducedDipoleReciprocalEnergy();
2055 }
2056
2057 pmeTimings.gkEnergyTotal -= System.nanoTime();
2058 realSpaceEnergyRegion.init(atoms, crystal,
2059 extendedSystem, esvTerm,
2060 coordinates, frame, axisAtom,
2061 globalMultipole, dMultipoledTirationESV, dMultipoledTautomerESV,
2062 vacuumInducedDipole, vacuumInducedDipoleCR,
2063 use, molecule,
2064 ip11, mask12, mask13, mask14, mask15, isSoft, ipdamp, thole,
2065 realSpaceNeighborParameters,
2066 gradient, false, nnTerm, lambdaMode,
2067 polarization, ewaldParameters, scaleParameters, alchemicalParametersGK,
2068 pmeTimings.realSpaceEnergyTime,
2069
2070 gradGK,
2071 torqueGK,
2072 null,
2073 null,
2074 shareddEdLambda,
2075 sharedd2EdLambda2);
2076 realSpaceEnergyRegion.executeWith(parallelTeam);
2077 eGK += realSpaceEnergyRegion.getPolarizationEnergy();
2078
2079
2080 alchemicalParametersGK.polarizationScale = 1.0;
2081
2082 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
2083 reciprocalEnergyRegion.init(
2084 atoms, crystal, gradient, false, esvTerm, use,
2085 globalMultipole, fractionalMultipole,
2086 dMultipoledTirationESV, dMultipoledTautomerESV,
2087 cartesianMultipolePhi, fracMultipolePhi,
2088 polarization, inducedDipole, inducedDipoleCR,
2089 cartesianInducedDipolePhi, cartesianInducedDipolePhiCR,
2090 fractionalInducedDipolePhi, fractionalInducedDipolePhiCR,
2091 reciprocalSpace, alchemicalParametersGK, extendedSystem,
2092 gradGK, torqueGK, null, null,
2093 shareddEdLambda, sharedd2EdLambda2);
2094 reciprocalEnergyRegion.executeWith(parallelTeam);
2095 eGK += reciprocalEnergyRegion.getInducedDipoleSelfEnergy()
2096 + reciprocalEnergyRegion.getInducedDipoleReciprocalEnergy();
2097 }
2098
2099 pmeTimings.gkEnergyTotal -= System.nanoTime();
2100 realSpaceEnergyRegion.init(
2101 atoms,
2102 crystal,
2103 extendedSystem,
2104 esvTerm,
2105 coordinates,
2106 frame,
2107 axisAtom,
2108 globalMultipole,
2109 dMultipoledTirationESV,
2110 dMultipoledTautomerESV,
2111 inducedDipole,
2112 inducedDipoleCR,
2113 use,
2114 molecule,
2115 ip11,
2116 mask12,
2117 mask13,
2118 mask14,
2119 mask15,
2120 isSoft,
2121 ipdamp,
2122 thole,
2123 realSpaceNeighborParameters,
2124 gradient,
2125 false,
2126 nnTerm,
2127 lambdaMode,
2128 polarization,
2129 ewaldParameters,
2130 scaleParameters,
2131 alchemicalParametersGK,
2132 pmeTimings.realSpaceEnergyTime,
2133
2134 gradGK,
2135 torqueGK,
2136 null,
2137 null,
2138 shareddEdLambda,
2139 sharedd2EdLambda2);
2140 realSpaceEnergyRegion.executeWith(parallelTeam);
2141 eGK += realSpaceEnergyRegion.getPolarizationEnergy();
2142
2143
2144 solvationEnergy += generalizedKirkwood.solvationEnergy(eGK, gradient, print);
2145 if (gradient) {
2146
2147 generalizedKirkwood.reduce(grad, torque, lambdaGrad, lambdaTorque);
2148 }
2149 gkInteractions += generalizedKirkwood.getInteractions();
2150 pmeTimings.gkEnergyTotal += System.nanoTime();
2151 }
2152
2153
2154 permanentRealSpaceEnergy += ereal;
2155 permanentSelfEnergy += eself;
2156 permanentReciprocalEnergy += erecip;
2157 permanentChargeCorrectionEnergy += ecorrect;
2158 inducedRealSpaceEnergy += ereali;
2159 inducedSelfEnergy += eselfi;
2160 inducedReciprocalEnergy += erecipi;
2161 permanentMultipoleEnergy += eself + erecip + ereal + ecorrect;
2162 polarizationEnergy += eselfi + erecipi + ereali;
2163 totalMultipoleEnergy += ereal + eself + erecip + ereali + eselfi + erecipi;
2164
2165
2166 if (logger.isLoggable(Level.FINE)) {
2167 StringBuilder sb = new StringBuilder();
2168 sb.append(format("\n Global Cartesian PME, lambdaMode=%s\n", lambdaMode.toString()));
2169 sb.append(format(" Multipole Reciprocal: %16.8f\n", erecip));
2170 sb.append(format(" Multipole Self-Energy: %16.8f\n", eself));
2171 sb.append(format(" Multipole Correction: %16.8f\n", ecorrect));
2172 sb.append(format(" Multipole Real Space: %16.8f\n", ereal));
2173 sb.append(format(" Polarization Reciprocal: %16.8f\n", erecipi));
2174 sb.append(format(" Polarization Self-Energy:%16.8f\n", eselfi));
2175 sb.append(format(" Polarization Real Space: %16.8f\n", ereali));
2176 if (generalizedKirkwoodTerm) {
2177 sb.append(format(" Generalized Kirkwood: %16.8f\n", solvationEnergy));
2178 }
2179 logger.info(sb.toString());
2180 }
2181
2182 return permanentMultipoleEnergy + polarizationEnergy + solvationEnergy;
2183 }
2184
2185
2186
2187
2188 private void permanentMultipoleField() {
2189 try {
2190
2191 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
2192 reciprocalSpace.computeBSplines();
2193 reciprocalSpace.splinePermanentMultipoles(globalMultipole, fractionalMultipole, use);
2194 }
2195
2196 field.reset(parallelTeam);
2197 fieldCR.reset(parallelTeam);
2198 permanentFieldRegion.init(atoms, crystal, coordinates,
2199 globalMultipole, inducedDipole, inducedDipoleCR, neighborLists,
2200 scaleParameters, use, molecule, ipdamp, thole, ip11,
2201 mask12, mask13, mask14, lambdaMode, reciprocalSpaceTerm, reciprocalSpace,
2202 ewaldParameters, pcgSolver, permanentSchedule, realSpaceNeighborParameters,
2203 field, fieldCR);
2204
2205
2206 sectionTeam.execute(permanentFieldRegion);
2207
2208
2209 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
2210 reciprocalSpace.computePermanentPhi(cartesianMultipolePhi, fracMultipolePhi);
2211 }
2212 } catch (RuntimeException e) {
2213 String message = "Fatal exception computing the permanent multipole field.\n";
2214 logger.log(Level.WARNING, message, e);
2215 throw e;
2216 } catch (Exception e) {
2217 String message = "Fatal exception computing the permanent multipole field.\n";
2218 logger.log(Level.SEVERE, message, e);
2219 }
2220 }
2221
2222 private void saveInducedDipolesToVacuumDipoles() {
2223 for (int i = 0; i < nAtoms; i++) {
2224 arraycopy(directDipole[i], 0, vacuumDirectDipole[i], 0, 3);
2225 arraycopy(directDipoleCR[i], 0, vacuumDirectDipoleCR[i], 0, 3);
2226 arraycopy(inducedDipole[0][i], 0, vacuumInducedDipole[0][i], 0, 3);
2227 arraycopy(inducedDipoleCR[0][i], 0, vacuumInducedDipoleCR[0][i], 0, 3);
2228 arraycopy(cartesianInducedDipolePhi[i], 0, cartesianVacuumDipolePhi[i], 0, tensorCount);
2229 arraycopy(cartesianInducedDipolePhiCR[i], 0, cartesianVacuumDipolePhiCR[i], 0,
2230 tensorCount);
2231 arraycopy(fractionalInducedDipolePhi[i], 0, fractionalVacuumDipolePhi[i], 0,
2232 tensorCount);
2233 arraycopy(fractionalInducedDipolePhiCR[i], 0, fractionalVacuumDipolePhiCR[i], 0,
2234 tensorCount);
2235 if (nSymm > 1) {
2236 for (int s = 1; s < nSymm; s++) {
2237 arraycopy(inducedDipole[s][i], 0, vacuumInducedDipole[s][i], 0, 3);
2238 arraycopy(inducedDipoleCR[s][i], 0, vacuumInducedDipoleCR[s][i], 0, 3);
2239 }
2240 }
2241 }
2242 }
2243
2244
2245
2246
2247 private int selfConsistentField(boolean print) {
2248 if (polarization == Polarization.NONE) {
2249 return -1;
2250 }
2251 long startTime = System.nanoTime();
2252
2253
2254 if (generalizedKirkwoodTerm) {
2255 pmeTimings.gkEnergyTotal = -System.nanoTime();
2256 generalizedKirkwood.computePermanentGKField();
2257 pmeTimings.gkEnergyTotal += System.nanoTime();
2258 logger.fine(format(" Computed GK permanent field %8.3f (sec)", pmeTimings.gkEnergyTotal * 1.0e-9));
2259 }
2260 directRegion.init(atoms, polarizability, globalMultipole, cartesianMultipolePhi,
2261 field, fieldCR, generalizedKirkwoodTerm, generalizedKirkwood, ewaldParameters,
2262 soluteDielectric, inducedDipole, inducedDipoleCR, directDipole, directDipoleCR,
2263 directField, directFieldCR);
2264 directRegion.executeWith(parallelTeam);
2265
2266
2267 if (polarization != Polarization.MUTUAL) {
2268 expandInducedDipoles();
2269 return 0;
2270 }
2271
2272
2273 if (scfPredictorParameters.scfPredictor != SCFPredictor.NONE) {
2274 switch (scfPredictorParameters.scfPredictor) {
2275 case ASPC -> scfPredictorParameters.aspcPredictor(lambdaMode, inducedDipole, inducedDipoleCR);
2276 case LS -> scfPredictorParameters.leastSquaresPredictor(lambdaMode, inducedDipole, inducedDipoleCR);
2277 case POLY -> scfPredictorParameters.polynomialPredictor(lambdaMode, inducedDipole, inducedDipoleCR);
2278 default -> {
2279 }
2280 }
2281 }
2282
2283
2284 expandInducedDipoles();
2285
2286
2287 try {
2288 return switch (scfAlgorithm) {
2289 case SOR -> scfBySOR(print, startTime);
2290 case EPT -> scfByEPT(print, startTime);
2291 default -> {
2292
2293 pcgSolver.init(atoms, coordinates, polarizability, ipdamp, thole,
2294 use, crystal, inducedDipole, inducedDipoleCR, directDipole, directDipoleCR,
2295 field, fieldCR, ewaldParameters, soluteDielectric, parallelTeam,
2296 realSpaceNeighborParameters.realSpaceSchedule, pmeTimings);
2297 yield pcgSolver.scfByPCG(print, startTime, this);
2298 }
2299 };
2300 } catch (EnergyException ex) {
2301 logger.warning(ex.toString());
2302 if (sorFallback) {
2303
2304 SCFAlgorithm currentSCFAlgorithm = scfAlgorithm;
2305 scfAlgorithm = SOR;
2306 double currentPolarSOR = sorRegion.getPolarSOR();
2307
2308 sorRegion.setPolarSOR(sorFallBackPolarSOR);
2309 int ret = selfConsistentField(print);
2310
2311 scfAlgorithm = currentSCFAlgorithm;
2312 sorRegion.setPolarSOR(currentPolarSOR);
2313 logger.info(" SOR induced dipoles computed due to SCF failure.");
2314 return ret;
2315 } else if (directFallback) {
2316
2317 polarization = Polarization.DIRECT;
2318 int ret = selfConsistentField(print);
2319 polarization = Polarization.MUTUAL;
2320 logger.info(" Direct induced dipoles computed due to SCF failure.");
2321 return ret;
2322 } else {
2323 throw ex;
2324 }
2325 }
2326 }
2327
2328
2329
2330
2331 private int scfBySOR(boolean print, long startTime) {
2332 long directTime = System.nanoTime() - startTime;
2333
2334
2335 StringBuilder sb = null;
2336 if (print) {
2337 sb = new StringBuilder("\n Self-Consistent Field\n Iter RMS Change (Debye) Time\n");
2338 }
2339 int completedSCFCycles = 0;
2340 int maxSCFCycles = 1000;
2341 double eps = 100.0;
2342 double previousEps;
2343 boolean done = false;
2344 while (!done) {
2345 long cycleTime = -System.nanoTime();
2346 try {
2347 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
2348 reciprocalSpace.splineInducedDipoles(inducedDipole, inducedDipoleCR, use);
2349 }
2350 field.reset(parallelTeam);
2351 fieldCR.reset(parallelTeam);
2352 inducedDipoleFieldRegion.init(
2353 atoms, crystal, use, molecule, ipdamp, thole, coordinates,
2354 realSpaceNeighborParameters, inducedDipole, inducedDipoleCR,
2355 reciprocalSpaceTerm, reciprocalSpace, lambdaMode, ewaldParameters,
2356 field, fieldCR, pmeTimings);
2357 inducedDipoleFieldRegion.executeWith(sectionTeam);
2358 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
2359 reciprocalSpace.computeInducedPhi(
2360 cartesianInducedDipolePhi, cartesianInducedDipolePhiCR,
2361 fractionalInducedDipolePhi, fractionalInducedDipolePhiCR);
2362 }
2363
2364 if (generalizedKirkwoodTerm) {
2365
2366 pmeTimings.gkEnergyTotal = -System.nanoTime();
2367 generalizedKirkwood.computeInducedGKField();
2368 pmeTimings.gkEnergyTotal += System.nanoTime();
2369 logger.fine(
2370 format(" Computed GK induced field %8.3f (sec)", pmeTimings.gkEnergyTotal * 1.0e-9));
2371 }
2372
2373 sorRegion.init(atoms, polarizability, inducedDipole, inducedDipoleCR,
2374 directDipole, directDipoleCR, cartesianInducedDipolePhi, cartesianInducedDipolePhiCR,
2375 field, fieldCR, generalizedKirkwoodTerm, generalizedKirkwood, ewaldParameters);
2376 parallelTeam.execute(sorRegion);
2377
2378 expandInducedDipoles();
2379 } catch (Exception e) {
2380 String message = "Exception computing mutual induced dipoles.";
2381 logger.log(Level.SEVERE, message, e);
2382 }
2383 completedSCFCycles++;
2384 previousEps = eps;
2385 eps = sorRegion.getEps();
2386 eps = Constants.ELEC_ANG_TO_DEBYE * sqrt(eps / (double) nAtoms);
2387 cycleTime += System.nanoTime();
2388 if (print) {
2389 sb.append(format(" %4d %15.10f %7.4f\n", completedSCFCycles, eps, cycleTime * NS2SEC));
2390 }
2391
2392
2393 if (eps > previousEps) {
2394 if (sb != null) {
2395 logger.warning(sb.toString());
2396 }
2397 String message = format(" SCF convergence failure: (%10.5f > %10.5f)\n", eps, previousEps);
2398 throw new EnergyException(message);
2399 }
2400
2401
2402
2403 if (completedSCFCycles >= maxSCFCycles) {
2404 if (sb != null) {
2405 logger.warning(sb.toString());
2406 }
2407 String message = format(" Maximum SCF iterations reached: (%d)\n", completedSCFCycles);
2408 throw new EnergyException(message);
2409 }
2410
2411
2412 if (eps < poleps) {
2413 done = true;
2414 }
2415 }
2416 if (print) {
2417 sb.append(format(" Direct: %7.4f\n", NS2SEC * directTime));
2418 startTime = System.nanoTime() - startTime;
2419 sb.append(format(" Total: %7.4f", startTime * NS2SEC));
2420 logger.info(sb.toString());
2421 }
2422 return completedSCFCycles;
2423 }
2424
2425
2426
2427
2428 private int scfByEPT(boolean print, long startTime) {
2429
2430
2431 for (int i = 0; i < nAtoms; i++) {
2432 for (int j = 0; j < 3; j++) {
2433 optRegion.optDipole[0][i][j] = directDipole[i][j];
2434 optRegion.optDipoleCR[0][i][j] = directDipoleCR[i][j];
2435 }
2436 }
2437
2438 int optOrder = optRegion.optOrder;
2439
2440 for (int currentOptOrder = 1; currentOptOrder <= optOrder; currentOptOrder++) {
2441 try {
2442 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
2443 reciprocalSpace.splineInducedDipoles(inducedDipole, inducedDipoleCR, use);
2444 }
2445 field.reset(parallelTeam);
2446 fieldCR.reset(parallelTeam);
2447 inducedDipoleFieldRegion.init(
2448 atoms, crystal, use, molecule, ipdamp, thole, coordinates, realSpaceNeighborParameters,
2449 inducedDipole, inducedDipoleCR, reciprocalSpaceTerm, reciprocalSpace, lambdaMode,
2450 ewaldParameters, field, fieldCR, pmeTimings);
2451 inducedDipoleFieldRegion.executeWith(sectionTeam);
2452
2453 if (reciprocalSpaceTerm && ewaldParameters.aewald > 0.0) {
2454 reciprocalSpace.computeInducedPhi(
2455 cartesianInducedDipolePhi, cartesianInducedDipolePhiCR,
2456 fractionalInducedDipolePhi, fractionalInducedDipolePhiCR);
2457 }
2458
2459 if (generalizedKirkwoodTerm) {
2460
2461 pmeTimings.gkEnergyTotal = -System.nanoTime();
2462 generalizedKirkwood.computeInducedGKField();
2463 pmeTimings.gkEnergyTotal += System.nanoTime();
2464 logger.fine(
2465 format(" Computed GK induced field %8.3f (sec)", pmeTimings.gkEnergyTotal * 1.0e-9));
2466 }
2467
2468 optRegion.init(
2469 currentOptOrder,
2470 atoms,
2471 polarizability,
2472 inducedDipole,
2473 inducedDipoleCR,
2474 cartesianInducedDipolePhi,
2475 cartesianInducedDipolePhiCR,
2476 field,
2477 fieldCR,
2478 generalizedKirkwoodTerm,
2479 generalizedKirkwood,
2480 ewaldParameters,
2481 soluteDielectric);
2482 parallelTeam.execute(optRegion);
2483
2484 expandInducedDipoles();
2485
2486 } catch (Exception e) {
2487 String message = "Exception computing opt induced dipoles.";
2488 logger.log(Level.SEVERE, message, e);
2489 }
2490 }
2491
2492
2493 for (int i = 0; i < nAtoms; i++) {
2494 for (int j = 0; j < 3; j++) {
2495 inducedDipole[0][i][j] = 0.0;
2496 inducedDipoleCR[0][i][j] = 0.0;
2497 double sum = 0.0;
2498 double sump = 0.0;
2499 for (int k = 0; k <= optOrder; k++) {
2500 sum += optRegion.optDipole[k][i][j];
2501 sump += optRegion.optDipoleCR[k][i][j];
2502 inducedDipole[0][i][j] += optRegion.optCoefficients[k] * sum;
2503 inducedDipoleCR[0][i][j] += optRegion.optCoefficients[k] * sump;
2504 }
2505 }
2506 }
2507
2508
2509 expandInducedDipoles();
2510
2511 return optOrder;
2512 }
2513
2514
2515
2516
2517
2518
2519 public double getTotalMultipoleEnergy() {
2520 return permanentMultipoleEnergy + polarizationEnergy;
2521 }
2522
2523
2524
2525
2526
2527
2528 public void setPolarization(Polarization set) {
2529 this.polarization = set;
2530 }
2531
2532
2533
2534
2535 private void assignMultipoles() {
2536 if (forceField == null) {
2537 String message = "No force field is defined.\n";
2538 logger.log(Level.SEVERE, message);
2539 return;
2540 }
2541
2542 if (forceField.getForceFieldTypeCount(ForceFieldType.MULTIPOLE) < 1
2543 && forceField.getForceFieldTypeCount(ForceFieldType.CHARGE) < 1) {
2544 String message = "Force field has no permanent electrostatic types.\n";
2545 logger.log(Level.SEVERE, message);
2546 return;
2547 }
2548 if (nAtoms < 1) {
2549 String message = "No atoms are defined.\n";
2550 logger.log(Level.SEVERE, message);
2551 return;
2552 }
2553 for (int i = 0; i < nAtoms; i++) {
2554 Atom atom = atoms[i];
2555
2556 if (!assignMultipole(elecForm, atom, forceField, localMultipole[i], i, axisAtom, frame)) {
2557 logger.info(
2558 format("No MultipoleType could be assigned:\n %s --> %s", atom, atom.getAtomType()));
2559 StringBuilder sb = new StringBuilder();
2560 List<Bond> bonds = atom.getBonds();
2561 for (Bond bond12 : bonds) {
2562 Atom a2 = bond12.get1_2(atom);
2563 AtomType aType2 = a2.getAtomType();
2564 sb.append(format("\n 1-2 %s --> %s", a2, aType2));
2565 }
2566 for (Bond bond12 : bonds) {
2567 Atom atom2 = bond12.get1_2(atom);
2568 bonds = atom2.getBonds();
2569 for (Bond bond23 : bonds) {
2570 Atom a2 = bond23.get1_2(atom2);
2571 AtomType aType2 = a2.getAtomType();
2572 sb.append(format("\n 1-3 %s --> %s", a2, aType2));
2573 }
2574 }
2575
2576 List<MultipoleType> multipoleTypes = forceField.getMultipoleTypes(atom.getAtomType().getKey());
2577 if (multipoleTypes != null && !multipoleTypes.isEmpty()) {
2578 sb.append("\n Similar Multipole types:");
2579 for (MultipoleType multipoleType : multipoleTypes) {
2580 sb.append(format("\n %s", multipoleType));
2581 }
2582 }
2583
2584 logger.log(Level.SEVERE, sb.toString());
2585 }
2586 }
2587
2588
2589 StringBuilder sb = new StringBuilder();
2590 for (int i = 0; i < nAtoms; i++) {
2591 boolean flag = false;
2592 for (int j = 0; j < 10; j++) {
2593 if (Double.isNaN(localMultipole[i][j])) {
2594 flag = true;
2595 break;
2596 }
2597 }
2598 if (flag) {
2599 sb.append("\n").append(atoms[i].toString()).append("\n");
2600 sb.append(format("%d", i + 1));
2601 for (int j = 0; j < 10; j++) {
2602 sb.append(format(" %8.3f", localMultipole[i][j]));
2603 }
2604 sb.append("\n");
2605 }
2606 }
2607 if (!sb.isEmpty()) {
2608 String message = "Fatal exception: Error assigning multipoles. " + sb;
2609 logger.log(Level.SEVERE, message);
2610 }
2611 }
2612
2613
2614
2615
2616 protected void assignPolarizationGroups() {
2617
2618 List<Integer> group = new ArrayList<>();
2619 for (int i = 0; i < nAtoms; i++) {
2620 Atom a = atoms[i];
2621 if (a.getIndex() - 1 != i) {
2622 logger.info(format(" PME Index: %d Atom Index: %d\n Atom: %s",
2623 i, a.getIndex() - 1, a));
2624 logger.severe(" Atom indexing is not consistent in PME.");
2625 }
2626 }
2627 for (Atom ai : atoms) {
2628 group.clear();
2629 int index = ai.getIndex() - 1;
2630 group.add(index);
2631 PolarizeType polarizeType = ai.getPolarizeType();
2632 if (polarizeType != null) {
2633 if (polarizeType.polarizationGroup != null) {
2634 PolarizeType.growGroup(group, ai);
2635 sort(group);
2636 ip11[index] = new int[group.size()];
2637 int j = 0;
2638 for (int k : group) {
2639 ip11[index][j++] = k;
2640 }
2641 } else {
2642 ip11[index] = new int[group.size()];
2643 int j = 0;
2644 for (int k : group) {
2645 ip11[index][j++] = k;
2646 }
2647 }
2648 } else {
2649 String message = "The polarize keyword was not found for atom "
2650 + (index + 1) + " with type " + ai.getType();
2651 logger.severe(message);
2652 }
2653 }
2654
2655 int[] mask = new int[nAtoms];
2656 List<Integer> list = new ArrayList<>();
2657 List<Integer> keep = new ArrayList<>();
2658 for (int i = 0; i < nAtoms; i++) {
2659 mask[i] = -1;
2660 }
2661 for (int i = 0; i < nAtoms; i++) {
2662 list.clear();
2663 for (int j : ip11[i]) {
2664 list.add(j);
2665 mask[j] = i;
2666 }
2667 keep.clear();
2668 for (int j : list) {
2669 Atom aj = atoms[j];
2670 List<Bond> bonds = aj.getBonds();
2671 for (Bond b : bonds) {
2672 Atom ak = b.get1_2(aj);
2673 int k = ak.getIndex() - 1;
2674 if (mask[k] != i) {
2675 keep.add(k);
2676 }
2677 }
2678 }
2679 list.clear();
2680 for (int j : keep) {
2681 for (int k : ip11[j]) {
2682 list.add(k);
2683 }
2684 }
2685 sort(list);
2686 ip12[i] = new int[list.size()];
2687 int j = 0;
2688 for (int k : list) {
2689 ip12[i][j++] = k;
2690 }
2691 }
2692
2693
2694 for (int i = 0; i < nAtoms; i++) {
2695 mask[i] = -1;
2696 }
2697 for (int i = 0; i < nAtoms; i++) {
2698 for (int j : ip11[i]) {
2699 mask[j] = i;
2700 }
2701 for (int j : ip12[i]) {
2702 mask[j] = i;
2703 }
2704 list.clear();
2705 for (int j : ip12[i]) {
2706 for (int k : ip12[j]) {
2707 if (mask[k] != i) {
2708 if (!list.contains(k)) {
2709 list.add(k);
2710 }
2711 }
2712 }
2713 }
2714 ip13[i] = new int[list.size()];
2715 sort(list);
2716 int j = 0;
2717 for (int k : list) {
2718 ip13[i][j++] = k;
2719 }
2720 }
2721 }
2722
2723
2724
2725
2726
2727
2728
2729
2730 public double[] computeMoments(Atom[] activeAtoms, boolean forceEnergy) {
2731
2732 var netchg = 0.0;
2733 var netdpl = 0.0;
2734 var xdpl = 0.0;
2735 var ydpl = 0.0;
2736 var zdpl = 0.0;
2737 var xxqdp = 0.0;
2738 var xyqdp = 0.0;
2739 var xzqdp = 0.0;
2740 var yxqdp = 0.0;
2741 var yyqdp = 0.0;
2742 var yzqdp = 0.0;
2743 var zxqdp = 0.0;
2744 var zyqdp = 0.0;
2745 var zzqdp = 0.0;
2746
2747
2748 double xmid = 0.0;
2749 double ymid = 0.0;
2750 double zmid = 0.0;
2751 double totalMass = 0;
2752 for (Atom atom : activeAtoms) {
2753 var m = atom.getMass();
2754 totalMass += m;
2755 xmid = xmid + atom.getX() * m;
2756 ymid = ymid + atom.getY() * m;
2757 zmid = zmid + atom.getZ() * m;
2758 }
2759 if (totalMass > 0) {
2760 xmid /= totalMass;
2761 ymid /= totalMass;
2762 zmid /= totalMass;
2763 }
2764 int n = activeAtoms.length;
2765 double[] xcm = new double[n];
2766 double[] ycm = new double[n];
2767 double[] zcm = new double[n];
2768 int k = 0;
2769 for (Atom atom : activeAtoms) {
2770 xcm[k] = atom.getX() - xmid;
2771 ycm[k] = atom.getY() - ymid;
2772 zcm[k] = atom.getZ() - zmid;
2773 k++;
2774 }
2775
2776 if (forceEnergy) {
2777 energy(false, false);
2778 }
2779
2780
2781 k = 0;
2782 for (Atom atom : activeAtoms) {
2783 int i = atom.getIndex() - 1;
2784 double[] globalMultipolei = globalMultipole[0][i];
2785 double[] inducedDipolei = inducedDipole[0][i];
2786
2787 var ci = globalMultipolei[t000];
2788 var dix = globalMultipolei[t100];
2789 var diy = globalMultipolei[t010];
2790 var diz = globalMultipolei[t001];
2791 var uix = inducedDipolei[0];
2792 var uiy = inducedDipolei[1];
2793 var uiz = inducedDipolei[2];
2794
2795 netchg += ci;
2796 xdpl += xcm[k] * ci + dix + uix;
2797 ydpl += ycm[k] * ci + diy + uiy;
2798 zdpl += zcm[k] * ci + diz + uiz;
2799 xxqdp += xcm[k] * xcm[k] * ci + 2.0 * xcm[k] * (dix + uix);
2800 xyqdp += xcm[k] * ycm[k] * ci + xcm[k] * (diy + uiy) + ycm[k] * (dix + uix);
2801 xzqdp += xcm[k] * zcm[k] * ci + xcm[k] * (diz + uiz) + zcm[k] * (dix + uix);
2802 yxqdp += ycm[k] * xcm[k] * ci + ycm[k] * (dix + uix) + xcm[k] * (diy + uiy);
2803 yyqdp += ycm[k] * ycm[k] * ci + 2.0 * ycm[k] * (diy + uiy);
2804 yzqdp += ycm[k] * zcm[k] * ci + ycm[k] * (diz + uiz) + zcm[k] * (diy + uiy);
2805 zxqdp += zcm[k] * xcm[k] * ci + zcm[k] * (dix + uix) + xcm[k] * (diz + uiz);
2806 zyqdp += zcm[k] * ycm[k] * ci + zcm[k] * (diy + uiy) + ycm[k] * (diz + uiz);
2807 zzqdp += zcm[k] * zcm[k] * ci + 2.0 * zcm[k] * (diz + uiz);
2808 k++;
2809 }
2810
2811
2812 var qave = (xxqdp + yyqdp + zzqdp) / 3.0;
2813 xxqdp = 1.5 * (xxqdp - qave);
2814 xyqdp = 1.5 * xyqdp;
2815 xzqdp = 1.5 * xzqdp;
2816 yxqdp = 1.5 * yxqdp;
2817 yyqdp = 1.5 * (yyqdp - qave);
2818 yzqdp = 1.5 * yzqdp;
2819 zxqdp = 1.5 * zxqdp;
2820 zyqdp = 1.5 * zyqdp;
2821 zzqdp = 1.5 * (zzqdp - qave);
2822
2823
2824 for (Atom atom : activeAtoms) {
2825 int i = atom.getIndex() - 1;
2826 double[] globalMultipolei = globalMultipole[0][i];
2827 var qixx = globalMultipolei[t200];
2828 var qiyy = globalMultipolei[t020];
2829 var qizz = globalMultipolei[t002];
2830 var qixy = globalMultipolei[t110];
2831 var qixz = globalMultipolei[t101];
2832 var qiyz = globalMultipolei[t011];
2833 xxqdp += qixx;
2834 xyqdp += qixy;
2835 xzqdp += qixz;
2836 yxqdp += qixy;
2837 yyqdp += qiyy;
2838 yzqdp += qiyz;
2839 zxqdp += qixz;
2840 zyqdp += qiyz;
2841 zzqdp += qizz;
2842 }
2843
2844
2845 xdpl = xdpl * ELEC_ANG_TO_DEBYE;
2846 ydpl = ydpl * ELEC_ANG_TO_DEBYE;
2847 zdpl = zdpl * ELEC_ANG_TO_DEBYE;
2848 xxqdp = xxqdp * ELEC_ANG_TO_DEBYE;
2849 xyqdp = xyqdp * ELEC_ANG_TO_DEBYE;
2850 xzqdp = xzqdp * ELEC_ANG_TO_DEBYE;
2851 yxqdp = yxqdp * ELEC_ANG_TO_DEBYE;
2852 yyqdp = yyqdp * ELEC_ANG_TO_DEBYE;
2853 yzqdp = yzqdp * ELEC_ANG_TO_DEBYE;
2854 zxqdp = zxqdp * ELEC_ANG_TO_DEBYE;
2855 zyqdp = zyqdp * ELEC_ANG_TO_DEBYE;
2856 zzqdp = zzqdp * ELEC_ANG_TO_DEBYE;
2857
2858
2859 netdpl = sqrt(xdpl * xdpl + ydpl * ydpl + zdpl * zdpl);
2860 double[][] a = new double[3][3];
2861 a[0][0] = xxqdp;
2862 a[0][1] = xyqdp;
2863 a[0][2] = xzqdp;
2864 a[1][0] = yxqdp;
2865 a[1][1] = yyqdp;
2866 a[1][2] = yzqdp;
2867 a[2][0] = zxqdp;
2868 a[2][1] = zyqdp;
2869 a[2][2] = zzqdp;
2870 EigenDecomposition e = new EigenDecomposition(new Array2DRowRealMatrix(a));
2871
2872 var netqdp = e.getRealEigenvalues();
2873
2874
2875 logger.info("\n Electric Moments\n");
2876 logger.info(format(" Total Electric Charge: %13.5f Electrons\n", netchg));
2877 logger.info(format(" Dipole Moment Magnitude: %13.5f Debye\n", netdpl));
2878 logger.info(format(" Dipole X,Y,Z-Components: %13.5f %13.5f %13.5f\n", xdpl, ydpl, zdpl));
2879 logger.info(format(" Quadrupole Moment Tensor: %13.5f %13.5f %13.5f", xxqdp, xyqdp, xzqdp));
2880 logger.info(format(" (Buckinghams) %13.5f %13.5f %13.5f", yxqdp, yyqdp, yzqdp));
2881 logger.info(format(" %13.5f %13.5f %13.5f\n", zxqdp, zyqdp, zzqdp));
2882 logger.info(
2883 format(
2884 " Principal Axes Quadrupole %13.5f %13.5f %13.5f\n", netqdp[2], netqdp[1], netqdp[0]));
2885
2886
2887 return new double[] {
2888 netchg, netdpl, xdpl, ydpl, zdpl, xxqdp, xyqdp, xzqdp, yxqdp, yyqdp, yzqdp, zxqdp, zyqdp, zzqdp
2889 };
2890 }
2891
2892
2893
2894
2895
2896
2897
2898
2899
2900
2901 private void log(int i, int k, double r, double eij) {
2902 logger.info(
2903 format(
2904 "%s %6d-%s %6d-%s %10.4f %10.4f",
2905 "ELEC",
2906 atoms[i].getIndex(),
2907 atoms[i].getAtomType().name,
2908 atoms[k].getIndex(),
2909 atoms[k].getAtomType().name,
2910 r,
2911 eij));
2912 }
2913
2914
2915
2916
2917
2918 private boolean esvTerm = false;
2919
2920
2921
2922
2923 private ExtendedSystem extendedSystem = null;
2924
2925
2926
2927
2928
2929
2930 public void attachExtendedSystem(ExtendedSystem system) {
2931
2932 esvTerm = true;
2933 extendedSystem = system;
2934
2935 setAtoms(extendedSystem.getExtendedAtoms(), extendedSystem.getExtendedMolecule());
2936
2937 if (dMultipoledTirationESV == null || dMultipoledTirationESV.length != nSymm
2938 || dMultipoledTirationESV[0].length != nAtoms) {
2939 dMultipoledTirationESV = new double[nSymm][nAtoms][10];
2940 dMultipoledTautomerESV = new double[nSymm][nAtoms][10];
2941 }
2942
2943 if (dPolardTitrationESV == null || dPolardTitrationESV.length != nAtoms) {
2944 dPolardTitrationESV = new double[nAtoms];
2945 dPolardTautomerESV = new double[nAtoms];
2946 }
2947
2948 logger.info(format(" Attached extended system (%d variables) to PME.\n", extendedSystem.getNumberOfVariables()));
2949 }
2950
2951
2952
2953
2954 double[] perAtomTitrationESV = null;
2955
2956
2957
2958
2959
2960 double[][][] dMultipoledTirationESV = null;
2961 double[][][] dMultipoledTautomerESV = null;
2962
2963
2964
2965
2966
2967 double[] dPolardTitrationESV = null;
2968 double[] dPolardTautomerESV = null;
2969
2970
2971
2972 LambdaFactors[] lambdaFactors = null;
2973
2974
2975
2976
2977
2978
2979
2980
2981
2982
2983
2984
2985 public class LambdaFactors {
2986
2987
2988
2989
2990 protected double lambdaProd = 1.0;
2991
2992
2993
2994 protected double lfAlpha = 0.0;
2995
2996
2997
2998 protected double dlfAlpha = 0.0;
2999
3000
3001
3002 protected double d2lfAlpha = 0.0;
3003
3004
3005
3006 protected double lfPowPerm = 1.0;
3007
3008
3009
3010 protected double dlfPowPerm = 0.0;
3011
3012
3013
3014 protected double d2lfPowPerm = 0.0;
3015
3016
3017
3018 protected double lfPowPol = 1.0;
3019
3020
3021
3022 protected double dlfPowPol = 0.0;
3023
3024
3025
3026 protected double d2lfPowPol = 0.0;
3027
3028
3029
3030 protected double dLpdL = 1.0;
3031
3032
3033
3034 protected double dLpdLi = 1.0;
3035
3036
3037
3038 protected double dLpdLk = 1.0;
3039
3040
3041
3042 protected int[] ik = new int[2];
3043
3044 public void print() {
3045 StringBuilder sb = new StringBuilder();
3046 if (this instanceof LambdaFactorsESV) {
3047 sb.append(format(" (QI-ESV) i,k:%d,%d", ik[0], ik[1]));
3048 } else {
3049 sb.append(" (QI-OST)");
3050 }
3051 sb.append(format(
3052 " lambda:%.2f lAlpha:%.2f,%.2f,%.2f lPowPerm:%.2f,%.2f,%.2f lPowPol:%.2f,%.2f,%.2f",
3053 lambdaProd,
3054 lfAlpha,
3055 dlfAlpha,
3056 d2lfAlpha,
3057 lfPowPerm,
3058 dlfPowPerm,
3059 d2lfPowPerm,
3060 lfPowPol,
3061 dlfPowPol,
3062 d2lfPowPol));
3063
3064 sb.append(format(
3065 "\n permExp:%.2f permAlpha:%.2f permWindow:%.2f,%.2f polExp:%.2f polWindow:%.2f,%.2f",
3066 alchemicalParameters.permLambdaExponent,
3067 alchemicalParameters.permLambdaAlpha,
3068 alchemicalParameters.permLambdaStart,
3069 alchemicalParameters.permLambdaEnd,
3070 alchemicalParameters.polLambdaExponent,
3071 alchemicalParameters.polLambdaStart,
3072 alchemicalParameters.polLambdaEnd));
3073 logger.info(sb.toString());
3074 }
3075
3076
3077
3078
3079
3080
3081
3082
3083 public void setFactors(int i, int k, LambdaMode mode) {
3084
3085 }
3086
3087
3088
3089
3090 public void setFactors() {
3091
3092 }
3093 }
3094
3095 public class LambdaFactorsOST extends LambdaFactors {
3096
3097 @Override
3098 public void setFactors() {
3099 lambdaProd = lambda;
3100 lfAlpha = alchemicalParameters.lAlpha;
3101 dlfAlpha = alchemicalParameters.dlAlpha;
3102 d2lfAlpha = alchemicalParameters.d2lAlpha;
3103 lfPowPerm = alchemicalParameters.permanentScale;
3104 dlfPowPerm =
3105 alchemicalParameters.dlPowPerm * alchemicalParameters.dEdLSign;
3106 d2lfPowPerm = alchemicalParameters.d2lPowPerm
3107 * alchemicalParameters.dEdLSign;
3108 lfPowPol = alchemicalParameters.polarizationScale;
3109 dlfPowPol =
3110 alchemicalParameters.dlPowPol * alchemicalParameters.dEdLSign;
3111 d2lfPowPol =
3112 alchemicalParameters.d2lPowPol * alchemicalParameters.dEdLSign;
3113 }
3114 }
3115
3116 public class LambdaFactorsESV extends LambdaFactors {
3117
3118 @Override
3119 public void setFactors(int i, int k, LambdaMode mode) {
3120 logger.info(format("Invoked Qi setFactors() method with i,k=%d,%d", i, k));
3121 ik[0] = i;
3122 ik[1] = k;
3123 final double L = lambda;
3124 lambdaProd = L * perAtomTitrationESV[i] * perAtomTitrationESV[k];
3125 final double esvli = perAtomTitrationESV[i];
3126 final double esvlk = perAtomTitrationESV[k];
3127 dLpdL = esvli * esvlk;
3128 dLpdLi = L * esvlk;
3129 dLpdLk = L * esvli;
3130
3131
3132 double permLambdaExponent = alchemicalParameters.permLambdaExponent;
3133 double permLambdaStart = alchemicalParameters.permLambdaStart;
3134 double permLambdaEnd = alchemicalParameters.permLambdaEnd;
3135
3136 double permWindow = 1.0 / (permLambdaEnd - permLambdaStart);
3137 double permLambda = (lambdaProd - permLambdaStart) * permWindow;
3138 lfPowPerm = pow(permLambda, alchemicalParameters.permLambdaExponent);
3139 dlfPowPerm = (permLambdaExponent < 1)
3140 ? 0.0 : permLambdaExponent * pow(permLambda, permLambdaExponent - 1) * permWindow;
3141 d2lfPowPerm = (permLambdaExponent < 2)
3142 ? 0.0 : permLambdaExponent
3143 * (permLambdaExponent - 1)
3144 * pow(permLambda, permLambdaExponent - 2)
3145 * permWindow
3146 * permWindow;
3147
3148
3149 double polLambdaExponent = alchemicalParameters.polLambdaExponent;
3150 double polLambdaStart = alchemicalParameters.polLambdaStart;
3151 double polLambdaEnd = alchemicalParameters.polLambdaEnd;
3152
3153 double polWindow = 1.0 / (polLambdaEnd - polLambdaStart);
3154 double polLambda = (lambdaProd - polLambdaStart) * polWindow;
3155 lfPowPol = pow(polLambda, polLambdaExponent);
3156 dlfPowPol = (polLambdaExponent < 1)
3157 ? 0.0 : polLambdaExponent * pow(polLambda, polLambdaExponent - 1) * polWindow;
3158 d2lfPowPol = (polLambdaExponent < 2)
3159 ? 0.0 : polLambdaExponent
3160 * (polLambdaExponent - 1)
3161 * pow(polLambda, polLambdaExponent - 2)
3162 * polWindow
3163 * polWindow;
3164
3165
3166 double permLambdaAlpha = alchemicalParameters.permLambdaAlpha;
3167 lfAlpha = permLambdaAlpha * (1.0 - permLambda) * (1.0 - permLambda);
3168 dlfAlpha = permLambdaAlpha * (1.0 - permLambda) * permWindow;
3169 d2lfAlpha = -permLambdaAlpha * permWindow * permWindow;
3170
3171
3172
3173
3174
3175
3176
3177
3178
3179
3180
3181 switch (mode) {
3182 case CONDENSED -> {
3183 lfPowPerm = 1.0 - lfPowPerm;
3184 dlfPowPerm = -dlfPowPerm;
3185 d2lfPowPerm = -d2lfPowPerm;
3186 if (polarization == Polarization.NONE || lambda < polLambdaStart) {
3187 lfPowPol = 0.0;
3188 dlfPowPol = 0.0;
3189 d2lfPowPol = 0.0;
3190 } else if (lambda <= polLambdaEnd) {
3191
3192 } else {
3193 lfPowPol = 1.0;
3194 dlfPowPol = 0.0;
3195 d2lfPowPol = 0.0;
3196 }
3197 }
3198 case CONDENSED_NO_LIGAND -> {
3199
3200 lfPowPerm = 1.0 - lfPowPerm;
3201 dlfPowPerm = -dlfPowPerm;
3202 d2lfPowPerm = -d2lfPowPerm;
3203 if (polarization == Polarization.NONE) {
3204 lfPowPol = 0.0;
3205 dlfPowPol = 0.0;
3206 d2lfPowPol = 0.0;
3207 } else if (lambda <= polLambdaEnd) {
3208 lfPowPol = 1.0 - lfPowPol;
3209 dlfPowPol = -dlfPowPol;
3210 d2lfPowPol = -d2lfPowPol;
3211 } else {
3212 lfPowPol = 0.0;
3213 dlfPowPol = 0.0;
3214 d2lfPowPol = 0.0;
3215 }
3216 }
3217 case VAPOR -> {
3218
3219 lfPowPerm = 1.0 - lfPowPerm;
3220 dlfPowPerm = -dlfPowPerm;
3221 d2lfPowPerm = -d2lfPowPerm;
3222 lfAlpha = 0.0;
3223 dlfAlpha = 0.0;
3224 d2lfAlpha = 0.0;
3225 if (polarization == Polarization.NONE || lambdaProd > polLambdaEnd) {
3226 lfPowPol = 0.0;
3227 dlfPowPol = 0.0;
3228 d2lfPowPol = 0.0;
3229 } else if (lambdaProd <= polLambdaEnd) {
3230 lfPowPol = 1.0 - lfPowPol;
3231 dlfPowPol = -dlfPowPol;
3232 d2lfPowPol = -d2lfPowPol;
3233 }
3234 }
3235 default -> {
3236 }
3237 }
3238 }
3239 }
3240
3241
3242
3243
3244
3245 public final LambdaFactors LambdaDefaults = new LambdaFactors();
3246 }