2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014, 2015 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.
42 #include "tng/tng_io.h"
45 #include "gromacs/fileio/gmxfio.h"
46 #include "gromacs/legacyheaders/copyrite.h"
47 #include "gromacs/legacyheaders/types/ifunc.h"
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/basedefinitions.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/programcontext.h"
55 #include "gromacs/utility/sysinfo.h"
57 static const char *modeToVerb(char mode)
72 gmx_fatal(FARGS, "Invalid file opening mode %c", mode);
79 void gmx_tng_open(const char *filename,
81 tng_trajectory_t *tng)
84 /* First check whether we have to make a backup,
85 * only for writing, not for read or append.
89 make_backup(filename);
92 /* tng must not be pointing at already allocated memory.
93 * Memory will be allocated by tng_util_trajectory_open() and must
94 * later on be freed by tng_util_trajectory_close(). */
95 if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
97 /* TNG does return more than one degree of error, but there is
98 no use case for GROMACS handling the non-fatal errors
101 "%s while opening %s for %s",
102 gmx_strerror("file"),
107 if (mode == 'w' || mode == 'a')
110 gmx_gethostname(hostname, 256);
113 tng_first_computer_name_set(*tng, hostname);
117 tng_last_computer_name_set(*tng, hostname);
120 char programInfo[256];
121 const char *precisionString = "";
123 precisionString = " (double precision)";
125 sprintf(programInfo, "%.100s, %.128s%.24s",
126 gmx::getProgramContext().displayName(),
127 GromacsVersion(), precisionString);
130 tng_first_program_name_set(*tng, programInfo);
134 tng_last_program_name_set(*tng, programInfo);
138 if (!gmx_getusername(username, 256))
142 tng_first_user_name_set(*tng, username);
146 tng_last_user_name_set(*tng, username);
147 tng_file_headers_write(*tng, TNG_USE_HASH);
152 gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
153 GMX_UNUSED_VALUE(filename);
154 GMX_UNUSED_VALUE(mode);
155 GMX_UNUSED_VALUE(tng);
159 void gmx_tng_close(tng_trajectory_t *tng)
161 /* We have to check that tng is set because
162 * tng_util_trajectory_close wants to return a NULL in it, and
163 * gives a fatal error if it is NULL. */
167 tng_util_trajectory_close(tng);
170 GMX_UNUSED_VALUE(tng);
175 static void addTngMoleculeFromTopology(tng_trajectory_t tng,
176 const char *moleculeName,
177 const t_atoms *atoms,
178 gmx_int64_t numMolecules,
179 tng_molecule_t *tngMol)
181 tng_chain_t tngChain = NULL;
182 tng_residue_t tngRes = NULL;
184 if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
186 gmx_file("Cannot add molecule to TNG molecular system.");
189 /* FIXME: The TNG atoms should contain mass and atomB info (for free
190 * energy calculations), i.e. in when it's available in TNG (2.0). */
191 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++)
193 const t_atom *at = &atoms->atom[atomIndex];
194 /* FIXME: Currently the TNG API can only add atoms belonging to a
195 * residue and chain. Wait for TNG 2.0*/
198 const t_resinfo *resInfo = &atoms->resinfo[at->resind];
199 char chainName[2] = {resInfo->chainid, 0};
200 tng_atom_t tngAtom = NULL;
205 prevAtom = &atoms->atom[atomIndex - 1];
212 /* If this is the first atom or if the residue changed add the
213 * residue to the TNG molecular system. */
214 if (!prevAtom || resInfo != &atoms->resinfo[prevAtom->resind])
216 /* If this is the first atom or if the chain changed add
217 * the chain to the TNG molecular system. */
218 if (!prevAtom || resInfo->chainid !=
219 atoms->resinfo[prevAtom->resind].chainid)
221 tng_molecule_chain_add(tng, *tngMol, chainName,
224 /* FIXME: When TNG supports both residue index and residue
225 * number the latter should be used. Wait for TNG 2.0*/
226 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
228 tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIndex]), *(atoms->atomtype[atomIndex]), &tngAtom);
231 tng_molecule_cnt_set(tng, *tngMol, numMolecules);
234 void gmx_tng_add_mtop(tng_trajectory_t tng,
235 const gmx_mtop_t *mtop)
238 const t_ilist *ilist;
243 /* No topology information available to add. */
247 for (int molIndex = 0; molIndex < mtop->nmolblock; molIndex++)
249 tng_molecule_t tngMol = NULL;
250 const gmx_moltype_t *molType =
251 &mtop->moltype[mtop->molblock[molIndex].type];
253 /* Add a molecule to the TNG trajectory with the same name as the
254 * current molecule. */
255 addTngMoleculeFromTopology(tng,
258 mtop->molblock[molIndex].nmol,
261 /* Bonds have to be deduced from interactions (constraints etc). Different
262 * interactions have different sets of parameters. */
263 /* Constraints are specified using two atoms */
264 for (i = 0; i < F_NRE; i++)
268 ilist = &molType->ilist[i];
272 while (j < ilist->nr)
274 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
280 /* Settle is described using three atoms */
281 ilist = &molType->ilist[F_SETTLE];
285 while (j < ilist->nr)
287 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
288 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+2], &tngBond);
295 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
296 * if they are positive.
298 * If only one of n1 and n2 is positive, then return it.
299 * If neither n1 or n2 is positive, then return -1. */
301 greatest_common_divisor_if_positive(int n1, int n2)
305 return (0 >= n2) ? -1 : n2;
312 /* We have a non-trivial greatest common divisor to compute. */
313 return gmx_greatest_common_divisor(n1, n2);
316 /* By default try to write 100 frames (of actual output) in each frame set.
317 * This number is the number of outputs of the most frequently written data
318 * type per frame set.
319 * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
320 * setups regarding compression efficiency and compression time. Make this
321 * a hidden command-line option? */
322 const int defaultFramesPerFrameSet = 100;
324 /*! \libinternal \brief Set the number of frames per frame
325 * set according to output intervals.
326 * The default is that 100 frames are written of the data
327 * that is written most often. */
328 static void tng_set_frames_per_frame_set(tng_trajectory_t tng,
329 const gmx_bool bUseLossyCompression,
330 const t_inputrec *ir)
334 /* Set the number of frames per frame set to contain at least
335 * defaultFramesPerFrameSet of the lowest common denominator of
336 * the writing interval of positions and velocities. */
337 /* FIXME after 5.0: consider nstenergy also? */
338 if (bUseLossyCompression)
340 gcd = ir->nstxout_compressed;
344 gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
345 gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
352 tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
355 /*! \libinternal \brief Set the data-writing intervals, and number of
356 * frames per frame set */
357 static void set_writing_intervals(tng_trajectory_t tng,
358 const gmx_bool bUseLossyCompression,
359 const t_inputrec *ir)
361 /* Define pointers to specific writing functions depending on if we
362 * write float or double data */
363 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
371 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
373 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
375 int xout, vout, fout;
376 // int gcd = -1, lowest = -1;
379 tng_set_frames_per_frame_set(tng, bUseLossyCompression, ir);
381 if (bUseLossyCompression)
383 xout = ir->nstxout_compressed;
386 compression = TNG_TNG_COMPRESSION;
393 compression = TNG_GZIP_COMPRESSION;
397 set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS,
398 "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
400 /* The design of TNG makes it awkward to try to write a box
401 * with multiple periodicities, which might be co-prime. Since
402 * the use cases for the box with a frame consisting only of
403 * velocities seem low, for now we associate box writing with
404 * position writing. */
405 set_writing_interval(tng, xout, 9, TNG_TRAJ_BOX_SHAPE,
406 "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA,
407 TNG_GZIP_COMPRESSION);
408 /* TODO: if/when we write energies to TNG also, reconsider how
409 * and when box information is written, because GROMACS
410 * behaviour pre-5.0 was to write the box with every
411 * trajectory frame and every energy frame, and probably
412 * people depend on this. */
414 /* TODO: If we need to write lambda values at steps when
415 * positions (or other data) are not also being written, then
416 * code in mdoutf.c will need to match however that is
418 set_writing_interval(tng, xout, 1, TNG_GMX_LAMBDA,
419 "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
420 TNG_GZIP_COMPRESSION);
422 /* FIXME: gcd and lowest currently not used. */
423 // gcd = greatest_common_divisor_if_positive(gcd, xout);
424 // if (lowest < 0 || xout < lowest)
431 set_writing_interval(tng, ir->nstvout, 3, TNG_TRAJ_VELOCITIES,
432 "VELOCITIES", TNG_PARTICLE_BLOCK_DATA,
435 /* FIXME: gcd and lowest currently not used. */
436 // gcd = greatest_common_divisor_if_positive(gcd, vout);
437 // if (lowest < 0 || vout < lowest)
444 set_writing_interval(tng, ir->nstfout, 3, TNG_TRAJ_FORCES,
445 "FORCES", TNG_PARTICLE_BLOCK_DATA,
446 TNG_GZIP_COMPRESSION);
448 /* FIXME: gcd and lowest currently not used. */
449 // gcd = greatest_common_divisor_if_positive(gcd, fout);
450 // if (lowest < 0 || fout < lowest)
455 /* FIXME: See above. gcd interval for lambdas is disabled. */
458 // /* Lambdas written at an interval of the lowest common denominator
459 // * of other output */
460 // set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA,
461 // "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
462 // TNG_GZIP_COMPRESSION);
464 // if (gcd < lowest / 10)
466 // gmx_warning("The lowest common denominator of trajectory output is "
467 // "every %d step(s), whereas the shortest output interval "
468 // "is every %d steps.", gcd, lowest);
474 void gmx_tng_prepare_md_writing(tng_trajectory_t tng,
475 const gmx_mtop_t *mtop,
476 const t_inputrec *ir)
479 gmx_tng_add_mtop(tng, mtop);
480 set_writing_intervals(tng, FALSE, ir);
481 tng_time_per_frame_set(tng, ir->delta_t * PICO);
483 GMX_UNUSED_VALUE(tng);
484 GMX_UNUSED_VALUE(mtop);
485 GMX_UNUSED_VALUE(ir);
490 /* Check if all atoms in the molecule system are selected
491 * by a selection group of type specified by gtype. */
492 static gmx_bool all_atoms_selected(const gmx_mtop_t *mtop,
495 const gmx_moltype_t *molType;
496 const t_atoms *atoms;
498 /* Iterate through all atoms in the molecule system and
499 * check if they belong to a selection group of the
501 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
503 molType = &mtop->moltype[mtop->molblock[molIndex].type];
505 atoms = &molType->atoms;
507 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
509 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
511 if (ggrpnr(&mtop->groups, gtype, i) != 0)
521 /* Create TNG molecules which will represent each of the selection
522 * groups for writing. But do that only if there is actually a
523 * specified selection group and it is not the whole system.
524 * TODO: Currently the only selection that is taken into account
525 * is egcCompressedX, but other selections should be added when
526 * e.g. writing energies is implemented.
528 static void add_selection_groups(tng_trajectory_t tng,
529 const gmx_mtop_t *mtop)
531 const gmx_moltype_t *molType;
532 const t_atoms *atoms;
534 const t_resinfo *resInfo;
535 const t_ilist *ilist;
538 tng_molecule_t mol, iterMol;
546 /* TODO: When the TNG molecules block is more flexible TNG selection
547 * groups should not need all atoms specified. It should be possible
548 * just to specify what molecules are selected (and/or which atoms in
549 * the molecule) and how many (if applicable). */
551 /* If no atoms are selected we do not need to create a
552 * TNG selection group molecule. */
553 if (mtop->groups.ngrpnr[egcCompressedX] == 0)
558 /* If all atoms are selected we do not have to create a selection
559 * group molecule in the TNG molecule system. */
560 if (all_atoms_selected(mtop, egcCompressedX))
565 /* The name of the TNG molecule containing the selection group is the
566 * same as the name of the selection group. */
567 nameIndex = *mtop->groups.grps[egcCompressedX].nm_ind;
568 groupName = *mtop->groups.grpname[nameIndex];
570 tng_molecule_alloc(tng, &mol);
571 tng_molecule_name_set(tng, mol, groupName);
572 tng_molecule_chain_add(tng, mol, "", &chain);
573 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
575 molType = &mtop->moltype[mtop->molblock[molIndex].type];
577 atoms = &molType->atoms;
579 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
581 bool bAtomsAdded = FALSE;
582 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
587 if (ggrpnr(&mtop->groups, egcCompressedX, i) != 0)
591 at = &atoms->atom[atomIndex];
594 resInfo = &atoms->resinfo[at->resind];
595 /* FIXME: When TNG supports both residue index and residue
596 * number the latter should be used. */
597 res_name = *resInfo->name;
598 res_id = at->resind + 1;
602 res_name = (char *)"";
605 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res)
608 /* Since there is ONE chain for selection groups do not keep the
609 * original residue IDs - otherwise there might be conflicts. */
610 tng_chain_residue_add(tng, chain, res_name, &res);
612 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
613 *(atoms->atomtype[atomIndex]),
614 atom_offset + atomIndex, &atom);
620 for (int k = 0; k < F_NRE; k++)
624 ilist = &molType->ilist[k];
628 while (l < ilist->nr)
631 atom1 = ilist->iatoms[l] + atom_offset;
632 atom2 = ilist->iatoms[l+1] + atom_offset;
633 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0 &&
634 ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
636 tng_molecule_bond_add(tng, mol, ilist->iatoms[l],
637 ilist->iatoms[l+1], &tngBond);
644 /* Settle is described using three atoms */
645 ilist = &molType->ilist[F_SETTLE];
649 while (l < ilist->nr)
651 int atom1, atom2, atom3;
652 atom1 = ilist->iatoms[l] + atom_offset;
653 atom2 = ilist->iatoms[l+1] + atom_offset;
654 atom3 = ilist->iatoms[l+2] + atom_offset;
655 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0)
657 if (ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
659 tng_molecule_bond_add(tng, mol, atom1,
662 if (ggrpnr(&mtop->groups, egcCompressedX, atom3) == 0)
664 tng_molecule_bond_add(tng, mol, atom1,
672 atom_offset += atoms->nr;
675 tng_molecule_existing_add(tng, &mol);
676 tng_molecule_cnt_set(tng, mol, 1);
677 tng_num_molecule_types_get(tng, &nMols);
678 for (gmx_int64_t k = 0; k < nMols; k++)
680 tng_molecule_of_index_get(tng, k, &iterMol);
685 tng_molecule_cnt_set(tng, iterMol, 0);
690 void gmx_tng_set_compression_precision(tng_trajectory_t tng,
694 tng_compression_precision_set(tng, prec);
696 GMX_UNUSED_VALUE(tng);
697 GMX_UNUSED_VALUE(prec);
701 void gmx_tng_prepare_low_prec_writing(tng_trajectory_t tng,
702 const gmx_mtop_t *mtop,
703 const t_inputrec *ir)
706 gmx_tng_add_mtop(tng, mtop);
707 add_selection_groups(tng, mtop);
708 set_writing_intervals(tng, TRUE, ir);
709 tng_time_per_frame_set(tng, ir->delta_t * PICO);
710 gmx_tng_set_compression_precision(tng, ir->x_compression_precision);
712 GMX_UNUSED_VALUE(tng);
713 GMX_UNUSED_VALUE(mtop);
714 GMX_UNUSED_VALUE(ir);
718 void gmx_fwrite_tng(tng_trajectory_t tng,
719 const gmx_bool bUseLossyCompression,
721 real elapsedPicoSeconds,
730 typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
740 static write_data_func_pointer write_data = tng_util_generic_with_time_double_write;
742 static write_data_func_pointer write_data = tng_util_generic_with_time_write;
744 double elapsedSeconds = elapsedPicoSeconds * PICO;
745 gmx_int64_t nParticles;
751 /* This function might get called when the type of the
752 compressed trajectory is actually XTC. So we exit and move
757 tng_num_particles_get(tng, &nParticles);
758 if (nAtoms != (int)nParticles)
760 tng_implicit_num_particles_set(tng, nAtoms);
763 if (bUseLossyCompression)
765 compression = TNG_TNG_COMPRESSION;
769 compression = TNG_GZIP_COMPRESSION;
772 /* The writing is done using write_data, which writes float or double
773 * depending on the GROMACS compilation. */
776 GMX_ASSERT(box, "Need a non-NULL box if positions are written");
778 if (write_data(tng, step, elapsedSeconds,
779 reinterpret_cast<const real *>(x),
780 3, TNG_TRAJ_POSITIONS, "POSITIONS",
781 TNG_PARTICLE_BLOCK_DATA,
782 compression) != TNG_SUCCESS)
784 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
786 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
787 * compression for box shape regardless of output mode */
788 if (write_data(tng, step, elapsedSeconds,
789 reinterpret_cast<const real *>(box),
790 9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
791 TNG_NON_PARTICLE_BLOCK_DATA,
792 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
794 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
800 if (write_data(tng, step, elapsedSeconds,
801 reinterpret_cast<const real *>(v),
802 3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
803 TNG_PARTICLE_BLOCK_DATA,
804 compression) != TNG_SUCCESS)
806 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
812 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
813 * compression for forces regardless of output mode */
814 if (write_data(tng, step, elapsedSeconds,
815 reinterpret_cast<const real *>(f),
816 3, TNG_TRAJ_FORCES, "FORCES",
817 TNG_PARTICLE_BLOCK_DATA,
818 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
820 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
824 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
825 * compression for lambdas regardless of output mode */
826 if (write_data(tng, step, elapsedSeconds,
827 reinterpret_cast<const real *>(&lambda),
828 1, TNG_GMX_LAMBDA, "LAMBDAS",
829 TNG_NON_PARTICLE_BLOCK_DATA,
830 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
832 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
835 GMX_UNUSED_VALUE(tng);
836 GMX_UNUSED_VALUE(bUseLossyCompression);
837 GMX_UNUSED_VALUE(step);
838 GMX_UNUSED_VALUE(elapsedPicoSeconds);
839 GMX_UNUSED_VALUE(lambda);
840 GMX_UNUSED_VALUE(box);
841 GMX_UNUSED_VALUE(nAtoms);
848 void fflush_tng(tng_trajectory_t tng)
855 tng_frame_set_premature_write(tng, TNG_USE_HASH);
857 GMX_UNUSED_VALUE(tng);
861 float gmx_tng_get_time_of_final_frame(tng_trajectory_t tng)
868 tng_num_frames_get(tng, &nFrames);
869 tng_util_time_of_frame_get(tng, nFrames - 1, &time);
874 GMX_UNUSED_VALUE(tng);