View Javadoc
1   // ******************************************************************************
2   //
3   // Title:       Force Field X.
4   // Description: Force Field X - Software for Molecular Biophysics.
5   // Copyright:   Copyright (c) Michael J. Schnieders 2001-2024.
6   //
7   // This file is part of Force Field X.
8   //
9   // Force Field X is free software; you can redistribute it and/or modify it
10  // under the terms of the GNU General Public License version 3 as published by
11  // the Free Software Foundation.
12  //
13  // Force Field X is distributed in the hope that it will be useful, but WITHOUT
14  // ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
15  // FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
16  // details.
17  //
18  // You should have received a copy of the GNU General Public License along with
19  // Force Field X; if not, write to the Free Software Foundation, Inc., 59 Temple
20  // Place, Suite 330, Boston, MA 02111-1307 USA
21  //
22  // Linking this library statically or dynamically with other modules is making a
23  // combined work based on this library. Thus, the terms and conditions of the
24  // GNU General Public License cover the whole combination.
25  //
26  // As a special exception, the copyright holders of this library give you
27  // permission to link this library with independent modules to produce an
28  // executable, regardless of the license terms of these independent modules, and
29  // to copy and distribute the resulting executable under terms of your choice,
30  // provided that you also meet, for each linked independent module, the terms
31  // and conditions of the license of that module. An independent module is a
32  // module which is not derived from or based on this library. If you modify this
33  // library, you may extend this exception to your version of the library, but
34  // you are not obligated to do so. If you do not wish to do so, delete this
35  // exception statement from your version.
36  //
37  // ******************************************************************************
38  package ffx.algorithms.thermodynamics;
39  
40  import ffx.utilities.FileUtils;
41  import ffx.utilities.HistogramXmlAdapter;
42  import jakarta.xml.bind.JAXBContext;
43  import jakarta.xml.bind.Marshaller;
44  import jakarta.xml.bind.Unmarshaller;
45  import jakarta.xml.bind.annotation.XmlAccessOrder;
46  import jakarta.xml.bind.annotation.XmlAccessType;
47  import jakarta.xml.bind.annotation.XmlAccessorOrder;
48  import jakarta.xml.bind.annotation.XmlAccessorType;
49  import jakarta.xml.bind.annotation.XmlElement;
50  import jakarta.xml.bind.annotation.XmlRootElement;
51  import jakarta.xml.bind.annotation.adapters.XmlJavaTypeAdapter;
52  import org.apache.commons.configuration2.CompositeConfiguration;
53  
54  import javax.annotation.Nullable;
55  import java.io.File;
56  import java.util.logging.Logger;
57  
58  import static java.lang.Math.pow;
59  import static java.lang.Math.round;
60  import static java.lang.String.format;
61  
62  @XmlRootElement(name = "HistogramData")
63  @XmlAccessorOrder(XmlAccessOrder.ALPHABETICAL)
64  @XmlAccessorType(XmlAccessType.NONE)
65  public class HistogramData {
66  
67    private static final Logger logger = Logger.getLogger(HistogramData.class.getName());
68  
69    /**
70     * The default Gaussian bias magnitude is 0.05 kcal/mol.
71     */
72    private static final double DEFAULT_BIAS_MAG = 0.05;
73  
74    /**
75     * The default dU/dL bias cutoff is 5 bins.
76     */
77    private static final int DEFAULT_DUDL_BIAS_CUTOFF = 5;
78  
79    /**
80     * The default lambda bias cutoff is 5 bins for continuous lambda and 0 bins for discrete lambda.
81     */
82    private static final int DEFAULT_LAMBDA_BIAS_CUTOFF = 5;
83  
84    /**
85     * The default count interval is 10 MD steps.
86     */
87    private static final int DEFAULT_COUNT_INTERVAL = 10;
88  
89    /**
90     * The default lambda bin width is 0.005.
91     * There are (1.0 / width - 1) full size bins (e.g., 199 for width = 0.005).
92     * The first and last bins are half size (e.g., 0.0025 for width = 0.005).
93     */
94    private static final double DEFAULT_LAMBDA_BIN_WIDTH = 0.005;
95  
96    /**
97     * The default number of lambda bins is 201 based on DEFAULT_LAMBDA_BIN_WIDTH = 0.005.
98     * The first and last bin are of size DEFAULT_LAMBDA_BIN_WIDTH / 2.
99     * Bin 0:   0.0 to DEFAULT_LAMBDA_BIN_WIDTH / 2.
100    * Bin 200: 1.0 - DEFAULT_LAMBDA_BIN_WIDTH / 2 to 1.0.
101    * Bins 1 to 199 have size DEFAULT_LAMBDA_BIN_WIDTH.
102    */
103   private static final int DEFAULT_LAMBDA_BINS = 201;
104 
105   /**
106    * The default dU/dL bin width is 2.0 kcal/mol.
107    */
108   private static final double DEFAULT_DUDL_BIN_WIDTH = 2.0;
109 
110   /**
111    * Default number of dU/dL bins is 101.
112    */
113   private static final int DEFAULT_DUDL_BINS = 101;
114 
115   /**
116    * The default dU/dL range is from -DEFAULT_DUDL_BINS to DEFAULT_DUDL_BINS.
117    * This depends on the 1) default bin width of 2.0 kcal/mol and 2) the central bin range of -1.0 to 1.0.
118    * The default dU/dL minimum is -101 kcal/mol and the maximum is 101 kcal/mol.
119    */
120   private static final double DEFAULT_DUDL_MINIMUM = -DEFAULT_DUDL_BINS;
121 
122   /**
123    * Default is that meta-dynamics is not used.
124    */
125   private static final boolean DEFAULT_METADYNAMICS = false;
126 
127   /**
128    * The default tempering factor is 2.0 (in units of kBT).
129    */
130   private static final double DEFAULT_TEMPERING_FACTOR = 2.0;
131 
132   /**
133    * The default lambda reset value is at lambda = 0.99.
134    */
135   private static final double DEFAULT_RESET_HISTOGRAM_AT_LAMBDA = 0.99;
136 
137   /**
138    * The default reset statistics flag is false.
139    */
140   private static final boolean DEFAULT_RESET_HISTOGRAM = false;
141 
142   /**
143    * The default tempering offset is DEFAULT_BIAS_TO_TEMPERING_OFFSET * DEFAULT_BIAS_MAG = 1.0 kcal/mol.
144    */
145   private final double DEFAULT_BIAS_TO_TEMPERING_OFFSET = 20.0;
146 
147   /**
148    * By default, each walker contributes to the same histogram.
149    */
150   private final boolean DEFAULT_INDEPENDENT_WALKERS = false;
151 
152   /**
153    * By default, only the Rank 0 process writes the histogram.
154    * <p>
155    * Note that if independent walkers are used, each walker will write its own histogram.
156    */
157   private final boolean DEFAULT_INDEPENDENT_WRITE = false;
158 
159   /**
160    * By default, the lambda state variable is continuous. If discreteLambda is true, the lambda state is discrete.
161    */
162   private static final boolean DEFAULT_DISCRETE_LAMBDA = false;
163 
164   /**
165    * Magnitude of each hill (not including tempering). The default biasMag = 0.05 (kcal/mol).
166    */
167   @XmlElement(name = "BiasMag", defaultValue = "" + DEFAULT_BIAS_MAG)
168   double biasMag = DEFAULT_BIAS_MAG;
169 
170   /**
171    * When evaluating the biasing potential, contributions from a Gaussian centered on a bin more
172    * than lambdaBiasCutoff away will be neglected.
173    * The continuous lambda simulations, the default lambdaBiasCutoff = 5.
174    * For discrete lambda simulations, the default lambdaBiasCutoff = 0.
175    */
176   @XmlElement(name = "LambdaBiasCutoff", defaultValue = "" + DEFAULT_LAMBDA_BIAS_CUTOFF)
177   int lambdaBiasCutoff = DEFAULT_LAMBDA_BIAS_CUTOFF;
178 
179   /**
180    * When evaluating the biasing potential, contributions from a Gaussian centered on a bin more
181    * than dUdLBiasCutoff" away will be neglected. The default dUdLBiasCutoff = 5.
182    */
183   @XmlElement(name = "dUdLBiasCutoff", defaultValue = "" + DEFAULT_DUDL_BIAS_CUTOFF)
184   int dUdLBiasCutoff = DEFAULT_DUDL_BIAS_CUTOFF;
185 
186   /**
187    * Interval between adding a count to the Recursion kernel in MD steps. The default countInterval = 10.
188    */
189   @XmlElement(name = "CountInterval", defaultValue = "" + DEFAULT_COUNT_INTERVAL)
190   int countInterval = DEFAULT_COUNT_INTERVAL;
191 
192   /**
193    * Width of a lambda bin, or the distance between discrete lambda values.
194    * The default is 1.0 / (lambdaBins - 1).
195    */
196   @XmlElement(name = "LambdaBinWidth", defaultValue = "" + DEFAULT_LAMBDA_BIN_WIDTH)
197   double lambdaBinWidth = DEFAULT_LAMBDA_BIN_WIDTH;
198 
199   /**
200    * If true, use discrete lambda values instead of continuous lambda values.
201    */
202   @XmlElement(name = "DiscreteLambda", defaultValue = "" + DEFAULT_DISCRETE_LAMBDA)
203   boolean discreteLambda = DEFAULT_DISCRETE_LAMBDA;
204 
205   /**
206    * The width of the FLambda bin. The default is 2.0 (kcal/mol).
207    */
208   @XmlElement(name = "dUdLBinWidth", defaultValue = "" + DEFAULT_DUDL_BIN_WIDTH)
209   double dUdLBinWidth = DEFAULT_DUDL_BIN_WIDTH;
210 
211   /**
212    * The minimum value of the first dU/dL bin.
213    * Initially this is set to: mindUdL = -(dUdLBinWidth * dUdLBins) / 2.0.
214    * However, mindUdL may be updated as the count matrix is dynamically resized.
215    */
216   @XmlElement(name = "dUdLMinimum", defaultValue = "" + DEFAULT_DUDL_MINIMUM)
217   double dUdLMinimum = DEFAULT_DUDL_MINIMUM;
218 
219   /**
220    * Flag to indicate use of 1D meta-dynamics instead of OST.
221    */
222   @XmlElement(name = "MetaDynamics", defaultValue = "" + DEFAULT_METADYNAMICS)
223   boolean metaDynamics = DEFAULT_METADYNAMICS;
224 
225   /**
226    * The Dama et al. transition-tempering rate parameter. A reasonable value is about 2 to 8 kT,
227    * with larger values being resulting in slower decay. The default temperingFactor = 2.0.
228    */
229   @XmlElement(name = "TemperingFactor", defaultValue = "" + DEFAULT_TEMPERING_FACTOR)
230   double temperingFactor = DEFAULT_TEMPERING_FACTOR;
231 
232   /**
233    * An offset applied before calculating tempering weight.
234    * <p>
235    * First, for every Lambda, we calculate the maximum bias at that lambda by searching all
236    * populated dU/dL bins: maxdUdL(L) = max[ G(L,F_L) ] where the max operator searches over all F_L
237    * bins.
238    * <p>
239    * Then, the minimum bias coverage is determined by searching the maxdUdL(L) array over Lambda.
240    * minBias = min[ maxdUdL(L) ] where the min operator searches all entries in the array.
241    * <p>
242    * Then the temperOffset is applied to the minBias: biasHeight = max[minBias - temperingOffset, 0]
243    * <p>
244    * The default temperingOffset = 20x the Gaussian bias magnitude.
245    */
246   @XmlElement(name = "TemperingOffset", defaultValue = "-1.0")
247   private double temperingOffset = -1.0;
248 
249   /**
250    * Once the lambda reset value is reached, OST statistics are reset.
251    */
252   @XmlElement(name = "ResetHistogramAtLambda", defaultValue = "" + DEFAULT_RESET_HISTOGRAM_AT_LAMBDA)
253   double resetHistogramAtLambda = DEFAULT_RESET_HISTOGRAM_AT_LAMBDA;
254 
255   /**
256    * Flag set to false once OST statistics are reset at lambdaResetValue.
257    */
258   @XmlElement(name = "ResetHistogram", defaultValue = "" + DEFAULT_RESET_HISTOGRAM)
259   boolean resetHistogram = DEFAULT_RESET_HISTOGRAM;
260 
261   @XmlElement(name = "IndependentWalkers", defaultValue = "" + DEFAULT_INDEPENDENT_WALKERS)
262   boolean independentWalkers = DEFAULT_INDEPENDENT_WALKERS;
263 
264   @XmlElement(name = "IndependentWrite", defaultValue = "" + DEFAULT_INDEPENDENT_WRITE)
265   boolean independentWrite = DEFAULT_INDEPENDENT_WRITE;
266 
267   /**
268    * Flag to indicate if OST should send and receive counts between processes synchronously or
269    * asynchronously. The latter can be faster by ~40% because simulation with Lambda &gt; 0.75 must
270    * compute two condensed phase self-consistent fields to interpolate polarization.
271    */
272   @XmlElement(name = "Asynchronous", defaultValue = "" + true)
273   boolean asynchronous = true;
274 
275   /**
276    * It is useful to have an odd number of bins, so that there is a bin from FL=-dFL/2 to dFL/2 so
277    * that as FL approaches zero its contribution to thermodynamic integration goes to zero.
278    * <p>
279    * Otherwise a contribution of zero from a L bin can only result from equal sampling of the
280    * ranges -dFL to 0 and 0 to dFL.
281    * <p>
282    * The default FLambdaBins = 101.
283    */
284   @XmlElement(name = "dUdLBins", defaultValue = "" + DEFAULT_DUDL_BINS)
285   public int dUdLBins = DEFAULT_DUDL_BINS;
286 
287   /**
288    * The number of hills added to the recursion kernel.
289    */
290   @XmlElement(name = "Counts", defaultValue = "0")
291   public int counts = 0;
292 
293   /**
294    * The recursion kernel stores the weight of each [lambda][dU/dL] bin.
295    * Note that the variable name begins with Z so that its last in the histogram XML file.
296    */
297   @XmlElement(name = "Data", required = true)
298   @XmlJavaTypeAdapter(HistogramXmlAdapter.class)
299   double[][] zHistogram = new double[DEFAULT_LAMBDA_BINS][DEFAULT_DUDL_BINS];
300 
301   /**
302    * The histogram file to write to / read from.
303    */
304   private File histogramFile = null;
305 
306   /**
307    * The name of the histogram file.
308    */
309   private String histogramFileName = null;
310 
311   /**
312    * Flag indicating if this histogram data came from reading a file.
313    */
314   private boolean histogramRead = false;
315 
316   /**
317    * For continuous lambda: The first Lambda bin is centered on 0.0 (-0.0025 to 0.0025). The final
318    * Lambda bin is centered on 1.0 ( 0.9975 to 1.0025).
319    * <p>
320    * With this scheme, the maximum of biasing Gaussian hills is at exactly 0.0 and 1.0.
321    * <p>
322    * For discrete lambda: The first value of lambda is 0.0 and last value is 1.0.
323    * <p>
324    * The default lambdaBins = 201.
325    */
326   public int lambdaBins = DEFAULT_LAMBDA_BINS;
327 
328   /**
329    * Half the width of a lambda bin, or zero for discrete lambda values.
330    */
331   public double lambdaBinWidth_2 = DEFAULT_LAMBDA_BIN_WIDTH / 2.0;
332 
333   /**
334    * The minimum value of the first lambda bin.
335    * <p>
336    * minLambda = -dL_2 for continuous lambda.
337    * <p>
338    * minLambda = 0 for discrete lambda.
339    */
340   public double minLambda = -DEFAULT_LAMBDA_BIN_WIDTH / 2.0;
341 
342   /**
343    * Either the discrete lambda values used, or null (continuous lambda).
344    */
345   public double[] lambdaLadder = null;
346 
347   /**
348    * The variance for the Gaussian bias in the lambda dimension.
349    * lambdaVariance = 2.0 * dL * 2.0 * dL;
350    */
351   public double lambdaVariance = pow(2.0 * DEFAULT_LAMBDA_BIN_WIDTH, 2);
352 
353 
354   /**
355    * The maximum value of the last dUdL bin.
356    * <p>
357    * maxdUdL = mindUdL + dUdLBins * dUdLBinWidth.
358    * <p>
359    * default = -101 + 101 * 2.0 = 101 = DEFAULT_DUDL_BINS.
360    */
361   public double dUdLMaximum = DEFAULT_DUDL_BINS;
362 
363   /**
364    * Half the width of the F_lambda bin.
365    */
366   public double dUdLBinWidth_2 = DEFAULT_DUDL_BIN_WIDTH / 2.0;
367 
368   /**
369    * The variance for the Gaussian bias in the dU/dL dimension.
370    * dUdLVariance = 2.0 * dFL * 2.0 * dFL;
371    */
372   public double dUdLVariance = pow(2.0 * dUdLBinWidth, 2);
373 
374   /**
375    * Marshall the histogram data to a file.
376    *
377    * @param histogramData The Histogram data.
378    * @param file          The file to write to.
379    */
380   public static void writeHistogram(HistogramData histogramData, File file) {
381     try {
382       JAXBContext jaxbContext = JAXBContext.newInstance(HistogramData.class);
383       Marshaller jaxbMarshaller = jaxbContext.createMarshaller();
384       jaxbMarshaller.setProperty(Marshaller.JAXB_FORMATTED_OUTPUT, true);
385       jaxbMarshaller.marshal(histogramData, file);
386     } catch (Exception e) {
387       logger.warning(" Exception writing histogram:\n " + e);
388     }
389   }
390 
391   /**
392    * Unmarshall the histogram data.
393    *
394    * @param file The file to read from.
395    * @return The Histogram data.
396    */
397   public static HistogramData readHistogram(@Nullable File file) {
398     if (file != null && file.exists() && file.canRead()) {
399       try {
400         JAXBContext jaxbContext = JAXBContext.newInstance(HistogramData.class);
401         Unmarshaller unmarshaller = jaxbContext.createUnmarshaller();
402         HistogramData histogramData = (HistogramData) unmarshaller.unmarshal(file);
403         histogramData.setHistogramFile(file);
404         histogramData.setHistogramRead(true);
405         histogramData.rectify();
406         return histogramData;
407       } catch (Exception e) {
408         logger.warning(" Exception reading histogram:\n " + e);
409       }
410     }
411     HistogramData histogramData = new HistogramData();
412     histogramData.setHistogramFile(file);
413     histogramData.rectify();
414     return histogramData;
415   }
416 
417   private void rectify() {
418     // Check that histogram size matches the number of lambda and dU/dL bins.
419     int hisLambdaBins = zHistogram.length;
420     if (lambdaBins != hisLambdaBins) {
421       logger.fine(" Current LambdaBins: " + lambdaBins);
422       logger.fine(" Updated to:         " + hisLambdaBins);
423       lambdaBins = hisLambdaBins;
424       lambdaBinWidth = 1.0 / (lambdaBins - 1);
425     }
426     int hisDUDLBins = zHistogram[0].length;
427     if (dUdLBins != hisDUDLBins) {
428       logger.fine(" Current dUdLBins:   " + dUdLBins);
429       logger.fine(" Updated to:         " + hisDUDLBins);
430       dUdLBins = hisDUDLBins;
431     }
432     rectifyLambdaVariables();
433     rectifyDUDLVariables();
434   }
435 
436   /**
437    * All transient variables that depend on dUdLBinWidth are set here.
438    */
439   private void rectifyDUDLVariables() {
440     // The center of the central bin is at 0.
441     dUdLBinWidth_2 = dUdLBinWidth / 2.0;
442     dUdLMaximum = dUdLMinimum + (dUdLBins * dUdLBinWidth);
443     dUdLVariance = pow(2.0 * dUdLBinWidth, 2);
444   }
445 
446   /**
447    * All transient variables that depend on lambdaBinWidth or discreteLambda are set here.
448    */
449   private void rectifyLambdaVariables() {
450     // Get the closest number of full width bins based on the current lambdaBinWidth.
451     int n = (int) round(1.0 / lambdaBinWidth);
452     // Now set the lambdaBinWidth exactly.
453     lambdaBinWidth = 1.0 / n;
454     // Divide a bin into two half size bins (the first and last bins are centered on 0 and 1, respectively).
455     lambdaBins = n + 1;
456     if (discreteLambda) {
457       lambdaLadder = new double[lambdaBins];
458       lambdaLadder[0] = 0.0;
459       lambdaLadder[lambdaBins - 1] = 1.0;
460       for (int i = 1; i < lambdaBins - 1; i++) {
461         lambdaLadder[i] = lambdaBinWidth * i;
462       }
463       lambdaBinWidth_2 = 0.0;
464       minLambda = 0.0;
465       lambdaBiasCutoff = 0;
466     } else {
467       lambdaLadder = null;
468       lambdaBinWidth_2 = lambdaBinWidth / 2.0;
469       minLambda = -lambdaBinWidth_2;
470       // Leave the lambda bias cutoff at it current value, unless its zero.
471       if (lambdaBiasCutoff == 0) {
472         // In this case, set it back to the default value.
473         lambdaBiasCutoff = DEFAULT_LAMBDA_BIAS_CUTOFF;
474       }
475     }
476     lambdaVariance = pow(2.0 * lambdaBinWidth, 2);
477   }
478 
479   public void setCountInterval(int countInterval) {
480     this.countInterval = countInterval;
481   }
482 
483   /**
484    * Marshall this histogram data to a file.
485    */
486   public void writeHistogram() {
487     writeHistogram(this, histogramFile);
488   }
489 
490   public String getHistogramFileName() {
491     return histogramFileName;
492   }
493 
494   public File getHistogramFile() {
495     return histogramFile;
496   }
497 
498   public void setHistogramFile(@Nullable File histogramFile) {
499     this.histogramFile = histogramFile;
500     if (histogramFile != null) {
501       histogramFileName = FileUtils.relativePathTo(histogramFile).toString();
502     }
503   }
504 
505   public boolean wasHistogramRead() {
506     return histogramRead;
507   }
508 
509   public void setHistogramRead(boolean histogramRead) {
510     this.histogramRead = histogramRead;
511   }
512 
513   /**
514    * Sets the value of independentWalkers; if true, it also sets writeIndependent to true.
515    *
516    * @param independentWalkers Value to set independentWalkers to.
517    */
518   public void setIndependentWalkers(boolean independentWalkers) {
519     this.independentWalkers = independentWalkers;
520     if (independentWalkers) {
521       setWriteIndependent(true);
522     }
523   }
524 
525   /**
526    * Sets the value of writeIndependent.
527    *
528    * @param independentWrite Value to set writeIndependent to.
529    * @throws IllegalArgumentException If independentWrite is false and independentWalkers is true.
530    */
531   public void setWriteIndependent(boolean independentWrite) throws IllegalArgumentException {
532     if (!independentWrite && independentWalkers) {
533       throw new IllegalArgumentException(" Independent walkers implies independent writing.");
534     }
535     this.independentWrite = independentWrite;
536   }
537 
538   /**
539    * Returns the value of independentWrite.
540    *
541    * @return independentWrite.
542    */
543   public boolean independentWrite() {
544     return independentWrite;
545   }
546 
547   /**
548    * Gets the Gaussian bias magnitude in kcal/mol.
549    *
550    * @return Gaussian bias magnitude in kcal/mol.
551    */
552   public double getBiasMag() {
553     return biasMag;
554   }
555 
556   /**
557    * Set the bias magnitude. If temperingOffset has not been explicitly set,
558    * the temperOffset is set to bias magnitude * 20.
559    *
560    * @param biasMag Gaussian bias magnitude in kcal/mol
561    */
562   public void setBiasMag(double biasMag) {
563     this.biasMag = biasMag;
564   }
565 
566   public void setTemperingFactor(double temperingFactor) {
567     this.temperingFactor = temperingFactor;
568   }
569 
570   public void setMetaDynamics(boolean metaDynamics) {
571     this.metaDynamics = metaDynamics;
572   }
573 
574   public double getLambdaBinWidth() {
575     return lambdaBinWidth;
576   }
577 
578   public int getLambdaBins() {
579     return lambdaBins;
580   }
581 
582   public int getDUDLBins() {
583     return dUdLBins;
584   }
585 
586   /**
587    * Sets lambdaBinWidth; if an invalid value is provided (not 0-1), resets it to default 0.005.
588    *
589    * @param lambdaBinWidth Lambda bin width.
590    */
591   public void setLambdaBinWidth(double lambdaBinWidth) {
592     if (lambdaBinWidth <= 0.0 || lambdaBinWidth > 1.0) {
593       this.lambdaBinWidth = DEFAULT_LAMBDA_BIN_WIDTH;
594     } else {
595       this.lambdaBinWidth = lambdaBinWidth;
596     }
597     rectifyLambdaVariables();
598   }
599 
600   /**
601    * Gets the tempering offset.
602    *
603    * @return Gets the tempering offset in kcal/mol.
604    */
605   public double getTemperingOffset() {
606     if (temperingOffset >= 0.0) {
607       return temperingOffset;
608     }
609     return biasMag * DEFAULT_BIAS_TO_TEMPERING_OFFSET;
610   }
611 
612   /**
613    * Sets an explicit tempering offset in kcal/mol.
614    *
615    * @param temperingOffset Tempering offset in kcal/mol.
616    */
617   public void setTemperingOffset(double temperingOffset) {
618     this.temperingOffset = temperingOffset;
619     this.temperingOffset = temperingOffset;
620   }
621 
622   public void setAsynchronous(boolean asynchronous) {
623     this.asynchronous = asynchronous;
624   }
625 
626   public double resetHistogramAtLambda() {
627     return resetHistogramAtLambda;
628   }
629 
630   /**
631    * Set the lambda value at which to reset the histogram.
632    *
633    * @param resetHistogramAtLambda The lambda value at which to reset the histogram.
634    */
635   public void setResetHistogramAtLambda(double resetHistogramAtLambda) {
636     if (resetHistogramAtLambda < 0.0 || resetHistogramAtLambda > 1.0) {
637       this.resetHistogramAtLambda = DEFAULT_RESET_HISTOGRAM_AT_LAMBDA;
638       this.resetHistogram = false;
639     } else {
640       this.resetHistogramAtLambda = resetHistogramAtLambda;
641       this.resetHistogram = true;
642     }
643   }
644 
645   public void setDiscreteLambda(boolean discreteLambda) {
646     this.discreteLambda = discreteLambda;
647     rectifyLambdaVariables();
648   }
649 
650   public void applyProperties(CompositeConfiguration properties) {
651     if (properties.containsKey("lambda-bias-cutoff")) {
652       lambdaBiasCutoff = properties.getInt("lambda-bias-cutoff");
653     }
654     if (properties.containsKey("dudl-bias-cutoff")) {
655       dUdLBiasCutoff = properties.getInt("dudl-bias-cutoff");
656     }
657     if (properties.containsKey("ost-bias-mag")) {
658       double bias = properties.getDouble("ost-bias-mag");
659       setBiasMag(bias);
660     }
661     if (properties.containsKey("ost-temperOffset")) {
662       double offset = properties.getDouble("ost-temperOffset");
663       setTemperingOffset(offset);
664     }
665     if (properties.containsKey("lambda-bin-width")) {
666       double width = properties.getDouble("lambda-bin-width");
667       setLambdaBinWidth(width);
668     }
669     if (properties.containsKey("flambda-bin-width")) {
670       dUdLBinWidth = properties.getDouble("flambda-bin-width");
671     }
672     if (properties.containsKey("discrete-lambda")) {
673       boolean discreteLambda = properties.getBoolean("discrete-lambda");
674       setDiscreteLambda(discreteLambda);
675     }
676   }
677 
678   public String toString() {
679     return format("  Lambda bins:      %6d\n", lambdaBins)
680         + format("  Lambda bin width: %6.3f\n", lambdaBinWidth)
681         + format("  dU/dL bins:       %6d\n", dUdLBins)
682         + format("  dU/dL bin width:  %6.3f (kcal/mol)\n", dUdLBinWidth)
683         + format("  Bias magnitude:   %6.3f (kcal/mol)\n", biasMag)
684         + format("  Tempering offset: %6.3f (kcal/mol)\n", getTemperingOffset())
685         + format("  Tempering rate:   %6.3f\n", temperingFactor)
686         + format("  Discrete lambda:  %6B\n", discreteLambda)
687         + format("  Meta-dynamics:    %6B\n\n", metaDynamics);
688   }
689 
690 }