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.
46 #include "tng/tng_io.h"
49 #include "gromacs/legacyheaders/copyrite.h"
50 #include "gromacs/legacyheaders/types/ifunc.h"
52 #include "gromacs/fileio/gmxfio.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/basenetwork.h"
57 #include "gromacs/utility/common.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/programcontext.h"
62 static const char *modeToVerb(char mode)
77 gmx_fatal(FARGS, "Invalid file opening mode %c", mode);
84 void gmx_tng_open(const char *filename,
86 tng_trajectory_t *tng)
89 /* First check whether we have to make a backup,
90 * only for writing, not for read or append.
95 /* only make backups for normal gromacs */
96 make_backup(filename);
100 /* tng must not be pointing at already allocated memory.
101 * Memory will be allocated by tng_util_trajectory_open() and must
102 * later on be freed by tng_util_trajectory_close(). */
103 if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
105 /* TNG does return more than one degree of error, but there is
106 no use case for GROMACS handling the non-fatal errors
109 "%s while opening %s for %s",
110 gmx_strerror("file"),
115 if (mode == 'w' || mode == 'a')
117 /* FIXME in TNG: When adding data to the header, subsequent blocks might get
118 * overwritten. This could be solved by moving the first trajectory
119 * frame set(s) to the end of the file. Could that cause other problems,
120 * e.g. when continuing a simulation? */
122 gmx_gethostname(hostname, 256);
125 tng_first_computer_name_set(*tng, hostname);
127 /* TODO: This should be implemented when the above fixme is done (adding data to
131 // tng_last_computer_name_set(*tng, hostname);
134 char programInfo[256];
135 const char *precisionString = "";
137 precisionString = " (double precision)";
139 sprintf(programInfo, "%.100s, %.128s%.24s",
140 gmx::getProgramContext().displayName(),
141 GromacsVersion(), precisionString);
144 tng_first_program_name_set(*tng, programInfo);
146 /* TODO: This should be implemented when the above fixme is done (adding data to
150 // tng_last_program_name_set(*tng, programInfo);
155 getlogin_r(username, 256);
158 tng_first_user_name_set(*tng, username);
160 /* TODO: This should be implemented when the above fixme is done (adding data to
164 // tng_last_user_name_set(*tng, username);
169 gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
170 GMX_UNUSED_VALUE(filename);
171 GMX_UNUSED_VALUE(mode);
172 GMX_UNUSED_VALUE(tng);
176 void gmx_tng_close(tng_trajectory_t *tng)
178 /* We have to check that tng is set because
179 * tng_util_trajectory_close wants to return a NULL in it, and
180 * gives a fatal error if it is NULL. */
184 tng_util_trajectory_close(tng);
187 GMX_UNUSED_VALUE(tng);
192 static void addTngMoleculeFromTopology(tng_trajectory_t tng,
193 const char *moleculeName,
194 const t_atoms *atoms,
195 gmx_int64_t numMolecules,
196 tng_molecule_t *tngMol)
198 tng_chain_t tngChain = NULL;
199 tng_residue_t tngRes = NULL;
201 if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
203 gmx_file("Cannot add molecule to TNG molecular system.");
206 /* FIXME: The TNG atoms should contain mass and atomB info (for free
207 * energy calculations), i.e. in when it's available in TNG (2.0). */
208 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++)
210 const t_atom *at = &atoms->atom[atomIndex];
211 /* FIXME: Currently the TNG API can only add atoms belonging to a
212 * residue and chain. Wait for TNG 2.0*/
215 const t_resinfo *resInfo = &atoms->resinfo[at->resind];
216 char chainName[2] = {resInfo->chainid, 0};
217 tng_atom_t tngAtom = NULL;
222 prevAtom = &atoms->atom[atomIndex - 1];
229 /* If this is the first atom or if the residue changed add the
230 * residue to the TNG molecular system. */
231 if (!prevAtom || resInfo != &atoms->resinfo[prevAtom->resind])
233 /* If this is the first atom or if the chain changed add
234 * the chain to the TNG molecular system. */
235 if (!prevAtom || resInfo->chainid !=
236 atoms->resinfo[prevAtom->resind].chainid)
238 tng_molecule_chain_add(tng, *tngMol, chainName,
241 /* FIXME: When TNG supports both residue index and residue
242 * number the latter should be used. Wait for TNG 2.0*/
243 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
245 tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIndex]), *(atoms->atomtype[atomIndex]), &tngAtom);
248 tng_molecule_cnt_set(tng, *tngMol, numMolecules);
251 void gmx_tng_add_mtop(tng_trajectory_t tng,
252 const gmx_mtop_t *mtop)
255 const t_ilist *ilist;
260 /* No topology information available to add. */
264 for (int molIndex = 0; molIndex < mtop->nmolblock; molIndex++)
266 tng_molecule_t tngMol = NULL;
267 const gmx_moltype_t *molType =
268 &mtop->moltype[mtop->molblock[molIndex].type];
270 /* Add a molecule to the TNG trajectory with the same name as the
271 * current molecule. */
272 addTngMoleculeFromTopology(tng,
275 mtop->molblock[molIndex].nmol,
278 /* Bonds have to be deduced from interactions (constraints etc). Different
279 * interactions have different sets of parameters. */
280 /* Constraints are specified using two atoms */
281 for (i = 0; i < F_NRE; i++)
285 ilist = &molType->ilist[i];
289 while (j < ilist->nr)
291 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
297 /* Settle is described using three atoms */
298 ilist = &molType->ilist[F_SETTLE];
302 while (j < ilist->nr)
304 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
305 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+2], &tngBond);
312 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
313 * if they are positive.
315 * If only one of n1 and n2 is positive, then return it.
316 * If neither n1 or n2 is positive, then return -1. */
318 greatest_common_divisor_if_positive(int n1, int n2)
322 return (0 >= n2) ? -1 : n2;
329 /* We have a non-trivial greatest common divisor to compute. */
330 return gmx_greatest_common_divisor(n1, n2);
333 /* By default try to write 100 frames (of actual output) in each frame set.
334 * This number is the number of outputs of the most frequently written data
335 * type per frame set.
336 * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
337 * setups regarding compression efficiency and compression time. Make this
338 * a hidden command-line option? */
339 const int defaultFramesPerFrameSet = 100;
341 /*! \libinternal \brief Set the number of frames per frame
342 * set according to output intervals.
343 * The default is that 100 frames are written of the data
344 * that is written most often. */
345 static void tng_set_frames_per_frame_set(tng_trajectory_t tng,
346 const gmx_bool bUseLossyCompression,
347 const t_inputrec *ir)
351 /* Set the number of frames per frame set to contain at least
352 * defaultFramesPerFrameSet of the lowest common denominator of
353 * the writing interval of positions and velocities. */
354 /* FIXME after 5.0: consider nstenergy also? */
355 if (bUseLossyCompression)
357 gcd = ir->nstxout_compressed;
361 gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
362 gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
369 tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
372 /*! \libinternal \brief Set the data-writing intervals, and number of
373 * frames per frame set */
374 static void set_writing_intervals(tng_trajectory_t tng,
375 const gmx_bool bUseLossyCompression,
376 const t_inputrec *ir)
378 /* Define pointers to specific writing functions depending on if we
379 * write float or double data */
380 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
388 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
390 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
392 int xout, vout, fout;
393 // int gcd = -1, lowest = -1;
396 tng_set_frames_per_frame_set(tng, bUseLossyCompression, ir);
398 if (bUseLossyCompression)
400 xout = ir->nstxout_compressed;
403 compression = TNG_TNG_COMPRESSION;
410 compression = TNG_GZIP_COMPRESSION;
414 set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS,
415 "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
417 /* The design of TNG makes it awkward to try to write a box
418 * with multiple periodicities, which might be co-prime. Since
419 * the use cases for the box with a frame consisting only of
420 * velocities seem low, for now we associate box writing with
421 * position writing. */
422 set_writing_interval(tng, xout, 9, TNG_TRAJ_BOX_SHAPE,
423 "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA,
424 TNG_GZIP_COMPRESSION);
425 /* TODO: if/when we write energies to TNG also, reconsider how
426 * and when box information is written, because GROMACS
427 * behaviour pre-5.0 was to write the box with every
428 * trajectory frame and every energy frame, and probably
429 * people depend on this. */
431 /* TODO: If we need to write lambda values at steps when
432 * positions (or other data) are not also being written, then
433 * code in mdoutf.c will need to match however that is
435 set_writing_interval(tng, xout, 1, TNG_GMX_LAMBDA,
436 "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
437 TNG_GZIP_COMPRESSION);
439 /* FIXME: gcd and lowest currently not used. */
440 // gcd = greatest_common_divisor_if_positive(gcd, xout);
441 // if (lowest < 0 || xout < lowest)
448 set_writing_interval(tng, ir->nstvout, 3, TNG_TRAJ_VELOCITIES,
449 "VELOCITIES", TNG_PARTICLE_BLOCK_DATA,
452 /* FIXME: gcd and lowest currently not used. */
453 // gcd = greatest_common_divisor_if_positive(gcd, vout);
454 // if (lowest < 0 || vout < lowest)
461 set_writing_interval(tng, ir->nstfout, 3, TNG_TRAJ_FORCES,
462 "FORCES", TNG_PARTICLE_BLOCK_DATA,
463 TNG_GZIP_COMPRESSION);
465 /* FIXME: gcd and lowest currently not used. */
466 // gcd = greatest_common_divisor_if_positive(gcd, fout);
467 // if (lowest < 0 || fout < lowest)
472 /* FIXME: See above. gcd interval for lambdas is disabled. */
475 // /* Lambdas written at an interval of the lowest common denominator
476 // * of other output */
477 // set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA,
478 // "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
479 // TNG_GZIP_COMPRESSION);
481 // if (gcd < lowest / 10)
483 // gmx_warning("The lowest common denominator of trajectory output is "
484 // "every %d step(s), whereas the shortest output interval "
485 // "is every %d steps.", gcd, lowest);
491 void gmx_tng_prepare_md_writing(tng_trajectory_t tng,
492 const gmx_mtop_t *mtop,
493 const t_inputrec *ir)
496 gmx_tng_add_mtop(tng, mtop);
497 set_writing_intervals(tng, FALSE, ir);
498 tng_time_per_frame_set(tng, ir->delta_t * PICO);
500 GMX_UNUSED_VALUE(tng);
501 GMX_UNUSED_VALUE(mtop);
502 GMX_UNUSED_VALUE(ir);
507 /* Check if all atoms in the molecule system are selected
508 * by a selection group of type specified by gtype. */
509 static gmx_bool all_atoms_selected(const gmx_mtop_t *mtop,
512 const gmx_moltype_t *molType;
513 const t_atoms *atoms;
515 /* Iterate through all atoms in the molecule system and
516 * check if they belong to a selection group of the
518 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
520 molType = &mtop->moltype[mtop->molblock[molIndex].type];
522 atoms = &molType->atoms;
524 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
526 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
528 if (ggrpnr(&mtop->groups, gtype, i) != 0)
538 /* Create TNG molecules which will represent each of the selection
539 * groups for writing. But do that only if there is actually a
540 * specified selection group and it is not the whole system.
541 * TODO: Currently the only selection that is taken into account
542 * is egcCompressedX, but other selections should be added when
543 * e.g. writing energies is implemented.
545 static void add_selection_groups(tng_trajectory_t tng,
546 const gmx_mtop_t *mtop)
548 const gmx_moltype_t *molType;
549 const t_atoms *atoms;
551 const t_resinfo *resInfo;
552 const t_ilist *ilist;
555 tng_molecule_t mol, iterMol;
563 /* TODO: When the TNG molecules block is more flexible TNG selection
564 * groups should not need all atoms specified. It should be possible
565 * just to specify what molecules are selected (and/or which atoms in
566 * the molecule) and how many (if applicable). */
568 /* If no atoms are selected we do not need to create a
569 * TNG selection group molecule. */
570 if (mtop->groups.ngrpnr[egcCompressedX] == 0)
575 /* If all atoms are selected we do not have to create a selection
576 * group molecule in the TNG molecule system. */
577 if (all_atoms_selected(mtop, egcCompressedX))
582 /* The name of the TNG molecule containing the selection group is the
583 * same as the name of the selection group. */
584 nameIndex = *mtop->groups.grps[egcCompressedX].nm_ind;
585 groupName = *mtop->groups.grpname[nameIndex];
587 tng_molecule_alloc(tng, &mol);
588 tng_molecule_name_set(tng, mol, groupName);
589 tng_molecule_chain_add(tng, mol, "", &chain);
590 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
592 molType = &mtop->moltype[mtop->molblock[molIndex].type];
594 atoms = &molType->atoms;
596 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
598 bool bAtomsAdded = FALSE;
599 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
604 if (ggrpnr(&mtop->groups, egcCompressedX, i) != 0)
608 at = &atoms->atom[atomIndex];
611 resInfo = &atoms->resinfo[at->resind];
612 /* FIXME: When TNG supports both residue index and residue
613 * number the latter should be used. */
614 res_name = *resInfo->name;
615 res_id = at->resind + 1;
619 res_name = (char *)"";
622 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res)
625 /* Since there is ONE chain for selection groups do not keep the
626 * original residue IDs - otherwise there might be conflicts. */
627 tng_chain_residue_add(tng, chain, res_name, &res);
629 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
630 *(atoms->atomtype[atomIndex]),
631 atom_offset + atomIndex, &atom);
637 for (int k = 0; k < F_NRE; k++)
641 ilist = &molType->ilist[k];
645 while (l < ilist->nr)
648 atom1 = ilist->iatoms[l] + atom_offset;
649 atom2 = ilist->iatoms[l+1] + atom_offset;
650 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0 &&
651 ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
653 tng_molecule_bond_add(tng, mol, ilist->iatoms[l],
654 ilist->iatoms[l+1], &tngBond);
661 /* Settle is described using three atoms */
662 ilist = &molType->ilist[F_SETTLE];
666 while (l < ilist->nr)
668 int atom1, atom2, atom3;
669 atom1 = ilist->iatoms[l] + atom_offset;
670 atom2 = ilist->iatoms[l+1] + atom_offset;
671 atom3 = ilist->iatoms[l+2] + atom_offset;
672 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0)
674 if (ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
676 tng_molecule_bond_add(tng, mol, atom1,
679 if (ggrpnr(&mtop->groups, egcCompressedX, atom3) == 0)
681 tng_molecule_bond_add(tng, mol, atom1,
689 atom_offset += atoms->nr;
692 tng_molecule_existing_add(tng, &mol);
693 tng_molecule_cnt_set(tng, mol, 1);
694 tng_num_molecule_types_get(tng, &nMols);
695 for (gmx_int64_t k = 0; k < nMols; k++)
697 tng_molecule_of_index_get(tng, k, &iterMol);
702 tng_molecule_cnt_set(tng, iterMol, 0);
707 void gmx_tng_set_compression_precision(tng_trajectory_t tng,
711 tng_compression_precision_set(tng, prec);
713 GMX_UNUSED_VALUE(tng);
714 GMX_UNUSED_VALUE(prec);
718 void gmx_tng_prepare_low_prec_writing(tng_trajectory_t tng,
719 const gmx_mtop_t *mtop,
720 const t_inputrec *ir)
723 gmx_tng_add_mtop(tng, mtop);
724 add_selection_groups(tng, mtop);
725 set_writing_intervals(tng, TRUE, ir);
726 tng_time_per_frame_set(tng, ir->delta_t * PICO);
727 gmx_tng_set_compression_precision(tng, ir->x_compression_precision);
729 GMX_UNUSED_VALUE(tng);
730 GMX_UNUSED_VALUE(mtop);
731 GMX_UNUSED_VALUE(ir);
735 void gmx_fwrite_tng(tng_trajectory_t tng,
736 const gmx_bool bUseLossyCompression,
738 real elapsedPicoSeconds,
747 typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
757 static write_data_func_pointer write_data = tng_util_generic_with_time_double_write;
759 static write_data_func_pointer write_data = tng_util_generic_with_time_write;
761 double elapsedSeconds = elapsedPicoSeconds * PICO;
762 gmx_int64_t nParticles;
768 /* This function might get called when the type of the
769 compressed trajectory is actually XTC. So we exit and move
774 tng_num_particles_get(tng, &nParticles);
775 if (nAtoms != (int)nParticles)
777 tng_implicit_num_particles_set(tng, nAtoms);
780 if (bUseLossyCompression)
782 compression = TNG_TNG_COMPRESSION;
786 compression = TNG_GZIP_COMPRESSION;
789 /* The writing is done using write_data, which writes float or double
790 * depending on the GROMACS compilation. */
793 GMX_ASSERT(box, "Need a non-NULL box if positions are written");
795 if (write_data(tng, step, elapsedSeconds,
796 reinterpret_cast<const real *>(x),
797 3, TNG_TRAJ_POSITIONS, "POSITIONS",
798 TNG_PARTICLE_BLOCK_DATA,
799 compression) != TNG_SUCCESS)
801 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
803 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
804 * compression for box shape regardless of output mode */
805 if (write_data(tng, step, elapsedSeconds,
806 reinterpret_cast<const real *>(box),
807 9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
808 TNG_NON_PARTICLE_BLOCK_DATA,
809 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
811 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
817 if (write_data(tng, step, elapsedSeconds,
818 reinterpret_cast<const real *>(v),
819 3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
820 TNG_PARTICLE_BLOCK_DATA,
821 compression) != TNG_SUCCESS)
823 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
829 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
830 * compression for forces regardless of output mode */
831 if (write_data(tng, step, elapsedSeconds,
832 reinterpret_cast<const real *>(f),
833 3, TNG_TRAJ_FORCES, "FORCES",
834 TNG_PARTICLE_BLOCK_DATA,
835 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
837 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
841 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
842 * compression for lambdas regardless of output mode */
843 if (write_data(tng, step, elapsedSeconds,
844 reinterpret_cast<const real *>(&lambda),
845 1, TNG_GMX_LAMBDA, "LAMBDAS",
846 TNG_NON_PARTICLE_BLOCK_DATA,
847 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
849 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
852 GMX_UNUSED_VALUE(tng);
853 GMX_UNUSED_VALUE(bUseLossyCompression);
854 GMX_UNUSED_VALUE(step);
855 GMX_UNUSED_VALUE(elapsedPicoSeconds);
856 GMX_UNUSED_VALUE(lambda);
857 GMX_UNUSED_VALUE(box);
858 GMX_UNUSED_VALUE(nAtoms);
865 void fflush_tng(tng_trajectory_t tng)
872 tng_frame_set_premature_write(tng, TNG_USE_HASH);
874 GMX_UNUSED_VALUE(tng);
878 float gmx_tng_get_time_of_final_frame(tng_trajectory_t tng)
885 tng_num_frames_get(tng, &nFrames);
886 tng_util_time_of_frame_get(tng, nFrames - 1, &time);
891 GMX_UNUSED_VALUE(tng);