2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #ifndef GMX_MDLIB_MDEBIN_H
38 #define GMX_MDLIB_MDEBIN_H
42 #include "gromacs/fileio/enxio.h"
43 #include "gromacs/mdlib/ebin.h"
44 #include "gromacs/mdtypes/enerdata.h"
46 class energyhistory_t
;
47 struct gmx_ekindata_t
;
48 struct gmx_enerdata_t
;
50 struct gmx_output_env_t
;
63 extern const char *egrp_nm
[egNR
+1];
65 /* The functions & data structures here determine the content for outputting
66 the .edr file; the file format and actual writing is done with functions
69 /* forward declaration */
70 typedef struct t_mde_delta_h_coll t_mde_delta_h_coll
;
73 /* This is the collection of energy averages collected during mdrun, and to
74 be written out to the .edr file. */
75 typedef struct t_mdebin
{
78 int ie
, iconrmsd
, ib
, ivol
, idens
, ipv
, ienthalpy
;
79 int isvir
, ifvir
, ipres
, ivir
, isurft
, ipc
, itemp
, itc
, itcb
, iu
, imu
;
81 int nE
, nEg
, nEc
, nTC
, nTCP
, nU
, nNHC
;
90 gmx_bool bNHC_trotter
;
91 gmx_bool bPrintNHChains
;
93 gmx_bool bMu
; /* true if dipole is calculated */
101 gmx_bool bEner
[F_NRE
];
102 gmx_bool bEInd
[egNR
];
105 FILE *fp_dhdl
; /* the dhdl.xvg output file */
106 double *dE
; /* energy components for dhdl.xvg output */
107 t_mde_delta_h_coll
*dhc
; /* the delta U components (raw data + histogram) */
112 /* delta_h block type enum: the kinds of energies written out. */
115 dhbtDH
= 0, /* delta H BAR energy difference*/
116 dhbtDHDL
= 1, /* dH/dlambda derivative */
117 dhbtEN
, /* System energy */
118 dhbtPV
, /* pV term */
119 dhbtEXPANDED
, /* expanded ensemble statistics */
125 t_mdebin
*init_mdebin(ener_file_t fp_ene
,
126 const gmx_mtop_t
*mtop
,
127 const t_inputrec
*ir
,
129 bool isRerun
= false);
130 /* Initiate MD energy bin and write header to energy file. */
133 void done_mdebin(t_mdebin
*mdebin
);
135 FILE *open_dhdl(const char *filename
, const t_inputrec
*ir
,
136 const gmx_output_env_t
*oenv
);
137 /* Open the dhdl file for output */
139 /* update the averaging structures. Called every time
140 the energies are evaluated. */
141 void upd_mdebin(t_mdebin
*md
,
146 gmx_enerdata_t
*enerd
,
155 gmx_ekindata_t
*ekind
,
157 const gmx::Constraints
*constr
);
159 void upd_mdebin_step(t_mdebin
*md
);
160 /* Updates only the step count in md */
162 void print_ebin_header(FILE *log
, int64_t steps
, double time
);
164 void print_ebin(ener_file_t fp_ene
, gmx_bool bEne
, gmx_bool bDR
, gmx_bool bOR
,
166 int64_t step
, double time
,
168 t_mdebin
*md
, t_fcdata
*fcd
,
169 gmx_groups_t
*groups
, t_grpopts
*opts
,
174 /* Between .edr writes, the averages are history dependent,
175 and that history needs to be retained in checkpoints.
176 These functions set/read the energyhistory_t class
177 that is written to checkpoints in checkpoint.c */
179 /* Set the energyhistory_t data from a mdebin structure */
180 void update_energyhistory(energyhistory_t
* enerhist
, const t_mdebin
* mdebin
);
182 /* Read the energyhistory_t data to a mdebin structure*/
183 void restore_energyhistory_from_state(t_mdebin
* mdebin
,
184 const energyhistory_t
* enerhist
);