Merge branch release-4-6
[alexxy/gromacs.git] / src / gromacs / fileio / trajectory_writing.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013, 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 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "typedefs.h"
40 #include "smalloc.h"
41 #include "sysstuff.h"
42 #include "vec.h"
43 #include "sim_util.h"
44 #include "mdrun.h"
45 #include "confio.h"
46 #include "trajectory_writing.h"
47 #include "mdoutf.h"
48
49 #include "gromacs/timing/wallcycle.h"
50
51 void
52 do_trajectory_writing(FILE           *fplog,
53                       t_commrec      *cr,
54                       int             nfile,
55                       const t_filenm  fnm[],
56                       gmx_large_int_t step,
57                       gmx_large_int_t step_rel,
58                       double          t,
59                       t_inputrec     *ir,
60                       t_state        *state,
61                       t_state        *state_global,
62                       gmx_mtop_t     *top_global,
63                       t_forcerec     *fr,
64                       gmx_update_t    upd,
65                       gmx_mdoutf_t   *outf,
66                       t_mdebin       *mdebin,
67                       gmx_ekindata_t *ekind,
68                       rvec           *f,
69                       rvec           *f_global,
70                       gmx_wallcycle_t wcycle,
71                       gmx_rng_t       mcrng,
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     int   n_xtc    = -1;
82     rvec *x_xtc    = NULL;
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->nstxtcout))
98     {
99         mdof_flags |= MDOF_XTC;
100     }
101     if (bCPT)
102     {
103         mdof_flags |= MDOF_CPT;
104     }
105     ;
106
107 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
108     if (bLastStep)
109     {
110         /* Enforce writing positions and velocities at end of run */
111         mdof_flags |= (MDOF_X | MDOF_V);
112     }
113 #endif
114 #ifdef GMX_FAHCORE
115     if (MASTER(cr))
116     {
117         fcReportProgress( ir->nsteps, step );
118     }
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(wcycle, ewcTRAJ);
130         if (bCPT)
131         {
132             if (state->flags & (1<<estLD_RNG))
133             {
134                 get_stochd_state(upd, state);
135             }
136             if (state->flags  & (1<<estMC_RNG))
137             {
138                 get_mc_state(mcrng, state);
139             }
140             if (MASTER(cr))
141             {
142                 if (bSumEkinhOld)
143                 {
144                     state_global->ekinstate.bUpToDate = FALSE;
145                 }
146                 else
147                 {
148                     update_ekinstate(&state_global->ekinstate, ekind);
149                     state_global->ekinstate.bUpToDate = TRUE;
150                 }
151                 update_energyhistory(&state_global->enerhist, mdebin);
152             }
153         }
154         write_traj(fplog, cr, outf, mdof_flags, top_global,
155                    step, t, state, state_global, f, f_global, &n_xtc, &x_xtc);
156         if (bCPT)
157         {
158             (*nchkpt)++;
159             bCPT = FALSE;
160         }
161         debug_gmx();
162         if (bLastStep && step_rel == ir->nsteps &&
163             bDoConfOut && MASTER(cr) &&
164             !bRerunMD)
165         {
166             /* x and v have been collected in write_traj,
167              * because a checkpoint file will always be written
168              * at the last step.
169              */
170             fprintf(stderr, "\nWriting final coordinates.\n");
171             if (fr->bMolPBC)
172             {
173                 /* Make molecules whole only for confout writing */
174                 do_pbc_mtop(fplog, ir->ePBC, state->box, top_global, state_global->x);
175             }
176             write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm),
177                                 *top_global->name, top_global,
178                                 state_global->x, state_global->v,
179                                 ir->ePBC, state->box);
180             debug_gmx();
181         }
182         wallcycle_stop(wcycle, ewcTRAJ);
183     }
184 }