2257451305faf4c9d917256f94b53011af51b424
[alexxy/gromacs.git] / src / gromacs / mdlib / mdebin.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifndef GMX_MDLIB_MDEBIN_H
38 #define GMX_MDLIB_MDEBIN_H
39
40 #include <stdio.h>
41
42 #include "gromacs/fileio/enxio.h"
43 #include "gromacs/mdlib/ebin.h"
44 #include "gromacs/mdtypes/forcerec.h"
45
46 class energyhistory_t;
47 struct gmx_ekindata_t;
48 struct gmx_mtop_t;
49 struct gmx_output_env_t;
50 struct t_expanded;
51 struct t_fcdata;
52 struct t_grpopts;
53 struct t_lambda;
54 class t_state;
55
56 namespace gmx
57 {
58 class Awh;
59 class Constraints;
60 }
61
62 /* The functions & data structures here determine the content for outputting
63    the .edr file; the file format and actual writing is done with functions
64    defined in enxio.h */
65
66 /* forward declaration */
67 typedef struct t_mde_delta_h_coll t_mde_delta_h_coll;
68
69
70 /* This is the collection of energy averages collected during mdrun, and to
71    be written out to the .edr file. */
72 typedef struct t_mdebin {
73     double              delta_t;
74     t_ebin             *ebin;
75     int                 ie, iconrmsd, ib, ivol, idens, ipv, ienthalpy;
76     int                 isvir, ifvir, ipres, ivir, isurft, ipc, itemp, itc, itcb, iu, imu;
77     int                 ivcos, ivisc;
78     int                 nE, nEg, nEc, nTC, nTCP, nU, nNHC;
79     int                *igrp;
80     int                 mde_n, mdeb_n;
81     real               *tmp_r;
82     rvec               *tmp_v;
83     gmx_bool            bConstr;
84     gmx_bool            bConstrVir;
85     gmx_bool            bTricl;
86     gmx_bool            bDynBox;
87     gmx_bool            bNHC_trotter;
88     gmx_bool            bPrintNHChains;
89     gmx_bool            bMTTK;
90     gmx_bool            bMu; /* true if dipole is calculated */
91     gmx_bool            bDiagPres;
92     int                 f_nre;
93     int                 epc;
94     real                ref_p;
95     int                 etc;
96     int                 nCrmsd;
97     gmx_bool            bEner[F_NRE];
98     gmx_bool            bEInd[egNR];
99     char              **print_grpnms;
100
101     FILE               *fp_dhdl; /* the dhdl.xvg output file */
102     double             *dE;      /* energy components for dhdl.xvg output */
103     t_mde_delta_h_coll *dhc;     /* the delta U components (raw data + histogram) */
104     real               *temperatures;
105 } t_mdebin;
106
107
108 /* delta_h block type enum: the kinds of energies written out. */
109 enum
110 {
111     dhbtDH   = 0, /* delta H BAR energy difference*/
112     dhbtDHDL = 1, /* dH/dlambda derivative */
113     dhbtEN,       /* System energy */
114     dhbtPV,       /* pV term */
115     dhbtEXPANDED, /* expanded ensemble statistics */
116     dhbtNR
117 };
118
119
120
121 t_mdebin *init_mdebin(ener_file_t       fp_ene,
122                       const gmx_mtop_t *mtop,
123                       const t_inputrec *ir,
124                       FILE             *fp_dhdl);
125 /* Initiate MD energy bin and write header to energy file. */
126
127 //! Destroy mdebin
128 void done_mdebin(t_mdebin *mdebin);
129
130 FILE *open_dhdl(const char *filename, const t_inputrec *ir,
131                 const gmx_output_env_t *oenv);
132 /* Open the dhdl file for output */
133
134 /* update the averaging structures. Called every time
135    the energies are evaluated. */
136 void upd_mdebin(t_mdebin                 *md,
137                 gmx_bool                  bDoDHDL,
138                 gmx_bool                  bSum,
139                 double                    time,
140                 real                      tmass,
141                 gmx_enerdata_t           *enerd,
142                 t_state                  *state,
143                 t_lambda                 *fep,
144                 t_expanded               *expand,
145                 matrix                    lastbox,
146                 tensor                    svir,
147                 tensor                    fvir,
148                 tensor                    vir,
149                 tensor                    pres,
150                 gmx_ekindata_t           *ekind,
151                 rvec                      mu_tot,
152                 const gmx::Constraints   *constr);
153
154 void upd_mdebin_step(t_mdebin *md);
155 /* Updates only the step count in md */
156
157 void print_ebin_header(FILE *log, gmx_int64_t steps, double time);
158
159 void print_ebin(ener_file_t fp_ene, gmx_bool bEne, gmx_bool bDR, gmx_bool bOR,
160                 FILE *log,
161                 gmx_int64_t step, double time,
162                 int mode,
163                 t_mdebin *md, t_fcdata *fcd,
164                 gmx_groups_t *groups, t_grpopts *opts,
165                 gmx::Awh *awh);
166
167
168
169 /* Between .edr writes, the averages are history dependent,
170    and that history needs to be retained in checkpoints.
171    These functions set/read the energyhistory_t class
172    that is written to checkpoints in checkpoint.c */
173
174 /* Set the energyhistory_t data from a mdebin structure */
175 void update_energyhistory(energyhistory_t * enerhist, const t_mdebin * mdebin);
176
177 /* Read the energyhistory_t data to a mdebin structure*/
178 void restore_energyhistory_from_state(t_mdebin              * mdebin,
179                                       const energyhistory_t * enerhist);
180
181 #endif