2 * This file is part of the GROMACS molecular simulation package.
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.
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.
39 #include "gromacs/commandline/filenm.h"
40 #include "gromacs/domdec/domdec.h"
41 #include "gromacs/domdec/domdec_struct.h"
42 #include "gromacs/fileio/checkpoint.h"
43 #include "gromacs/fileio/gmxfio.h"
44 #include "gromacs/fileio/tngio.h"
45 #include "gromacs/fileio/trrio.h"
46 #include "gromacs/fileio/xtcio.h"
47 #include "gromacs/fileio/xvgr.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdlib/mdrun.h"
50 #include "gromacs/mdlib/trajectory_writing.h"
51 #include "gromacs/mdtypes/commrec.h"
52 #include "gromacs/mdtypes/imdoutputprovider.h"
53 #include "gromacs/mdtypes/inputrec.h"
54 #include "gromacs/mdtypes/md_enums.h"
55 #include "gromacs/mdtypes/state.h"
56 #include "gromacs/timing/wallcycle.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/pleasecite.h"
59 #include "gromacs/utility/smalloc.h"
65 tng_trajectory_t tng_low_prec
;
66 int x_compression_precision
; /* only used by XTC output */
69 gmx_bool bKeepAndNumCPT
;
76 int natoms_x_compressed
;
77 gmx_groups_t
*groups
; /* for compressed position writing */
78 gmx_wallcycle_t wcycle
;
80 gmx::IMDOutputProvider
*outputProvider
;
84 gmx_mdoutf_t
init_mdoutf(FILE *fplog
, int nfile
, const t_filenm fnm
[],
85 const MdrunOptions
&mdrunOptions
,
87 gmx::IMDOutputProvider
*outputProvider
,
88 const t_inputrec
*ir
, gmx_mtop_t
*top_global
,
89 const gmx_output_env_t
*oenv
, gmx_wallcycle_t wcycle
)
92 const char *appendMode
= "a+", *writeMode
= "w+", *filemode
;
93 gmx_bool bAppendFiles
, bCiteTng
= FALSE
;
100 of
->fp_xtc
= nullptr;
102 of
->tng_low_prec
= nullptr;
103 of
->fp_dhdl
= nullptr;
105 of
->eIntegrator
= ir
->eI
;
106 of
->bExpanded
= ir
->bExpanded
;
107 of
->elamstats
= ir
->expandedvals
->elamstats
;
108 of
->simulation_part
= ir
->simulation_part
;
109 of
->x_compression_precision
= static_cast<int>(ir
->x_compression_precision
);
111 of
->f_global
= nullptr;
112 of
->outputProvider
= outputProvider
;
116 bAppendFiles
= mdrunOptions
.continuationOptions
.appendFiles
;
118 of
->bKeepAndNumCPT
= mdrunOptions
.checkpointOptions
.keepAndNumberCheckpointFiles
;
120 filemode
= bAppendFiles
? appendMode
: writeMode
;
122 if (EI_DYNAMICS(ir
->eI
) &&
123 ir
->nstxout_compressed
> 0)
125 const char *filename
;
126 filename
= ftp2fn(efCOMPRESSED
, nfile
, fnm
);
127 switch (fn2ftp(filename
))
130 of
->fp_xtc
= open_xtc(filename
, filemode
);
133 gmx_tng_open(filename
, filemode
[0], &of
->tng_low_prec
);
134 if (filemode
[0] == 'w')
136 gmx_tng_prepare_low_prec_writing(of
->tng_low_prec
, top_global
, ir
);
141 gmx_incons("Invalid reduced precision file format");
144 if ((EI_DYNAMICS(ir
->eI
) || EI_ENERGY_MINIMIZATION(ir
->eI
))
147 !(EI_DYNAMICS(ir
->eI
) &&
154 const char *filename
;
155 filename
= ftp2fn(efTRN
, nfile
, fnm
);
156 switch (fn2ftp(filename
))
160 /* If there is no uncompressed coordinate output and
161 there is compressed TNG output write forces
162 and/or velocities to the TNG file instead. */
163 if (ir
->nstxout
!= 0 || ir
->nstxout_compressed
== 0 ||
166 of
->fp_trn
= gmx_trr_open(filename
, filemode
);
170 gmx_tng_open(filename
, filemode
[0], &of
->tng
);
171 if (filemode
[0] == 'w')
173 gmx_tng_prepare_md_writing(of
->tng
, top_global
, ir
);
178 gmx_incons("Invalid full precision file format");
181 if (EI_DYNAMICS(ir
->eI
) || EI_ENERGY_MINIMIZATION(ir
->eI
))
183 of
->fp_ene
= open_enx(ftp2fn(efEDR
, nfile
, fnm
), filemode
);
185 of
->fn_cpt
= opt2fn("-cpo", nfile
, fnm
);
187 if ((ir
->efep
!= efepNO
|| ir
->bSimTemp
) && ir
->fepvals
->nstdhdl
> 0 &&
188 (ir
->fepvals
->separate_dhdl_file
== esepdhdlfileYES
) &&
193 of
->fp_dhdl
= gmx_fio_fopen(opt2fn("-dhdl", nfile
, fnm
), filemode
);
197 of
->fp_dhdl
= open_dhdl(opt2fn("-dhdl", nfile
, fnm
), ir
, oenv
);
201 outputProvider
->initOutput(fplog
, nfile
, fnm
, bAppendFiles
, oenv
);
203 /* Set up atom counts so they can be passed to actual
204 trajectory-writing routines later. Also, XTC writing needs
205 to know what (and how many) atoms might be in the XTC
206 groups, and how to look up later which ones they are. */
207 of
->natoms_global
= top_global
->natoms
;
208 of
->groups
= &top_global
->groups
;
209 of
->natoms_x_compressed
= 0;
210 for (i
= 0; (i
< top_global
->natoms
); i
++)
212 if (ggrpnr(of
->groups
, egcCompressedX
, i
) == 0)
214 of
->natoms_x_compressed
++;
218 if (ir
->nstfout
&& DOMAINDECOMP(cr
))
220 snew(of
->f_global
, top_global
->natoms
);
226 please_cite(fplog
, "Lundborg2014");
232 ener_file_t
mdoutf_get_fp_ene(gmx_mdoutf_t of
)
237 FILE *mdoutf_get_fp_dhdl(gmx_mdoutf_t of
)
242 gmx_wallcycle_t
mdoutf_get_wcycle(gmx_mdoutf_t of
)
247 void mdoutf_write_to_trajectory_files(FILE *fplog
, t_commrec
*cr
,
250 gmx_mtop_t
*top_global
,
251 gmx_int64_t step
, double t
,
252 t_state
*state_local
, t_state
*state_global
,
253 ObservablesHistory
*observablesHistory
,
254 gmx::ArrayRef
<gmx::RVec
> f_local
)
258 if (DOMAINDECOMP(cr
))
260 if (mdof_flags
& MDOF_CPT
)
262 dd_collect_state(cr
->dd
, state_local
, state_global
);
266 if (mdof_flags
& (MDOF_X
| MDOF_X_COMPRESSED
))
268 gmx::ArrayRef
<gmx::RVec
> globalXRef
= MASTER(cr
) ? gmx::makeArrayRef(state_global
->x
) : gmx::EmptyArrayRef();
269 dd_collect_vec(cr
->dd
, state_local
, state_local
->x
, globalXRef
);
271 if (mdof_flags
& MDOF_V
)
273 gmx::ArrayRef
<gmx::RVec
> globalVRef
= MASTER(cr
) ? gmx::makeArrayRef(state_global
->v
) : gmx::EmptyArrayRef();
274 dd_collect_vec(cr
->dd
, state_local
, state_local
->v
, globalVRef
);
277 f_global
= of
->f_global
;
278 if (mdof_flags
& MDOF_F
)
280 dd_collect_vec(cr
->dd
, state_local
, f_local
, gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec
*>(f_global
), f_local
.size()));
285 /* We have the whole state locally: copy the local state pointer */
286 state_global
= state_local
;
288 f_global
= as_rvec_array(f_local
.data());
293 if (mdof_flags
& MDOF_CPT
)
296 fflush_tng(of
->tng_low_prec
);
297 ivec one_ivec
= { 1, 1, 1 };
298 write_checkpoint(of
->fn_cpt
, of
->bKeepAndNumCPT
,
300 DOMAINDECOMP(cr
) ? cr
->dd
->nc
: one_ivec
,
301 DOMAINDECOMP(cr
) ? cr
->dd
->nnodes
: cr
->nnodes
,
302 of
->eIntegrator
, of
->simulation_part
,
303 of
->bExpanded
, of
->elamstats
, step
, t
,
304 state_global
, observablesHistory
);
307 if (mdof_flags
& (MDOF_X
| MDOF_V
| MDOF_F
))
309 const rvec
*x
= (mdof_flags
& MDOF_X
) ? as_rvec_array(state_global
->x
.data()) : nullptr;
310 const rvec
*v
= (mdof_flags
& MDOF_V
) ? as_rvec_array(state_global
->v
.data()) : nullptr;
311 const rvec
*f
= (mdof_flags
& MDOF_F
) ? f_global
: nullptr;
315 gmx_trr_write_frame(of
->fp_trn
, step
, t
, state_local
->lambda
[efptFEP
],
316 state_local
->box
, top_global
->natoms
,
318 if (gmx_fio_flush(of
->fp_trn
) != 0)
320 gmx_file("Cannot write trajectory; maybe you are out of disk space?");
324 /* If a TNG file is open for uncompressed coordinate output also write
325 velocities and forces to it. */
328 gmx_fwrite_tng(of
->tng
, FALSE
, step
, t
, state_local
->lambda
[efptFEP
],
333 /* If only a TNG file is open for compressed coordinate output (no uncompressed
334 coordinate output) also write forces and velocities to it. */
335 else if (of
->tng_low_prec
)
337 gmx_fwrite_tng(of
->tng_low_prec
, FALSE
, step
, t
, state_local
->lambda
[efptFEP
],
343 if (mdof_flags
& MDOF_X_COMPRESSED
)
345 rvec
*xxtc
= nullptr;
347 if (of
->natoms_x_compressed
== of
->natoms_global
)
349 /* We are writing the positions of all of the atoms to
350 the compressed output */
351 xxtc
= as_rvec_array(state_global
->x
.data());
355 /* We are writing the positions of only a subset of
356 the atoms to the compressed output, so we have to
357 make a copy of the subset of coordinates. */
360 snew(xxtc
, of
->natoms_x_compressed
);
361 for (i
= 0, j
= 0; (i
< of
->natoms_global
); i
++)
363 if (ggrpnr(of
->groups
, egcCompressedX
, i
) == 0)
365 copy_rvec(state_global
->x
[i
], xxtc
[j
++]);
369 if (write_xtc(of
->fp_xtc
, of
->natoms_x_compressed
, step
, t
,
370 state_local
->box
, xxtc
, of
->x_compression_precision
) == 0)
373 "XTC error. This indicates you are out of disk space, or a "
374 "simulation with major instabilities resulting in coordinates "
375 "that are NaN or too large to be represented in the XTC format.\n");
377 gmx_fwrite_tng(of
->tng_low_prec
,
381 state_local
->lambda
[efptFEP
],
383 of
->natoms_x_compressed
,
387 if (of
->natoms_x_compressed
!= of
->natoms_global
)
395 void mdoutf_tng_close(gmx_mdoutf_t of
)
397 if (of
->tng
|| of
->tng_low_prec
)
399 wallcycle_start(of
->wcycle
, ewcTRAJ
);
400 gmx_tng_close(&of
->tng
);
401 gmx_tng_close(&of
->tng_low_prec
);
402 wallcycle_stop(of
->wcycle
, ewcTRAJ
);
406 void done_mdoutf(gmx_mdoutf_t of
)
408 if (of
->fp_ene
!= nullptr)
410 close_enx(of
->fp_ene
);
414 close_xtc(of
->fp_xtc
);
418 gmx_trr_close(of
->fp_trn
);
420 if (of
->fp_dhdl
!= nullptr)
422 gmx_fio_fclose(of
->fp_dhdl
);
424 of
->outputProvider
->finishOutput();
425 if (of
->f_global
!= nullptr)
430 gmx_tng_close(&of
->tng
);
431 gmx_tng_close(&of
->tng_low_prec
);