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.algorithms.thermodynamics;
39
40 import edu.rit.mp.DoubleBuf;
41 import edu.rit.pj.Comm;
42 import ffx.algorithms.AlgorithmListener;
43 import ffx.algorithms.cli.DynamicsOptions;
44 import ffx.algorithms.cli.LambdaParticleOptions;
45 import ffx.algorithms.dynamics.Barostat;
46 import ffx.algorithms.dynamics.integrators.Stochastic;
47 import ffx.algorithms.optimize.Minimize;
48 import ffx.crystal.Crystal;
49 import ffx.crystal.CrystalPotential;
50 import ffx.numerics.Potential;
51 import ffx.numerics.integrate.DataSet;
52 import ffx.numerics.integrate.DoublesDataSet;
53 import ffx.numerics.integrate.Integrate1DNumeric;
54 import ffx.numerics.integrate.Integrate1DNumeric.IntegrationType;
55 import ffx.potential.MolecularAssembly;
56 import ffx.potential.SystemState;
57 import ffx.potential.Utilities;
58 import ffx.potential.bonded.LambdaInterface;
59 import ffx.potential.parsers.PDBFilter;
60 import ffx.potential.parsers.SystemFilter;
61 import ffx.potential.parsers.XYZFilter;
62 import org.apache.commons.configuration2.CompositeConfiguration;
63 import org.apache.commons.io.FilenameUtils;
64
65 import javax.annotation.Nullable;
66 import java.io.File;
67 import java.io.PrintWriter;
68 import java.util.ArrayList;
69 import java.util.Arrays;
70 import java.util.List;
71 import java.util.Optional;
72 import java.util.logging.Level;
73 import java.util.logging.Logger;
74
75 import static ffx.numerics.integrate.Integrate1DNumeric.IntegrationType.SIMPSONS;
76 import static ffx.utilities.Constants.R;
77 import static java.lang.String.format;
78 import static java.lang.System.arraycopy;
79 import static java.util.Arrays.fill;
80 import static org.apache.commons.math3.util.FastMath.PI;
81 import static org.apache.commons.math3.util.FastMath.abs;
82 import static org.apache.commons.math3.util.FastMath.asin;
83 import static org.apache.commons.math3.util.FastMath.exp;
84 import static org.apache.commons.math3.util.FastMath.floor;
85 import static org.apache.commons.math3.util.FastMath.sin;
86 import static org.apache.commons.math3.util.FastMath.sqrt;
87
88
89
90
91
92
93
94
95
96
97
98 public class OrthogonalSpaceTempering implements CrystalPotential, LambdaInterface {
99
100 private static final Logger logger = Logger.getLogger(OrthogonalSpaceTempering.class.getName());
101
102
103
104 protected final CrystalPotential potential;
105
106
107
108 protected final Barostat barostat;
109
110
111
112 protected final AlgorithmListener algorithmListener;
113
114
115
116 protected final boolean print = false;
117
118
119
120 protected final int nVariables;
121
122
123
124 private final LambdaInterface lambdaInterface;
125
126
127
128 private final List<Histogram> allHistograms = new ArrayList<>();
129
130
131
132 private final OptimizationParameters optimizationParameters;
133
134
135
136 private final CompositeConfiguration properties;
137
138
139
140 protected MolecularAssembly molecularAssembly;
141
142
143
144
145 protected Potential.STATE state = Potential.STATE.BOTH;
146
147
148
149 protected double forceFieldEnergy;
150
151
152
153 private Histogram histogram;
154
155
156
157 private int histogramIndex;
158
159
160
161 private boolean propagateLambda = true;
162
163
164
165 private final double[] dUdXdL;
166
167
168
169
170 private final double[] tempGradient;
171
172
173
174
175 private double dForceFieldEnergydL;
176
177
178
179 private double gLdEdL = 0.0;
180
181
182
183 private double biasEnergy;
184
185
186
187 private double totalEnergy;
188
189
190
191 private double dUdLambda;
192
193
194
195 private double d2UdL2;
196
197
198
199 private boolean hardWallConstraint = false;
200
201 private final DynamicsOptions dynamicsOptions;
202
203 private final LambdaParticleOptions lambdaParticleOptions;
204
205
206
207
208
209
210
211
212
213
214
215
216
217 public OrthogonalSpaceTempering(LambdaInterface lambdaInterface, CrystalPotential potential,
218 HistogramData histogramData, LambdaData lambdaData, CompositeConfiguration properties,
219 DynamicsOptions dynamicsOptions, LambdaParticleOptions lambdaParticleOptions,
220 AlgorithmListener algorithmListener) {
221 this.lambdaInterface = lambdaInterface;
222 this.potential = potential;
223 this.properties = properties;
224 this.dynamicsOptions = dynamicsOptions;
225 this.lambdaParticleOptions = lambdaParticleOptions;
226 this.algorithmListener = algorithmListener;
227 nVariables = potential.getNumberOfVariables();
228
229 if (potential instanceof Barostat) {
230 barostat = (Barostat) potential;
231 } else {
232 barostat = null;
233 }
234
235 dUdXdL = new double[nVariables];
236 tempGradient = new double[nVariables];
237
238
239 histogram = new Histogram(properties, histogramData, lambdaData);
240 histogramIndex = 0;
241 allHistograms.add(histogram);
242
243
244 optimizationParameters = new OptimizationParameters(properties);
245 }
246
247
248
249
250
251
252 public void addHistogram(HistogramData histogramData, LambdaData lambdaData) {
253 Histogram newHisto = new Histogram(properties, histogramData, lambdaData);
254 histogramData.asynchronous = this.histogram.hd.asynchronous;
255 allHistograms.add(newHisto);
256 }
257
258
259
260
261 @Override
262 public boolean dEdLZeroAtEnds() {
263 return false;
264 }
265
266
267
268
269 @Override
270 public boolean destroy() {
271
272 histogram.destroy();
273 return potential.destroy();
274 }
275
276
277
278
279 public double energy(double[] x) {
280
281
282 if (state == Potential.STATE.FAST) {
283 forceFieldEnergy = potential.energy(x);
284 return forceFieldEnergy;
285 }
286
287
288 if (propagateLambda) {
289 histogram.stochasticPreForce();
290 }
291
292
293 fill(tempGradient, 0.0);
294 forceFieldEnergy = potential.energyAndGradient(x, tempGradient);
295
296 dForceFieldEnergydL = lambdaInterface.getdEdL();
297 d2UdL2 = lambdaInterface.getd2EdL2();
298 int lambdaBin = histogram.indexForLambda();
299 dUdLambda = dForceFieldEnergydL;
300
301 gLdEdL = 0.0;
302 double bias1D;
303
304 histogram.updateFreeEnergyDifference(false, false);
305 if (histogram.hd.metaDynamics) {
306 bias1D = histogram.energyAndGradientMeta(true);
307 } else {
308
309 if (histogram.hd.biasMag > 0.0) {
310 double[] chainRule = new double[2];
311 gLdEdL = histogram.energyAndGradient2D(dUdLambda, chainRule);
312 double dGdLambda = chainRule[0];
313 double dGdFLambda = chainRule[1];
314 dUdLambda += dGdLambda + dGdFLambda * d2UdL2;
315 }
316
317
318 bias1D = histogram.energyAndGradient1D(true);
319 }
320
321
322 biasEnergy = bias1D + gLdEdL;
323
324 if (print) {
325 logger.info(format(" Bias Energy %16.8f", biasEnergy));
326 logger.info(format(" %s %16.8f (Kcal/mole)", "OST Potential ", forceFieldEnergy + biasEnergy));
327 }
328
329 if (propagateLambda) {
330 histogram.stochasticPostForce();
331
332 histogram.ld.stepsTaken++;
333 long energyCount = histogram.ld.stepsTaken;
334
335
336 int printFrequency = dynamicsOptions.getReportFrequency(100);
337 if (energyCount % printFrequency == 0) {
338 double dBdL = dUdLambda - dForceFieldEnergydL;
339 double lambda = histogram.ld.lambda;
340 int lambdaBins = histogram.hd.getLambdaBins();
341 if (lambdaBins < 1000) {
342 logger.info(format(" L=%6.4f (%3d) F_LU=%10.4f F_LB=%10.4f F_L=%10.4f V_L=%10.4f", lambda,
343 lambdaBin, dForceFieldEnergydL, dBdL, dUdLambda, histogram.ld.thetaVelocity));
344 } else {
345 logger.info(format(" L=%6.4f (%4d) F_LU=%10.4f F_LB=%10.4f F_L=%10.4f V_L=%10.4f", lambda,
346 lambdaBin, dForceFieldEnergydL, dBdL, dUdLambda, histogram.ld.thetaVelocity));
347 }
348 }
349
350
351 if (energyCount % histogram.hd.countInterval == 0) {
352 histogram.addBias(dForceFieldEnergydL);
353
354
355 if (optimizationParameters.doOptimization) {
356 optimizationParameters.optimize(forceFieldEnergy, x, tempGradient);
357 }
358
359 if (algorithmListener != null) {
360 algorithmListener.algorithmUpdate(molecularAssembly);
361 }
362 }
363 }
364
365 totalEnergy = forceFieldEnergy + biasEnergy;
366
367 return totalEnergy;
368 }
369
370
371
372
373 @Override
374 public double energyAndGradient(double[] x, double[] gradient) {
375
376
377 if (state != STATE.FAST && propagateLambda) {
378 histogram.stochasticPreForce();
379 }
380
381
382 forceFieldEnergy = potential.energyAndGradient(x, gradient);
383
384
385 if (state == STATE.FAST) {
386 return forceFieldEnergy;
387 }
388
389
390 dForceFieldEnergydL = lambdaInterface.getdEdL();
391 dUdLambda = dForceFieldEnergydL;
392 d2UdL2 = lambdaInterface.getd2EdL2();
393
394 gLdEdL = 0.0;
395 double bias1D;
396
397 histogram.updateFreeEnergyDifference(false, false);
398 if (histogram.hd.metaDynamics) {
399 bias1D = histogram.energyAndGradientMeta(true);
400 } else {
401 if (histogram.hd.biasMag > 0) {
402 double[] chainRule = new double[2];
403 gLdEdL = histogram.energyAndGradient2D(dUdLambda, chainRule);
404 double dGdLambda = chainRule[0];
405 double dGdFLambda = chainRule[1];
406
407
408 dUdLambda += dGdLambda + dGdFLambda * d2UdL2;
409
410
411 fill(dUdXdL, 0.0);
412 lambdaInterface.getdEdXdL(dUdXdL);
413 for (int i = 0; i < nVariables; i++) {
414 gradient[i] += dGdFLambda * dUdXdL[i];
415 }
416 }
417
418
419 bias1D = histogram.energyAndGradient1D(true);
420 }
421
422
423 biasEnergy = bias1D + gLdEdL;
424
425 if (print) {
426 logger.info(format(" %s %16.8f", "Bias Energy ", biasEnergy));
427 logger.info(format(" %s %16.8f %s", "OST Potential ", forceFieldEnergy + biasEnergy, "(Kcal/mole)"));
428 }
429
430 if (propagateLambda) {
431 histogram.stochasticPostForce();
432
433 histogram.ld.stepsTaken++;
434 long energyCount = histogram.ld.stepsTaken;
435
436
437 int printFrequency = dynamicsOptions.getReportFrequency(100);
438 if (energyCount % printFrequency == 0) {
439 double dBdL = dUdLambda - dForceFieldEnergydL;
440 int lambdaBin = histogram.indexForLambda();
441 double lambda = histogram.ld.lambda;
442 int lambdaBins = histogram.hd.getLambdaBins();
443 if (lambdaBins < 1000) {
444 logger.info(format(" L=%6.4f (%3d) F_LU=%10.4f F_LB=%10.4f F_L=%10.4f V_L=%10.4f", lambda,
445 lambdaBin, dForceFieldEnergydL, dBdL, dUdLambda, histogram.ld.thetaVelocity));
446 } else {
447 logger.info(format(" L=%6.4f (%4d) F_LU=%10.4f F_LB=%10.4f F_L=%10.4f V_L=%10.4f", lambda,
448 lambdaBin, dForceFieldEnergydL, dBdL, dUdLambda, histogram.ld.thetaVelocity));
449 }
450 }
451
452
453 if (energyCount % histogram.hd.countInterval == 0) {
454 histogram.addBias(dForceFieldEnergydL);
455
456
457 if (optimizationParameters.doOptimization) {
458 optimizationParameters.optimize(forceFieldEnergy, x, gradient);
459 }
460
461 if (algorithmListener != null) {
462 algorithmListener.algorithmUpdate(molecularAssembly);
463 }
464 }
465 }
466
467 totalEnergy = forceFieldEnergy + biasEnergy;
468
469 return totalEnergy;
470 }
471
472
473
474
475 @Override
476 public double[] getAcceleration(double[] acceleration) {
477 return potential.getAcceleration(acceleration);
478 }
479
480 public Histogram[] getAllHistograms() {
481 int nHisto = allHistograms.size();
482 Histogram[] ret = new Histogram[nHisto];
483 ret = allHistograms.toArray(ret);
484 return ret;
485 }
486
487
488
489
490 @Override
491 public double[] getCoordinates(double[] doubles) {
492 return potential.getCoordinates(doubles);
493 }
494
495
496
497
498 @Override
499 public void setCoordinates(double[] doubles) {
500 potential.setCoordinates(doubles);
501 }
502
503
504
505
506 @Override
507 public Crystal getCrystal() {
508 return potential.getCrystal();
509 }
510
511
512
513
514 @Override
515 public void setCrystal(Crystal crystal) {
516 potential.setCrystal(crystal);
517 }
518
519
520
521
522
523
524
525 public long getEnergyCount() {
526 return histogram.ld.stepsTaken;
527 }
528
529
530
531
532 @Override
533 public Potential.STATE getEnergyTermState() {
534 return state;
535 }
536
537
538
539
540 @Override
541 public void setEnergyTermState(Potential.STATE state) {
542 this.state = state;
543 potential.setEnergyTermState(state);
544 }
545
546
547
548
549
550
551 public double getForceFieldEnergy() {
552 return forceFieldEnergy;
553 }
554
555
556
557
558
559
560 public Histogram getHistogram() {
561 return histogram;
562 }
563
564
565
566
567
568
569 public double getLambda() {
570 return histogram.ld.lambda;
571 }
572
573
574
575
576
577
578 public void setLambda(double lambda) {
579 if (histogram != null) {
580 lambda = histogram.mapLambda(lambda);
581 } else {
582 logger.warning(" OrthogonalSpaceTempering.setLambda was called before histogram constructed!");
583 logger.info(Utilities.stackTraceToString(new RuntimeException()));
584 }
585 lambdaInterface.setLambda(lambda);
586 histogram.ld.lambda = lambda;
587 histogram.ld.theta = asin(sqrt(lambda));
588 }
589
590
591
592
593 @Override
594 public double[] getMass() {
595 return potential.getMass();
596 }
597
598
599
600
601 @Override
602 public int getNumberOfVariables() {
603 return potential.getNumberOfVariables();
604 }
605
606
607
608
609
610
611 public OptimizationParameters getOptimizationParameters() {
612 return optimizationParameters;
613 }
614
615
616
617
618
619
620 public Potential getPotentialEnergy() {
621 return potential;
622 }
623
624
625
626
627 @Override
628 public double[] getPreviousAcceleration(double[] previousAcceleration) {
629 return potential.getPreviousAcceleration(previousAcceleration);
630 }
631
632
633
634
635 @Override
636 public double[] getScaling() {
637 return potential.getScaling();
638 }
639
640
641
642
643 @Override
644 public void setScaling(double[] scaling) {
645 potential.setScaling(scaling);
646 }
647
648
649
650
651 @Override
652 public double getTotalEnergy() {
653 return totalEnergy;
654 }
655
656
657
658
659
660
661 public double getTotaldEdLambda() {
662 return dUdLambda;
663 }
664
665
666
667
668
669
670 @Override
671 public Potential.VARIABLE_TYPE[] getVariableTypes() {
672 return potential.getVariableTypes();
673 }
674
675
676
677
678 @Override
679 public double[] getVelocity(double[] velocity) {
680 return potential.getVelocity(velocity);
681 }
682
683
684
685
686 @Override
687 public double getd2EdL2() {
688 throw new UnsupportedOperationException(
689 " Second derivatives of the bias are not implemented, as they require third derivatives of the potential.");
690 }
691
692
693
694
695 @Override
696 public double getdEdL() {
697 return getTotaldEdLambda();
698 }
699
700
701
702
703 @Override
704 public void getdEdXdL(double[] gradient) {
705 throw new UnsupportedOperationException(
706 " Second derivatives of the bias are not implemented, as they require third derivatives of the potential.");
707 }
708
709
710 public void logOutputFiles(int index) {
711 logger.info(format(" OST: Lambda file %s, histogram %s", histogram.ld.getLambdaFile(),
712 allHistograms.get(index).hd.getHistogramFile()));
713 }
714
715
716
717
718 @Override
719 public void setAcceleration(double[] acceleration) {
720 potential.setAcceleration(acceleration);
721 }
722
723
724
725
726
727
728
729 public void setHardWallConstraint(boolean hardWallConstraint) {
730 this.hardWallConstraint = hardWallConstraint;
731 }
732
733 public void setMolecularAssembly(MolecularAssembly molecularAssembly) {
734 this.molecularAssembly = molecularAssembly;
735 }
736
737
738
739
740 @Override
741 public void setPreviousAcceleration(double[] previousAcceleration) {
742 potential.setPreviousAcceleration(previousAcceleration);
743 }
744
745
746
747
748
749
750 public void setPropagateLambda(boolean propagateLambda) {
751 this.propagateLambda = propagateLambda;
752 }
753
754
755
756
757 public boolean getPropagateLambda() {
758 return propagateLambda;
759 }
760
761
762
763
764 @Override
765 public void setVelocity(double[] velocity) {
766 potential.setVelocity(velocity);
767 }
768
769
770
771
772
773
774 public void switchHistogram(int index) {
775 histogramIndex = index;
776 histogram = allHistograms.get(histogramIndex);
777 logger.info(" OST switching to histogram " + histogramIndex);
778 }
779
780 @Override
781 public void writeAdditionalRestartInfo(boolean recursive) {
782 histogram.writeRestart();
783 if (recursive) {
784 potential.writeAdditionalRestartInfo(recursive);
785 }
786 }
787
788 double getLambdaWriteOut() {
789 return histogram.hd.resetHistogramAtLambda;
790 }
791
792
793
794
795
796
797 double getForceFielddEdL() {
798 return dForceFieldEnergydL;
799 }
800
801
802
803
804
805
806 double getBiasEnergy() {
807 return biasEnergy;
808 }
809
810
811
812
813
814
815
816
817
818
819 boolean insideHardWallConstraint(double lambda, double dUdL) {
820 if (hardWallConstraint) {
821 double weight = histogram.getRecursionKernelValue(lambda, dUdL);
822 return weight > 0.0;
823 }
824 return true;
825 }
826
827
828
829
830 public class OptimizationParameters {
831
832
833
834
835
836
837 private boolean doOptimization = false;
838
839
840
841 private final boolean doUnitCellReset;
842
843
844
845
846
847 private final double lambdaCutoff;
848
849
850
851
852
853 private double optimumEnergy = Double.MAX_VALUE;
854
855
856
857
858
859 private final int frequency;
860
861
862
863
864
865 private final double eps;
866
867
868
869
870
871 private final double tolerance;
872
873
874
875
876
877 private final double energyWindow;
878
879
880
881 private double[] optimumCoords;
882
883
884
885 private File optimizationFile;
886
887
888
889 private SystemFilter optimizationFilter;
890
891
892
893
894 OptimizationParameters(CompositeConfiguration properties) {
895 energyWindow = properties.getDouble("ost-opt-energy-window", 10.0);
896 eps = properties.getDouble("ost-opt-eps", 0.1);
897 tolerance = properties.getDouble("ost-opt-tolerance", 1.0e-8);
898 frequency = properties.getInt("ost-opt-frequency", 10000);
899 lambdaCutoff = properties.getDouble("ost-opt-lambda-cutoff", 0.8);
900 doUnitCellReset = properties.getBoolean("ost-opt-unitcell-reset", false);
901 }
902
903 private void log() {
904 logger.info("\n Optimization Parameters");
905 logger.info(format(" Energy Window: %6.3f (kcal/mol)", energyWindow));
906 logger.info(format(" EPS: %6.4f RMS (kcal/mol/Å)", eps));
907 logger.info(format(" Tolerance: %6.4f (kcal/mol)", tolerance));
908 logger.info(format(" Frequency: %6d (steps)", frequency));
909 logger.info(format(" Lambda Cutoff: %6.4f", lambdaCutoff));
910 logger.info(format(" Unit Cell Reset: %6B", doUnitCellReset));
911 logger.info(format(" File: %s", optimizationFile.getName()));
912 }
913
914
915
916
917
918
919 public double[] getOptimumCoordinates() {
920 if (optimumEnergy < Double.MAX_VALUE) {
921 return optimumCoords;
922 } else {
923 logger.info(
924 "Lambda optimization cutoff was not reached. Try increasing the number of timesteps.");
925 return null;
926 }
927 }
928
929
930
931
932
933
934 public double getOptimumEnergy() {
935 if (optimumEnergy == Double.MAX_VALUE) {
936 logger.info(
937 "Lambda optimization cutoff was not reached. Try increasing the number of timesteps.");
938 }
939 return optimumEnergy;
940 }
941
942
943
944
945
946
947
948
949 public void optimize(double e, double[] x, @Nullable double[] gradient) {
950
951
952 double lambda = histogram.ld.lambda;
953 long stepsTaken = histogram.ld.stepsTaken;
954 if (doOptimization && lambda > lambdaCutoff && stepsTaken % frequency == 0) {
955 if (gradient == null) {
956 gradient = new double[x.length];
957 }
958 } else {
959 return;
960 }
961
962 logger.info(format("\n OST Minimization (Step %d)", stepsTaken));
963
964
965 lambdaInterface.setLambda(1.0);
966
967
968 potential.setEnergyTermState(Potential.STATE.BOTH);
969
970
971 boolean origBaroActive = true;
972 if (barostat != null) {
973 origBaroActive = barostat.isActive();
974 barostat.setActive(false);
975 }
976
977
978 try {
979 double startingEnergy = potential.energy(x);
980 Minimize minimize = new Minimize(molecularAssembly, potential, algorithmListener);
981 minimize.minimize(eps);
982
983 double minEnergy = potential.getTotalEnergy();
984
985 if (minEnergy < optimumEnergy + energyWindow) {
986 if (minEnergy < optimumEnergy) {
987 optimumEnergy = minEnergy;
988 }
989 int n = potential.getNumberOfVariables();
990 optimumCoords = new double[n];
991 optimumCoords = potential.getCoordinates(optimumCoords);
992 double mass = molecularAssembly.getMass();
993 Crystal crystal = molecularAssembly.getCrystal();
994 double density = crystal.getDensity(mass);
995 optimizationFilter.writeFile(optimizationFile, true);
996 Crystal uc = crystal.getUnitCell();
997 logger.info(format(" Minimum: %12.6f %s (%12.6f g/cc) optimized from %12.6f at step %d.",
998 minEnergy, uc.toShortString(), density, startingEnergy, stepsTaken));
999 }
1000 } catch (Exception ex) {
1001 String message = ex.getMessage();
1002 logger.info(
1003 format(" Exception minimizing coordinates at lambda=%8.6f\n %s.", lambda, message));
1004 logger.info(" Sampling will continue.");
1005 }
1006
1007
1008 lambdaInterface.setLambda(lambda);
1009
1010
1011 potential.setEnergyTermState(state);
1012
1013
1014 if (barostat != null) {
1015 barostat.setActive(origBaroActive);
1016 }
1017
1018 if (doUnitCellReset) {
1019 logger.info("\n Resetting Unit Cell");
1020 double mass = molecularAssembly.getMass();
1021 double density = molecularAssembly.getCrystal().getDensity(mass);
1022 molecularAssembly.applyRandomDensity(density);
1023 molecularAssembly.applyRandomSymOp(0.0);
1024 lambda = 0.0;
1025 lambdaInterface.setLambda(lambda);
1026 } else {
1027
1028 double eCheck = potential.energyAndGradient(x, gradient);
1029
1030 if (abs(eCheck - e) > tolerance) {
1031 logger.warning(
1032 format(" Optimization could not revert coordinates %16.8f vs. %16.8f.", e, eCheck));
1033 }
1034 }
1035 }
1036
1037
1038
1039
1040
1041
1042
1043 public void setOptimization(boolean doOptimization, MolecularAssembly molecularAssembly) {
1044 this.doOptimization = doOptimization;
1045 OrthogonalSpaceTempering.this.molecularAssembly = molecularAssembly;
1046 File file = molecularAssembly.getFile();
1047 String fileName = FilenameUtils.removeExtension(file.getAbsolutePath());
1048 String ext = FilenameUtils.getExtension(file.getAbsolutePath());
1049 if (optimizationFilter == null) {
1050 if (ext.toUpperCase().contains("XYZ")) {
1051 optimizationFile = new File(fileName + "_opt.arc");
1052 optimizationFilter = new XYZFilter(optimizationFile, molecularAssembly,
1053 molecularAssembly.getForceField(), molecularAssembly.getProperties());
1054 } else {
1055 optimizationFile = new File(fileName + "_opt.pdb");
1056 PDBFilter pdbFilter = new PDBFilter(optimizationFile, molecularAssembly,
1057 molecularAssembly.getForceField(), molecularAssembly.getProperties());
1058 int models = pdbFilter.countNumModels();
1059 pdbFilter.setModelNumbering(models);
1060 optimizationFilter = pdbFilter;
1061 }
1062 }
1063
1064 if (this.doOptimization) {
1065 log();
1066 }
1067 }
1068 }
1069
1070
1071
1072
1073
1074
1075
1076 public class Histogram {
1077
1078
1079
1080
1081 private static final double TWO_D_NORMALIZATION = 1.0;
1082
1083
1084
1085
1086 private static final double ONE_D_NORMALIZATION = Math.sqrt(2.0 * Math.PI);
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096 private static final double PSEUDO_BIAS_MAGNITUDE = 1.0;
1097
1098
1099
1100 protected final Comm world;
1101
1102
1103
1104 protected final int rank;
1105
1106
1107
1108 final double[] ensembleAveragedUdL;
1109
1110
1111
1112
1113
1114
1115 private final double deltaT;
1116
1117
1118
1119 private final IntegrationType integrationType;
1120
1121
1122
1123 private final SendAsynchronous sendAsynchronous;
1124
1125
1126
1127 private final SendSynchronous sendSynchronous;
1128
1129
1130
1131
1132 private double temperingWeight = 1.0;
1133
1134
1135
1136
1137
1138 private double[] kernelValues;
1139
1140
1141
1142 private double lastReceivedLambda;
1143
1144
1145
1146 private double lastReceiveddUdL;
1147 private final boolean spreadBias;
1148
1149
1150
1151 final HistogramData hd;
1152
1153
1154
1155 final LambdaData ld;
1156
1157
1158
1159 private final Stochastic stochastic;
1160
1161
1162
1163 private final SystemState lambdaState;
1164
1165
1166
1167 private int currentNumberOfHills = 0;
1168
1169
1170
1171 private double currentFreeEnergyDifference = 0.0;
1172
1173
1174
1175
1176
1177
1178
1179 Histogram(CompositeConfiguration properties, HistogramData histogramData, LambdaData lambdaData) {
1180 hd = histogramData;
1181 ld = lambdaData;
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195 deltaT = hd.temperingFactor * R * dynamicsOptions.getTemperature();
1196
1197
1198 ensembleAveragedUdL = new double[hd.getLambdaBins()];
1199
1200
1201 kernelValues = new double[hd.getDUDLBins()];
1202
1203 String propString = properties.getString("ost-integrationType", "SIMPSONS");
1204 IntegrationType testType;
1205 try {
1206 testType = IntegrationType.valueOf(propString.toUpperCase());
1207 } catch (Exception ex) {
1208 logger.warning(format(" Invalid argument %s to ost-integrationType; resetting to SIMPSONS", propString));
1209 testType = SIMPSONS;
1210 }
1211 integrationType = testType;
1212
1213 spreadBias = properties.getBoolean("ost-spread-bias", false);
1214
1215
1216
1217
1218
1219 world = Comm.world();
1220 int numProc = world.size();
1221 rank = world.rank();
1222 if (hd.asynchronous) {
1223
1224 sendAsynchronous = new SendAsynchronous(this);
1225 sendAsynchronous.start();
1226 sendSynchronous = null;
1227 } else {
1228 Histogram[] histograms = new Histogram[numProc];
1229 int[] rankToHistogramMap = new int[numProc];
1230 for (int i = 0; i < numProc; i++) {
1231 histograms[i] = this;
1232 rankToHistogramMap[i] = 0;
1233 }
1234 sendSynchronous = new SendSynchronous(histograms, rankToHistogramMap);
1235 sendAsynchronous = null;
1236 }
1237
1238 lastReceivedLambda = ld.lambda;
1239 if (hd.discreteLambda) {
1240 lastReceivedLambda = mapLambda(lastReceivedLambda);
1241 logger.info(format(" Discrete lambda: initializing lambda to nearest bin %.5f", lastReceivedLambda));
1242 ld.lambda = lastReceivedLambda;
1243 ld.theta = asin(sqrt(lastReceivedLambda));
1244 lambdaInterface.setLambda(lastReceivedLambda);
1245 lambdaState = null;
1246 stochastic = null;
1247 } else {
1248
1249 lambdaState = new SystemState(1);
1250 stochastic = new Stochastic(lambdaParticleOptions.getLambdaFriction(), lambdaState);
1251 stochastic.setTemperature(dynamicsOptions.getTemperature());
1252 stochastic.setTimeStep(dynamicsOptions.getDtPsec());
1253 double[] mass = lambdaState.getMass();
1254 double[] thetaPosition = lambdaState.x();
1255 double[] thetaVelocity = lambdaState.v();
1256 double[] thetaAccel = lambdaState.a();
1257 mass[0] = lambdaParticleOptions.getLambdaMass();
1258 thetaPosition[0] = ld.theta;
1259 thetaVelocity[0] = ld.thetaVelocity;
1260 thetaAccel[0] = ld.thetaAcceleration;
1261 }
1262
1263 lastReceiveddUdL = getdEdL();
1264
1265 if (hd.wasHistogramRead()) {
1266 updateFreeEnergyDifference(true, false);
1267 }
1268 }
1269
1270 public void disableResetStatistics() {
1271 hd.resetHistogram = false;
1272 }
1273
1274
1275
1276
1277
1278
1279
1280 public StringBuffer evaluateTotalOSTBias(boolean bias) {
1281 StringBuffer sb = new StringBuffer();
1282 double[] chainRule = new double[2];
1283 for (int dUdLBin = 0; dUdLBin < hd.getDUDLBins(); dUdLBin++) {
1284 double currentdUdL = dUdLforBin(dUdLBin);
1285 sb.append(format(" %16.8f", currentdUdL));
1286 for (int lambdaBin = 0; lambdaBin < hd.getLambdaBins(); lambdaBin++) {
1287 double currentLambda = lambdaForIndex(lambdaBin);
1288 double bias1D = energyAndGradient1D(currentLambda, false);
1289 double bias2D = energyAndGradient2D(currentLambda, currentdUdL, chainRule, hd.biasMag);
1290 double totalBias = bias1D + bias2D;
1291
1292 if (!bias) {
1293 totalBias = -totalBias;
1294 }
1295 sb.append(format(" %16.8f", totalBias));
1296 }
1297 sb.append("\n");
1298 }
1299 return sb;
1300 }
1301
1302
1303
1304
1305
1306
1307
1308 public StringBuffer evaluate2DOSTBias(boolean bias) {
1309 StringBuffer sb = new StringBuffer();
1310 double[] chainRule = new double[2];
1311 for (int dUdLBin = 0; dUdLBin < hd.getDUDLBins(); dUdLBin++) {
1312 double currentdUdL = dUdLforBin(dUdLBin);
1313 sb.append(format(" %16.8f", currentdUdL));
1314 for (int lambdaBin = 0; lambdaBin < hd.getLambdaBins(); lambdaBin++) {
1315 double currentLambda = lambdaForIndex(lambdaBin);
1316 double bias2D = energyAndGradient2D(currentLambda, currentdUdL, chainRule, hd.biasMag);
1317 if (!bias) {
1318 bias2D = -bias2D;
1319 }
1320 sb.append(format(" %16.8f", bias2D));
1321 }
1322 sb.append("\n");
1323 }
1324 return sb;
1325 }
1326
1327
1328
1329
1330
1331
1332
1333 public boolean getIndependentWalkers() {
1334 return hd.independentWalkers;
1335 }
1336
1337 public double getLambdaResetValue() {
1338 return hd.resetHistogramAtLambda;
1339 }
1340
1341
1342
1343
1344
1345
1346 public int getRank() {
1347 return rank;
1348 }
1349
1350 public boolean getResetStatistics() {
1351 return hd.resetHistogram;
1352 }
1353
1354
1355
1356
1357
1358
1359 public Optional<SendSynchronous> getSynchronousSend() {
1360 return Optional.ofNullable(sendSynchronous);
1361 }
1362
1363 public void setLastReceiveddUdL(double lastReceiveddUdL) {
1364 this.lastReceiveddUdL = lastReceiveddUdL;
1365 }
1366
1367 public double updateFreeEnergyDifference(boolean print, boolean save) {
1368 if (hd.metaDynamics) {
1369 return updateMetaDynamicsFreeEnergyDifference(print, save);
1370 } else {
1371 return updateOSTFreeEnergyDifference(print, save);
1372 }
1373 }
1374
1375
1376
1377
1378
1379
1380
1381
1382 public double updateMetaDynamicsFreeEnergyDifference(boolean print, boolean save) {
1383 synchronized (this) {
1384 if (!print && currentNumberOfHills == hd.counts) {
1385 return currentFreeEnergyDifference;
1386 }
1387 double freeEnergyDifferenceMeta = 0.0;
1388 double freeEnergyDifferenceTI = 0.0;
1389
1390 double minFL = Double.MAX_VALUE;
1391
1392 int lambdaBins = hd.getLambdaBins();
1393 double[] metaFreeEnergy = new double[lambdaBins];
1394
1395
1396 double totalWeight = 0;
1397 StringBuilder stringBuilder = new StringBuilder();
1398
1399
1400 for (int lambdaBin = 0; lambdaBin < lambdaBins; lambdaBin++) {
1401 int firstdUdLBin = firstdUdLBin(lambdaBin);
1402 int lastdUdLBin = lastdUdLBin(lambdaBin);
1403 double lambdaCount = 0;
1404
1405 double mindUdLForLambda = 0.0;
1406 double maxdUdLforLambda = 0.0;
1407 double maxBias = 0;
1408 if (firstdUdLBin == -1 || lastdUdLBin == -1) {
1409 ensembleAveragedUdL[lambdaBin] = 0.0;
1410 minFL = 0.0;
1411 } else {
1412 double ensembleAverageFLambda = 0.0;
1413 double partitionFunction = 0.0;
1414
1415
1416 double offset = evaluateKernelForLambda(lambdaBin, firstdUdLBin, lastdUdLBin);
1417
1418 for (int dUdLBin = firstdUdLBin; dUdLBin <= lastdUdLBin; dUdLBin++) {
1419 double kernel = kernelValues[dUdLBin];
1420
1421 if (kernel - offset > maxBias) {
1422 maxBias = kernel - offset;
1423 }
1424
1425
1426 partitionFunction += kernel;
1427
1428 double currentdUdL = dUdLforBin(dUdLBin);
1429 ensembleAverageFLambda += currentdUdL * kernel;
1430 lambdaCount += getRecursionKernelValue(lambdaBin, dUdLBin);
1431 }
1432 if (minFL > maxBias) {
1433 minFL = maxBias;
1434 }
1435 ensembleAveragedUdL[lambdaBin] = (partitionFunction == 0) ? 0 : ensembleAverageFLambda / partitionFunction;
1436 mindUdLForLambda = hd.dUdLMinimum + firstdUdLBin * hd.dUdLBinWidth;
1437 maxdUdLforLambda = hd.dUdLMinimum + (lastdUdLBin + 1) * hd.dUdLBinWidth;
1438 }
1439
1440 double deltaFreeEnergy = ensembleAveragedUdL[lambdaBin] * deltaForLambdaBin(lambdaBin);
1441 freeEnergyDifferenceTI += deltaFreeEnergy;
1442 totalWeight += lambdaCount;
1443
1444
1445 double lambdaBinWidth = hd.lambdaBinWidth;
1446 double llL = lambdaBin * lambdaBinWidth - hd.lambdaBinWidth_2;
1447 double ulL = llL + lambdaBinWidth;
1448 if (llL < 0.0) {
1449 llL = 0.0;
1450 }
1451 if (ulL > 1.0) {
1452 ulL = 1.0;
1453 }
1454
1455 double midLambda = llL;
1456 if (!hd.discreteLambda) {
1457 midLambda = (llL + ulL) / 2.0;
1458 }
1459 double bias1D = energyAndGradient1D(midLambda, false);
1460 metaFreeEnergy[lambdaBin] = energyAndGradientMeta(midLambda, false);
1461 freeEnergyDifferenceMeta = -(metaFreeEnergy[lambdaBin] - metaFreeEnergy[0]);
1462
1463 if (print || save) {
1464 stringBuilder.append(format(" %6.2e %7.5f %7.1f %7.1f %9.2f %9.2f %9.2f %9.2f %9.2f\n", lambdaCount,
1465 midLambda, mindUdLForLambda, maxdUdLforLambda, ensembleAveragedUdL[lambdaBin],
1466 bias1D, freeEnergyDifferenceTI, metaFreeEnergy[lambdaBin], freeEnergyDifferenceMeta));
1467 }
1468 }
1469
1470 double temperingOffset = hd.getTemperingOffset();
1471 double temperEnergy = (minFL > temperingOffset) ? temperingOffset - minFL : 0;
1472 temperingWeight = exp(temperEnergy / deltaT);
1473
1474 if (print) {
1475 logger.info(" Weight Lambda dU/dL Bins <dU/dL> g(L) dG_OST(L) Meta(L) dG_Meta(L)");
1476 logger.info(stringBuilder.toString());
1477 logger.info(" Histogram Evaluation");
1478 logger.info(format(" Free Energy Difference: %12.4f kcal/mol", freeEnergyDifferenceMeta));
1479 logger.info(format(" Number of Hills: %12d", hd.counts));
1480 logger.info(format(" Total Bias Added: %12.4f kcal/mol", totalWeight * hd.biasMag));
1481 logger.info(format(" Minimum Bias: %12.4f kcal/mol", minFL));
1482 logger.info(format(" Tempering Percentage: %12.4f %%\n", temperingWeight * 100));
1483 }
1484
1485 if (save) {
1486 String modelFilename = molecularAssembly.getFile().getAbsolutePath();
1487 File saveDir = new File(FilenameUtils.getFullPath(modelFilename));
1488 String dirName = saveDir + File.separator;
1489 String fileName = dirName + "histogram.txt";
1490 try {
1491 logger.info(" Writing " + fileName);
1492 PrintWriter printWriter = new PrintWriter(fileName);
1493 printWriter.write(stringBuilder.toString());
1494 printWriter.close();
1495 } catch (Exception e) {
1496 logger.info(format(" Failed to write %s.", fileName));
1497 }
1498 }
1499
1500 currentNumberOfHills = hd.counts;
1501 currentFreeEnergyDifference = freeEnergyDifferenceMeta;
1502 return currentFreeEnergyDifference;
1503 }
1504 }
1505
1506
1507
1508
1509
1510
1511
1512
1513 public double updateOSTFreeEnergyDifference(boolean print, boolean save) {
1514 synchronized (this) {
1515 if (!print && currentNumberOfHills == hd.counts) {
1516 return currentFreeEnergyDifference;
1517 }
1518
1519 double freeEnergyDifferenceOST = 0.0;
1520 double minFL = Double.MAX_VALUE;
1521
1522
1523
1524 boolean biasMagZero = hd.biasMag <= 0.0;
1525
1526
1527 double totalWeight = 0;
1528 double beta = 1.0 / (R * dynamicsOptions.getTemperature());
1529 StringBuilder stringBuilder = new StringBuilder();
1530
1531
1532 for (int lambdaBin = 0; lambdaBin < hd.lambdaBins; lambdaBin++) {
1533 int firstdUdLBin = firstdUdLBin(lambdaBin);
1534 int lastdUdLBin = lastdUdLBin(lambdaBin);
1535 double lambdaCount = 0;
1536
1537 double mindUdLforLambda = 0.0;
1538 double maxdUdLforLambda = 0.0;
1539 double maxBias = 0;
1540 if (firstdUdLBin == -1 || lastdUdLBin == -1) {
1541 ensembleAveragedUdL[lambdaBin] = 0.0;
1542 minFL = 0.0;
1543 } else {
1544 double ensembleAverage = 0.0;
1545 double partitionFunction = 0.0;
1546
1547
1548 double offset = evaluateKernelForLambda(lambdaBin, firstdUdLBin, lastdUdLBin);
1549
1550 for (int dUdLBin = firstdUdLBin; dUdLBin <= lastdUdLBin; dUdLBin++) {
1551 double kernel = kernelValues[dUdLBin];
1552
1553 if (kernel - offset > maxBias) {
1554 maxBias = kernel - offset;
1555 }
1556
1557 double weight;
1558 if (biasMagZero) {
1559
1560 weight = kernel;
1561 } else {
1562
1563 weight = exp(kernel * beta);
1564 }
1565 partitionFunction += weight;
1566
1567 double currentdUdL = dUdLforBin(dUdLBin);
1568 ensembleAverage += currentdUdL * weight;
1569 lambdaCount += getRecursionKernelValue(lambdaBin, dUdLBin);
1570 }
1571 if (minFL > maxBias) {
1572 minFL = maxBias;
1573 }
1574 ensembleAveragedUdL[lambdaBin] = (partitionFunction == 0) ? 0 : ensembleAverage / partitionFunction;
1575 mindUdLforLambda = hd.dUdLMinimum + firstdUdLBin * hd.dUdLBinWidth;
1576 maxdUdLforLambda = hd.dUdLMinimum + (lastdUdLBin + 1) * hd.dUdLBinWidth;
1577 }
1578
1579 double deltaFreeEnergy = ensembleAveragedUdL[lambdaBin] * deltaForLambdaBin(lambdaBin);
1580 freeEnergyDifferenceOST += deltaFreeEnergy;
1581 totalWeight += lambdaCount;
1582
1583 if (print || save) {
1584 double llL = lambdaBin * hd.lambdaBinWidth - hd.lambdaBinWidth_2;
1585 double ulL = llL + hd.lambdaBinWidth;
1586 if (llL < 0.0) {
1587 llL = 0.0;
1588 }
1589 if (ulL > 1.0) {
1590 ulL = 1.0;
1591 }
1592
1593 double midLambda = llL;
1594 if (!hd.discreteLambda) {
1595 midLambda = (llL + ulL) / 2.0;
1596 }
1597 double bias1D = energyAndGradient1D(midLambda, false);
1598
1599 double bias2D = 0.0;
1600 if (!biasMagZero) {
1601 bias2D = computeBiasEnergy(midLambda, ensembleAveragedUdL[lambdaBin]) - bias1D;
1602 }
1603
1604 stringBuilder.append(
1605 format(" %6.2e %7.5f %7.1f %7.1f %8.2f %8.2f %8.2f %8.2f %8.2f %8.2f\n", lambdaCount,
1606 midLambda, mindUdLforLambda, maxdUdLforLambda, ensembleAveragedUdL[lambdaBin],
1607 bias1D, bias2D, bias1D + bias2D, freeEnergyDifferenceOST, bias1D + bias2D + freeEnergyDifferenceOST));
1608 }
1609 }
1610
1611 if (!biasMagZero) {
1612 double temperingOffset = hd.getTemperingOffset();
1613 double temperEnergy = (minFL > temperingOffset) ? temperingOffset - minFL : 0;
1614 temperingWeight = exp(temperEnergy / deltaT);
1615 }
1616
1617 freeEnergyDifferenceOST = integrateNumeric(ensembleAveragedUdL, integrationType);
1618
1619 if (print) {
1620 logger.info(" Weight Lambda dU/dL Bins <dU/dL> g(L) f(L,<dU/dL>) Bias dG(L) Bias+dG(L)");
1621 logger.info(stringBuilder.toString());
1622 logger.info(" Histogram Evaluation");
1623 logger.info(format(" Free Energy Difference: %12.4f kcal/mol", freeEnergyDifferenceOST));
1624 logger.info(format(" Number of Hills: %12d", hd.counts));
1625 logger.info(format(" Total Bias Added: %12.4f kcal/mol", totalWeight * hd.biasMag));
1626 logger.info(format(" Minimum Bias: %12.4f kcal/mol", minFL));
1627 logger.info(format(" Tempering Percentage: %12.4f %%\n", temperingWeight * 100));
1628 }
1629
1630 if (save) {
1631 String modelFilename = molecularAssembly.getFile().getAbsolutePath();
1632 File saveDir = new File(FilenameUtils.getFullPath(modelFilename));
1633 String dirName = saveDir + File.separator;
1634 String fileName = dirName + "histogram.txt";
1635 try {
1636 logger.info("\n Writing " + fileName);
1637 PrintWriter printWriter = new PrintWriter(fileName);
1638 printWriter.write(stringBuilder.toString());
1639 printWriter.close();
1640 } catch (Exception e) {
1641 logger.info(format(" Failed to write %s.", fileName));
1642 }
1643 }
1644
1645 currentNumberOfHills = hd.counts;
1646 currentFreeEnergyDifference = freeEnergyDifferenceOST;
1647 return currentFreeEnergyDifference;
1648 }
1649 }
1650
1651
1652
1653
1654
1655
1656
1657 private double deltaForLambdaBin(int lambdaBin) {
1658 if (!hd.discreteLambda && (lambdaBin == 0 || lambdaBin == hd.lambdaBins - 1)) {
1659
1660 return hd.lambdaBinWidth_2;
1661 } else if (hd.discreteLambda && lambdaBin == 0) {
1662
1663 return 0.0;
1664 }
1665
1666 return hd.lambdaBinWidth;
1667 }
1668
1669
1670
1671
1672
1673
1674
1675
1676 private int firstdUdLBin(int lambdaBin) {
1677
1678 synchronized (this) {
1679
1680 for (int jFL = 0; jFL < hd.dUdLBins; jFL++) {
1681 double count = hd.zHistogram[lambdaBin][jFL];
1682 if (count > 0) {
1683 return jFL;
1684 }
1685 }
1686 return -1;
1687 }
1688 }
1689
1690
1691
1692
1693
1694
1695
1696
1697 private int lastdUdLBin(int lambdaBin) {
1698
1699 synchronized (this) {
1700
1701 for (int jFL = hd.dUdLBins - 1; jFL >= 0; jFL--) {
1702 double count = hd.zHistogram[lambdaBin][jFL];
1703 if (count > 0) {
1704 return jFL;
1705 }
1706 }
1707 return -1;
1708 }
1709 }
1710
1711
1712
1713
1714
1715
1716 void writeRestart() {
1717 if (rank == 0 || hd.independentWrite) {
1718 updateFreeEnergyDifference(true, false);
1719 try {
1720 hd.writeHistogram();
1721 logger.info(format(" Wrote histogram restart to: %s.", hd.getHistogramFileName()));
1722 } catch (Exception ex) {
1723 String message = format(" Exception writing histogram restart file %s.", hd.getHistogramFileName());
1724 logger.log(Level.INFO, Utilities.stackTraceToString(ex));
1725 logger.log(Level.SEVERE, message, ex);
1726 }
1727 }
1728
1729
1730 try {
1731 ld.writeLambdaData();
1732 logger.info(format(" Wrote lambda restart to: %s.", ld.getLambdaFileName()));
1733 } catch (Exception ex) {
1734 String message = format(" Exception writing lambda restart file %s.", ld.getLambdaFileName());
1735 logger.log(Level.INFO, Utilities.stackTraceToString(ex));
1736 logger.log(Level.SEVERE, message, ex);
1737 }
1738 }
1739
1740 private double mapLambda(double lambda) {
1741 if (hd.discreteLambda) {
1742 return hd.lambdaLadder[indexForDiscreteLambda(lambda)];
1743 } else {
1744 return lambda;
1745 }
1746 }
1747
1748
1749
1750
1751
1752
1753
1754
1755 int indexForLambda(double lambda) {
1756 if (hd.discreteLambda) {
1757 return indexForDiscreteLambda(lambda);
1758 } else {
1759 return indexForContinuousLambda(lambda);
1760 }
1761 }
1762
1763 int indexForLambda() {
1764 return indexForLambda(ld.lambda);
1765 }
1766
1767 private double lambdaForIndex(int bin) {
1768 if (hd.discreteLambda) {
1769 return hd.lambdaLadder[bin];
1770 } else {
1771 return bin * hd.lambdaBinWidth;
1772 }
1773 }
1774
1775 private int indexForContinuousLambda(double lambda) {
1776 int lambdaBin = (int) floor((lambda - hd.minLambda) / hd.lambdaBinWidth);
1777 if (lambdaBin < 0) {
1778 lambdaBin = 0;
1779 }
1780 if (lambdaBin >= hd.lambdaBins) {
1781 lambdaBin = hd.lambdaBins - 1;
1782 }
1783 return lambdaBin;
1784 }
1785
1786 private int indexForDiscreteLambda(double lambda) {
1787 assert hd.discreteLambda && hd.lambdaLadder != null && hd.lambdaLadder.length > 0;
1788
1789 int initialGuess = indexForContinuousLambda(lambda);
1790 double minErr = Double.MAX_VALUE;
1791 int minErrBin = -1;
1792 for (int i = -1; i < 2; i++) {
1793 int guessBin = i + initialGuess;
1794 if (guessBin < 0 || guessBin >= hd.lambdaBins) {
1795 continue;
1796 }
1797 double guessLam = hd.lambdaLadder[guessBin];
1798 double guessErr = Math.abs(guessLam - lambda);
1799 if (guessErr < minErr) {
1800 minErr = guessErr;
1801 minErrBin = guessBin;
1802 }
1803 }
1804
1805 assert minErr < 1.0E-6;
1806 return minErrBin;
1807 }
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817 int binFordUdL(double dUdL) {
1818
1819
1820 if (dUdL < hd.dUdLMinimum) {
1821 return -1;
1822 }
1823
1824
1825 if (dUdL > hd.dUdLMaximum) {
1826 return -1;
1827 }
1828
1829 int bin = (int) floor((dUdL - hd.dUdLMinimum) / hd.dUdLBinWidth);
1830
1831 if (bin == hd.dUdLBins) {
1832 bin = hd.dUdLBins - 1;
1833 }
1834
1835 if (bin < 0) {
1836 bin = 0;
1837 }
1838
1839 return bin;
1840 }
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850 double getRecursionKernelValue(int lambdaBin, int dUdLBin) {
1851
1852 synchronized (this) {
1853
1854 lambdaBin = lambdaMirror(lambdaBin);
1855
1856
1857 if (dUdLBin < 0 || dUdLBin >= hd.dUdLBins) {
1858 return 0.0;
1859 }
1860
1861 return hd.zHistogram[lambdaBin][dUdLBin];
1862 }
1863 }
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873 double getRecursionKernelValue(double lambda, double dUdL) {
1874 return getRecursionKernelValue(indexForLambda(lambda), binFordUdL(dUdL));
1875 }
1876
1877
1878
1879
1880
1881
1882
1883
1884 void addToRecursionKernelValue(double currentLambda, double currentdUdL, double value) {
1885 synchronized (this) {
1886 if (spreadBias) {
1887
1888 int dUdLbiasCutoff = hd.dUdLBiasCutoff;
1889 checkRecursionKernelSize(dUdLforBin(binFordUdL(currentdUdL) - dUdLbiasCutoff));
1890 checkRecursionKernelSize(dUdLforBin(binFordUdL(currentdUdL) + dUdLbiasCutoff));
1891
1892 int currentLambdaBin = indexForLambda(currentLambda);
1893 int currentdUdLBin = binFordUdL(currentdUdL);
1894
1895
1896 double invLs2 = 0.5 / hd.lambdaVariance;
1897 double invFLs2 = 0.5 / hd.dUdLVariance;
1898
1899
1900 double normalize = 0.0;
1901 int lambdaBiasCutoff = hd.lambdaBiasCutoff;
1902 for (int iL = -lambdaBiasCutoff; iL <= lambdaBiasCutoff; iL++) {
1903 int lambdaBin = currentLambdaBin + iL;
1904 double deltaL = currentLambda - lambdaBin * hd.lambdaBinWidth;
1905 double deltaL2 = deltaL * deltaL;
1906
1907 double L2exp = exp(-deltaL2 * invLs2);
1908 for (int jFL = -dUdLbiasCutoff; jFL <= dUdLbiasCutoff; jFL++) {
1909 int dUdLBin = currentdUdLBin + jFL;
1910 double deltaFL = currentdUdL - dUdLforBin(dUdLBin);
1911 double deltaFL2 = deltaFL * deltaFL;
1912 normalize += L2exp * exp(-deltaFL2 * invFLs2);
1913 }
1914 }
1915
1916 for (int iL = -lambdaBiasCutoff; iL <= lambdaBiasCutoff; iL++) {
1917 int lambdaBin = currentLambdaBin + iL;
1918 double deltaL = currentLambda - lambdaBin * hd.lambdaBinWidth;
1919 double deltaL2 = deltaL * deltaL;
1920
1921 double L2exp = exp(-deltaL2 * invLs2);
1922 lambdaBin = lambdaMirror(lambdaBin);
1923 for (int jFL = -dUdLbiasCutoff; jFL <= dUdLbiasCutoff; jFL++) {
1924 int dUdLBin = currentdUdLBin + jFL;
1925 double deltaFL = currentdUdL - dUdLforBin(dUdLBin);
1926 double deltaFL2 = deltaFL * deltaFL;
1927 double weight = value / normalize * L2exp * exp(-deltaFL2 * invFLs2);
1928 hd.zHistogram[lambdaBin][dUdLBin] += weight;
1929 }
1930 }
1931 hd.counts++;
1932 } else {
1933
1934 checkRecursionKernelSize(currentdUdL);
1935 int lambdaBin = indexForLambda(currentLambda);
1936 int dUdLBin = binFordUdL(currentdUdL);
1937 try {
1938 hd.zHistogram[lambdaBin][dUdLBin] += value;
1939 hd.counts++;
1940 } catch (IndexOutOfBoundsException e) {
1941 logger.info(format(" Histogram dimensions %d %d", hd.lambdaBins, hd.dUdLBins));
1942 logger.info(format(
1943 " Count skipped in addToRecursionKernelValue due to an index out of bounds exception.\n L=%10.8f (%d), dU/dL=%10.8f (%d) and count=%10.8f",
1944 currentLambda, lambdaBin, currentdUdL, dUdLBin, value));
1945
1946 }
1947 }
1948 }
1949 }
1950
1951
1952
1953
1954 void allocateRecursionKernel() {
1955
1956 synchronized (this) {
1957 hd.zHistogram = new double[hd.lambdaBins][hd.dUdLBins];
1958 kernelValues = new double[hd.dUdLBins];
1959 }
1960 }
1961
1962
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973 private double integrateNumeric(double[] dUdLs, IntegrationType type) {
1974 double val;
1975 if (hd.discreteLambda) {
1976
1977 double[] lams = Integrate1DNumeric.generateXPoints(0.0, 1.0, hd.lambdaBins, false);
1978 DataSet dSet = new DoublesDataSet(lams, dUdLs, false);
1979 val = Integrate1DNumeric.integrateData(dSet, Integrate1DNumeric.IntegrationSide.LEFT, type);
1980 } else {
1981
1982 double[] midLams = Integrate1DNumeric.generateXPoints(hd.lambdaBinWidth, 1.0 - hd.lambdaBinWidth,
1983 (hd.lambdaBins - 2), false);
1984 double[] midVals = Arrays.copyOfRange(dUdLs, 1, (hd.lambdaBins - 1));
1985 DataSet dSet = new DoublesDataSet(midLams, midVals, false);
1986
1987 val = Integrate1DNumeric.integrateData(dSet, Integrate1DNumeric.IntegrationSide.LEFT, type);
1988
1989
1990
1991 double dL_4 = hd.lambdaBinWidth_2 * 0.5;
1992
1993
1994
1995 double val0 = 0;
1996 double val1 = 0;
1997
1998
1999 if (!lambdaInterface.dEdLZeroAtEnds()) {
2000 double recipSlopeLen = 1.0 / (hd.lambdaBinWidth * 0.75);
2001
2002 double slope = dUdLs[0] - dUdLs[1];
2003 slope *= recipSlopeLen;
2004 val0 = dUdLs[0] + (slope * dL_4);
2005
2006 slope = dUdLs[hd.lambdaBins - 1] - dUdLs[hd.lambdaBins - 2];
2007 slope *= recipSlopeLen;
2008 val1 = dUdLs[hd.lambdaBins - 1] + (slope * dL_4);
2009 logger.fine(format(" Inferred dU/dL values at 0 and 1: %10.5g , %10.5g", val0, val1));
2010 }
2011
2012
2013
2014 val += trapezoid(0, dL_4, val0, dUdLs[0]);
2015 val += trapezoid(dL_4, hd.lambdaBinWidth, dUdLs[0], dUdLs[1]);
2016 val += trapezoid(1.0 - hd.lambdaBinWidth, 1.0 - dL_4, dUdLs[hd.lambdaBins - 2],
2017 dUdLs[hd.lambdaBins - 1]);
2018 val += trapezoid(1.0 - dL_4, 1.0, dUdLs[hd.lambdaBins - 1], val1);
2019 }
2020
2021 return val;
2022 }
2023
2024
2025
2026
2027
2028
2029
2030
2031
2032
2033 private double trapezoid(double x0, double x1, double fx0, double fx1) {
2034 double val = 0.5 * (fx0 + fx1);
2035 val *= (x1 - x0);
2036 return val;
2037 }
2038
2039
2040
2041
2042
2043
2044
2045
2046
2047
2048
2049
2050
2051
2052
2053 private double evaluateKernelForLambda(int lambda, int llFL, int ulFL) {
2054 double maxKernel = Double.MIN_VALUE;
2055
2056 double gaussianBiasMagnitude = hd.biasMag;
2057 if (hd.biasMag <= 0.0) {
2058 gaussianBiasMagnitude = PSEUDO_BIAS_MAGNITUDE;
2059 }
2060
2061 for (int jFL = llFL; jFL <= ulFL; jFL++) {
2062 double value = evaluateKernel(lambda, jFL, gaussianBiasMagnitude);
2063 kernelValues[jFL] = value;
2064 if (value > maxKernel) {
2065 maxKernel = value;
2066 }
2067 }
2068
2069
2070 double offset = 0.0;
2071 if (hd.biasMag > 0.0 && !hd.metaDynamics) {
2072 offset = -maxKernel;
2073 for (int jFL = llFL; jFL <= ulFL; jFL++) {
2074 kernelValues[jFL] += offset;
2075 }
2076 }
2077 return offset;
2078 }
2079
2080
2081
2082
2083
2084
2085
2086 private int lambdaMirror(int bin) {
2087 if (bin < 0) {
2088 return -bin;
2089 }
2090 int lastBin = hd.lambdaBins - 1;
2091 if (bin > lastBin) {
2092
2093 bin -= lastBin;
2094
2095 return lastBin - bin;
2096 }
2097
2098 return bin;
2099 }
2100
2101
2102
2103
2104
2105
2106
2107
2108 private double mirrorFactor(int bin) {
2109 if (hd.discreteLambda) {
2110 return 1.0;
2111 }
2112 if (bin == 0 || bin == hd.lambdaBins - 1) {
2113 return 2.0;
2114 }
2115 return 1.0;
2116 }
2117
2118
2119
2120
2121
2122
2123
2124 private double dUdLforBin(int dUdLBin) {
2125 return hd.dUdLMinimum + dUdLBin * hd.dUdLBinWidth + hd.dUdLBinWidth_2;
2126 }
2127
2128
2129
2130
2131 private double evaluateKernel(int currentLambdaBin, int currentdUdLBin, double gaussianBiasMagnitude) {
2132
2133 double currentLambda = currentLambdaBin * hd.lambdaBinWidth;
2134 double currentdUdL = dUdLforBin(currentdUdLBin);
2135
2136
2137 double invLs2 = 0.5 / hd.lambdaVariance;
2138 double invFLs2 = 0.5 / hd.dUdLVariance;
2139
2140 double sum = 0.0;
2141 int lambdaBiasCutoff = hd.lambdaBiasCutoff;
2142 for (int iL = -lambdaBiasCutoff; iL <= lambdaBiasCutoff; iL++) {
2143 int lambdaBin = currentLambdaBin + iL;
2144 double deltaL = currentLambda - lambdaBin * hd.lambdaBinWidth;
2145 double deltaL2 = deltaL * deltaL;
2146
2147
2148 double L2exp = exp(-deltaL2 * invLs2);
2149
2150
2151 lambdaBin = lambdaMirror(lambdaBin);
2152 double mirrorFactor = mirrorFactor(lambdaBin);
2153
2154 int dUdLbiasCutoff = hd.dUdLBiasCutoff;
2155 for (int jFL = -dUdLbiasCutoff; jFL <= dUdLbiasCutoff; jFL++) {
2156 int dUdLBin = currentdUdLBin + jFL;
2157 double weight = mirrorFactor * getRecursionKernelValue(lambdaBin, dUdLBin);
2158 if (weight <= 0.0) {
2159 continue;
2160 }
2161 double deltaFL = currentdUdL - dUdLforBin(dUdLBin);
2162 double deltaFL2 = deltaFL * deltaFL;
2163 double e = weight * L2exp * exp(-deltaFL2 * invFLs2);
2164 sum += e;
2165 }
2166 }
2167 return gaussianBiasMagnitude * sum;
2168 }
2169
2170
2171
2172
2173
2174
2175
2176
2177 double computeBiasEnergy(double currentLambda, double currentdUdL) {
2178 synchronized (this) {
2179 int currentLambdaBin = indexForLambda(currentLambda);
2180 int currentdUdLBin = binFordUdL(currentdUdL);
2181
2182 double bias1D;
2183 double bias2D = 0.0;
2184
2185 if (!hd.metaDynamics) {
2186 if (hd.biasMag > 0.0) {
2187 int lambdaBiasCutoff = hd.lambdaBiasCutoff;
2188 for (int iL = -lambdaBiasCutoff; iL <= lambdaBiasCutoff; iL++) {
2189
2190 int lambdaBin = currentLambdaBin + iL;
2191 double deltaL = currentLambda - (lambdaBin * hd.lambdaBinWidth);
2192 double deltaL2 = deltaL * deltaL;
2193 double expL2 = exp(-deltaL2 / (2.0 * hd.lambdaVariance));
2194
2195
2196 lambdaBin = lambdaMirror(lambdaBin);
2197 double mirrorFactor = mirrorFactor(lambdaBin);
2198
2199 int dUdLbiasCutoff = hd.dUdLBiasCutoff;
2200 for (int iFL = -dUdLbiasCutoff; iFL <= dUdLbiasCutoff; iFL++) {
2201 int dUdLBin = currentdUdLBin + iFL;
2202
2203 double weight = mirrorFactor * getRecursionKernelValue(lambdaBin, dUdLBin);
2204 if (weight <= 0.0) {
2205 continue;
2206 }
2207
2208 double deltaFL = currentdUdL - dUdLforBin(dUdLBin);
2209 double deltaFL2 = deltaFL * deltaFL;
2210 double bias = weight * hd.biasMag * expL2 * exp(-deltaFL2 / (2.0 * hd.dUdLVariance));
2211 bias2D += bias;
2212 }
2213 }
2214 }
2215
2216 bias1D = energyAndGradient1D(currentLambda, false);
2217 } else {
2218 bias1D = energyAndGradientMeta(currentLambda, false);
2219 }
2220
2221
2222 return bias1D + bias2D;
2223 }
2224 }
2225
2226
2227
2228
2229
2230
2231
2232
2233 double energyAndGradient2D(double currentdUdLambda, double[] chainRule) {
2234 return energyAndGradient2D(ld.lambda, currentdUdLambda, chainRule, hd.biasMag);
2235 }
2236
2237
2238
2239
2240
2241
2242
2243
2244
2245
2246 double energyAndGradient2D(double currentLambda, double currentdUdLambda, double[] chainRule,
2247 double gaussianBiasMagnitude) {
2248
2249 if (gaussianBiasMagnitude <= 0.0) {
2250 chainRule[0] = 0.0;
2251 chainRule[1] = 0.0;
2252 return 0;
2253 }
2254
2255
2256 double gLdEdL = 0.0;
2257 double dGdLambda = 0.0;
2258 double dGdFLambda = 0.0;
2259 int currentLambdaBin = indexForLambda(currentLambda);
2260 int currentdUdLBin = binFordUdL(currentdUdLambda);
2261
2262 int lambdaBiasCutoff = hd.lambdaBiasCutoff;
2263 for (int iL = -lambdaBiasCutoff; iL <= lambdaBiasCutoff; iL++) {
2264 int lambdaBin = currentLambdaBin + iL;
2265
2266
2267 double deltaL = currentLambda - (lambdaBin * hd.lambdaBinWidth);
2268 double deltaL2 = deltaL * deltaL;
2269 double expL2 = exp(-deltaL2 / (2.0 * hd.lambdaVariance));
2270
2271
2272 double mirrorFactor = mirrorFactor(lambdaBin);
2273 int dUdLbiasCutoff = hd.dUdLBiasCutoff;
2274 for (int iFL = -dUdLbiasCutoff; iFL <= dUdLbiasCutoff; iFL++) {
2275 int dUdLBin = currentdUdLBin + iFL;
2276
2277 double weight = mirrorFactor * getRecursionKernelValue(lambdaBin, dUdLBin);
2278 if (weight <= 0.0) {
2279 continue;
2280 }
2281
2282 double deltaFL = currentdUdLambda - dUdLforBin(dUdLBin);
2283 double deltaFL2 = deltaFL * deltaFL;
2284 double bias = weight * gaussianBiasMagnitude * expL2 * exp(-deltaFL2 / (2.0 * hd.dUdLVariance));
2285 gLdEdL += bias;
2286 dGdLambda -= deltaL / hd.lambdaVariance * bias;
2287 dGdFLambda -= deltaFL / hd.dUdLVariance * bias;
2288 }
2289 }
2290
2291 chainRule[0] = dGdLambda;
2292 chainRule[1] = dGdFLambda;
2293 return gLdEdL;
2294 }
2295
2296
2297
2298
2299
2300
2301
2302
2303
2304
2305
2306
2307
2308 private double energyAndGradient1D(boolean gradient) {
2309 return energyAndGradient1D(ld.lambda, gradient);
2310 }
2311
2312
2313
2314
2315
2316
2317
2318
2319
2320
2321
2322
2323
2324 private double energyAndGradient1D(double currentLambda, boolean gradient) {
2325 synchronized (this) {
2326 double biasEnergy = 0.0;
2327 for (int iL0 = 0; iL0 < hd.lambdaBins - 1; iL0++) {
2328 int iL1 = iL0 + 1;
2329
2330
2331 double L0 = iL0 * hd.lambdaBinWidth;
2332 double L1 = L0 + hd.lambdaBinWidth;
2333 double FL0 = ensembleAveragedUdL[iL0];
2334 double FL1 = ensembleAveragedUdL[iL1];
2335 double deltaFL = FL1 - FL0;
2336
2337
2338
2339
2340
2341 boolean done = false;
2342 if (currentLambda <= L1) {
2343 done = true;
2344 L1 = currentLambda;
2345 }
2346
2347
2348 biasEnergy += (FL0 * L1 + deltaFL * L1 * (0.5 * L1 - L0) / hd.lambdaBinWidth);
2349 biasEnergy -= (FL0 * L0 + deltaFL * L0 * (-0.5 * L0) / hd.lambdaBinWidth);
2350 if (done) {
2351
2352 if (gradient) {
2353 dUdLambda -= FL0 + (L1 - L0) * deltaFL / hd.lambdaBinWidth;
2354 }
2355 break;
2356 }
2357 }
2358 return -biasEnergy;
2359 }
2360 }
2361
2362
2363
2364
2365
2366
2367
2368 private double energyAndGradientMeta(boolean gradient) {
2369 return energyAndGradientMeta(ld.lambda, gradient);
2370 }
2371
2372
2373
2374
2375
2376
2377
2378
2379 private double energyAndGradientMeta(double currentLambda, boolean gradient) {
2380
2381 double biasEnergy = 0.0;
2382
2383
2384 int currentBin = indexForLambda(currentLambda);
2385
2386
2387 int lambdaBiasCutoff = hd.lambdaBiasCutoff;
2388 for (int iL = -lambdaBiasCutoff; iL <= lambdaBiasCutoff; iL++) {
2389 int lambdaBin = currentBin + iL;
2390 double deltaL = currentLambda - (lambdaBin * hd.lambdaBinWidth);
2391 double deltaL2 = deltaL * deltaL;
2392
2393 lambdaBin = lambdaMirror(lambdaBin);
2394 double mirrorFactor = mirrorFactor(lambdaBin);
2395 double weight = mirrorFactor * hd.biasMag * countsForLambda(lambdaBin);
2396 if (weight > 0) {
2397 double e = weight * exp(-deltaL2 / (2.0 * hd.lambdaVariance));
2398 biasEnergy += e;
2399
2400 if (gradient) {
2401 dUdLambda -= deltaL / hd.lambdaVariance * e;
2402 }
2403 }
2404 }
2405
2406 return biasEnergy;
2407 }
2408
2409
2410
2411
2412
2413
2414
2415 private double countsForLambda(int lambdaBin) {
2416 synchronized (this) {
2417 double count = 0.0;
2418 for (int i = 0; i < hd.dUdLBins; i++) {
2419 count += hd.zHistogram[lambdaBin][i];
2420 }
2421 return count;
2422 }
2423 }
2424
2425
2426
2427
2428 void checkRecursionKernelSize(double currentdUdL) {
2429
2430 synchronized (this) {
2431 if (currentdUdL > hd.dUdLMaximum) {
2432 logger.info(format(" Current F_lambda %8.2f > maximum histogram size %8.2f.", currentdUdL, hd.dUdLMaximum));
2433
2434 double origDeltaG = updateFreeEnergyDifference(false, false);
2435
2436 int newdUdLBins = hd.dUdLBins;
2437 while (hd.dUdLMinimum + newdUdLBins * hd.dUdLBinWidth < currentdUdL) {
2438 newdUdLBins += 100;
2439 }
2440 double[][] newRecursionKernel = new double[hd.lambdaBins][newdUdLBins];
2441
2442
2443 for (int i = 0; i < hd.lambdaBins; i++) {
2444 arraycopy(hd.zHistogram[i], 0, newRecursionKernel[i], 0, hd.dUdLBins);
2445 }
2446
2447 hd.zHistogram = newRecursionKernel;
2448 hd.dUdLBins = newdUdLBins;
2449 kernelValues = new double[hd.dUdLBins];
2450 hd.dUdLMaximum = hd.dUdLMinimum + hd.dUdLBinWidth * hd.dUdLBins;
2451 logger.info(format(" New histogram %8.2f to %8.2f with %d bins.\n", hd.dUdLMinimum, hd.dUdLMaximum, hd.dUdLBins));
2452
2453 double newFreeEnergy = updateFreeEnergyDifference(true, false);
2454 assert (origDeltaG == newFreeEnergy);
2455 }
2456 if (currentdUdL < hd.dUdLMinimum) {
2457 logger.info(format(" Current F_lambda %8.2f < minimum histogram size %8.2f.", currentdUdL, hd.dUdLMinimum));
2458
2459 double origDeltaG = updateFreeEnergyDifference(false, false);
2460
2461 int offset = 100;
2462 while (currentdUdL < hd.dUdLMinimum - offset * hd.dUdLBinWidth) {
2463 offset += 100;
2464 }
2465 int newdUdLBins = hd.dUdLBins + offset;
2466 double[][] newRecursionKernel = new double[hd.lambdaBins][newdUdLBins];
2467
2468
2469
2470 for (int i = 0; i < hd.lambdaBins; i++) {
2471 arraycopy(hd.zHistogram[i], 0, newRecursionKernel[i], offset, hd.dUdLBins);
2472 }
2473
2474 hd.zHistogram = newRecursionKernel;
2475 hd.dUdLMinimum = hd.dUdLMinimum - offset * hd.dUdLBinWidth;
2476 hd.dUdLBins = newdUdLBins;
2477 kernelValues = new double[hd.dUdLBins];
2478 logger.info(format(" New histogram %8.2f to %8.2f with %d bins.\n", hd.dUdLMinimum, hd.dUdLMaximum, hd.dUdLBins));
2479
2480 double newFreeEnergy = updateFreeEnergyDifference(true, false);
2481 assert (origDeltaG == newFreeEnergy);
2482 }
2483 }
2484 }
2485
2486
2487
2488
2489
2490
2491 void addBias(double dUdL) {
2492
2493 if (hd.asynchronous) {
2494 sendAsynchronous.send(ld.lambda, dUdL, temperingWeight);
2495 } else {
2496 sendSynchronous.send(ld.lambda, dUdL, temperingWeight);
2497 }
2498 }
2499
2500
2501
2502
2503 private void stochasticPreForce() {
2504
2505 double[] thetaPosition = lambdaState.x();
2506 thetaPosition[0] = ld.theta;
2507 stochastic.preForce(OrthogonalSpaceTempering.this);
2508 ld.theta = thetaPosition[0];
2509
2510
2511 if (ld.theta > PI) {
2512 ld.theta -= 2.0 * PI;
2513 } else if (ld.theta <= -PI) {
2514 ld.theta += 2.0 * PI;
2515 }
2516
2517
2518 double sinTheta = sin(ld.theta);
2519 ld.lambda = sinTheta * sinTheta;
2520 lambdaInterface.setLambda(ld.lambda);
2521 }
2522
2523
2524
2525
2526 private void stochasticPostForce() {
2527 double[] thetaPosition = lambdaState.x();
2528 double[] gradient = lambdaState.gradient();
2529
2530 gradient[0] = dUdLambda * sin(2 * ld.theta);
2531 thetaPosition[0] = ld.theta;
2532 stochastic.postForce(gradient);
2533
2534
2535 double[] thetaVelocity = lambdaState.v();
2536 double[] thetaAcceleration = lambdaState.a();
2537 ld.theta = thetaPosition[0];
2538 ld.thetaVelocity = thetaVelocity[0];
2539 ld.thetaAcceleration = thetaAcceleration[0];
2540
2541
2542 if (ld.theta > PI) {
2543 ld.theta -= 2.0 * PI;
2544 } else if (ld.theta <= -PI) {
2545 ld.theta += 2.0 * PI;
2546 }
2547
2548
2549 double sinTheta = sin(ld.theta);
2550 ld.lambda = sinTheta * sinTheta;
2551 lambdaInterface.setLambda(ld.lambda);
2552 }
2553
2554
2555
2556
2557
2558
2559
2560 double getLastReceivedLambda() {
2561 return lastReceivedLambda;
2562 }
2563
2564 public void setLastReceivedLambda(double lastReceivedLambda) {
2565 this.lastReceivedLambda = lastReceivedLambda;
2566 }
2567
2568
2569
2570
2571
2572
2573
2574 double getLastReceivedDUDL() {
2575 return lastReceiveddUdL;
2576 }
2577
2578 void destroy() {
2579 if (sendAsynchronous != null && sendAsynchronous.isAlive()) {
2580 double[] killMessage = new double[]{Double.NaN, Double.NaN, Double.NaN, Double.NaN};
2581 DoubleBuf killBuf = DoubleBuf.buffer(killMessage);
2582 try {
2583 logger.fine(" Sending the termination message.");
2584 world.send(rank, killBuf);
2585 logger.fine(" Termination message was sent successfully.");
2586 logger.fine(format(" Receive thread alive %b status %s", sendAsynchronous.isAlive(),
2587 sendAsynchronous.getState()));
2588 } catch (Exception ex) {
2589 String message = format(" Asynchronous Multiwalker OST termination signal "
2590 + "failed to be sent for process %d.", rank);
2591 logger.log(Level.SEVERE, message, ex);
2592 }
2593 } else {
2594 logger.fine(
2595 " CountReceiveThread was either not initialized, or is not alive. This is the case for the Histogram script.");
2596 }
2597 }
2598
2599 }
2600
2601 }