Merge commit d30f2cb6 from release-2020 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / trajectory_writing.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014,2015,2016,2017,2018,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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "trajectory_writing.h"
38
39 #include "gromacs/commandline/filenm.h"
40 #include "gromacs/fileio/confio.h"
41 #include "gromacs/fileio/tngio.h"
42 #include "gromacs/math/vec.h"
43 #include "gromacs/mdlib/mdoutf.h"
44 #include "gromacs/mdlib/stat.h"
45 #include "gromacs/mdlib/update.h"
46 #include "gromacs/mdtypes/commrec.h"
47 #include "gromacs/mdtypes/forcerec.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/mdtypes/observableshistory.h"
50 #include "gromacs/mdtypes/state.h"
51 #include "gromacs/pbcutil/pbc.h"
52 #include "gromacs/timing/wallcycle.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/smalloc.h"
55
56 void do_md_trajectory_writing(FILE*                    fplog,
57                               t_commrec*               cr,
58                               int                      nfile,
59                               const t_filenm           fnm[],
60                               int64_t                  step,
61                               int64_t                  step_rel,
62                               double                   t,
63                               t_inputrec*              ir,
64                               t_state*                 state,
65                               t_state*                 state_global,
66                               ObservablesHistory*      observablesHistory,
67                               const gmx_mtop_t*        top_global,
68                               t_forcerec*              fr,
69                               gmx_mdoutf_t             outf,
70                               const gmx::EnergyOutput& energyOutput,
71                               gmx_ekindata_t*          ekind,
72                               gmx::ArrayRef<gmx::RVec> f,
73                               gmx_bool                 bCPT,
74                               gmx_bool                 bRerunMD,
75                               gmx_bool                 bLastStep,
76                               gmx_bool                 bDoConfOut,
77                               gmx_bool                 bSumEkinhOld)
78 {
79     int   mdof_flags;
80     rvec* x_for_confout = nullptr;
81
82     mdof_flags = 0;
83     if (do_per_step(step, ir->nstxout))
84     {
85         mdof_flags |= MDOF_X;
86     }
87     if (do_per_step(step, ir->nstvout))
88     {
89         mdof_flags |= MDOF_V;
90     }
91     if (do_per_step(step, ir->nstfout))
92     {
93         mdof_flags |= MDOF_F;
94     }
95     if (do_per_step(step, ir->nstxout_compressed))
96     {
97         mdof_flags |= MDOF_X_COMPRESSED;
98     }
99     if (bCPT)
100     {
101         mdof_flags |= MDOF_CPT;
102     }
103     if (do_per_step(step, mdoutf_get_tng_box_output_interval(outf)))
104     {
105         mdof_flags |= MDOF_BOX;
106     }
107     if (do_per_step(step, mdoutf_get_tng_lambda_output_interval(outf)))
108     {
109         mdof_flags |= MDOF_LAMBDA;
110     }
111     if (do_per_step(step, mdoutf_get_tng_compressed_box_output_interval(outf)))
112     {
113         mdof_flags |= MDOF_BOX_COMPRESSED;
114     }
115     if (do_per_step(step, mdoutf_get_tng_compressed_lambda_output_interval(outf)))
116     {
117         mdof_flags |= MDOF_LAMBDA_COMPRESSED;
118     }
119
120 #if GMX_FAHCORE
121     if (bLastStep)
122     {
123         /* Enforce writing positions and velocities at end of run */
124         mdof_flags |= (MDOF_X | MDOF_V);
125     }
126     if (MASTER(cr))
127     {
128         fcReportProgress(ir->nsteps, step);
129     }
130
131 #    if defined(__native_client__)
132     fcCheckin(MASTER(cr));
133 #    endif
134
135     /* sync bCPT and fc record-keeping */
136     if (bCPT && MASTER(cr))
137     {
138         fcRequestCheckPoint();
139     }
140 #endif
141
142     if (mdof_flags != 0)
143     {
144         wallcycle_start(mdoutf_get_wcycle(outf), ewcTRAJ);
145         if (bCPT)
146         {
147             if (MASTER(cr))
148             {
149                 if (bSumEkinhOld)
150                 {
151                     state_global->ekinstate.bUpToDate = FALSE;
152                 }
153                 else
154                 {
155                     update_ekinstate(&state_global->ekinstate, ekind);
156                     state_global->ekinstate.bUpToDate = TRUE;
157                 }
158
159                 energyOutput.fillEnergyHistory(observablesHistory->energyHistory.get());
160             }
161         }
162         // Note that part of the following code is duplicated in StatePropagatorData::trajectoryWriterTeardown.
163         // This duplication is needed while both legacy and modular code paths are in use.
164         // TODO: Remove duplication asap, make sure to keep in sync in the meantime.
165         mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags, top_global->natoms, step, t,
166                                          state, state_global, observablesHistory, f);
167         if (bLastStep && step_rel == ir->nsteps && bDoConfOut && MASTER(cr) && !bRerunMD)
168         {
169             if (fr->bMolPBC && state == state_global)
170             {
171                 /* This (single-rank) run needs to allocate a
172                    temporary array of size natoms so that any
173                    periodicity removal for mdrun -confout does not
174                    perturb the update and thus the final .edr
175                    output. This makes .cpt restarts look binary
176                    identical, and makes .edr restarts binary
177                    identical. */
178                 snew(x_for_confout, state_global->natoms);
179                 copy_rvecn(state_global->x.rvec_array(), x_for_confout, 0, state_global->natoms);
180             }
181             else
182             {
183                 /* With DD, or no bMolPBC, it doesn't matter if
184                    we change state_global->x.rvec_array() */
185                 x_for_confout = state_global->x.rvec_array();
186             }
187
188             /* x and v have been collected in mdoutf_write_to_trajectory_files,
189              * because a checkpoint file will always be written
190              * at the last step.
191              */
192             fprintf(stderr, "\nWriting final coordinates.\n");
193             if (fr->bMolPBC && !ir->bPeriodicMols)
194             {
195                 /* Make molecules whole only for confout writing */
196                 do_pbc_mtop(ir->ePBC, state->box, top_global, x_for_confout);
197             }
198             write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm), *top_global->name, top_global,
199                                 x_for_confout, state_global->v.rvec_array(), ir->ePBC, state->box);
200             if (fr->bMolPBC && state == state_global)
201             {
202                 sfree(x_for_confout);
203             }
204         }
205         wallcycle_stop(mdoutf_get_wcycle(outf), ewcTRAJ);
206     }
207 }