2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 /*! \libinternal \file
36 * \brief Header for the code that writes energy-like quantities.
38 * \author Mark Abraham <mark.j.abraham@gmail.com>
39 * \author Paul Bauer <paul.bauer.q@gmail.com>
40 * \author Artem Zhmurov <zhmurov@gmail.com>
42 * \ingroup module_mdlib
45 #ifndef GMX_MDLIB_ENERGYOUTPUT_H
46 #define GMX_MDLIB_ENERGYOUTPUT_H
50 #include "gromacs/mdtypes/enerdata.h"
52 class energyhistory_t;
54 struct gmx_ekindata_t;
55 struct gmx_enerdata_t;
56 struct SimulationGroups;
58 struct gmx_output_env_t;
68 struct t_mde_delta_h_coll;
76 //! \brief Printed names for intergroup energies
77 extern const char *egrp_nm[egNR+1];
79 /* \brief delta_h block type enum: the kinds of energies written out. */
82 //! Delta H BAR energy difference
84 //! dH/dlambda derivative
90 //! Expanded ensemble statistics
92 //! Total number of energy types in this enum
99 /* Energy output class
101 * This is the collection of energy averages collected during mdrun, and to
102 * be written out to the .edr file.
104 * \todo Use more std containers.
105 * \todo Remove GMX_CONSTRAINTVIR
106 * \todo Write free-energy output also to energy file (after adding more tests)
112 /*! \brief Initiate MD energy bin
114 * \param[in] fp_ene Energy output file.
115 * \param[in] mtop Topology.
116 * \param[in] ir Input parameters.
117 * \param[in] pull_work Pulling simulations data
118 * \param[in] fp_dhdl FEP file.
119 * \param[in] isRerun Is this is a rerun instead of the simulations.
121 EnergyOutput(ener_file *fp_ene,
122 const gmx_mtop_t *mtop,
123 const t_inputrec *ir,
124 const pull_t *pull_work,
130 /*! \brief Update the averaging structures.
132 * Called every step on which the thermodynamic values are evaluated.
134 * \param[in] bDoDHDL Whether the FEP is enabled.
135 * \param[in] bSum If this stepshould be recorded to compute sums and averaes.
136 * \param[in] time Current simulation time.
137 * \param[in] tmass Total mass
138 * \param[in] enerd Energy data object.
139 * \param[in] state System state.
140 * \param[in] fep FEP data.
141 * \param[in] expand Expanded ensemble (for FEP).
142 * \param[in] lastbox PBC data.
143 * \param[in] svir Constraint virial.
144 * \param[in] fvir Force virial.
145 * \param[in] vir Total virial.
146 * \param[in] pres Pressure.
147 * \param[in] ekind Kinetic energy data.
148 * \param[in] mu_tot Total dipole.
149 * \param[in] constr Constraints object to get RMSD from (for LINCS only).
151 void addDataAtEnergyStep(bool bDoDHDL,
155 gmx_enerdata_t *enerd,
164 gmx_ekindata_t *ekind,
166 const gmx::Constraints *constr);
168 /*! \brief Update the data averaging structure counts.
170 * Updates the number of steps, the values have not being computed.
172 void recordNonEnergyStep();
174 /*! \brief Writes current quantites to log and energy files.
176 * Prints current values of energies, pressure, temperature, restraint
177 * data, etc. to energy output file and to the log file (if not nullptr).
179 * This function only does something useful when bEne || bDR || bOR || log.
181 * \todo Perhaps this responsibility should involve some other
182 * object visiting all the contributing objects.
184 * \param[in] fp_ene Energy file for the output.
185 * \param[in] bEne If it is a step for energy output or last step.
186 * \param[in] bDR If it is a step of writing distance restraints.
187 * \param[in] bOR If it is a step of writing orientation restraints.
188 * \param[in] log Pointer to the log file.
189 * \param[in] step Current step.
190 * \param[in] time Current simulation time.
191 * \param[in] fcd Bonded force computation data,
192 * including orientation and distance restraints.
193 * \param[in] awh AWH data.
195 void printStepToEnergyFile(ener_file *fp_ene, bool bEne, bool bDR, bool bOR,
197 int64_t step, double time,
201 /*! \brief Print reference temperatures for annealing groups.
203 * Nothing is done if log is nullptr.
205 * \param[in] log Log file to print to.
206 * \param[in] groups Information on atom groups.
207 * \param[in] opts Atom temperature coupling groups options
208 * (annealing is done by groups).
210 void printAnnealingTemperatures(FILE *log, SimulationGroups *groups, t_grpopts *opts);
212 /*! \brief Prints average values to log file.
214 * This is called at the end of the simulation run to print accumulated average values.
215 * Nothing it done if log is nullptr.
217 * \param[in] log Where to print.
218 * \param[in] groups Atom groups.
220 void printAverages(FILE *log, const SimulationGroups *groups);
222 /*! \brief Get the number of thermodynamic terms recorded.
224 * \todo Refactor this to return the expected output size,
225 * rather than exposing the implementation details about
226 * thermodynamic terms.
228 * \returns Number of thermodynamic terms.
230 int numEnergyTerms() const;
232 /*! \brief Fill the energyhistory_t data.
234 * Between .edr writes, the averages are history dependent,
235 * and that history needs to be retained in checkpoints.
236 * These functions set/read the energyhistory_t class
237 * that is written to checkpoints.
239 * \param[out] enerhist Energy history data structure.
241 void fillEnergyHistory(energyhistory_t * enerhist) const;
243 /*! \brief Restore from energyhistory_t data.
245 * \param[in] enerhist Energy history data structure.
247 void restoreFromEnergyHistory(const energyhistory_t &enerhist);
249 //! Print an output header to the log file.
250 void printHeader(FILE *log, int64_t steps, double time);
256 //! Structure to store energy components and their running averages
257 t_ebin *ebin_ = nullptr;
259 //! Is the periodic box triclinic
260 bool bTricl_ = false;
261 //! NHC trotter is used
262 bool bNHC_trotter_ = false;
263 //! If Nose-Hoover chains should be printed
264 bool bPrintNHChains_ = false;
265 //! If MTTK integrator was used
268 //! Temperature control scheme
271 //! Which of the main energy terms should be printed
272 bool bEner_[F_NRE] = {false};
273 //! Index for main energy terms
275 //! Number of energy terms from F_NRE list to be saved (i.e. number of 'true' in bEner)
278 //! Index for constraints RMSD
280 /* !\brief How many constraints RMSD terms there are.
281 * Usually 1 if LINCS is used and 0 otherwise)
282 * nCrmsd > 0 indicates when constraints RMSD is saves, hence no boolean
286 //! Is the periodic box dynamic
287 bool bDynBox_ = false;
288 //! Index for box dimensions
290 //! Index for box volume
292 //! Index for density
294 //! Triclinic box and not a rerun
295 bool bDiagPres_ = false;
296 //! Reference pressure, averaged over dimensions
298 //! Index for thermodynamic work (pV)
300 //! Index for entalpy (pV + total energy)
303 /*! \brief If the constraints virial should be printed.
304 * Can only be true if "GMX_CONSTRAINTVIR" environmental variable is set */
305 bool bConstrVir_ = false;
306 //! Index for constrains virial
308 //! Index for force virial
311 //! If we have pressure computed
313 //! Index for total virial
315 //! Index for pressure
317 /*! \brief Index for surface tension
318 * [(pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ]]*/
321 //! Pressure control scheme
323 //! Index for velocity of the box borders
326 //! If dipole was calculated
328 //! Index for the dipole
331 //! Index for coseine acceleration used for viscocity calculation
333 //! Index for viscocity
336 //! Which energy terms from egNR list should be printed in group-to-group block
337 bool bEInd_[egNR] = {false};
338 //! Number of energy terms to be printed (these, for which bEInd[] == true)
340 //! Number of energy output groups
342 //! Number of intergroup energy sets to be printed for each energy term (nE = (nEg*(nEg+1))/2)
344 //! Indexes for integroup energy sets (each set with nEc energies)
345 int *igrp_ = nullptr;
347 //! Number of temperature coupling groups
349 //! Index for temperature
352 //! Number of Nose-Hoover chains
354 //! Number of temperature coupling coefficients in case of NH Chains
356 //! Index for temperature coupling coefficient in case of NH chains
359 //! Number of temperature coupling terms if the temperature control is dealt by barostat (MTTK case)
361 //! Scalling factors for temperaturs control in MTTK
363 //! Index for scalling factor of MTTK
366 //! Number of acceleration groups
368 //! Index for group velocities
371 //! Array to accumulate values during update
372 real *tmp_r_ = nullptr;
373 //! Array to accumulate values during update
374 rvec *tmp_v_ = nullptr;
376 //! The dhdl.xvg output file
377 FILE *fp_dhdl_ = nullptr;
378 //! Energy components for dhdl.xvg output
379 double *dE_ = nullptr;
380 //! The delta U components (raw data + histogram)
381 t_mde_delta_h_coll *dhc_ = nullptr;
382 //! Temperatures for simulated tempering groups
383 real *temperatures_ = nullptr;
384 //! Number of temperatures actually saved
385 int numTemperatures_ = 0;
390 //! Open the dhdl file for output
391 FILE *open_dhdl(const char *filename, const t_inputrec *ir,
392 const gmx_output_env_t *oenv);