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