cb27b6a626e68dff30f2e73ddf67c2a27d713f32
[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, 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/mdrun.h"
45 #include "gromacs/mdlib/sim_util.h"
46 #include "gromacs/mdlib/update.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/mdtypes/observableshistory.h"
50 #include "gromacs/mdtypes/state.h"
51 #include "gromacs/timing/wallcycle.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/smalloc.h"
54
55 void
56 do_md_trajectory_writing(FILE                    *fplog,
57                          t_commrec               *cr,
58                          int                      nfile,
59                          const t_filenm           fnm[],
60                          gmx_int64_t              step,
61                          gmx_int64_t              step_rel,
62                          double                   t,
63                          t_inputrec              *ir,
64                          t_state                 *state,
65                          t_state                 *state_global,
66                          ObservablesHistory      *observablesHistory,
67                          gmx_mtop_t              *top_global,
68                          t_forcerec              *fr,
69                          gmx_mdoutf_t             outf,
70                          t_mdebin                *mdebin,
71                          gmx_ekindata_t          *ekind,
72                          gmx::ArrayRef<gmx::RVec> f,
73                          int                     *nchkpt,
74                          gmx_bool                 bCPT,
75                          gmx_bool                 bRerunMD,
76                          gmx_bool                 bLastStep,
77                          gmx_bool                 bDoConfOut,
78                          gmx_bool                 bSumEkinhOld
79                          )
80 {
81     int   mdof_flags;
82     rvec *x_for_confout = nullptr;
83
84     mdof_flags = 0;
85     if (do_per_step(step, ir->nstxout))
86     {
87         mdof_flags |= MDOF_X;
88     }
89     if (do_per_step(step, ir->nstvout))
90     {
91         mdof_flags |= MDOF_V;
92     }
93     if (do_per_step(step, ir->nstfout))
94     {
95         mdof_flags |= MDOF_F;
96     }
97     if (do_per_step(step, ir->nstxout_compressed))
98     {
99         mdof_flags |= MDOF_X_COMPRESSED;
100     }
101     if (bCPT)
102     {
103         mdof_flags |= MDOF_CPT;
104     }
105     if (do_per_step(step, mdoutf_get_tng_box_output_interval(outf)))
106     {
107         mdof_flags |= MDOF_BOX;
108     }
109     if (do_per_step(step, mdoutf_get_tng_lambda_output_interval(outf)))
110     {
111         mdof_flags |= MDOF_LAMBDA;
112     }
113     if (do_per_step(step, mdoutf_get_tng_compressed_box_output_interval(outf)))
114     {
115         mdof_flags |= MDOF_BOX_COMPRESSED;
116     }
117     if (do_per_step(step, mdoutf_get_tng_compressed_lambda_output_interval(outf)))
118     {
119         mdof_flags |= MDOF_LAMBDA_COMPRESSED;
120     }
121
122 #if defined(GMX_FAHCORE)
123     if (bLastStep)
124     {
125         /* Enforce writing positions and velocities at end of run */
126         mdof_flags |= (MDOF_X | MDOF_V);
127     }
128     if (MASTER(cr))
129     {
130         fcReportProgress( ir->nsteps, step );
131     }
132
133 #if defined(__native_client__)
134     fcCheckin(MASTER(cr));
135 #endif
136
137     /* sync bCPT and fc record-keeping */
138     if (bCPT && MASTER(cr))
139     {
140         fcRequestCheckPoint();
141     }
142 #endif
143
144     if (mdof_flags != 0)
145     {
146         wallcycle_start(mdoutf_get_wcycle(outf), ewcTRAJ);
147         if (bCPT)
148         {
149             if (MASTER(cr))
150             {
151                 if (bSumEkinhOld)
152                 {
153                     state_global->ekinstate.bUpToDate = FALSE;
154                 }
155                 else
156                 {
157                     update_ekinstate(&state_global->ekinstate, ekind);
158                     state_global->ekinstate.bUpToDate = TRUE;
159                 }
160
161                 update_energyhistory(observablesHistory->energyHistory.get(), mdebin);
162             }
163         }
164         mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags, top_global,
165                                          step, t, state, state_global, observablesHistory, f);
166         if (bCPT)
167         {
168             (*nchkpt)++;
169         }
170         if (bLastStep && step_rel == ir->nsteps &&
171             bDoConfOut && MASTER(cr) &&
172             !bRerunMD)
173         {
174             if (fr->bMolPBC && state == state_global)
175             {
176                 /* This (single-rank) run needs to allocate a
177                    temporary array of size natoms so that any
178                    periodicity removal for mdrun -confout does not
179                    perturb the update and thus the final .edr
180                    output. This makes .cpt restarts look binary
181                    identical, and makes .edr restarts binary
182                    identical. */
183                 snew(x_for_confout, state_global->natoms);
184                 copy_rvecn(as_rvec_array(state_global->x.data()), x_for_confout, 0, state_global->natoms);
185             }
186             else
187             {
188                 /* With DD, or no bMolPBC, it doesn't matter if
189                    we change as_rvec_array(state_global->x.data()) */
190                 x_for_confout = as_rvec_array(state_global->x.data());
191             }
192
193             /* x and v have been collected in mdoutf_write_to_trajectory_files,
194              * because a checkpoint file will always be written
195              * at the last step.
196              */
197             fprintf(stderr, "\nWriting final coordinates.\n");
198             if (fr->bMolPBC && !ir->bPeriodicMols)
199             {
200                 /* Make molecules whole only for confout writing */
201                 do_pbc_mtop(fplog, ir->ePBC, state->box, top_global, x_for_confout);
202             }
203             write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm),
204                                 *top_global->name, top_global,
205                                 x_for_confout, as_rvec_array(state_global->v.data()),
206                                 ir->ePBC, state->box);
207             if (fr->bMolPBC && state == state_global)
208             {
209                 sfree(x_for_confout);
210             }
211         }
212         wallcycle_stop(mdoutf_get_wcycle(outf), ewcTRAJ);
213     }
214 }