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.cli;
39  
40  import ffx.algorithms.AlgorithmListener;
41  import ffx.algorithms.dynamics.MDEngine;
42  import ffx.algorithms.dynamics.MolecularDynamics;
43  import ffx.algorithms.dynamics.integrators.Integrator;
44  import ffx.algorithms.dynamics.integrators.IntegratorEnum;
45  import ffx.algorithms.dynamics.thermostats.Thermostat;
46  import ffx.algorithms.dynamics.thermostats.ThermostatEnum;
47  import ffx.numerics.Potential;
48  import ffx.potential.MolecularAssembly;
49  import ffx.potential.cli.WriteoutOptions;
50  import ffx.utilities.Constants;
51  import picocli.CommandLine.ArgGroup;
52  import picocli.CommandLine.Option;
53  
54  import javax.annotation.Nullable;
55  import java.util.logging.Logger;
56  
57  /**
58   * Represents command line options for scripts that run molecular dynamics.
59   *
60   * @author Michael J. Schnieders
61   * @author Hernan V. Bernabe
62   * @since 1.0
63   */
64  public class DynamicsOptions {
65  
66    private static final Logger logger = Logger.getLogger(DynamicsOptions.class.getName());
67    /**
68     * Thermostat.
69     */
70    public ThermostatEnum thermostat;
71    /**
72     * Integrator.
73     */
74    public IntegratorEnum integrator;
75  
76    /**
77     * The ArgGroup keeps the DynamicsOptions together when printing help.
78     */
79    @ArgGroup(heading = "%n Dynamics Options%n", validate = false)
80    private final DynamicsOptionGroup group = new DynamicsOptionGroup();
81    private MDEngine engine = null;
82  
83    /**
84     * The restart save frequency in picoseconds (1.0 psec default).
85     *
86     * @return Restart file interval in picoseconds.
87     */
88    public double getCheckpoint() {
89      return group.checkpoint;
90    }
91  
92    /**
93     * The checkpoint frequency in steps.
94     *
95     * @param defaultFrequency The default frequency if the checkpoint interval is less than the time step.
96     * @return The checkpoint frequency in steps.
97     */
98    public int getCheckpointFrequency(int defaultFrequency) {
99      if (group.checkpoint > getDtPsec()) {
100       return (int) (group.checkpoint / getDtPsec());
101     }
102     return defaultFrequency;
103   }
104 
105   public void setCheckpoint(double checkpoint) {
106     group.checkpoint = checkpoint;
107   }
108 
109   /**
110    * The time step in femtoseconds (default of 1.0). A value of 2.0 is possible for the RESPA
111    * integrator.
112    *
113    * @return Time step in femtoseconds.
114    */
115   public double getDt() {
116     return group.dt;
117   }
118 
119   public double getDtPsec() {
120     return group.dt * Constants.FSEC_TO_PSEC;
121   }
122 
123   public void setDt(double dt) {
124     group.dt = dt;
125   }
126 
127   /**
128    * Initialize a MolecularDynamics from the parsed options.
129    *
130    * @param potential         a {@link ffx.numerics.Potential} object.
131    * @param activeAssembly    a {@link ffx.potential.MolecularAssembly} object.
132    * @param algorithmListener a {@link ffx.algorithms.AlgorithmListener} object.
133    * @param writeoutOptions   a {@link WriteoutOptions} object.
134    * @return a {@link MolecularDynamics} object.
135    */
136   public MolecularDynamics getDynamics(WriteoutOptions writeoutOptions, Potential potential,
137                                        MolecularAssembly activeAssembly, AlgorithmListener algorithmListener) {
138     return getDynamics(writeoutOptions, potential, activeAssembly, algorithmListener, engine);
139   }
140 
141   /**
142    * Initialize a MolecularDynamics from the parsed options.
143    *
144    * @param potential         a {@link ffx.numerics.Potential} object.
145    * @param activeAssembly    a {@link ffx.potential.MolecularAssembly} object.
146    * @param algorithmListener a {@link ffx.algorithms.AlgorithmListener} object.
147    * @param writeoutOptions   a {@link WriteoutOptions} object.
148    * @param requestedEngine   The requested engine (either FFX or OpenMM).
149    * @return a {@link MolecularDynamics} object.
150    */
151   public MolecularDynamics getDynamics(WriteoutOptions writeoutOptions, Potential potential,
152                                        MolecularAssembly activeAssembly, AlgorithmListener algorithmListener,
153                                        @Nullable MDEngine requestedEngine) {
154     MolecularDynamics molDyn;
155 
156     if (requestedEngine == null) {
157       molDyn = MolecularDynamics.dynamicsFactory(activeAssembly, potential, algorithmListener,
158           thermostat, integrator);
159     } else {
160       molDyn = MolecularDynamics.dynamicsFactory(activeAssembly, potential, algorithmListener,
161           thermostat, integrator, requestedEngine);
162     }
163     molDyn.setFileType(writeoutOptions.getFileType());
164     molDyn.setRestartFrequency(group.checkpoint);
165     molDyn.setIntervalSteps(group.trajSteps);
166 
167     return molDyn;
168   }
169 
170   public long getNumSteps() {
171     return group.steps;
172   }
173 
174   /**
175    * Getter for the field <code>optimize</code>.
176    *
177    * @return Whether to optimize structures.
178    */
179   public boolean getOptimize() {
180     return group.optimize;
181   }
182 
183   /**
184    * The thermodynamics reporting frequency in picoseconds (0.1 psec default).
185    *
186    * @return Thermodynamics logging interval in picoseconds.
187    */
188   public double getReport() {
189     return group.report;
190   }
191 
192   /**
193    * The molecular dynamics reporting frequency in steps.
194    *
195    * @param defaultFrequency The default frequency if the report interval is less than the time step.
196    * @return The reporting frequency in steps.
197    */
198   public int getReportFrequency(int defaultFrequency) {
199     if (group.report > getDtPsec()) {
200       return (int) (group.report / getDtPsec());
201     }
202     return defaultFrequency;
203   }
204 
205   public void setReport(double report) {
206     group.report = report;
207   }
208 
209   /**
210    * Write/snapshot appending interval.
211    *
212    * @return Interval between appending snapshots in psec.
213    */
214   public double getSnapshotInterval() {
215     return group.write;
216   }
217 
218   /**
219    * The simulation temperature (Kelvin).
220    *
221    * @return Temperature in Kelvins.
222    */
223   public double getTemperature() {
224     return group.temperature;
225   }
226 
227   public void setTemperature(double temperature) {
228     group.temperature = temperature;
229   }
230 
231   /**
232    * Parse the thermostat and integrator.
233    */
234   public void init() {
235     thermostat = Thermostat.parseThermostat(group.thermostatString);
236     integrator = Integrator.parseIntegrator(group.integratorString);
237     if (group.engineString != null) {
238       try {
239         engine = MDEngine.valueOf(group.engineString.toUpperCase());
240       } catch (Exception ex) {
241         logger.warning(String.format(
242             " Could not parse %s as a valid dynamics engine! Defaulting to the Platform-recommended engine.",
243             group.engineString));
244         engine = null;
245       }
246     }
247   }
248 
249   /**
250    * The desired thermostat: current choices are Adiabatic, Berendsen, or Bussi.
251    *
252    * @return Returns a String for the requested thermostat.
253    */
254   public String getThermostatString() {
255     return group.thermostatString;
256   }
257 
258   /**
259    * Set the thermostat.
260    *
261    * @param thermostat Thermostat to replace the requested one with.
262    */
263   public void setThermostat(ThermostatEnum thermostat) {
264     this.thermostat = thermostat;
265   }
266 
267   /**
268    * Set the integrator.
269    *
270    * @param integrator Integrator to replace the requested one with.
271    */
272   public void setIntegrator(IntegratorEnum integrator) {
273     this.integrator = integrator;
274   }
275 
276   public void setThermostatString(String thermostatString) {
277     group.thermostatString = thermostatString;
278   }
279 
280   /**
281    * The integrator: current choices are Beeman, RESPA, Stochastic (Langevin) or Verlet.
282    *
283    * @return Returns a String for the requested integrator.
284    */
285   public String getIntegratorString() {
286     return group.integratorString;
287   }
288 
289   public void setIntegratorString(String integratorString) {
290     group.integratorString = integratorString;
291   }
292 
293   /**
294    * Snapshot save frequency in picoseconds (1.0 psec default).
295    *
296    * @return Returns the frequency to save snapshots.
297    */
298   public double getWrite() {
299     return group.write;
300   }
301 
302   public void setWrite(double write) {
303     group.write = write;
304   }
305 
306   /**
307    * The number of molecular dynamics steps (default is 1,000,000).
308    *
309    * @return Returns the number of MD time steps.
310    */
311   public long getSteps() {
312     return group.steps;
313   }
314 
315   public void setSteps(long steps) {
316     group.steps = steps;
317   }
318 
319   /**
320    * Number of steps for each OpenMM MD cycle.
321    *
322    * @return Returns the number of steps for OpenMM MD cycles.
323    */
324   public int getTrajSteps() {
325     return group.trajSteps;
326   }
327 
328   public void setTrajSteps(int trajSteps) {
329     group.trajSteps = trajSteps;
330   }
331 
332   /**
333    * Saves low-energy snapshots discovered (only for single topology simulations).
334    *
335    * @return Returns true if low-energy snapshots should be saved.
336    */
337   public boolean isOptimize() {
338     return group.optimize;
339   }
340 
341   public void setOptimize(boolean optimize) {
342     group.optimize = optimize;
343   }
344 
345   /**
346    * The default engine choice for integrating the equations of motion
347    *
348    * @return Returns a String for the requested engine.
349    */
350   public String getEngineString() {
351     return group.engineString;
352   }
353 
354   public void setEngineString(String engineString) {
355     group.engineString = engineString;
356   }
357 
358   /**
359    * Collection of Dynamics Options.
360    */
361   private static class DynamicsOptionGroup {
362 
363     /**
364      * -d or --dt sets the time step in femtoseconds (default of 1.0). A value of 2.0 is possible for
365      * the RESPA integrator.
366      */
367     @Option(names = {"-d", "--dt"}, paramLabel = "1.0", defaultValue = "1.0",
368         description = "Time discretization step in femtoseconds.")
369     private double dt = 1.0;
370 
371     /**
372      * -b or --thermostat sets the desired thermostat: current choices are Adiabatic, Berendsen, or
373      * Bussi.
374      */
375     @Option(names = {"-b", "--thermostat"}, paramLabel = "Bussi", defaultValue = "Bussi",
376         description = "Thermostat: [Adiabatic / Berendsen / Bussi].")
377     private String thermostatString = "Bussi";
378 
379     /**
380      * -i or --integrator sets the desired integrator: current choices are Beeman, RESPA, Stochastic
381      * (i.e., Langevin dynamics) or Verlet.
382      */
383     @Option(names = {"-i", "--integrator"}, paramLabel = "Verlet", defaultValue = "Verlet",
384         description = "Integrator: [Beeman / Respa / Stochastic / Verlet].")
385     private String integratorString = "Verlet";
386 
387     /**
388      * -r or --report sets the thermodynamics reporting frequency in picoseconds (0.1 psec default).
389      */
390     @Option(names = {"-r", "--report"}, paramLabel = "0.25", defaultValue = "0.25",
391         description = "Interval in psec to report thermodynamics (psec).")
392     private double report = 0.25;
393 
394     /**
395      * -w or --write sets snapshot save frequency in picoseconds (1.0 psec default).
396      */
397     @Option(names = {"-w", "--write"}, paramLabel = "10.0", defaultValue = "10.0",
398         description = "Interval in psec to write out coordinates (psec).")
399     private double write = 10.0;
400 
401     /**
402      * -t or --temperature sets the simulation temperature (Kelvin).
403      */
404     @Option(names = {"-t", "--temperature"}, paramLabel = "298.15", defaultValue = "298.15",
405         description = "Temperature (Kelvin).")
406     private double temperature = 298.15;
407 
408     /**
409      * -n or --steps sets the number of molecular dynamics steps (default is 1 nsec).
410      */
411     @Option(names = {"-n", "--numberOfSteps"}, paramLabel = "1000000", defaultValue = "1000000",
412         description = "Number of molecular dynamics steps.")
413     private long steps = 1000000;
414 
415     /**
416      * -z or --trajSteps Number of steps for each OpenMM MD cycle.
417      */
418     @Option(names = {"-z", "--trajSteps"}, paramLabel = "100", defaultValue = "100",
419         description = "Number of steps per MD cycle (--mdE = OpenMM only).")
420     private int trajSteps = 100;
421 
422     /**
423      * -o or --optimize saves low-energy snapshots discovered (only for single topology simulations).
424      */
425     @Option(names = {"-o", "--optimize"}, defaultValue = "false",
426         description = "Optimize and save low-energy snapshots.")
427     private boolean optimize = false;
428 
429     /**
430      * -k or --checkpoint sets the restart save frequency in picoseconds (1.0 psec default).
431      */
432     @Option(names = {"-k", "--checkpoint"}, paramLabel = "1.0", defaultValue = "1.0",
433         description = "Interval in psec to write out restart files (.dyn, .his, etc).")
434     private double checkpoint = 1.0;
435 
436     /**
437      * --mdE or --molecularDynamicsEngine over-rides the default engine choice for integrating the
438      * equations of motion
439      */
440     @Option(names = {"--mdE", "--molecularDynamicsEngine"}, paramLabel = "FFX",
441         description = "Use FFX or OpenMM to integrate dynamics.")
442     private String engineString = "FFX";
443   }
444 }