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