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