2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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;
74 struct MdModulesNotifier;
75 enum class StartingBehavior;
78 //! \brief Printed names for intergroup energies
79 extern const char* egrp_nm[egNR + 1];
81 /* \brief delta_h block type enum: the kinds of energies written out. */
84 //! Delta H BAR energy difference
86 //! dH/dlambda derivative
92 //! Expanded ensemble statistics
94 //! Total number of energy types in this enum
102 * \brief Arrays connected to Pressure and Temperature coupling
104 struct PTCouplingArrays
106 //! Box velocities for Parrinello-Rahman P-coupling.
108 //! Nose-Hoover coordinates.
109 ArrayRef<const double> nosehoover_xi;
110 //! Nose-Hoover velocities.
111 ArrayRef<const double> nosehoover_vxi;
112 //! Pressure Nose-Hoover coordinates.
113 ArrayRef<const double> nhpres_xi;
114 //! Pressure Nose-Hoover velocities.
115 ArrayRef<const double> nhpres_vxi;
118 /* Energy output class
120 * This is the collection of energy averages collected during mdrun, and to
121 * be written out to the .edr file.
123 * \todo Use more std containers.
124 * \todo Remove GMX_CONSTRAINTVIR
125 * \todo Write free-energy output also to energy file (after adding more tests)
130 /*! \brief Initiate MD energy bin
132 * \param[in] fp_ene Energy output file.
133 * \param[in] mtop Topology.
134 * \param[in] ir Input parameters.
135 * \param[in] pull_work Pulling simulations data
136 * \param[in] fp_dhdl FEP file.
137 * \param[in] isRerun Is this is a rerun instead of the simulations.
138 * \param[in] startingBehavior Run starting behavior.
139 * \param[in] mdModulesNotifier Notifications to MD modules.
141 EnergyOutput(ener_file* fp_ene,
142 const gmx_mtop_t* mtop,
143 const t_inputrec* ir,
144 const pull_t* pull_work,
147 StartingBehavior startingBehavior,
148 const MdModulesNotifier& mdModulesNotifier);
152 /*! \brief Update the averaging structures.
154 * Called every step on which the thermodynamic values are evaluated.
156 * \param[in] bDoDHDL Whether the FEP is enabled.
157 * \param[in] bSum If this stepshould be recorded to compute sums and averages.
158 * \param[in] time Current simulation time.
159 * \param[in] tmass Total mass
160 * \param[in] enerd Energy data object.
161 * \param[in] fep FEP data.
162 * \param[in] expand Expanded ensemble (for FEP).
163 * \param[in] lastbox PBC data.
164 * \param[in] ptCouplingArrays Arrays connected to pressure and temperature coupling.
165 * \param[in] fep_state The current alchemical state we are in.
166 * \param[in] svir Constraint virial.
167 * \param[in] fvir Force virial.
168 * \param[in] vir Total virial.
169 * \param[in] pres Pressure.
170 * \param[in] ekind Kinetic energy data.
171 * \param[in] mu_tot Total dipole.
172 * \param[in] constr Constraints object to get RMSD from (for LINCS only).
174 void addDataAtEnergyStep(bool bDoDHDL,
178 const gmx_enerdata_t* enerd,
180 const t_expanded* expand,
181 const matrix lastbox,
182 PTCouplingArrays ptCouplingArrays,
188 const gmx_ekindata_t* ekind,
190 const gmx::Constraints* constr);
192 /*! \brief Update the data averaging structure counts.
194 * Updates the number of steps, the values have not being computed.
196 void recordNonEnergyStep();
198 /*! \brief Writes current quantites to log and energy files.
200 * Prints current values of energies, pressure, temperature, restraint
201 * data, etc. to energy output file and to the log file (if not nullptr).
203 * This function only does something useful when bEne || bDR || bOR || log.
205 * \todo Perhaps this responsibility should involve some other
206 * object visiting all the contributing objects.
208 * \param[in] fp_ene Energy file for the output.
209 * \param[in] bEne If it is a step for energy output or last step.
210 * \param[in] bDR If it is a step of writing distance restraints.
211 * \param[in] bOR If it is a step of writing orientation restraints.
212 * \param[in] log Pointer to the log file.
213 * \param[in] step Current step.
214 * \param[in] time Current simulation time.
215 * \param[in] fcd Bonded force computation data,
216 * including orientation and distance restraints.
217 * \param[in] awh AWH data.
219 void printStepToEnergyFile(ener_file* fp_ene,
229 /*! \brief Print reference temperatures for annealing groups.
231 * Nothing is done if log is nullptr.
233 * \param[in] log Log file to print to.
234 * \param[in] groups Information on atom groups.
235 * \param[in] opts Atom temperature coupling groups options
236 * (annealing is done by groups).
238 static void printAnnealingTemperatures(FILE* log, const SimulationGroups* groups, t_grpopts* opts);
240 /*! \brief Prints average values to log file.
242 * This is called at the end of the simulation run to print accumulated average values.
243 * Nothing it done if log is nullptr.
245 * \param[in] log Where to print.
246 * \param[in] groups Atom groups.
248 void printAverages(FILE* log, const SimulationGroups* groups);
250 /*! \brief Get the number of thermodynamic terms recorded.
252 * \todo Refactor this to return the expected output size,
253 * rather than exposing the implementation details about
254 * thermodynamic terms.
256 * \returns Number of thermodynamic terms.
258 int numEnergyTerms() const;
260 /*! \brief Fill the energyhistory_t data.
262 * Between .edr writes, the averages are history dependent,
263 * and that history needs to be retained in checkpoints.
264 * These functions set/read the energyhistory_t class
265 * that is written to checkpoints.
267 * \param[out] enerhist Energy history data structure.
269 void fillEnergyHistory(energyhistory_t* enerhist) const;
271 /*! \brief Restore from energyhistory_t data.
273 * \param[in] enerhist Energy history data structure.
275 void restoreFromEnergyHistory(const energyhistory_t& enerhist);
277 //! Print an output header to the log file.
278 static void printHeader(FILE* log, int64_t steps, double time);
284 //! Structure to store energy components and their running averages
285 t_ebin* ebin_ = nullptr;
287 //! Is the periodic box triclinic
288 bool bTricl_ = false;
289 //! NHC trotter is used
290 bool bNHC_trotter_ = false;
291 //! If Nose-Hoover chains should be printed
292 bool bPrintNHChains_ = false;
293 //! If MTTK integrator was used
296 //! Temperature control scheme
299 //! Which of the main energy terms should be printed
300 bool bEner_[F_NRE] = { false };
301 //! Index for main energy terms
303 //! Number of energy terms from F_NRE list to be saved (i.e. number of 'true' in bEner)
306 //! Index for constraints RMSD
308 /* !\brief How many constraints RMSD terms there are.
309 * Usually 1 if LINCS is used and 0 otherwise)
310 * nCrmsd > 0 indicates when constraints RMSD is saves, hence no boolean
314 //! Is the periodic box dynamic
315 bool bDynBox_ = false;
316 //! Index for box dimensions
318 //! Index for box volume
320 //! Index for density
322 //! Triclinic box and not a rerun
323 bool bDiagPres_ = false;
324 //! Reference pressure, averaged over dimensions
326 //! Index for thermodynamic work (pV)
328 //! Index for entalpy (pV + total energy)
331 /*! \brief If the constraints virial should be printed.
332 * Can only be true if "GMX_CONSTRAINTVIR" environmental variable is set */
333 bool bConstrVir_ = false;
334 //! Index for constrains virial
336 //! Index for force virial
339 //! If we have pressure computed
341 //! Index for total virial
343 //! Index for pressure
345 /*! \brief Index for surface tension
346 * [(pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ]]*/
349 //! Pressure control scheme
351 //! Index for velocity of the box borders
354 //! If dipole was calculated
356 //! Index for the dipole
359 //! Index for coseine acceleration used for viscocity calculation
361 //! Index for viscocity
364 //! Which energy terms from egNR list should be printed in group-to-group block
365 bool bEInd_[egNR] = { false };
366 //! Number of energy terms to be printed (these, for which bEInd[] == true)
368 //! Number of energy output groups
370 //! Number of intergroup energy sets to be printed for each energy term (nE = (nEg*(nEg+1))/2)
372 //! Indexes for integroup energy sets (each set with nEc energies)
373 int* igrp_ = nullptr;
375 //! Number of temperature coupling groups
377 //! Index for temperature
380 //! Number of Nose-Hoover chains
382 //! Number of temperature coupling coefficients in case of NH Chains
384 //! Index for temperature coupling coefficient in case of NH chains
387 //! Number of temperature coupling terms if the temperature control is dealt by barostat (MTTK case)
389 //! Scalling factors for temperaturs control in MTTK
391 //! Index for scalling factor of MTTK
394 //! Number of acceleration groups
396 //! Index for group velocities
399 //! Array to accumulate values during update
400 real* tmp_r_ = nullptr;
401 //! Array to accumulate values during update
402 rvec* tmp_v_ = nullptr;
404 //! The dhdl.xvg output file
405 FILE* fp_dhdl_ = nullptr;
406 //! Energy components for dhdl.xvg output
407 double* dE_ = nullptr;
408 //! The delta U components (raw data + histogram)
409 t_mde_delta_h_coll* dhc_ = nullptr;
410 //! Temperatures for simulated tempering groups
411 real* temperatures_ = nullptr;
412 //! Number of temperatures actually saved
413 int numTemperatures_ = 0;
418 //! Open the dhdl file for output
419 FILE* open_dhdl(const char* filename, const t_inputrec* ir, const gmx_output_env_t* oenv);