Update instructions in containers.rst
[gromacs.git] / src / gromacs / mdlib / mdoutf.cpp
blob6ba96c89503b6d69924eb544949906a3898e968d
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2016,2017 The GROMACS development team.
5 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
36 #include "gmxpre.h"
38 #include "mdoutf.h"
40 #include "config.h"
42 #include "gromacs/commandline/filenm.h"
43 #include "gromacs/domdec/collect.h"
44 #include "gromacs/domdec/domdec_struct.h"
45 #include "gromacs/fileio/checkpoint.h"
46 #include "gromacs/fileio/gmxfio.h"
47 #include "gromacs/fileio/tngio.h"
48 #include "gromacs/fileio/trrio.h"
49 #include "gromacs/fileio/xtcio.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdlib/trajectory_writing.h"
53 #include "gromacs/mdrunutility/handlerestart.h"
54 #include "gromacs/mdrunutility/multisim.h"
55 #include "gromacs/mdtypes/awh_history.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/df_history.h"
58 #include "gromacs/mdtypes/edsamhistory.h"
59 #include "gromacs/mdtypes/energyhistory.h"
60 #include "gromacs/mdtypes/imdoutputprovider.h"
61 #include "gromacs/mdtypes/inputrec.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/mdtypes/mdrunoptions.h"
64 #include "gromacs/mdtypes/observableshistory.h"
65 #include "gromacs/mdtypes/state.h"
66 #include "gromacs/mdtypes/swaphistory.h"
67 #include "gromacs/timing/wallcycle.h"
68 #include "gromacs/topology/topology.h"
69 #include "gromacs/utility/baseversion.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/pleasecite.h"
72 #include "gromacs/utility/programcontext.h"
73 #include "gromacs/utility/smalloc.h"
74 #include "gromacs/utility/sysinfo.h"
76 struct gmx_mdoutf
78 t_fileio* fp_trn;
79 t_fileio* fp_xtc;
80 gmx_tng_trajectory_t tng;
81 gmx_tng_trajectory_t tng_low_prec;
82 int x_compression_precision; /* only used by XTC output */
83 ener_file_t fp_ene;
84 const char* fn_cpt;
85 gmx_bool bKeepAndNumCPT;
86 int eIntegrator;
87 gmx_bool bExpanded;
88 int elamstats;
89 int simulation_part;
90 FILE* fp_dhdl;
91 int natoms_global;
92 int natoms_x_compressed;
93 const SimulationGroups* groups; /* for compressed position writing */
94 gmx_wallcycle_t wcycle;
95 rvec* f_global;
96 gmx::IMDOutputProvider* outputProvider;
97 const gmx::MdModulesNotifier* mdModulesNotifier;
98 bool simulationsShareState;
99 MPI_Comm mastersComm;
103 gmx_mdoutf_t init_mdoutf(FILE* fplog,
104 int nfile,
105 const t_filenm fnm[],
106 const gmx::MdrunOptions& mdrunOptions,
107 const t_commrec* cr,
108 gmx::IMDOutputProvider* outputProvider,
109 const gmx::MdModulesNotifier& mdModulesNotifier,
110 const t_inputrec* ir,
111 const gmx_mtop_t* top_global,
112 const gmx_output_env_t* oenv,
113 gmx_wallcycle_t wcycle,
114 const gmx::StartingBehavior startingBehavior,
115 bool simulationsShareState,
116 const gmx_multisim_t* ms)
118 gmx_mdoutf_t of;
119 const char * appendMode = "a+", *writeMode = "w+", *filemode;
120 gmx_bool bCiteTng = FALSE;
121 int i;
122 bool restartWithAppending = (startingBehavior == gmx::StartingBehavior::RestartWithAppending);
124 snew(of, 1);
126 of->fp_trn = nullptr;
127 of->fp_ene = nullptr;
128 of->fp_xtc = nullptr;
129 of->tng = nullptr;
130 of->tng_low_prec = nullptr;
131 of->fp_dhdl = nullptr;
133 of->eIntegrator = ir->eI;
134 of->bExpanded = ir->bExpanded;
135 of->elamstats = ir->expandedvals->elamstats;
136 of->simulation_part = ir->simulation_part;
137 of->x_compression_precision = static_cast<int>(ir->x_compression_precision);
138 of->wcycle = wcycle;
139 of->f_global = nullptr;
140 of->outputProvider = outputProvider;
142 GMX_RELEASE_ASSERT(!simulationsShareState || ms != nullptr,
143 "Need valid multisim object when simulations share state");
144 of->simulationsShareState = simulationsShareState;
145 if (of->simulationsShareState)
147 of->mastersComm = ms->mastersComm_;
150 if (MASTER(cr))
152 of->bKeepAndNumCPT = mdrunOptions.checkpointOptions.keepAndNumberCheckpointFiles;
154 filemode = restartWithAppending ? appendMode : writeMode;
156 if (EI_DYNAMICS(ir->eI) && ir->nstxout_compressed > 0)
158 const char* filename;
159 filename = ftp2fn(efCOMPRESSED, nfile, fnm);
160 switch (fn2ftp(filename))
162 case efXTC: of->fp_xtc = open_xtc(filename, filemode); break;
163 case efTNG:
164 gmx_tng_open(filename, filemode[0], &of->tng_low_prec);
165 if (filemode[0] == 'w')
167 gmx_tng_prepare_low_prec_writing(of->tng_low_prec, top_global, ir);
169 bCiteTng = TRUE;
170 break;
171 default: gmx_incons("Invalid reduced precision file format");
174 if ((EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
175 && (!GMX_FAHCORE
176 && !(EI_DYNAMICS(ir->eI) && ir->nstxout == 0 && ir->nstvout == 0 && ir->nstfout == 0)))
178 const char* filename;
179 filename = ftp2fn(efTRN, nfile, fnm);
180 switch (fn2ftp(filename))
182 case efTRR:
183 case efTRN:
184 /* If there is no uncompressed coordinate output and
185 there is compressed TNG output write forces
186 and/or velocities to the TNG file instead. */
187 if (ir->nstxout != 0 || ir->nstxout_compressed == 0 || !of->tng_low_prec)
189 of->fp_trn = gmx_trr_open(filename, filemode);
191 break;
192 case efTNG:
193 gmx_tng_open(filename, filemode[0], &of->tng);
194 if (filemode[0] == 'w')
196 gmx_tng_prepare_md_writing(of->tng, top_global, ir);
198 bCiteTng = TRUE;
199 break;
200 default: gmx_incons("Invalid full precision file format");
203 if (EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
205 of->fp_ene = open_enx(ftp2fn(efEDR, nfile, fnm), filemode);
207 of->fn_cpt = opt2fn("-cpo", nfile, fnm);
209 if ((ir->efep != efepNO || ir->bSimTemp) && ir->fepvals->nstdhdl > 0
210 && (ir->fepvals->separate_dhdl_file == esepdhdlfileYES) && EI_DYNAMICS(ir->eI))
212 if (restartWithAppending)
214 of->fp_dhdl = gmx_fio_fopen(opt2fn("-dhdl", nfile, fnm), filemode);
216 else
218 of->fp_dhdl = open_dhdl(opt2fn("-dhdl", nfile, fnm), ir, oenv);
222 outputProvider->initOutput(fplog, nfile, fnm, restartWithAppending, oenv);
223 of->mdModulesNotifier = &mdModulesNotifier;
225 /* Set up atom counts so they can be passed to actual
226 trajectory-writing routines later. Also, XTC writing needs
227 to know what (and how many) atoms might be in the XTC
228 groups, and how to look up later which ones they are. */
229 of->natoms_global = top_global->natoms;
230 of->groups = &top_global->groups;
231 of->natoms_x_compressed = 0;
232 for (i = 0; (i < top_global->natoms); i++)
234 if (getGroupType(*of->groups, SimulationAtomGroupType::CompressedPositionOutput, i) == 0)
236 of->natoms_x_compressed++;
240 if (ir->nstfout && DOMAINDECOMP(cr))
242 snew(of->f_global, top_global->natoms);
246 if (bCiteTng)
248 please_cite(fplog, "Lundborg2014");
251 return of;
254 ener_file_t mdoutf_get_fp_ene(gmx_mdoutf_t of)
256 return of->fp_ene;
259 FILE* mdoutf_get_fp_dhdl(gmx_mdoutf_t of)
261 return of->fp_dhdl;
264 gmx_wallcycle_t mdoutf_get_wcycle(gmx_mdoutf_t of)
266 return of->wcycle;
269 static void mpiBarrierBeforeRename(const bool applyMpiBarrierBeforeRename, MPI_Comm mpiBarrierCommunicator)
271 if (applyMpiBarrierBeforeRename)
273 #if GMX_MPI
274 MPI_Barrier(mpiBarrierCommunicator);
275 #else
276 GMX_RELEASE_ASSERT(false, "Should not request a barrier without MPI");
277 GMX_UNUSED_VALUE(mpiBarrierCommunicator);
278 #endif
281 /*! \brief Write a checkpoint to the filename
283 * Appends the _step<step>.cpt with bNumberAndKeep, otherwise moves
284 * the previous checkpoint filename with suffix _prev.cpt.
286 static void write_checkpoint(const char* fn,
287 gmx_bool bNumberAndKeep,
288 FILE* fplog,
289 const t_commrec* cr,
290 ivec domdecCells,
291 int nppnodes,
292 int eIntegrator,
293 int simulation_part,
294 gmx_bool bExpanded,
295 int elamstats,
296 int64_t step,
297 double t,
298 t_state* state,
299 ObservablesHistory* observablesHistory,
300 const gmx::MdModulesNotifier& mdModulesNotifier,
301 gmx::WriteCheckpointDataHolder* modularSimulatorCheckpointData,
302 bool applyMpiBarrierBeforeRename,
303 MPI_Comm mpiBarrierCommunicator)
305 t_fileio* fp;
306 char* fntemp; /* the temporary checkpoint file name */
307 int npmenodes;
308 char buf[1024], suffix[5 + STEPSTRSIZE], sbuf[STEPSTRSIZE];
309 t_fileio* ret;
311 if (DOMAINDECOMP(cr))
313 npmenodes = cr->npmenodes;
315 else
317 npmenodes = 0;
320 #if !GMX_NO_RENAME
321 /* make the new temporary filename */
322 snew(fntemp, std::strlen(fn) + 5 + STEPSTRSIZE);
323 std::strcpy(fntemp, fn);
324 fntemp[std::strlen(fn) - std::strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
325 sprintf(suffix, "_%s%s", "step", gmx_step_str(step, sbuf));
326 std::strcat(fntemp, suffix);
327 std::strcat(fntemp, fn + std::strlen(fn) - std::strlen(ftp2ext(fn2ftp(fn))) - 1);
328 #else
329 /* if we can't rename, we just overwrite the cpt file.
330 * dangerous if interrupted.
332 snew(fntemp, std::strlen(fn));
333 std::strcpy(fntemp, fn);
334 #endif
335 std::string timebuf = gmx_format_current_time();
337 if (fplog)
339 fprintf(fplog, "Writing checkpoint, step %s at %s\n\n", gmx_step_str(step, buf), timebuf.c_str());
342 /* Get offsets for open files */
343 auto outputfiles = gmx_fio_get_output_file_positions();
345 fp = gmx_fio_open(fntemp, "w");
347 /* We can check many more things now (CPU, acceleration, etc), but
348 * it is highly unlikely to have two separate builds with exactly
349 * the same version, user, time, and build host!
352 int nlambda = (state->dfhist ? state->dfhist->nlambda : 0);
354 edsamhistory_t* edsamhist = observablesHistory->edsamHistory.get();
355 int nED = (edsamhist ? edsamhist->nED : 0);
357 swaphistory_t* swaphist = observablesHistory->swapHistory.get();
358 int eSwapCoords = (swaphist ? swaphist->eSwapCoords : eswapNO);
360 CheckpointHeaderContents headerContents = { 0,
361 { 0 },
362 { 0 },
363 { 0 },
364 { 0 },
365 GMX_DOUBLE,
366 { 0 },
367 { 0 },
368 eIntegrator,
369 simulation_part,
370 step,
372 nppnodes,
373 { 0 },
374 npmenodes,
375 state->natoms,
376 state->ngtc,
377 state->nnhpres,
378 state->nhchainlength,
379 nlambda,
380 state->flags,
386 nED,
387 eSwapCoords,
388 false };
389 std::strcpy(headerContents.version, gmx_version());
390 std::strcpy(headerContents.fprog, gmx::getProgramContext().fullBinaryPath());
391 std::strcpy(headerContents.ftime, timebuf.c_str());
392 if (DOMAINDECOMP(cr))
394 copy_ivec(domdecCells, headerContents.dd_nc);
397 write_checkpoint_data(fp, headerContents, bExpanded, elamstats, state, observablesHistory,
398 mdModulesNotifier, &outputfiles, modularSimulatorCheckpointData);
400 /* we really, REALLY, want to make sure to physically write the checkpoint,
401 and all the files it depends on, out to disk. Because we've
402 opened the checkpoint with gmx_fio_open(), it's in our list
403 of open files. */
404 ret = gmx_fio_all_output_fsync();
406 if (ret)
408 char buf[STRLEN];
409 sprintf(buf, "Cannot fsync '%s'; maybe you are out of disk space?", gmx_fio_getname(ret));
411 if (getenv(GMX_IGNORE_FSYNC_FAILURE_ENV) == nullptr)
413 gmx_file(buf);
415 else
417 gmx_warning("%s", buf);
421 if (gmx_fio_close(fp) != 0)
423 gmx_file("Cannot read/write checkpoint; corrupt file, or maybe you are out of disk space?");
426 /* we don't move the checkpoint if the user specified they didn't want it,
427 or if the fsyncs failed */
428 #if !GMX_NO_RENAME
429 if (!bNumberAndKeep && !ret)
431 if (gmx_fexist(fn))
433 /* Rename the previous checkpoint file */
434 mpiBarrierBeforeRename(applyMpiBarrierBeforeRename, mpiBarrierCommunicator);
436 std::strcpy(buf, fn);
437 buf[std::strlen(fn) - std::strlen(ftp2ext(fn2ftp(fn))) - 1] = '\0';
438 std::strcat(buf, "_prev");
439 std::strcat(buf, fn + std::strlen(fn) - std::strlen(ftp2ext(fn2ftp(fn))) - 1);
440 if (!GMX_FAHCORE)
442 /* we copy here so that if something goes wrong between now and
443 * the rename below, there's always a state.cpt.
444 * If renames are atomic (such as in POSIX systems),
445 * this copying should be unneccesary.
447 gmx_file_copy(fn, buf, FALSE);
448 /* We don't really care if this fails:
449 * there's already a new checkpoint.
452 else
454 gmx_file_rename(fn, buf);
458 /* Rename the checkpoint file from the temporary to the final name */
459 mpiBarrierBeforeRename(applyMpiBarrierBeforeRename, mpiBarrierCommunicator);
461 if (gmx_file_rename(fntemp, fn) != 0)
463 gmx_file("Cannot rename checkpoint file; maybe you are out of disk space?");
466 #endif /* GMX_NO_RENAME */
468 sfree(fntemp);
470 #if GMX_FAHCORE
471 /*code for alternate checkpointing scheme. moved from top of loop over
472 steps */
473 fcRequestCheckPoint();
474 if (fcCheckPointParallel(cr->nodeid, NULL, 0) == 0)
476 gmx_fatal(3, __FILE__, __LINE__, "Checkpoint error on step %d\n", step);
478 #endif /* end GMX_FAHCORE block */
481 void mdoutf_write_checkpoint(gmx_mdoutf_t of,
482 FILE* fplog,
483 const t_commrec* cr,
484 int64_t step,
485 double t,
486 t_state* state_global,
487 ObservablesHistory* observablesHistory,
488 gmx::WriteCheckpointDataHolder* modularSimulatorCheckpointData)
490 fflush_tng(of->tng);
491 fflush_tng(of->tng_low_prec);
492 /* Write the checkpoint file.
493 * When simulations share the state, an MPI barrier is applied before
494 * renaming old and new checkpoint files to minimize the risk of
495 * checkpoint files getting out of sync.
497 ivec one_ivec = { 1, 1, 1 };
498 write_checkpoint(of->fn_cpt, of->bKeepAndNumCPT, fplog, cr,
499 DOMAINDECOMP(cr) ? cr->dd->numCells : one_ivec,
500 DOMAINDECOMP(cr) ? cr->dd->nnodes : cr->nnodes, of->eIntegrator,
501 of->simulation_part, of->bExpanded, of->elamstats, step, t, state_global,
502 observablesHistory, *(of->mdModulesNotifier), modularSimulatorCheckpointData,
503 of->simulationsShareState, of->mastersComm);
506 void mdoutf_write_to_trajectory_files(FILE* fplog,
507 const t_commrec* cr,
508 gmx_mdoutf_t of,
509 int mdof_flags,
510 int natoms,
511 int64_t step,
512 double t,
513 t_state* state_local,
514 t_state* state_global,
515 ObservablesHistory* observablesHistory,
516 gmx::ArrayRef<const gmx::RVec> f_local,
517 gmx::WriteCheckpointDataHolder* modularSimulatorCheckpointData)
519 const rvec* f_global;
521 if (DOMAINDECOMP(cr))
523 if (mdof_flags & MDOF_CPT)
525 dd_collect_state(cr->dd, state_local, state_global);
527 else
529 if (mdof_flags & (MDOF_X | MDOF_X_COMPRESSED))
531 auto globalXRef = MASTER(cr) ? state_global->x : gmx::ArrayRef<gmx::RVec>();
532 dd_collect_vec(cr->dd, state_local->ddp_count, state_local->ddp_count_cg_gl,
533 state_local->cg_gl, state_local->x, globalXRef);
535 if (mdof_flags & MDOF_V)
537 auto globalVRef = MASTER(cr) ? state_global->v : gmx::ArrayRef<gmx::RVec>();
538 dd_collect_vec(cr->dd, state_local->ddp_count, state_local->ddp_count_cg_gl,
539 state_local->cg_gl, state_local->v, globalVRef);
542 f_global = of->f_global;
543 if (mdof_flags & MDOF_F)
545 dd_collect_vec(
546 cr->dd, state_local->ddp_count, state_local->ddp_count_cg_gl, state_local->cg_gl, f_local,
547 gmx::arrayRefFromArray(reinterpret_cast<gmx::RVec*>(of->f_global), f_local.size()));
550 else
552 /* We have the whole state locally: copy the local state pointer */
553 state_global = state_local;
555 f_global = as_rvec_array(f_local.data());
558 if (MASTER(cr))
560 if (mdof_flags & MDOF_CPT)
562 mdoutf_write_checkpoint(of, fplog, cr, step, t, state_global, observablesHistory,
563 modularSimulatorCheckpointData);
566 if (mdof_flags & (MDOF_X | MDOF_V | MDOF_F))
568 const rvec* x = (mdof_flags & MDOF_X) ? state_global->x.rvec_array() : nullptr;
569 const rvec* v = (mdof_flags & MDOF_V) ? state_global->v.rvec_array() : nullptr;
570 const rvec* f = (mdof_flags & MDOF_F) ? f_global : nullptr;
572 if (of->fp_trn)
574 gmx_trr_write_frame(of->fp_trn, step, t, state_local->lambda[efptFEP],
575 state_local->box, natoms, x, v, f);
576 if (gmx_fio_flush(of->fp_trn) != 0)
578 gmx_file("Cannot write trajectory; maybe you are out of disk space?");
582 /* If a TNG file is open for uncompressed coordinate output also write
583 velocities and forces to it. */
584 else if (of->tng)
586 gmx_fwrite_tng(of->tng, FALSE, step, t, state_local->lambda[efptFEP],
587 state_local->box, natoms, x, v, f);
589 /* If only a TNG file is open for compressed coordinate output (no uncompressed
590 coordinate output) also write forces and velocities to it. */
591 else if (of->tng_low_prec)
593 gmx_fwrite_tng(of->tng_low_prec, FALSE, step, t, state_local->lambda[efptFEP],
594 state_local->box, natoms, x, v, f);
597 if (mdof_flags & MDOF_X_COMPRESSED)
599 rvec* xxtc = nullptr;
601 if (of->natoms_x_compressed == of->natoms_global)
603 /* We are writing the positions of all of the atoms to
604 the compressed output */
605 xxtc = state_global->x.rvec_array();
607 else
609 /* We are writing the positions of only a subset of
610 the atoms to the compressed output, so we have to
611 make a copy of the subset of coordinates. */
612 int i, j;
614 snew(xxtc, of->natoms_x_compressed);
615 auto x = makeArrayRef(state_global->x);
616 for (i = 0, j = 0; (i < of->natoms_global); i++)
618 if (getGroupType(*of->groups, SimulationAtomGroupType::CompressedPositionOutput, i) == 0)
620 copy_rvec(x[i], xxtc[j++]);
624 if (write_xtc(of->fp_xtc, of->natoms_x_compressed, step, t, state_local->box, xxtc,
625 of->x_compression_precision)
626 == 0)
628 gmx_fatal(FARGS,
629 "XTC error. This indicates you are out of disk space, or a "
630 "simulation with major instabilities resulting in coordinates "
631 "that are NaN or too large to be represented in the XTC format.\n");
633 gmx_fwrite_tng(of->tng_low_prec, TRUE, step, t, state_local->lambda[efptFEP],
634 state_local->box, of->natoms_x_compressed, xxtc, nullptr, nullptr);
635 if (of->natoms_x_compressed != of->natoms_global)
637 sfree(xxtc);
640 if (mdof_flags & (MDOF_BOX | MDOF_LAMBDA) && !(mdof_flags & (MDOF_X | MDOF_V | MDOF_F)))
642 if (of->tng)
644 real lambda = -1;
645 rvec* box = nullptr;
646 if (mdof_flags & MDOF_BOX)
648 box = state_local->box;
650 if (mdof_flags & MDOF_LAMBDA)
652 lambda = state_local->lambda[efptFEP];
654 gmx_fwrite_tng(of->tng, FALSE, step, t, lambda, box, natoms, nullptr, nullptr, nullptr);
657 if (mdof_flags & (MDOF_BOX_COMPRESSED | MDOF_LAMBDA_COMPRESSED)
658 && !(mdof_flags & (MDOF_X_COMPRESSED)))
660 if (of->tng_low_prec)
662 real lambda = -1;
663 rvec* box = nullptr;
664 if (mdof_flags & MDOF_BOX_COMPRESSED)
666 box = state_local->box;
668 if (mdof_flags & MDOF_LAMBDA_COMPRESSED)
670 lambda = state_local->lambda[efptFEP];
672 gmx_fwrite_tng(of->tng_low_prec, FALSE, step, t, lambda, box, natoms, nullptr,
673 nullptr, nullptr);
677 #if GMX_FAHCORE
678 /* Write a FAH checkpoint after writing any other data. We may end up
679 checkpointing twice but it's fast so it's ok. */
680 if ((mdof_flags & ~MDOF_CPT))
682 fcCheckpoint();
684 #endif
688 void mdoutf_tng_close(gmx_mdoutf_t of)
690 if (of->tng || of->tng_low_prec)
692 wallcycle_start(of->wcycle, ewcTRAJ);
693 gmx_tng_close(&of->tng);
694 gmx_tng_close(&of->tng_low_prec);
695 wallcycle_stop(of->wcycle, ewcTRAJ);
699 void done_mdoutf(gmx_mdoutf_t of)
701 if (of->fp_ene != nullptr)
703 done_ener_file(of->fp_ene);
705 if (of->fp_xtc)
707 close_xtc(of->fp_xtc);
709 if (of->fp_trn)
711 gmx_trr_close(of->fp_trn);
713 if (of->fp_dhdl != nullptr)
715 gmx_fio_fclose(of->fp_dhdl);
717 of->outputProvider->finishOutput();
718 if (of->f_global != nullptr)
720 sfree(of->f_global);
723 gmx_tng_close(&of->tng);
724 gmx_tng_close(&of->tng_low_prec);
726 sfree(of);
729 int mdoutf_get_tng_box_output_interval(gmx_mdoutf_t of)
731 if (of->tng)
733 return gmx_tng_get_box_output_interval(of->tng);
735 return 0;
738 int mdoutf_get_tng_lambda_output_interval(gmx_mdoutf_t of)
740 if (of->tng)
742 return gmx_tng_get_lambda_output_interval(of->tng);
744 return 0;
747 int mdoutf_get_tng_compressed_box_output_interval(gmx_mdoutf_t of)
749 if (of->tng_low_prec)
751 return gmx_tng_get_box_output_interval(of->tng_low_prec);
753 return 0;
756 int mdoutf_get_tng_compressed_lambda_output_interval(gmx_mdoutf_t of)
758 if (of->tng_low_prec)
760 return gmx_tng_get_lambda_output_interval(of->tng_low_prec);
762 return 0;