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