Remove the rest of the device coordinates management from PME
[gromacs.git] / src / gromacs / mdlib / trajectory_writing.cpp
blob9c87f46c4c237a84cbf8c13ef33ae17b91915c51
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
5 * Copyright (c) 2018,2019,2020, 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.
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.
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.
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.
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.
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.
36 #include "gmxpre.h"
38 #include "trajectory_writing.h"
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/mdoutf.h"
45 #include "gromacs/mdlib/stat.h"
46 #include "gromacs/mdlib/update.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/mdtypes/forcerec.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/observableshistory.h"
51 #include "gromacs/mdtypes/state.h"
52 #include "gromacs/pbcutil/pbc.h"
53 #include "gromacs/timing/wallcycle.h"
54 #include "gromacs/topology/topology.h"
55 #include "gromacs/utility/smalloc.h"
57 void do_md_trajectory_writing(FILE* fplog,
58 t_commrec* cr,
59 int nfile,
60 const t_filenm fnm[],
61 int64_t step,
62 int64_t step_rel,
63 double t,
64 t_inputrec* ir,
65 t_state* state,
66 t_state* state_global,
67 ObservablesHistory* observablesHistory,
68 const gmx_mtop_t* top_global,
69 t_forcerec* fr,
70 gmx_mdoutf_t outf,
71 const gmx::EnergyOutput& energyOutput,
72 gmx_ekindata_t* ekind,
73 gmx::ArrayRef<gmx::RVec> f,
74 gmx_bool bCPT,
75 gmx_bool bRerunMD,
76 gmx_bool bLastStep,
77 gmx_bool bDoConfOut,
78 gmx_bool bSumEkinhOld)
80 int mdof_flags;
81 rvec* x_for_confout = nullptr;
83 mdof_flags = 0;
84 if (do_per_step(step, ir->nstxout))
86 mdof_flags |= MDOF_X;
88 if (do_per_step(step, ir->nstvout))
90 mdof_flags |= MDOF_V;
92 if (do_per_step(step, ir->nstfout))
94 mdof_flags |= MDOF_F;
96 if (do_per_step(step, ir->nstxout_compressed))
98 mdof_flags |= MDOF_X_COMPRESSED;
100 if (bCPT)
102 mdof_flags |= MDOF_CPT;
104 if (do_per_step(step, mdoutf_get_tng_box_output_interval(outf)))
106 mdof_flags |= MDOF_BOX;
108 if (do_per_step(step, mdoutf_get_tng_lambda_output_interval(outf)))
110 mdof_flags |= MDOF_LAMBDA;
112 if (do_per_step(step, mdoutf_get_tng_compressed_box_output_interval(outf)))
114 mdof_flags |= MDOF_BOX_COMPRESSED;
116 if (do_per_step(step, mdoutf_get_tng_compressed_lambda_output_interval(outf)))
118 mdof_flags |= MDOF_LAMBDA_COMPRESSED;
121 #if GMX_FAHCORE
122 if (bLastStep)
124 /* Enforce writing positions and velocities at end of run */
125 mdof_flags |= (MDOF_X | MDOF_V);
127 if (MASTER(cr))
129 fcReportProgress(ir->nsteps, step);
132 # if defined(__native_client__)
133 fcCheckin(MASTER(cr));
134 # endif
136 /* sync bCPT and fc record-keeping */
137 if (bCPT && MASTER(cr))
139 fcRequestCheckPoint();
141 #endif
143 if (mdof_flags != 0)
145 wallcycle_start(mdoutf_get_wcycle(outf), ewcTRAJ);
146 if (bCPT)
148 if (MASTER(cr))
150 if (bSumEkinhOld)
152 state_global->ekinstate.bUpToDate = FALSE;
154 else
156 update_ekinstate(&state_global->ekinstate, ekind);
157 state_global->ekinstate.bUpToDate = TRUE;
160 energyOutput.fillEnergyHistory(observablesHistory->energyHistory.get());
163 // Note that part of the following code is duplicated in StatePropagatorData::trajectoryWriterTeardown.
164 // This duplication is needed while both legacy and modular code paths are in use.
165 // TODO: Remove duplication asap, make sure to keep in sync in the meantime.
166 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags, top_global->natoms, step, t,
167 state, state_global, observablesHistory, f);
168 if (bLastStep && step_rel == ir->nsteps && bDoConfOut && MASTER(cr) && !bRerunMD)
170 if (fr->bMolPBC && state == state_global)
172 /* This (single-rank) run needs to allocate a
173 temporary array of size natoms so that any
174 periodicity removal for mdrun -confout does not
175 perturb the update and thus the final .edr
176 output. This makes .cpt restarts look binary
177 identical, and makes .edr restarts binary
178 identical. */
179 snew(x_for_confout, state_global->natoms);
180 copy_rvecn(state_global->x.rvec_array(), x_for_confout, 0, state_global->natoms);
182 else
184 /* With DD, or no bMolPBC, it doesn't matter if
185 we change state_global->x.rvec_array() */
186 x_for_confout = state_global->x.rvec_array();
189 /* x and v have been collected in mdoutf_write_to_trajectory_files,
190 * because a checkpoint file will always be written
191 * at the last step.
193 fprintf(stderr, "\nWriting final coordinates.\n");
194 if (fr->bMolPBC && !ir->bPeriodicMols)
196 /* Make molecules whole only for confout writing */
197 do_pbc_mtop(ir->pbcType, state->box, top_global, x_for_confout);
199 write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm), *top_global->name, top_global,
200 x_for_confout, state_global->v.rvec_array(), ir->pbcType, state->box);
201 if (fr->bMolPBC && state == state_global)
203 sfree(x_for_confout);
206 wallcycle_stop(mdoutf_get_wcycle(outf), ewcTRAJ);