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.
44 #include "tng/tng_io.h"
47 #include "gromacs/legacyheaders/copyrite.h"
48 #include "gromacs/legacyheaders/types/ifunc.h"
50 #include "gromacs/fileio/gmxfio.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/utilities.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/basenetwork.h"
55 #include "gromacs/utility/common.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/gmxassert.h"
58 #include "gromacs/utility/programcontext.h"
60 static const char *modeToVerb(char mode)
74 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.
90 /* only make backups for normal gromacs */
91 make_backup(filename);
95 /* tng must not be pointing at already allocated memory.
96 * Memory will be allocated by tng_util_trajectory_open() and must
97 * later on be freed by tng_util_trajectory_close(). */
98 if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
100 /* TNG does return more than one degree of error, but there is
101 no use case for GROMACS handling the non-fatal errors
104 "%s while opening %s for %s",
105 gmx_strerror("file"),
110 if (mode == 'w' || mode == 'a')
112 /* FIXME in TNG: When adding data to the header, subsequent blocks might get
113 * overwritten. This could be solved by moving the first trajectory
114 * frame set(s) to the end of the file. Could that cause other problems,
115 * e.g. when continuing a simulation? */
117 gmx_gethostname(hostname, 256);
120 tng_first_computer_name_set(*tng, hostname);
122 /* TODO: This should be implemented when the above fixme is done (adding data to
126 // tng_last_computer_name_set(*tng, hostname);
129 char programInfo[256];
130 const char *precisionString = "";
132 precisionString = " (double precision)";
134 sprintf(programInfo, "%.100s, %.128s%.24s",
135 gmx::getProgramContext().displayName(),
136 GromacsVersion(), precisionString);
139 tng_first_program_name_set(*tng, programInfo);
141 /* TODO: This should be implemented when the above fixme is done (adding data to
145 // tng_last_program_name_set(*tng, programInfo);
150 getlogin_r(username, 256);
153 tng_first_user_name_set(*tng, username);
155 /* TODO: This should be implemented when the above fixme is done (adding data to
159 // tng_last_user_name_set(*tng, username);
164 gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
165 GMX_UNUSED_VALUE(filename);
166 GMX_UNUSED_VALUE(mode);
167 GMX_UNUSED_VALUE(tng);
171 void gmx_tng_close(tng_trajectory_t *tng)
173 /* We have to check that tng is set because
174 * tng_util_trajectory_close wants to return a NULL in it, and
175 * gives a fatal error if it is NULL. */
179 tng_util_trajectory_close(tng);
182 GMX_UNUSED_VALUE(tng);
187 static void addTngMoleculeFromTopology(tng_trajectory_t tng,
188 const char *moleculeName,
189 const t_atoms *atoms,
190 gmx_int64_t numMolecules,
191 tng_molecule_t *tngMol)
193 tng_chain_t tngChain = NULL;
194 tng_residue_t tngRes = NULL;
196 if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
198 gmx_file("Cannot add molecule to TNG molecular system.");
201 /* FIXME: The TNG atoms should contain mass and atomB info (for free
202 * energy calculations), i.e. in when it's available in TNG (2.0). */
203 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++)
205 const t_atom *at = &atoms->atom[atomIndex];
206 /* FIXME: Currently the TNG API can only add atoms belonging to a
207 * residue and chain. Wait for TNG 2.0*/
210 const t_resinfo *resInfo = &atoms->resinfo[at->resind];
211 char chainName[2] = {resInfo->chainid, 0};
212 tng_atom_t tngAtom = NULL;
217 prevAtom = &atoms->atom[atomIndex - 1];
224 /* If this is the first atom or if the residue changed add the
225 * residue to the TNG molecular system. */
226 if (!prevAtom || resInfo != &atoms->resinfo[prevAtom->resind])
228 /* If this is the first atom or if the chain changed add
229 * the chain to the TNG molecular system. */
230 if (!prevAtom || resInfo->chainid !=
231 atoms->resinfo[prevAtom->resind].chainid)
233 tng_molecule_chain_add(tng, *tngMol, chainName,
236 /* FIXME: When TNG supports both residue index and residue
237 * number the latter should be used. Wait for TNG 2.0*/
238 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
240 tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIndex]), *(atoms->atomtype[atomIndex]), &tngAtom);
243 tng_molecule_cnt_set(tng, *tngMol, numMolecules);
246 void gmx_tng_add_mtop(tng_trajectory_t tng,
247 const gmx_mtop_t *mtop)
250 const t_ilist *ilist;
255 /* No topology information available to add. */
259 for (int molIndex = 0; molIndex < mtop->nmolblock; molIndex++)
261 tng_molecule_t tngMol = NULL;
262 const gmx_moltype_t *molType =
263 &mtop->moltype[mtop->molblock[molIndex].type];
265 /* Add a molecule to the TNG trajectory with the same name as the
266 * current molecule. */
267 addTngMoleculeFromTopology(tng,
270 mtop->molblock[molIndex].nmol,
273 /* Bonds have to be deduced from interactions (constraints etc). Different
274 * interactions have different sets of parameters. */
275 /* Constraints are specified using two atoms */
276 for (i = 0; i < F_NRE; i++)
280 ilist = &molType->ilist[i];
284 while (j < ilist->nr)
286 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
292 /* Settle is described using three atoms */
293 ilist = &molType->ilist[F_SETTLE];
297 while (j < ilist->nr)
299 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
300 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+2], &tngBond);
307 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
308 * if they are positive.
310 * If only one of n1 and n2 is positive, then return it.
311 * If neither n1 or n2 is positive, then return -1. */
313 greatest_common_divisor_if_positive(int n1, int n2)
317 return (0 >= n2) ? -1 : n2;
324 /* We have a non-trivial greatest common divisor to compute. */
325 return gmx_greatest_common_divisor(n1, n2);
328 /* By default try to write 100 frames (of actual output) in each frame set.
329 * This number is the number of outputs of the most frequently written data
330 * type per frame set.
331 * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
332 * setups regarding compression efficiency and compression time. Make this
333 * a hidden command-line option? */
334 const int defaultFramesPerFrameSet = 100;
336 /*! \libinternal \brief Set the number of frames per frame
337 * set according to output intervals.
338 * The default is that 100 frames are written of the data
339 * that is written most often. */
340 static void tng_set_frames_per_frame_set(tng_trajectory_t tng,
341 const gmx_bool bUseLossyCompression,
342 const t_inputrec *ir)
346 /* Set the number of frames per frame set to contain at least
347 * defaultFramesPerFrameSet of the lowest common denominator of
348 * the writing interval of positions and velocities. */
349 /* FIXME after 5.0: consider nstenergy also? */
350 if (bUseLossyCompression)
352 gcd = ir->nstxout_compressed;
356 gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
357 gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
364 tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
367 /*! \libinternal \brief Set the data-writing intervals, and number of
368 * frames per frame set */
369 static void set_writing_intervals(tng_trajectory_t tng,
370 const gmx_bool bUseLossyCompression,
371 const t_inputrec *ir)
373 /* Define pointers to specific writing functions depending on if we
374 * write float or double data */
375 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
383 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
385 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
387 int xout, vout, fout;
388 // int gcd = -1, lowest = -1;
391 tng_set_frames_per_frame_set(tng, bUseLossyCompression, ir);
393 if (bUseLossyCompression)
395 xout = ir->nstxout_compressed;
398 compression = TNG_TNG_COMPRESSION;
405 compression = TNG_GZIP_COMPRESSION;
409 set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS,
410 "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
412 /* The design of TNG makes it awkward to try to write a box
413 * with multiple periodicities, which might be co-prime. Since
414 * the use cases for the box with a frame consisting only of
415 * velocities seem low, for now we associate box writing with
416 * position writing. */
417 set_writing_interval(tng, xout, 9, TNG_TRAJ_BOX_SHAPE,
418 "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA,
419 TNG_GZIP_COMPRESSION);
420 /* TODO: if/when we write energies to TNG also, reconsider how
421 * and when box information is written, because GROMACS
422 * behaviour pre-5.0 was to write the box with every
423 * trajectory frame and every energy frame, and probably
424 * people depend on this. */
426 /* TODO: If we need to write lambda values at steps when
427 * positions (or other data) are not also being written, then
428 * code in mdoutf.c will need to match however that is
430 set_writing_interval(tng, xout, 1, TNG_GMX_LAMBDA,
431 "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
432 TNG_GZIP_COMPRESSION);
434 /* FIXME: gcd and lowest currently not used. */
435 // gcd = greatest_common_divisor_if_positive(gcd, xout);
436 // if (lowest < 0 || xout < lowest)
443 set_writing_interval(tng, ir->nstvout, 3, TNG_TRAJ_VELOCITIES,
444 "VELOCITIES", TNG_PARTICLE_BLOCK_DATA,
447 /* FIXME: gcd and lowest currently not used. */
448 // gcd = greatest_common_divisor_if_positive(gcd, vout);
449 // if (lowest < 0 || vout < lowest)
456 set_writing_interval(tng, ir->nstfout, 3, TNG_TRAJ_FORCES,
457 "FORCES", TNG_PARTICLE_BLOCK_DATA,
458 TNG_GZIP_COMPRESSION);
460 /* FIXME: gcd and lowest currently not used. */
461 // gcd = greatest_common_divisor_if_positive(gcd, fout);
462 // if (lowest < 0 || fout < lowest)
467 /* FIXME: See above. gcd interval for lambdas is disabled. */
470 // /* Lambdas written at an interval of the lowest common denominator
471 // * of other output */
472 // set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA,
473 // "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
474 // TNG_GZIP_COMPRESSION);
476 // if (gcd < lowest / 10)
478 // gmx_warning("The lowest common denominator of trajectory output is "
479 // "every %d step(s), whereas the shortest output interval "
480 // "is every %d steps.", gcd, lowest);
486 void gmx_tng_prepare_md_writing(tng_trajectory_t tng,
487 const gmx_mtop_t *mtop,
488 const t_inputrec *ir)
491 gmx_tng_add_mtop(tng, mtop);
492 set_writing_intervals(tng, FALSE, ir);
493 tng_time_per_frame_set(tng, ir->delta_t * PICO);
495 GMX_UNUSED_VALUE(tng);
496 GMX_UNUSED_VALUE(mtop);
497 GMX_UNUSED_VALUE(ir);
502 /* Check if all atoms in the molecule system are selected
503 * by a selection group of type specified by gtype. */
504 static gmx_bool all_atoms_selected(const gmx_mtop_t *mtop,
507 const gmx_moltype_t *molType;
508 const t_atoms *atoms;
510 /* Iterate through all atoms in the molecule system and
511 * check if they belong to a selection group of the
513 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
515 molType = &mtop->moltype[mtop->molblock[molIndex].type];
517 atoms = &molType->atoms;
519 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
521 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
523 if (ggrpnr(&mtop->groups, gtype, i) != 0)
533 /* Create TNG molecules which will represent each of the selection
534 * groups for writing. But do that only if there is actually a
535 * specified selection group and it is not the whole system.
536 * TODO: Currently the only selection that is taken into account
537 * is egcCompressedX, but other selections should be added when
538 * e.g. writing energies is implemented.
540 static void add_selection_groups(tng_trajectory_t tng,
541 const gmx_mtop_t *mtop)
543 const gmx_moltype_t *molType;
544 const t_atoms *atoms;
546 const t_resinfo *resInfo;
547 const t_ilist *ilist;
550 tng_molecule_t mol, iterMol;
558 /* TODO: When the TNG molecules block is more flexible TNG selection
559 * groups should not need all atoms specified. It should be possible
560 * just to specify what molecules are selected (and/or which atoms in
561 * the molecule) and how many (if applicable). */
563 /* If no atoms are selected we do not need to create a
564 * TNG selection group molecule. */
565 if (mtop->groups.ngrpnr[egcCompressedX] == 0)
570 /* If all atoms are selected we do not have to create a selection
571 * group molecule in the TNG molecule system. */
572 if (all_atoms_selected(mtop, egcCompressedX))
577 /* The name of the TNG molecule containing the selection group is the
578 * same as the name of the selection group. */
579 nameIndex = *mtop->groups.grps[egcCompressedX].nm_ind;
580 groupName = *mtop->groups.grpname[nameIndex];
582 tng_molecule_alloc(tng, &mol);
583 tng_molecule_name_set(tng, mol, groupName);
584 tng_molecule_chain_add(tng, mol, "", &chain);
585 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
587 molType = &mtop->moltype[mtop->molblock[molIndex].type];
589 atoms = &molType->atoms;
591 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
593 bool bAtomsAdded = FALSE;
594 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
599 if (ggrpnr(&mtop->groups, egcCompressedX, i) != 0)
603 at = &atoms->atom[atomIndex];
606 resInfo = &atoms->resinfo[at->resind];
607 /* FIXME: When TNG supports both residue index and residue
608 * number the latter should be used. */
609 res_name = *resInfo->name;
610 res_id = at->resind + 1;
614 res_name = (char *)"";
617 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res)
620 /* Since there is ONE chain for selection groups do not keep the
621 * original residue IDs - otherwise there might be conflicts. */
622 tng_chain_residue_add(tng, chain, res_name, &res);
624 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
625 *(atoms->atomtype[atomIndex]),
626 atom_offset + atomIndex, &atom);
632 for (int k = 0; k < F_NRE; k++)
636 ilist = &molType->ilist[k];
640 while (l < ilist->nr)
643 atom1 = ilist->iatoms[l] + atom_offset;
644 atom2 = ilist->iatoms[l+1] + atom_offset;
645 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0 &&
646 ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
648 tng_molecule_bond_add(tng, mol, ilist->iatoms[l],
649 ilist->iatoms[l+1], &tngBond);
656 /* Settle is described using three atoms */
657 ilist = &molType->ilist[F_SETTLE];
661 while (l < ilist->nr)
663 int atom1, atom2, atom3;
664 atom1 = ilist->iatoms[l] + atom_offset;
665 atom2 = ilist->iatoms[l+1] + atom_offset;
666 atom3 = ilist->iatoms[l+2] + atom_offset;
667 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0)
669 if (ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
671 tng_molecule_bond_add(tng, mol, atom1,
674 if (ggrpnr(&mtop->groups, egcCompressedX, atom3) == 0)
676 tng_molecule_bond_add(tng, mol, atom1,
684 atom_offset += atoms->nr;
687 tng_molecule_existing_add(tng, &mol);
688 tng_molecule_cnt_set(tng, mol, 1);
689 tng_num_molecule_types_get(tng, &nMols);
690 for (gmx_int64_t k = 0; k < nMols; k++)
692 tng_molecule_of_index_get(tng, k, &iterMol);
697 tng_molecule_cnt_set(tng, iterMol, 0);
702 void gmx_tng_set_compression_precision(tng_trajectory_t tng,
706 tng_compression_precision_set(tng, prec);
708 GMX_UNUSED_VALUE(tng);
709 GMX_UNUSED_VALUE(prec);
713 void gmx_tng_prepare_low_prec_writing(tng_trajectory_t tng,
714 const gmx_mtop_t *mtop,
715 const t_inputrec *ir)
718 gmx_tng_add_mtop(tng, mtop);
719 add_selection_groups(tng, mtop);
720 set_writing_intervals(tng, TRUE, ir);
721 tng_time_per_frame_set(tng, ir->delta_t * PICO);
722 gmx_tng_set_compression_precision(tng, ir->x_compression_precision);
724 GMX_UNUSED_VALUE(tng);
725 GMX_UNUSED_VALUE(mtop);
726 GMX_UNUSED_VALUE(ir);
730 void gmx_fwrite_tng(tng_trajectory_t tng,
731 const gmx_bool bUseLossyCompression,
733 real elapsedPicoSeconds,
742 typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
752 static write_data_func_pointer write_data = tng_util_generic_with_time_double_write;
754 static write_data_func_pointer write_data = tng_util_generic_with_time_write;
756 double elapsedSeconds = elapsedPicoSeconds * PICO;
757 gmx_int64_t nParticles;
763 /* This function might get called when the type of the
764 compressed trajectory is actually XTC. So we exit and move
769 tng_num_particles_get(tng, &nParticles);
770 if (nAtoms != (int)nParticles)
772 tng_implicit_num_particles_set(tng, nAtoms);
775 if (bUseLossyCompression)
777 compression = TNG_TNG_COMPRESSION;
781 compression = TNG_GZIP_COMPRESSION;
784 /* The writing is done using write_data, which writes float or double
785 * depending on the GROMACS compilation. */
788 GMX_ASSERT(box, "Need a non-NULL box if positions are written");
790 if (write_data(tng, step, elapsedSeconds,
791 reinterpret_cast<const real *>(x),
792 3, TNG_TRAJ_POSITIONS, "POSITIONS",
793 TNG_PARTICLE_BLOCK_DATA,
794 compression) != TNG_SUCCESS)
796 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
798 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
799 * compression for box shape regardless of output mode */
800 if (write_data(tng, step, elapsedSeconds,
801 reinterpret_cast<const real *>(box),
802 9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
803 TNG_NON_PARTICLE_BLOCK_DATA,
804 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
806 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
812 if (write_data(tng, step, elapsedSeconds,
813 reinterpret_cast<const real *>(v),
814 3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
815 TNG_PARTICLE_BLOCK_DATA,
816 compression) != TNG_SUCCESS)
818 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 forces regardless of output mode */
826 if (write_data(tng, step, elapsedSeconds,
827 reinterpret_cast<const real *>(f),
828 3, TNG_TRAJ_FORCES, "FORCES",
829 TNG_PARTICLE_BLOCK_DATA,
830 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
832 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
836 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
837 * compression for lambdas regardless of output mode */
838 if (write_data(tng, step, elapsedSeconds,
839 reinterpret_cast<const real *>(&lambda),
840 1, TNG_GMX_LAMBDA, "LAMBDAS",
841 TNG_NON_PARTICLE_BLOCK_DATA,
842 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
844 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
847 GMX_UNUSED_VALUE(tng);
848 GMX_UNUSED_VALUE(bUseLossyCompression);
849 GMX_UNUSED_VALUE(step);
850 GMX_UNUSED_VALUE(elapsedPicoSeconds);
851 GMX_UNUSED_VALUE(lambda);
852 GMX_UNUSED_VALUE(box);
853 GMX_UNUSED_VALUE(nAtoms);
860 void fflush_tng(tng_trajectory_t tng)
867 tng_frame_set_premature_write(tng, TNG_USE_HASH);
869 GMX_UNUSED_VALUE(tng);
873 float gmx_tng_get_time_of_final_frame(tng_trajectory_t tng)
880 tng_num_frames_get(tng, &nFrames);
881 tng_util_time_of_frame_get(tng, nFrames - 1, &time);
886 GMX_UNUSED_VALUE(tng);