2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020,2021, 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
52 #include "gromacs/math/vectypes.h"
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/mdtypes/enerdata.h"
55 #include "gromacs/topology/ifunc.h"
56 #include "gromacs/utility/arrayref.h"
57 #include "gromacs/utility/real.h"
59 class energyhistory_t;
62 struct gmx_enerdata_t;
63 struct SimulationGroups;
65 struct gmx_output_env_t;
75 struct t_mde_delta_h_coll;
81 struct MDModulesNotifiers;
82 enum class StartingBehavior;
85 //! \brief Printed names for intergroup energies
86 const char* enumValueToString(NonBondedEnergyTerms enumValue);
88 /* \brief delta_h block type enum: the kinds of energies written out. */
91 //! Delta H BAR energy difference
93 //! dH/dlambda derivative
99 //! Expanded ensemble statistics
101 //! Total number of energy types in this enum
107 class EnergyDriftTracker;
110 * \brief Arrays connected to Pressure and Temperature coupling
112 struct PTCouplingArrays
114 //! Box velocities for Parrinello-Rahman P-coupling.
116 //! Nose-Hoover coordinates.
117 ArrayRef<const double> nosehoover_xi;
118 //! Nose-Hoover velocities.
119 ArrayRef<const double> nosehoover_vxi;
120 //! Pressure Nose-Hoover coordinates.
121 ArrayRef<const double> nhpres_xi;
122 //! Pressure Nose-Hoover velocities.
123 ArrayRef<const double> nhpres_vxi;
126 /* Energy output class
128 * This is the collection of energy averages collected during mdrun, and to
129 * be written out to the .edr file.
131 * \todo Use more std containers.
132 * \todo Write free-energy output also to energy file (after adding more tests)
137 /*! \brief Initiate MD energy bin
139 * \param[in] fp_ene Energy output file.
140 * \param[in] mtop Topology.
141 * \param[in] inputrec Input parameters.
142 * \param[in] pull_work Pulling simulations data
143 * \param[in] fp_dhdl FEP file.
144 * \param[in] isRerun Is this is a rerun instead of the simulations.
145 * \param[in] startingBehavior Run starting behavior.
146 * \param[in] simulationsShareState Tells whether the physical state is shared over simulations
147 * \param[in] mdModulesNotifiers Notifications to MD modules.
149 EnergyOutput(ener_file* fp_ene,
150 const gmx_mtop_t& mtop,
151 const t_inputrec& inputrec,
152 const pull_t* pull_work,
155 StartingBehavior startingBehavior,
156 bool simulationsShareState,
157 const MDModulesNotifiers& mdModulesNotifiers);
161 /*! \brief Update the averaging structures.
163 * Called every step on which the thermodynamic values are evaluated.
165 * \param[in] bDoDHDL Whether the FEP is enabled.
166 * \param[in] bSum If this stepshould be recorded to compute sums and averages.
167 * \param[in] time Current simulation time.
168 * \param[in] tmass Total mass
169 * \param[in] enerd Energy data object.
170 * \param[in] fep FEP data.
171 * \param[in] expand Expanded ensemble (for FEP).
172 * \param[in] lastbox PBC data.
173 * \param[in] ptCouplingArrays Arrays connected to pressure and temperature coupling.
174 * \param[in] fep_state The current alchemical state we are in.
175 * \param[in] vir Total virial.
176 * \param[in] pres Pressure.
177 * \param[in] ekind Kinetic energy data.
178 * \param[in] mu_tot Total dipole.
179 * \param[in] constr Constraints object to get RMSD from (for LINCS only).
181 void addDataAtEnergyStep(bool bDoDHDL,
185 const gmx_enerdata_t* enerd,
187 const t_expanded* expand,
188 const matrix lastbox,
189 PTCouplingArrays ptCouplingArrays,
193 const gmx_ekindata_t* ekind,
195 const gmx::Constraints* constr);
197 /*! \brief Update the data averaging structure counts.
199 * Updates the number of steps, the values have not being computed.
201 void recordNonEnergyStep();
203 /*! \brief Writes current quantites to log and energy files.
205 * Prints current values of energies, pressure, temperature, restraint
206 * data, etc. to energy output file and to the log file (if not nullptr).
208 * This function only does something useful when bEne || bDR || bOR || log.
210 * \todo Perhaps this responsibility should involve some other
211 * object visiting all the contributing objects.
213 * \param[in] fp_ene Energy file for the output.
214 * \param[in] bEne If it is a step for energy output or last step.
215 * \param[in] bDR If it is a step of writing distance restraints.
216 * \param[in] bOR If it is a step of writing orientation restraints.
217 * \param[in] log Pointer to the log file.
218 * \param[in] step Current step.
219 * \param[in] time Current simulation time.
220 * \param[in] fcd Bonded force computation data,
221 * including orientation and distance restraints.
222 * \param[in] awh AWH data.
224 void printStepToEnergyFile(ener_file* fp_ene,
234 /*! \brief Print reference temperatures for annealing groups.
236 * Nothing is done if log is nullptr.
238 * \param[in] log Log file to print to.
239 * \param[in] groups Information on atom groups.
240 * \param[in] opts Atom temperature coupling groups options
241 * (annealing is done by groups).
243 static void printAnnealingTemperatures(FILE* log, const SimulationGroups* groups, const t_grpopts* opts);
245 /*! \brief Prints average values to log file.
247 * This is called at the end of the simulation run to print accumulated average values.
248 * Nothing it done if log is nullptr.
250 * \param[in] log Where to print.
251 * \param[in] groups Atom groups.
253 void printAverages(FILE* log, const SimulationGroups* groups);
255 /*! \brief Get the number of thermodynamic terms recorded.
257 * \todo Refactor this to return the expected output size,
258 * rather than exposing the implementation details about
259 * thermodynamic terms.
261 * \returns Number of thermodynamic terms.
263 int numEnergyTerms() const;
265 /*! \brief Fill the energyhistory_t data.
267 * Between .edr writes, the averages are history dependent,
268 * and that history needs to be retained in checkpoints.
269 * These functions set/read the energyhistory_t class
270 * that is written to checkpoints.
272 * \param[out] enerhist Energy history data structure.
274 void fillEnergyHistory(energyhistory_t* enerhist) const;
276 /*! \brief Restore from energyhistory_t data.
278 * \param[in] enerhist Energy history data structure.
280 void restoreFromEnergyHistory(const energyhistory_t& enerhist);
282 //! Print an output header to the log file.
283 static void printHeader(FILE* log, int64_t steps, double time);
285 /*! \brief Print conserved energy drift message to \p fplog
287 * Note that this is only over the current run (times are printed),
288 * this is not from the original start time for runs with continuation.
289 * This has the advantage that one can find if conservation issues are
290 * from the current run with the current settings on the current hardware.
292 void printEnergyConservation(FILE* fplog, int simulationPart, bool usingMdIntegrator) const;
298 //! Structure to store energy components and their running averages
299 t_ebin* ebin_ = nullptr;
301 //! Is the periodic box triclinic
302 bool bTricl_ = false;
303 //! NHC trotter is used
304 bool bNHC_trotter_ = false;
305 //! If Nose-Hoover chains should be printed
306 bool bPrintNHChains_ = false;
307 //! If MTTK integrator was used
310 //! Temperature control scheme
311 TemperatureCoupling etc_ = TemperatureCoupling::No;
313 //! Which of the main energy terms should be printed
314 bool bEner_[F_NRE] = { false };
315 //! Index for main energy terms
317 //! Number of energy terms from F_NRE list to be saved (i.e. number of 'true' in bEner)
320 //! Index for constraints RMSD
322 /* !\brief How many constraints RMSD terms there are.
323 * Usually 1 if LINCS is used and 0 otherwise)
324 * nCrmsd > 0 indicates when constraints RMSD is saves, hence no boolean
328 //! Is the periodic box dynamic
329 bool bDynBox_ = false;
330 //! Index for box dimensions
332 //! Index for box volume
334 //! Index for density
336 //! Triclinic box and not a rerun
337 bool bDiagPres_ = false;
338 //! Reference pressure, averaged over dimensions
340 //! Index for thermodynamic work (pV)
342 //! Index for entalpy (pV + total energy)
345 //! If we have pressure computed
347 //! Index for total virial
349 //! Index for pressure
351 /*! \brief Index for surface tension
352 * [(pres[ZZ][ZZ]-(pres[XX][XX]+pres[YY][YY])*0.5)*box[ZZ][ZZ]]*/
355 //! Pressure control scheme
356 PressureCoupling epc_ = PressureCoupling::No;
357 //! Index for velocity of the box borders
360 //! If dipole was calculated
362 //! Index for the dipole
365 //! Index for coseine acceleration used for viscocity calculation
367 //! Index for viscocity
370 //! Which energy terms from NonBondedEnergyTerms list should be printed in group-to-group block
371 gmx::EnumerationArray<NonBondedEnergyTerms, bool> bEInd_;
372 //! Number of energy terms to be printed (these, for which bEInd[] == true)
374 //! Number of energy output groups
376 //! Number of intergroup energy sets to be printed for each energy term (nE = (nEg*(nEg+1))/2)
378 //! Indexes for integroup energy sets (each set with nEc energies)
379 std::vector<int> igrp_;
381 //! Number of temperature coupling groups
383 //! Index for temperature
386 //! Number of Nose-Hoover chains
388 //! Number of temperature coupling coefficients in case of NH Chains
390 //! Index for temperature coupling coefficient in case of NH chains
393 //! Number of temperature coupling terms if the temperature control is dealt by barostat (MTTK case)
395 //! Scalling factors for temperaturs control in MTTK
397 //! Index for scalling factor of MTTK
400 //! Array to accumulate values during update
401 std::vector<real> tmp_r_;
403 //! The dhdl.xvg output file
404 FILE* fp_dhdl_ = nullptr;
405 //! Energy components for dhdl.xvg output
406 std::vector<double> dE_;
407 //! The delta U components (raw data + histogram)
408 std::unique_ptr<t_mde_delta_h_coll> dhc_;
409 //! Temperatures for simulated tempering groups
410 std::vector<real> temperatures_;
412 //! For tracking the conserved or total energy
413 std::unique_ptr<EnergyDriftTracker> conservedEnergyTracker_;
418 //! Open the dhdl file for output
419 FILE* open_dhdl(const char* filename, const t_inputrec* ir, const gmx_output_env_t* oenv);