2 * This file is part of the GROMACS molecular simulation package.
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.
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.
46 #include "../../external/tng_io/include/tng_io.h"
49 #include "gromacs/legacyheaders/copyrite.h"
50 #include "gromacs/legacyheaders/gmx_fatal.h"
51 #include "gromacs/legacyheaders/main.h"
52 #include "gromacs/legacyheaders/physics.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/programinfo.h"
55 #include "gromacs/math/utilities.h"
58 static const char *modeToVerb(char mode)
72 gmx_fatal(FARGS, "Invalid file opening mode %c", mode);
77 void gmx_tng_open(const char *filename,
79 tng_trajectory_t *tng)
82 /* First check whether we have to make a backup,
83 * only for writing, not for read or append.
88 /* only make backups for normal gromacs */
89 make_backup(filename);
93 /* tng must not be pointing at already allocated memory.
94 * Memory will be allocated by tng_util_trajectory_open() and must
95 * later on be freed by tng_util_trajectory_close(). */
96 if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
98 /* TNG does return more than one degree of error, but there is
99 no use case for GROMACS handling the non-fatal errors
102 "%s while opening %s for %s",
103 gmx_strerror("file"),
108 if (mode == 'w' || mode == 'a')
110 /* FIXME in TNG: When adding data to the header, subsequent blocks might get
111 * overwritten. This could be solved by moving the first trajectory
112 * frame set(s) to the end of the file. Could that cause other problems,
113 * e.g. when continuing a simulation? */
115 gmx_gethostname(hostname, 256);
118 tng_first_computer_name_set(*tng, hostname);
120 /* TODO: This should be implemented when the above fixme is done (adding data to
124 // tng_last_computer_name_set(*tng, hostname);
127 char programInfo[256];
128 const char *precisionString = "";
130 precisionString = " (double precision)";
132 sprintf(programInfo, "%.100s, %.128s%.24s", gmx::ProgramInfo::getInstance().displayName().c_str(),
133 GromacsVersion(), precisionString);
136 tng_first_program_name_set(*tng, programInfo);
138 /* TODO: This should be implemented when the above fixme is done (adding data to
142 // tng_last_program_name_set(*tng, programInfo);
147 getlogin_r(username, 256);
150 tng_first_user_name_set(*tng, username);
152 /* TODO: This should be implemented when the above fixme is done (adding data to
156 // tng_last_user_name_set(*tng, username);
161 gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
162 GMX_UNUSED_VALUE(filename);
163 GMX_UNUSED_VALUE(mode);
164 GMX_UNUSED_VALUE(tng);
168 void gmx_tng_close(tng_trajectory_t *tng)
170 /* We have to check that tng is set because
171 * tng_util_trajectory_close wants to return a NULL in it, and
172 * gives a fatal error if it is NULL. */
176 tng_util_trajectory_close(tng);
179 GMX_UNUSED_VALUE(tng);
184 static void addTngMoleculeFromTopology(tng_trajectory_t tng,
185 const char *moleculeName,
186 const t_atoms *atoms,
187 gmx_int64_t numMolecules,
188 tng_molecule_t *tngMol)
190 if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
192 gmx_file("Cannot add molecule to TNG molecular system.");
195 /* FIXME: The TNG atoms should contain mass and atomB info (for free
196 * energy calculations), i.e. in when it's available in TNG (2.0). */
197 for (int atomIt = 0; atomIt < atoms->nr; atomIt++)
199 const t_atom *at = &atoms->atom[atomIt];
200 /* FIXME: Currently the TNG API can only add atoms belonging to a
201 * residue and chain. Wait for TNG 2.0*/
204 const t_resinfo *resInfo = &atoms->resinfo[at->resind];
205 char chainName[2] = {resInfo->chainid, 0};
206 tng_chain_t tngChain = NULL;
207 tng_residue_t tngRes = NULL;
208 tng_atom_t tngAtom = NULL;
210 if (tng_molecule_chain_find (tng, *tngMol, chainName,
211 (gmx_int64_t)-1, &tngChain) !=
214 tng_molecule_chain_add (tng, *tngMol, chainName,
218 /* FIXME: When TNG supports both residue index and residue
219 * number the latter should be used. Wait for TNG 2.0*/
220 if (tng_chain_residue_find(tng, tngChain, *resInfo->name,
221 at->resind + 1, &tngRes)
224 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
226 tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIt]), *(atoms->atomtype[atomIt]), &tngAtom);
229 tng_molecule_cnt_set(tng, *tngMol, numMolecules);
232 void gmx_tng_add_mtop(tng_trajectory_t tng,
233 const gmx_mtop_t *mtop)
236 const t_ilist *ilist;
241 /* No topology information available to add. */
245 for (int molIt = 0; molIt < mtop->nmolblock; molIt++)
247 tng_molecule_t tngMol = NULL;
248 const gmx_moltype_t *molType =
249 &mtop->moltype[mtop->molblock[molIt].type];
251 /* Add a molecule to the TNG trajectory with the same name as the
252 * current molecule. */
253 addTngMoleculeFromTopology(tng,
256 mtop->molblock[molIt].nmol,
259 /* Bonds have to be deduced from interactions (constraints etc). Different
260 * interactions have different sets of parameters. */
261 /* Constraints are specified using two atoms */
262 for (i = 0; i < F_NRE; i++)
266 ilist = &molType->ilist[i];
270 while (j < ilist->nr)
272 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
278 /* Settle is described using three atoms */
279 ilist = &molType->ilist[F_SETTLE];
283 while (j < ilist->nr)
285 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
286 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+2], &tngBond);
293 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
294 * if they are positive.
296 * If only one of n1 and n2 is positive, then return it.
297 * If neither n1 or n2 is positive, then return -1. */
299 greatest_common_divisor_if_positive(int n1, int n2)
303 return (0 >= n2) ? -1 : n2;
310 /* We have a non-trivial greatest common divisor to compute. */
311 return gmx_greatest_common_divisor(n1, n2);
314 /* By default try to write 100 frames (of actual output) in each frame set.
315 * This number is the number of outputs of the most frequently written data
316 * type per frame set.
317 * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
318 * setups regarding compression efficiency and compression time. Make this
319 * a hidden command-line option? */
320 const int defaultFramesPerFrameSet = 100;
322 /*! \libinternal \brief Set the number of frames per frame
323 * set according to output intervals.
324 * The default is that 100 frames are written of the data
325 * that is written most often. */
326 static void tng_set_frames_per_frame_set(tng_trajectory_t tng,
327 const gmx_bool bUseLossyCompression,
328 const t_inputrec *ir)
332 /* Set the number of frames per frame set to contain at least
333 * defaultFramesPerFrameSet of the lowest common denominator of
334 * the writing interval of positions and velocities. */
335 /* FIXME after 5.0: consider nstenergy also? */
336 if (bUseLossyCompression)
338 gcd = ir->nstxout_compressed;
342 gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
343 gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
350 tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
353 /*! \libinternal \brief Set the data-writing intervals, and number of
354 * frames per frame set */
355 static void set_writing_intervals(tng_trajectory_t tng,
356 const gmx_bool bUseLossyCompression,
357 const t_inputrec *ir)
359 /* Define pointers to specific writing functions depending on if we
360 * write float or double data */
361 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
369 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
371 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
373 int xout, vout, fout;
374 // int gcd = -1, lowest = -1;
377 tng_set_frames_per_frame_set(tng, bUseLossyCompression, ir);
379 if (bUseLossyCompression)
381 xout = ir->nstxout_compressed;
384 compression = TNG_TNG_COMPRESSION;
391 compression = TNG_GZIP_COMPRESSION;
395 set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS,
396 "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
398 /* The design of TNG makes it awkward to try to write a box
399 * with multiple periodicities, which might be co-prime. Since
400 * the use cases for the box with a frame consisting only of
401 * velocities seem low, for now we associate box writing with
402 * position writing. */
403 set_writing_interval(tng, xout, 9, TNG_TRAJ_BOX_SHAPE,
404 "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA,
405 TNG_GZIP_COMPRESSION);
406 /* TODO: if/when we write energies to TNG also, reconsider how
407 * and when box information is written, because GROMACS
408 * behaviour pre-5.0 was to write the box with every
409 * trajectory frame and every energy frame, and probably
410 * people depend on this. */
412 /* TODO: If we need to write lambda values at steps when
413 * positions (or other data) are not also being written, then
414 * code in mdoutf.c will need to match however that is
416 set_writing_interval(tng, xout, 1, TNG_GMX_LAMBDA,
417 "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
418 TNG_GZIP_COMPRESSION);
420 /* FIXME: gcd and lowest currently not used. */
421 // gcd = greatest_common_divisor_if_positive(gcd, xout);
422 // if (lowest < 0 || xout < lowest)
429 set_writing_interval(tng, ir->nstvout, 3, TNG_TRAJ_VELOCITIES,
430 "VELOCITIES", TNG_PARTICLE_BLOCK_DATA,
433 /* FIXME: gcd and lowest currently not used. */
434 // gcd = greatest_common_divisor_if_positive(gcd, vout);
435 // if (lowest < 0 || vout < lowest)
442 set_writing_interval(tng, ir->nstfout, 3, TNG_TRAJ_FORCES,
443 "FORCES", TNG_PARTICLE_BLOCK_DATA,
444 TNG_GZIP_COMPRESSION);
446 /* FIXME: gcd and lowest currently not used. */
447 // gcd = greatest_common_divisor_if_positive(gcd, fout);
448 // if (lowest < 0 || fout < lowest)
453 /* FIXME: See above. gcd interval for lambdas is disabled. */
456 // /* Lambdas written at an interval of the lowest common denominator
457 // * of other output */
458 // set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA,
459 // "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
460 // TNG_GZIP_COMPRESSION);
462 // if (gcd < lowest / 10)
464 // gmx_warning("The lowest common denominator of trajectory output is "
465 // "every %d step(s), whereas the shortest output interval "
466 // "is every %d steps.", gcd, lowest);
472 void gmx_tng_prepare_md_writing(tng_trajectory_t tng,
473 const gmx_mtop_t *mtop,
474 const t_inputrec *ir)
477 gmx_tng_add_mtop(tng, mtop);
478 set_writing_intervals(tng, FALSE, ir);
479 tng_time_per_frame_set(tng, ir->delta_t * PICO);
481 GMX_UNUSED_VALUE(tng);
482 GMX_UNUSED_VALUE(mtop);
483 GMX_UNUSED_VALUE(ir);
488 /* Create a TNG molecule representing the selection groups
490 static void add_selection_groups(tng_trajectory_t tng,
491 const gmx_mtop_t *mtop)
493 const gmx_moltype_t *molType;
494 const t_atoms *atoms;
496 const t_resinfo *resInfo;
497 const t_ilist *ilist;
498 int nAtoms = 0, i = 0, j, molIt, atomIt, nameIndex;
500 tng_molecule_t mol, iterMol;
508 /* The name of the TNG molecule containing the selection group is the
509 * same as the name of the selection group. */
510 nameIndex = *mtop->groups.grps[egcCompressedX].nm_ind;
511 groupName = *mtop->groups.grpname[nameIndex];
513 tng_molecule_alloc(tng, &mol);
514 tng_molecule_name_set(tng, mol, groupName);
515 tng_molecule_chain_add(tng, mol, "", &chain);
516 for (molIt = 0; molIt < mtop->nmoltype; molIt++)
518 molType = &mtop->moltype[mtop->molblock[molIt].type];
520 atoms = &molType->atoms;
522 for (j = 0; j < mtop->molblock[molIt].nmol; j++)
524 bool bAtomsAdded = FALSE;
525 for (atomIt = 0; atomIt < atoms->nr; atomIt++, i++)
530 if (ggrpnr(&mtop->groups, egcCompressedX, i) != 0)
534 at = &atoms->atom[atomIt];
537 resInfo = &atoms->resinfo[at->resind];
538 /* FIXME: When TNG supports both residue index and residue
539 * number the latter should be used. */
540 res_name = *resInfo->name;
541 res_id = at->resind + 1;
545 res_name = (char *)"";
548 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res)
551 /* Since there is ONE chain for selection groups do not keep the
552 * original residue IDs - otherwise there might be conflicts. */
553 tng_chain_residue_add(tng, chain, res_name, &res);
555 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIt]),
556 *(atoms->atomtype[atomIt]),
557 atom_offset + atomIt, &atom);
564 for (int k = 0; k < F_NRE; k++)
568 ilist = &molType->ilist[k];
572 while (l < ilist->nr)
575 atom1 = ilist->iatoms[l] + atom_offset;
576 atom2 = ilist->iatoms[l+1] + atom_offset;
577 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0 &&
578 ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
580 tng_molecule_bond_add(tng, mol, ilist->iatoms[l],
581 ilist->iatoms[l+1], &tngBond);
588 /* Settle is described using three atoms */
589 ilist = &molType->ilist[F_SETTLE];
593 while (l < ilist->nr)
595 int atom1, atom2, atom3;
596 atom1 = ilist->iatoms[l] + atom_offset;
597 atom2 = ilist->iatoms[l+1] + atom_offset;
598 atom3 = ilist->iatoms[l+2] + atom_offset;
599 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0)
601 if (ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
603 tng_molecule_bond_add(tng, mol, atom1,
606 if (ggrpnr(&mtop->groups, egcCompressedX, atom3) == 0)
608 tng_molecule_bond_add(tng, mol, atom1,
616 atom_offset += atoms->nr;
621 tng_molecule_existing_add(tng, &mol);
622 tng_molecule_cnt_set(tng, mol, 1);
623 tng_num_molecule_types_get(tng, &nMols);
624 for (gmx_int64_t k = 0; k < nMols; k++)
626 tng_molecule_of_index_get(tng, k, &iterMol);
631 tng_molecule_cnt_set(tng, iterMol, 0);
636 tng_molecule_free(tng, &mol);
641 void gmx_tng_set_compression_precision(tng_trajectory_t tng,
645 tng_compression_precision_set(tng, 1.0/prec);
647 GMX_UNUSED_VALUE(tng);
648 GMX_UNUSED_VALUE(prec);
652 void gmx_tng_prepare_low_prec_writing(tng_trajectory_t tng,
653 const gmx_mtop_t *mtop,
654 const t_inputrec *ir)
657 gmx_tng_add_mtop(tng, mtop);
658 add_selection_groups(tng, mtop);
659 set_writing_intervals(tng, TRUE, ir);
660 tng_time_per_frame_set(tng, ir->delta_t * PICO);
661 gmx_tng_set_compression_precision(tng, ir->x_compression_precision);
663 GMX_UNUSED_VALUE(tng);
664 GMX_UNUSED_VALUE(mtop);
665 GMX_UNUSED_VALUE(ir);
669 void gmx_fwrite_tng(tng_trajectory_t tng,
670 const gmx_bool bUseLossyCompression,
672 real elapsedPicoSeconds,
681 typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
691 static write_data_func_pointer write_data = tng_util_generic_with_time_double_write;
693 static write_data_func_pointer write_data = tng_util_generic_with_time_write;
695 double elapsedSeconds = elapsedPicoSeconds * PICO;
696 gmx_int64_t nParticles;
702 /* This function might get called when the type of the
703 compressed trajectory is actually XTC. So we exit and move
708 tng_num_particles_get(tng, &nParticles);
709 if (nAtoms != (int)nParticles)
711 tng_implicit_num_particles_set(tng, nAtoms);
714 if (bUseLossyCompression)
716 compression = TNG_TNG_COMPRESSION;
720 compression = TNG_GZIP_COMPRESSION;
723 /* The writing is done using write_data, which writes float or double
724 * depending on the GROMACS compilation. */
727 GMX_ASSERT(box, "Need a non-NULL box if positions are written");
729 if (write_data(tng, step, elapsedSeconds,
730 reinterpret_cast<const real *>(x),
731 3, TNG_TRAJ_POSITIONS, "POSITIONS",
732 TNG_PARTICLE_BLOCK_DATA,
733 compression) != TNG_SUCCESS)
735 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
737 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
738 * compression for box shape regardless of output mode */
739 if (write_data(tng, step, elapsedSeconds,
740 reinterpret_cast<const real *>(box),
741 9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
742 TNG_NON_PARTICLE_BLOCK_DATA,
743 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
745 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
751 if (write_data(tng, step, elapsedSeconds,
752 reinterpret_cast<const real *>(v),
753 3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
754 TNG_PARTICLE_BLOCK_DATA,
755 compression) != TNG_SUCCESS)
757 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
763 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
764 * compression for forces regardless of output mode */
765 if (write_data(tng, step, elapsedSeconds,
766 reinterpret_cast<const real *>(f),
767 3, TNG_TRAJ_FORCES, "FORCES",
768 TNG_PARTICLE_BLOCK_DATA,
769 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
771 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
775 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
776 * compression for lambdas regardless of output mode */
777 if (write_data(tng, step, elapsedSeconds,
778 reinterpret_cast<const real *>(&lambda),
779 1, TNG_GMX_LAMBDA, "LAMBDAS",
780 TNG_NON_PARTICLE_BLOCK_DATA,
781 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
783 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
786 GMX_UNUSED_VALUE(tng);
787 GMX_UNUSED_VALUE(bUseLossyCompression);
788 GMX_UNUSED_VALUE(step);
789 GMX_UNUSED_VALUE(elapsedPicoSeconds);
790 GMX_UNUSED_VALUE(lambda);
791 GMX_UNUSED_VALUE(box);
792 GMX_UNUSED_VALUE(nAtoms);
799 void fflush_tng(tng_trajectory_t tng)
806 tng_frame_set_premature_write(tng, TNG_USE_HASH);
808 GMX_UNUSED_VALUE(tng);
812 float gmx_tng_get_time_of_final_frame(tng_trajectory_t tng)
819 tng_num_frames_get(tng, &nFrames);
820 tng_util_time_of_frame_get(tng, nFrames - 1, &time);
825 GMX_UNUSED_VALUE(tng);