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)
76 gmx_fatal(FARGS, "Invalid file opening mode %c", mode);
81 void gmx_tng_open(const char *filename,
83 tng_trajectory_t *tng)
86 /* First check whether we have to make a backup,
87 * only for writing, not for read or append.
92 /* only make backups for normal gromacs */
93 make_backup(filename);
97 /* tng must not be pointing at already allocated memory.
98 * Memory will be allocated by tng_util_trajectory_open() and must
99 * later on be freed by tng_util_trajectory_close(). */
100 if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
102 /* TNG does return more than one degree of error, but there is
103 no use case for GROMACS handling the non-fatal errors
106 "%s while opening %s for %s",
107 gmx_strerror("file"),
112 if (mode == 'w' || mode == 'a')
114 /* FIXME in TNG: When adding data to the header, subsequent blocks might get
115 * overwritten. This could be solved by moving the first trajectory
116 * frame set(s) to the end of the file. Could that cause other problems,
117 * e.g. when continuing a simulation? */
119 gmx_gethostname(hostname, 256);
122 tng_first_computer_name_set(*tng, hostname);
124 /* TODO: This should be implemented when the above fixme is done (adding data to
128 // tng_last_computer_name_set(*tng, hostname);
131 char programInfo[256];
132 const char *precisionString = "";
134 precisionString = " (double precision)";
136 sprintf(programInfo, "%.100s, %.128s%.24s",
137 gmx::getProgramContext().displayName(),
138 GromacsVersion(), precisionString);
141 tng_first_program_name_set(*tng, programInfo);
143 /* TODO: This should be implemented when the above fixme is done (adding data to
147 // tng_last_program_name_set(*tng, programInfo);
152 getlogin_r(username, 256);
155 tng_first_user_name_set(*tng, username);
157 /* TODO: This should be implemented when the above fixme is done (adding data to
161 // tng_last_user_name_set(*tng, username);
166 gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
167 GMX_UNUSED_VALUE(filename);
168 GMX_UNUSED_VALUE(mode);
169 GMX_UNUSED_VALUE(tng);
173 void gmx_tng_close(tng_trajectory_t *tng)
175 /* We have to check that tng is set because
176 * tng_util_trajectory_close wants to return a NULL in it, and
177 * gives a fatal error if it is NULL. */
181 tng_util_trajectory_close(tng);
184 GMX_UNUSED_VALUE(tng);
189 static void addTngMoleculeFromTopology(tng_trajectory_t tng,
190 const char *moleculeName,
191 const t_atoms *atoms,
192 gmx_int64_t numMolecules,
193 tng_molecule_t *tngMol)
195 tng_chain_t tngChain = NULL;
196 tng_residue_t tngRes = NULL;
198 if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
200 gmx_file("Cannot add molecule to TNG molecular system.");
203 /* FIXME: The TNG atoms should contain mass and atomB info (for free
204 * energy calculations), i.e. in when it's available in TNG (2.0). */
205 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++)
207 const t_atom *at = &atoms->atom[atomIndex];
208 /* FIXME: Currently the TNG API can only add atoms belonging to a
209 * residue and chain. Wait for TNG 2.0*/
212 const t_resinfo *resInfo = &atoms->resinfo[at->resind];
213 char chainName[2] = {resInfo->chainid, 0};
214 tng_atom_t tngAtom = NULL;
219 prevAtom = &atoms->atom[atomIndex - 1];
226 /* If this is the first atom or if the residue changed add the
227 * residue to the TNG molecular system. */
228 if (!prevAtom || resInfo != &atoms->resinfo[prevAtom->resind])
230 /* If this is the first atom or if the chain changed add
231 * the chain to the TNG molecular system. */
232 if (!prevAtom || resInfo->chainid !=
233 atoms->resinfo[prevAtom->resind].chainid)
235 tng_molecule_chain_add(tng, *tngMol, chainName,
238 /* FIXME: When TNG supports both residue index and residue
239 * number the latter should be used. Wait for TNG 2.0*/
240 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
242 tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIndex]), *(atoms->atomtype[atomIndex]), &tngAtom);
245 tng_molecule_cnt_set(tng, *tngMol, numMolecules);
248 void gmx_tng_add_mtop(tng_trajectory_t tng,
249 const gmx_mtop_t *mtop)
252 const t_ilist *ilist;
257 /* No topology information available to add. */
261 for (int molIndex = 0; molIndex < mtop->nmolblock; molIndex++)
263 tng_molecule_t tngMol = NULL;
264 const gmx_moltype_t *molType =
265 &mtop->moltype[mtop->molblock[molIndex].type];
267 /* Add a molecule to the TNG trajectory with the same name as the
268 * current molecule. */
269 addTngMoleculeFromTopology(tng,
272 mtop->molblock[molIndex].nmol,
275 /* Bonds have to be deduced from interactions (constraints etc). Different
276 * interactions have different sets of parameters. */
277 /* Constraints are specified using two atoms */
278 for (i = 0; i < F_NRE; i++)
282 ilist = &molType->ilist[i];
286 while (j < ilist->nr)
288 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
294 /* Settle is described using three atoms */
295 ilist = &molType->ilist[F_SETTLE];
299 while (j < ilist->nr)
301 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
302 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+2], &tngBond);
309 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
310 * if they are positive.
312 * If only one of n1 and n2 is positive, then return it.
313 * If neither n1 or n2 is positive, then return -1. */
315 greatest_common_divisor_if_positive(int n1, int n2)
319 return (0 >= n2) ? -1 : n2;
326 /* We have a non-trivial greatest common divisor to compute. */
327 return gmx_greatest_common_divisor(n1, n2);
330 /* By default try to write 100 frames (of actual output) in each frame set.
331 * This number is the number of outputs of the most frequently written data
332 * type per frame set.
333 * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
334 * setups regarding compression efficiency and compression time. Make this
335 * a hidden command-line option? */
336 const int defaultFramesPerFrameSet = 100;
338 /*! \libinternal \brief Set the number of frames per frame
339 * set according to output intervals.
340 * The default is that 100 frames are written of the data
341 * that is written most often. */
342 static void tng_set_frames_per_frame_set(tng_trajectory_t tng,
343 const gmx_bool bUseLossyCompression,
344 const t_inputrec *ir)
348 /* Set the number of frames per frame set to contain at least
349 * defaultFramesPerFrameSet of the lowest common denominator of
350 * the writing interval of positions and velocities. */
351 /* FIXME after 5.0: consider nstenergy also? */
352 if (bUseLossyCompression)
354 gcd = ir->nstxout_compressed;
358 gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
359 gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
366 tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
369 /*! \libinternal \brief Set the data-writing intervals, and number of
370 * frames per frame set */
371 static void set_writing_intervals(tng_trajectory_t tng,
372 const gmx_bool bUseLossyCompression,
373 const t_inputrec *ir)
375 /* Define pointers to specific writing functions depending on if we
376 * write float or double data */
377 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
385 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
387 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
389 int xout, vout, fout;
390 // int gcd = -1, lowest = -1;
393 tng_set_frames_per_frame_set(tng, bUseLossyCompression, ir);
395 if (bUseLossyCompression)
397 xout = ir->nstxout_compressed;
400 compression = TNG_TNG_COMPRESSION;
407 compression = TNG_GZIP_COMPRESSION;
411 set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS,
412 "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
414 /* The design of TNG makes it awkward to try to write a box
415 * with multiple periodicities, which might be co-prime. Since
416 * the use cases for the box with a frame consisting only of
417 * velocities seem low, for now we associate box writing with
418 * position writing. */
419 set_writing_interval(tng, xout, 9, TNG_TRAJ_BOX_SHAPE,
420 "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA,
421 TNG_GZIP_COMPRESSION);
422 /* TODO: if/when we write energies to TNG also, reconsider how
423 * and when box information is written, because GROMACS
424 * behaviour pre-5.0 was to write the box with every
425 * trajectory frame and every energy frame, and probably
426 * people depend on this. */
428 /* TODO: If we need to write lambda values at steps when
429 * positions (or other data) are not also being written, then
430 * code in mdoutf.c will need to match however that is
432 set_writing_interval(tng, xout, 1, TNG_GMX_LAMBDA,
433 "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
434 TNG_GZIP_COMPRESSION);
436 /* FIXME: gcd and lowest currently not used. */
437 // gcd = greatest_common_divisor_if_positive(gcd, xout);
438 // if (lowest < 0 || xout < lowest)
445 set_writing_interval(tng, ir->nstvout, 3, TNG_TRAJ_VELOCITIES,
446 "VELOCITIES", TNG_PARTICLE_BLOCK_DATA,
449 /* FIXME: gcd and lowest currently not used. */
450 // gcd = greatest_common_divisor_if_positive(gcd, vout);
451 // if (lowest < 0 || vout < lowest)
458 set_writing_interval(tng, ir->nstfout, 3, TNG_TRAJ_FORCES,
459 "FORCES", TNG_PARTICLE_BLOCK_DATA,
460 TNG_GZIP_COMPRESSION);
462 /* FIXME: gcd and lowest currently not used. */
463 // gcd = greatest_common_divisor_if_positive(gcd, fout);
464 // if (lowest < 0 || fout < lowest)
469 /* FIXME: See above. gcd interval for lambdas is disabled. */
472 // /* Lambdas written at an interval of the lowest common denominator
473 // * of other output */
474 // set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA,
475 // "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
476 // TNG_GZIP_COMPRESSION);
478 // if (gcd < lowest / 10)
480 // gmx_warning("The lowest common denominator of trajectory output is "
481 // "every %d step(s), whereas the shortest output interval "
482 // "is every %d steps.", gcd, lowest);
488 void gmx_tng_prepare_md_writing(tng_trajectory_t tng,
489 const gmx_mtop_t *mtop,
490 const t_inputrec *ir)
493 gmx_tng_add_mtop(tng, mtop);
494 set_writing_intervals(tng, FALSE, ir);
495 tng_time_per_frame_set(tng, ir->delta_t * PICO);
497 GMX_UNUSED_VALUE(tng);
498 GMX_UNUSED_VALUE(mtop);
499 GMX_UNUSED_VALUE(ir);
504 /* Check if all atoms in the molecule system are selected
505 * by a selection group of type specified by gtype. */
506 static gmx_bool all_atoms_selected(const gmx_mtop_t *mtop,
509 const gmx_moltype_t *molType;
510 const t_atoms *atoms;
512 /* Iterate through all atoms in the molecule system and
513 * check if they belong to a selection group of the
515 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
517 molType = &mtop->moltype[mtop->molblock[molIndex].type];
519 atoms = &molType->atoms;
521 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
523 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
525 if (ggrpnr(&mtop->groups, gtype, i) != 0)
535 /* Create TNG molecules which will represent each of the selection
536 * groups for writing. But do that only if there is actually a
537 * specified selection group and it is not the whole system.
538 * TODO: Currently the only selection that is taken into account
539 * is egcCompressedX, but other selections should be added when
540 * e.g. writing energies is implemented.
542 static void add_selection_groups(tng_trajectory_t tng,
543 const gmx_mtop_t *mtop)
545 const gmx_moltype_t *molType;
546 const t_atoms *atoms;
548 const t_resinfo *resInfo;
549 const t_ilist *ilist;
552 tng_molecule_t mol, iterMol;
560 /* TODO: When the TNG molecules block is more flexible TNG selection
561 * groups should not need all atoms specified. It should be possible
562 * just to specify what molecules are selected (and/or which atoms in
563 * the molecule) and how many (if applicable). */
565 /* If no atoms are selected we do not need to create a
566 * TNG selection group molecule. */
567 if (mtop->groups.ngrpnr[egcCompressedX] == 0)
572 /* If all atoms are selected we do not have to create a selection
573 * group molecule in the TNG molecule system. */
574 if (all_atoms_selected(mtop, egcCompressedX))
579 /* The name of the TNG molecule containing the selection group is the
580 * same as the name of the selection group. */
581 nameIndex = *mtop->groups.grps[egcCompressedX].nm_ind;
582 groupName = *mtop->groups.grpname[nameIndex];
584 tng_molecule_alloc(tng, &mol);
585 tng_molecule_name_set(tng, mol, groupName);
586 tng_molecule_chain_add(tng, mol, "", &chain);
587 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
589 molType = &mtop->moltype[mtop->molblock[molIndex].type];
591 atoms = &molType->atoms;
593 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
595 bool bAtomsAdded = FALSE;
596 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
601 if (ggrpnr(&mtop->groups, egcCompressedX, i) != 0)
605 at = &atoms->atom[atomIndex];
608 resInfo = &atoms->resinfo[at->resind];
609 /* FIXME: When TNG supports both residue index and residue
610 * number the latter should be used. */
611 res_name = *resInfo->name;
612 res_id = at->resind + 1;
616 res_name = (char *)"";
619 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res)
622 /* Since there is ONE chain for selection groups do not keep the
623 * original residue IDs - otherwise there might be conflicts. */
624 tng_chain_residue_add(tng, chain, res_name, &res);
626 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
627 *(atoms->atomtype[atomIndex]),
628 atom_offset + atomIndex, &atom);
634 for (int k = 0; k < F_NRE; k++)
638 ilist = &molType->ilist[k];
642 while (l < ilist->nr)
645 atom1 = ilist->iatoms[l] + atom_offset;
646 atom2 = ilist->iatoms[l+1] + atom_offset;
647 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0 &&
648 ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
650 tng_molecule_bond_add(tng, mol, ilist->iatoms[l],
651 ilist->iatoms[l+1], &tngBond);
658 /* Settle is described using three atoms */
659 ilist = &molType->ilist[F_SETTLE];
663 while (l < ilist->nr)
665 int atom1, atom2, atom3;
666 atom1 = ilist->iatoms[l] + atom_offset;
667 atom2 = ilist->iatoms[l+1] + atom_offset;
668 atom3 = ilist->iatoms[l+2] + atom_offset;
669 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0)
671 if (ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
673 tng_molecule_bond_add(tng, mol, atom1,
676 if (ggrpnr(&mtop->groups, egcCompressedX, atom3) == 0)
678 tng_molecule_bond_add(tng, mol, atom1,
686 atom_offset += atoms->nr;
689 tng_molecule_existing_add(tng, &mol);
690 tng_molecule_cnt_set(tng, mol, 1);
691 tng_num_molecule_types_get(tng, &nMols);
692 for (gmx_int64_t k = 0; k < nMols; k++)
694 tng_molecule_of_index_get(tng, k, &iterMol);
699 tng_molecule_cnt_set(tng, iterMol, 0);
704 void gmx_tng_set_compression_precision(tng_trajectory_t tng,
708 tng_compression_precision_set(tng, prec);
710 GMX_UNUSED_VALUE(tng);
711 GMX_UNUSED_VALUE(prec);
715 void gmx_tng_prepare_low_prec_writing(tng_trajectory_t tng,
716 const gmx_mtop_t *mtop,
717 const t_inputrec *ir)
720 gmx_tng_add_mtop(tng, mtop);
721 add_selection_groups(tng, mtop);
722 set_writing_intervals(tng, TRUE, ir);
723 tng_time_per_frame_set(tng, ir->delta_t * PICO);
724 gmx_tng_set_compression_precision(tng, ir->x_compression_precision);
726 GMX_UNUSED_VALUE(tng);
727 GMX_UNUSED_VALUE(mtop);
728 GMX_UNUSED_VALUE(ir);
732 void gmx_fwrite_tng(tng_trajectory_t tng,
733 const gmx_bool bUseLossyCompression,
735 real elapsedPicoSeconds,
744 typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
754 static write_data_func_pointer write_data = tng_util_generic_with_time_double_write;
756 static write_data_func_pointer write_data = tng_util_generic_with_time_write;
758 double elapsedSeconds = elapsedPicoSeconds * PICO;
759 gmx_int64_t nParticles;
765 /* This function might get called when the type of the
766 compressed trajectory is actually XTC. So we exit and move
771 tng_num_particles_get(tng, &nParticles);
772 if (nAtoms != (int)nParticles)
774 tng_implicit_num_particles_set(tng, nAtoms);
777 if (bUseLossyCompression)
779 compression = TNG_TNG_COMPRESSION;
783 compression = TNG_GZIP_COMPRESSION;
786 /* The writing is done using write_data, which writes float or double
787 * depending on the GROMACS compilation. */
790 GMX_ASSERT(box, "Need a non-NULL box if positions are written");
792 if (write_data(tng, step, elapsedSeconds,
793 reinterpret_cast<const real *>(x),
794 3, TNG_TRAJ_POSITIONS, "POSITIONS",
795 TNG_PARTICLE_BLOCK_DATA,
796 compression) != TNG_SUCCESS)
798 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
800 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
801 * compression for box shape regardless of output mode */
802 if (write_data(tng, step, elapsedSeconds,
803 reinterpret_cast<const real *>(box),
804 9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
805 TNG_NON_PARTICLE_BLOCK_DATA,
806 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
808 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
814 if (write_data(tng, step, elapsedSeconds,
815 reinterpret_cast<const real *>(v),
816 3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
817 TNG_PARTICLE_BLOCK_DATA,
818 compression) != TNG_SUCCESS)
820 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
826 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
827 * compression for forces regardless of output mode */
828 if (write_data(tng, step, elapsedSeconds,
829 reinterpret_cast<const real *>(f),
830 3, TNG_TRAJ_FORCES, "FORCES",
831 TNG_PARTICLE_BLOCK_DATA,
832 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
834 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
838 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
839 * compression for lambdas regardless of output mode */
840 if (write_data(tng, step, elapsedSeconds,
841 reinterpret_cast<const real *>(&lambda),
842 1, TNG_GMX_LAMBDA, "LAMBDAS",
843 TNG_NON_PARTICLE_BLOCK_DATA,
844 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
846 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
849 GMX_UNUSED_VALUE(tng);
850 GMX_UNUSED_VALUE(bUseLossyCompression);
851 GMX_UNUSED_VALUE(step);
852 GMX_UNUSED_VALUE(elapsedPicoSeconds);
853 GMX_UNUSED_VALUE(lambda);
854 GMX_UNUSED_VALUE(box);
855 GMX_UNUSED_VALUE(nAtoms);
862 void fflush_tng(tng_trajectory_t tng)
869 tng_frame_set_premature_write(tng, TNG_USE_HASH);
871 GMX_UNUSED_VALUE(tng);
875 float gmx_tng_get_time_of_final_frame(tng_trajectory_t tng)
882 tng_num_frames_get(tng, &nFrames);
883 tng_util_time_of_frame_get(tng, nFrames - 1, &time);
888 GMX_UNUSED_VALUE(tng);