Make stepWorkload.useGpuXBufferOps flag consistent
[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, The GROMACS development team.
5  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 #include "gmxpre.h"
37
38 #include "trajectory_writing.h"
39
40 #include "gromacs/commandline/filenm.h"
41 #include "gromacs/fileio/confio.h"
42 #include "gromacs/fileio/tngio.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/mdlib/energyoutput.h"
45 #include "gromacs/mdlib/mdoutf.h"
46 #include "gromacs/mdlib/stat.h"
47 #include "gromacs/mdlib/update.h"
48 #include "gromacs/mdtypes/checkpointdata.h"
49 #include "gromacs/mdtypes/commrec.h"
50 #include "gromacs/mdtypes/forcerec.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/observableshistory.h"
53 #include "gromacs/mdtypes/state.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/timing/wallcycle.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/smalloc.h"
58
59 void do_md_trajectory_writing(FILE*                          fplog,
60                               t_commrec*                     cr,
61                               int                            nfile,
62                               const t_filenm                 fnm[],
63                               int64_t                        step,
64                               int64_t                        step_rel,
65                               double                         t,
66                               const t_inputrec*              ir,
67                               t_state*                       state,
68                               t_state*                       state_global,
69                               ObservablesHistory*            observablesHistory,
70                               const gmx_mtop_t&              top_global,
71                               t_forcerec*                    fr,
72                               gmx_mdoutf_t                   outf,
73                               const gmx::EnergyOutput&       energyOutput,
74                               gmx_ekindata_t*                ekind,
75                               gmx::ArrayRef<const gmx::RVec> f,
76                               gmx_bool                       bCPT,
77                               gmx_bool                       bRerunMD,
78                               gmx_bool                       bLastStep,
79                               gmx_bool                       bDoConfOut,
80                               gmx_bool                       bSumEkinhOld)
81 {
82     int   mdof_flags;
83     rvec* x_for_confout = nullptr;
84
85     mdof_flags = 0;
86     if (do_per_step(step, ir->nstxout))
87     {
88         mdof_flags |= MDOF_X;
89     }
90     if (do_per_step(step, ir->nstvout))
91     {
92         mdof_flags |= MDOF_V;
93     }
94     if (do_per_step(step, ir->nstfout))
95     {
96         mdof_flags |= MDOF_F;
97     }
98     if (do_per_step(step, ir->nstxout_compressed))
99     {
100         mdof_flags |= MDOF_X_COMPRESSED;
101     }
102     if (bCPT)
103     {
104         mdof_flags |= MDOF_CPT;
105     }
106     if (do_per_step(step, mdoutf_get_tng_box_output_interval(outf)))
107     {
108         mdof_flags |= MDOF_BOX;
109     }
110     if (do_per_step(step, mdoutf_get_tng_lambda_output_interval(outf)))
111     {
112         mdof_flags |= MDOF_LAMBDA;
113     }
114     if (do_per_step(step, mdoutf_get_tng_compressed_box_output_interval(outf)))
115     {
116         mdof_flags |= MDOF_BOX_COMPRESSED;
117     }
118     if (do_per_step(step, mdoutf_get_tng_compressed_lambda_output_interval(outf)))
119     {
120         mdof_flags |= MDOF_LAMBDA_COMPRESSED;
121     }
122
123     if (mdof_flags != 0)
124     {
125         wallcycle_start(mdoutf_get_wcycle(outf), WallCycleCounter::Traj);
126         if (bCPT)
127         {
128             if (MASTER(cr))
129             {
130                 if (bSumEkinhOld)
131                 {
132                     state_global->ekinstate.bUpToDate = FALSE;
133                 }
134                 else
135                 {
136                     update_ekinstate(&state_global->ekinstate, ekind);
137                     state_global->ekinstate.bUpToDate = TRUE;
138                 }
139
140                 energyOutput.fillEnergyHistory(observablesHistory->energyHistory.get());
141             }
142         }
143         // The current function is only called by legacy code, while
144         // mdoutf_write_to_trajectory_files is also called from modular simulator. Create a dummy
145         // modular simulator checkpointing object for compatibility.
146         gmx::WriteCheckpointDataHolder checkpointDataHolder;
147         // Note that part of the following code is duplicated in StatePropagatorData::trajectoryWriterTeardown.
148         // This duplication is needed while both legacy and modular code paths are in use.
149         // TODO: Remove duplication asap, make sure to keep in sync in the meantime.
150         mdoutf_write_to_trajectory_files(
151                 fplog, cr, outf, mdof_flags, top_global.natoms, step, t, state, state_global, observablesHistory, f, &checkpointDataHolder);
152         if (bLastStep && step_rel == ir->nsteps && bDoConfOut && MASTER(cr) && !bRerunMD)
153         {
154             if (fr->bMolPBC && state == state_global)
155             {
156                 /* This (single-rank) run needs to allocate a
157                    temporary array of size natoms so that any
158                    periodicity removal for mdrun -confout does not
159                    perturb the update and thus the final .edr
160                    output. This makes .cpt restarts look binary
161                    identical, and makes .edr restarts binary
162                    identical. */
163                 snew(x_for_confout, state_global->natoms);
164                 copy_rvecn(state_global->x.rvec_array(), x_for_confout, 0, state_global->natoms);
165             }
166             else
167             {
168                 /* With DD, or no bMolPBC, it doesn't matter if
169                    we change state_global->x.rvec_array() */
170                 x_for_confout = state_global->x.rvec_array();
171             }
172
173             /* x and v have been collected in mdoutf_write_to_trajectory_files,
174              * because a checkpoint file will always be written
175              * at the last step.
176              */
177             fprintf(stderr, "\nWriting final coordinates.\n");
178             if (fr->bMolPBC && !ir->bPeriodicMols)
179             {
180                 /* Make molecules whole only for confout writing */
181                 do_pbc_mtop(ir->pbcType, state->box, &top_global, x_for_confout);
182             }
183             write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm),
184                                 *top_global.name,
185                                 top_global,
186                                 x_for_confout,
187                                 state_global->v.rvec_array(),
188                                 ir->pbcType,
189                                 state->box);
190             if (fr->bMolPBC && state == state_global)
191             {
192                 sfree(x_for_confout);
193             }
194         }
195         wallcycle_stop(mdoutf_get_wcycle(outf), WallCycleCounter::Traj);
196     }
197 #if GMX_FAHCORE
198     if (MASTER(cr))
199     {
200         fcWriteVisFrame(ir->pbcType, state_global->box, top_global, state_global->x.rvec_array());
201     }
202 #endif
203 }