When writing TNG include file closing in wallcycle.
[gromacs/AngularHB.git] / src / gromacs / fileio / trajectory_writing.c
blob0d90e8ea5325ffbe41231c2baa538e13425bb9a2
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014, 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.
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.
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.
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.
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.
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.
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include "typedefs.h"
40 #include "gromacs/utility/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"
49 #include "gromacs/timing/wallcycle.h"
51 void
52 do_md_trajectory_writing(FILE *fplog,
53 t_commrec *cr,
54 int nfile,
55 const t_filenm fnm[],
56 gmx_int64_t step,
57 gmx_int64_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_mdoutf_t outf,
65 t_mdebin *mdebin,
66 gmx_ekindata_t *ekind,
67 rvec *f,
68 rvec *f_global,
69 int *nchkpt,
70 gmx_bool bCPT,
71 gmx_bool bRerunMD,
72 gmx_bool bLastStep,
73 gmx_bool bDoConfOut,
74 gmx_bool bSumEkinhOld
77 int mdof_flags;
79 mdof_flags = 0;
80 if (do_per_step(step, ir->nstxout))
82 mdof_flags |= MDOF_X;
84 if (do_per_step(step, ir->nstvout))
86 mdof_flags |= MDOF_V;
88 if (do_per_step(step, ir->nstfout))
90 mdof_flags |= MDOF_F;
92 if (do_per_step(step, ir->nstxout_compressed))
94 mdof_flags |= MDOF_X_COMPRESSED;
96 if (bCPT)
98 mdof_flags |= MDOF_CPT;
102 #if defined(GMX_FAHCORE) || defined(GMX_WRITELASTSTEP)
103 if (bLastStep)
105 /* Enforce writing positions and velocities at end of run */
106 mdof_flags |= (MDOF_X | MDOF_V);
108 #endif
109 #ifdef GMX_FAHCORE
110 if (MASTER(cr))
112 fcReportProgress( ir->nsteps, step );
115 #if defined(__native_client__)
116 fcCheckin(MASTER(cr));
117 #endif
119 /* sync bCPT and fc record-keeping */
120 if (bCPT && MASTER(cr))
122 fcRequestCheckPoint();
124 #endif
126 if (mdof_flags != 0)
128 wallcycle_start(mdoutf_get_wcycle(outf), ewcTRAJ);
129 if (bCPT)
131 if (MASTER(cr))
133 if (bSumEkinhOld)
135 state_global->ekinstate.bUpToDate = FALSE;
137 else
139 update_ekinstate(&state_global->ekinstate, ekind);
140 state_global->ekinstate.bUpToDate = TRUE;
142 update_energyhistory(&state_global->enerhist, mdebin);
145 mdoutf_write_to_trajectory_files(fplog, cr, outf, mdof_flags, top_global,
146 step, t, state, state_global, f, f_global);
147 if (bCPT)
149 (*nchkpt)++;
150 bCPT = FALSE;
152 debug_gmx();
153 if (bLastStep && step_rel == ir->nsteps &&
154 bDoConfOut && MASTER(cr) &&
155 !bRerunMD)
157 /* x and v have been collected in mdoutf_write_to_trajectory_files,
158 * because a checkpoint file will always be written
159 * at the last step.
161 fprintf(stderr, "\nWriting final coordinates.\n");
162 if (fr->bMolPBC)
164 /* Make molecules whole only for confout writing */
165 do_pbc_mtop(fplog, ir->ePBC, state->box, top_global, state_global->x);
167 write_sto_conf_mtop(ftp2fn(efSTO, nfile, fnm),
168 *top_global->name, top_global,
169 state_global->x, state_global->v,
170 ir->ePBC, state->box);
171 debug_gmx();
173 wallcycle_stop(mdoutf_get_wcycle(outf), ewcTRAJ);