OCL: Make variables const
[gromacs.git] / src / gromacs / mdlib / mdoutf.cpp
blob1ddf2c7ae050c2a4b3f53b4016e8f3baec5b18e5
1 /*
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.
35 #include "gmxpre.h"
37 #include "mdoutf.h"
39 #include "gromacs/commandline/filenm.h"
40 #include "gromacs/domdec/collect.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/topology/topology.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/pleasecite.h"
60 #include "gromacs/utility/smalloc.h"
62 struct gmx_mdoutf {
63 t_fileio *fp_trn;
64 t_fileio *fp_xtc;
65 gmx_tng_trajectory_t tng;
66 gmx_tng_trajectory_t tng_low_prec;
67 int x_compression_precision; /* only used by XTC output */
68 ener_file_t fp_ene;
69 const char *fn_cpt;
70 gmx_bool bKeepAndNumCPT;
71 int eIntegrator;
72 gmx_bool bExpanded;
73 int elamstats;
74 int simulation_part;
75 FILE *fp_dhdl;
76 int natoms_global;
77 int natoms_x_compressed;
78 gmx_groups_t *groups; /* for compressed position writing */
79 gmx_wallcycle_t wcycle;
80 rvec *f_global;
81 gmx::IMDOutputProvider *outputProvider;
85 gmx_mdoutf_t init_mdoutf(FILE *fplog, int nfile, const t_filenm fnm[],
86 const MdrunOptions &mdrunOptions,
87 const t_commrec *cr,
88 gmx::IMDOutputProvider *outputProvider,
89 const t_inputrec *ir, gmx_mtop_t *top_global,
90 const gmx_output_env_t *oenv, gmx_wallcycle_t wcycle)
92 gmx_mdoutf_t of;
93 const char *appendMode = "a+", *writeMode = "w+", *filemode;
94 gmx_bool bAppendFiles, bCiteTng = FALSE;
95 int i;
97 snew(of, 1);
99 of->fp_trn = nullptr;
100 of->fp_ene = nullptr;
101 of->fp_xtc = nullptr;
102 of->tng = nullptr;
103 of->tng_low_prec = nullptr;
104 of->fp_dhdl = nullptr;
106 of->eIntegrator = ir->eI;
107 of->bExpanded = ir->bExpanded;
108 of->elamstats = ir->expandedvals->elamstats;
109 of->simulation_part = ir->simulation_part;
110 of->x_compression_precision = static_cast<int>(ir->x_compression_precision);
111 of->wcycle = wcycle;
112 of->f_global = nullptr;
113 of->outputProvider = outputProvider;
115 if (MASTER(cr))
117 bAppendFiles = mdrunOptions.continuationOptions.appendFiles;
119 of->bKeepAndNumCPT = mdrunOptions.checkpointOptions.keepAndNumberCheckpointFiles;
121 filemode = bAppendFiles ? appendMode : writeMode;
123 if (EI_DYNAMICS(ir->eI) &&
124 ir->nstxout_compressed > 0)
126 const char *filename;
127 filename = ftp2fn(efCOMPRESSED, nfile, fnm);
128 switch (fn2ftp(filename))
130 case efXTC:
131 of->fp_xtc = open_xtc(filename, filemode);
132 break;
133 case efTNG:
134 gmx_tng_open(filename, filemode[0], &of->tng_low_prec);
135 if (filemode[0] == 'w')
137 gmx_tng_prepare_low_prec_writing(of->tng_low_prec, top_global, ir);
139 bCiteTng = TRUE;
140 break;
141 default:
142 gmx_incons("Invalid reduced precision file format");
145 if ((EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI)) &&
146 (!GMX_FAHCORE &&
147 !(EI_DYNAMICS(ir->eI) &&
148 ir->nstxout == 0 &&
149 ir->nstvout == 0 &&
150 ir->nstfout == 0)
154 const char *filename;
155 filename = ftp2fn(efTRN, nfile, fnm);
156 switch (fn2ftp(filename))
158 case efTRR:
159 case efTRN:
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 ||
164 !of->tng_low_prec)
166 of->fp_trn = gmx_trr_open(filename, filemode);
168 break;
169 case efTNG:
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);
175 bCiteTng = TRUE;
176 break;
177 default:
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 ) &&
189 EI_DYNAMICS(ir->eI))
191 if (bAppendFiles)
193 of->fp_dhdl = gmx_fio_fopen(opt2fn("-dhdl", nfile, fnm), filemode);
195 else
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 (getGroupType(of->groups, egcCompressedX, i) == 0)
214 of->natoms_x_compressed++;
218 if (ir->nstfout && DOMAINDECOMP(cr))
220 snew(of->f_global, top_global->natoms);
224 if (bCiteTng)
226 please_cite(fplog, "Lundborg2014");
229 return of;
232 ener_file_t mdoutf_get_fp_ene(gmx_mdoutf_t of)
234 return of->fp_ene;
237 FILE *mdoutf_get_fp_dhdl(gmx_mdoutf_t of)
239 return of->fp_dhdl;
242 gmx_wallcycle_t mdoutf_get_wcycle(gmx_mdoutf_t of)
244 return of->wcycle;
247 void mdoutf_write_to_trajectory_files(FILE *fplog, const t_commrec *cr,
248 gmx_mdoutf_t of,
249 int mdof_flags,
250 gmx_mtop_t *top_global,
251 int64_t step, double t,
252 t_state *state_local, t_state *state_global,
253 ObservablesHistory *observablesHistory,
254 gmx::ArrayRef<gmx::RVec> f_local)
256 rvec *f_global;
258 if (DOMAINDECOMP(cr))
260 if (mdof_flags & MDOF_CPT)
262 dd_collect_state(cr->dd, state_local, state_global);
264 else
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()));
283 else
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());
291 if (MASTER(cr))
293 if (mdof_flags & MDOF_CPT)
295 fflush_tng(of->tng);
296 fflush_tng(of->tng_low_prec);
297 ivec one_ivec = { 1, 1, 1 };
298 write_checkpoint(of->fn_cpt, of->bKeepAndNumCPT,
299 fplog, cr,
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;
313 if (of->fp_trn)
315 gmx_trr_write_frame(of->fp_trn, step, t, state_local->lambda[efptFEP],
316 state_local->box, top_global->natoms,
317 x, v, f);
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. */
326 else if (of->tng)
328 gmx_fwrite_tng(of->tng, FALSE, step, t, state_local->lambda[efptFEP],
329 state_local->box,
330 top_global->natoms,
331 x, v, f);
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],
338 state_local->box,
339 top_global->natoms,
340 x, v, f);
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());
353 else
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. */
358 int i, j;
360 snew(xxtc, of->natoms_x_compressed);
361 for (i = 0, j = 0; (i < of->natoms_global); i++)
363 if (getGroupType(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)
372 gmx_fatal(FARGS,
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,
378 TRUE,
379 step,
381 state_local->lambda[efptFEP],
382 state_local->box,
383 of->natoms_x_compressed,
384 xxtc,
385 nullptr,
386 nullptr);
387 if (of->natoms_x_compressed != of->natoms_global)
389 sfree(xxtc);
392 if (mdof_flags & (MDOF_BOX | MDOF_LAMBDA) && !(mdof_flags & (MDOF_X | MDOF_V | MDOF_F)) )
394 if (of->tng)
396 real lambda = -1;
397 rvec *box = nullptr;
398 if (mdof_flags & MDOF_BOX)
400 box = state_local->box;
402 if (mdof_flags & MDOF_LAMBDA)
404 lambda = state_local->lambda[efptFEP];
406 gmx_fwrite_tng(of->tng, FALSE, step, t, lambda,
407 box, top_global->natoms,
408 nullptr, nullptr, nullptr);
411 if (mdof_flags & (MDOF_BOX_COMPRESSED | MDOF_LAMBDA_COMPRESSED) && !(mdof_flags & (MDOF_X_COMPRESSED)) )
413 if (of->tng_low_prec)
415 real lambda = -1;
416 rvec *box = nullptr;
417 if (mdof_flags & MDOF_BOX_COMPRESSED)
419 box = state_local->box;
421 if (mdof_flags & MDOF_LAMBDA_COMPRESSED)
423 lambda = state_local->lambda[efptFEP];
425 gmx_fwrite_tng(of->tng_low_prec, FALSE, step, t, lambda,
426 box, top_global->natoms,
427 nullptr, nullptr, nullptr);
433 void mdoutf_tng_close(gmx_mdoutf_t of)
435 if (of->tng || of->tng_low_prec)
437 wallcycle_start(of->wcycle, ewcTRAJ);
438 gmx_tng_close(&of->tng);
439 gmx_tng_close(&of->tng_low_prec);
440 wallcycle_stop(of->wcycle, ewcTRAJ);
444 void done_mdoutf(gmx_mdoutf_t of)
446 if (of->fp_ene != nullptr)
448 close_enx(of->fp_ene);
450 if (of->fp_xtc)
452 close_xtc(of->fp_xtc);
454 if (of->fp_trn)
456 gmx_trr_close(of->fp_trn);
458 if (of->fp_dhdl != nullptr)
460 gmx_fio_fclose(of->fp_dhdl);
462 of->outputProvider->finishOutput();
463 if (of->f_global != nullptr)
465 sfree(of->f_global);
468 gmx_tng_close(&of->tng);
469 gmx_tng_close(&of->tng_low_prec);
471 sfree(of);
474 int mdoutf_get_tng_box_output_interval(gmx_mdoutf_t of)
476 if (of->tng)
478 return gmx_tng_get_box_output_interval(of->tng);
480 return 0;
483 int mdoutf_get_tng_lambda_output_interval(gmx_mdoutf_t of)
485 if (of->tng)
487 return gmx_tng_get_lambda_output_interval(of->tng);
489 return 0;
492 int mdoutf_get_tng_compressed_box_output_interval(gmx_mdoutf_t of)
494 if (of->tng_low_prec)
496 return gmx_tng_get_box_output_interval(of->tng_low_prec);
498 return 0;
501 int mdoutf_get_tng_compressed_lambda_output_interval(gmx_mdoutf_t of)
503 if (of->tng_low_prec)
505 return gmx_tng_get_lambda_output_interval(of->tng_low_prec);
507 return 0;