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.
49 #include "tng/tng_io.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/mdtypes/inputrec.h"
55 #include "gromacs/topology/ifunc.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/trajectory/trajectoryframe.h"
58 #include "gromacs/utility/basedefinitions.h"
59 #include "gromacs/utility/baseversion.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/gmxassert.h"
63 #include "gromacs/utility/programcontext.h"
64 #include "gromacs/utility/smalloc.h"
65 #include "gromacs/utility/sysinfo.h"
66 #include "gromacs/utility/unique_cptr.h"
68 static const char *modeToVerb(char mode)
83 gmx_fatal(FARGS, "Invalid file opening mode %c", mode);
90 void gmx_tng_open(const char *filename,
92 tng_trajectory_t *tng)
95 /* First check whether we have to make a backup,
96 * only for writing, not for read or append.
100 make_backup(filename);
103 /* tng must not be pointing at already allocated memory.
104 * Memory will be allocated by tng_util_trajectory_open() and must
105 * later on be freed by tng_util_trajectory_close(). */
106 if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
108 /* TNG does return more than one degree of error, but there is
109 no use case for GROMACS handling the non-fatal errors
112 "File I/O error while opening %s for %s",
117 if (mode == 'w' || mode == 'a')
120 gmx_gethostname(hostname, 256);
123 tng_first_computer_name_set(*tng, hostname);
127 tng_last_computer_name_set(*tng, hostname);
130 char programInfo[256];
131 const char *precisionString = "";
133 precisionString = " (double precision)";
135 sprintf(programInfo, "%.100s %.128s%.24s",
136 gmx::getProgramContext().displayName(),
137 gmx_version(), precisionString);
140 tng_first_program_name_set(*tng, programInfo);
144 tng_last_program_name_set(*tng, programInfo);
148 if (!gmx_getusername(username, 256))
152 tng_first_user_name_set(*tng, username);
156 tng_last_user_name_set(*tng, username);
157 tng_file_headers_write(*tng, TNG_USE_HASH);
162 gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
163 GMX_UNUSED_VALUE(filename);
164 GMX_UNUSED_VALUE(mode);
165 GMX_UNUSED_VALUE(tng);
169 void gmx_tng_close(tng_trajectory_t *tng)
171 /* We have to check that tng is set because
172 * tng_util_trajectory_close wants to return a NULL in it, and
173 * gives a fatal error if it is NULL. */
177 tng_util_trajectory_close(tng);
180 GMX_UNUSED_VALUE(tng);
185 static void addTngMoleculeFromTopology(tng_trajectory_t tng,
186 const char *moleculeName,
187 const t_atoms *atoms,
188 gmx_int64_t numMolecules,
189 tng_molecule_t *tngMol)
191 tng_chain_t tngChain = nullptr;
192 tng_residue_t tngRes = nullptr;
194 if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
196 gmx_file("Cannot add molecule to TNG molecular system.");
199 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++)
201 const t_atom *at = &atoms->atom[atomIndex];
202 /* FIXME: Currently the TNG API can only add atoms belonging to a
203 * residue and chain. Wait for TNG 2.0*/
206 const t_resinfo *resInfo = &atoms->resinfo[at->resind];
207 char chainName[2] = {resInfo->chainid, 0};
208 tng_atom_t tngAtom = nullptr;
213 prevAtom = &atoms->atom[atomIndex - 1];
220 /* If this is the first atom or if the residue changed add the
221 * residue to the TNG molecular system. */
222 if (!prevAtom || resInfo != &atoms->resinfo[prevAtom->resind])
224 /* If this is the first atom or if the chain changed add
225 * the chain to the TNG molecular system. */
226 if (!prevAtom || resInfo->chainid !=
227 atoms->resinfo[prevAtom->resind].chainid)
229 tng_molecule_chain_add(tng, *tngMol, chainName,
232 /* FIXME: When TNG supports both residue index and residue
233 * number the latter should be used. Wait for TNG 2.0*/
234 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
236 tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIndex]), *(atoms->atomtype[atomIndex]), &tngAtom);
239 tng_molecule_cnt_set(tng, *tngMol, numMolecules);
242 void gmx_tng_add_mtop(tng_trajectory_t tng,
243 const gmx_mtop_t *mtop)
247 std::vector<real> atomCharges;
248 std::vector<real> atomMasses;
249 const t_ilist *ilist;
255 /* No topology information available to add. */
260 datatype = TNG_DOUBLE_DATA;
262 datatype = TNG_FLOAT_DATA;
265 atomCharges.reserve(mtop->natoms);
266 atomMasses.reserve(mtop->natoms);
268 for (int molIndex = 0; molIndex < mtop->nmolblock; molIndex++)
270 tng_molecule_t tngMol = nullptr;
271 const gmx_moltype_t *molType = &mtop->moltype[mtop->molblock[molIndex].type];
273 /* Add a molecule to the TNG trajectory with the same name as the
274 * current molecule. */
275 addTngMoleculeFromTopology(tng,
278 mtop->molblock[molIndex].nmol,
281 /* Bonds have to be deduced from interactions (constraints etc). Different
282 * interactions have different sets of parameters. */
283 /* Constraints are specified using two atoms */
284 for (i = 0; i < F_NRE; i++)
288 ilist = &molType->ilist[i];
292 while (j < ilist->nr)
294 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
300 /* Settle is described using three atoms */
301 ilist = &molType->ilist[F_SETTLE];
305 while (j < ilist->nr)
307 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
308 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+2], &tngBond);
312 /* First copy atom charges and masses, first atom by atom and then copy the memory for the molecule instances.
313 * FIXME: Atom B state data should also be written to TNG (v 2.0?) */
314 for (int atomCounter = 0; atomCounter < molType->atoms.nr; atomCounter++)
316 atomCharges.push_back(molType->atoms.atom[atomCounter].q);
317 atomMasses.push_back(molType->atoms.atom[atomCounter].m);
319 for (int molCounter = 1; molCounter < mtop->molblock[molIndex].nmol; molCounter++)
321 std::copy_n(atomCharges.end() - molType->atoms.nr, molType->atoms.nr, std::back_inserter(atomCharges));
322 std::copy_n(atomMasses.end() - molType->atoms.nr, molType->atoms.nr, std::back_inserter(atomMasses));
325 /* Write the TNG data blocks. */
326 tng_particle_data_block_add(tng, TNG_TRAJ_PARTIAL_CHARGES, "PARTIAL CHARGES",
327 datatype, TNG_NON_TRAJECTORY_BLOCK,
328 1, 1, 1, 0, mtop->natoms,
329 TNG_GZIP_COMPRESSION, atomCharges.data());
330 tng_particle_data_block_add(tng, TNG_TRAJ_MASSES, "ATOM MASSES",
331 datatype, TNG_NON_TRAJECTORY_BLOCK,
332 1, 1, 1, 0, mtop->natoms,
333 TNG_GZIP_COMPRESSION, atomMasses.data());
336 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
337 * if they are positive.
339 * If only one of n1 and n2 is positive, then return it.
340 * If neither n1 or n2 is positive, then return -1. */
342 greatest_common_divisor_if_positive(int n1, int n2)
346 return (0 >= n2) ? -1 : n2;
353 /* We have a non-trivial greatest common divisor to compute. */
354 return gmx_greatest_common_divisor(n1, n2);
357 /* By default try to write 100 frames (of actual output) in each frame set.
358 * This number is the number of outputs of the most frequently written data
359 * type per frame set.
360 * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
361 * setups regarding compression efficiency and compression time. Make this
362 * a hidden command-line option? */
363 const int defaultFramesPerFrameSet = 100;
365 /*! \libinternal \brief Set the number of frames per frame
366 * set according to output intervals.
367 * The default is that 100 frames are written of the data
368 * that is written most often. */
369 static void tng_set_frames_per_frame_set(tng_trajectory_t tng,
370 const gmx_bool bUseLossyCompression,
371 const t_inputrec *ir)
375 /* Set the number of frames per frame set to contain at least
376 * defaultFramesPerFrameSet of the lowest common denominator of
377 * the writing interval of positions and velocities. */
378 /* FIXME after 5.0: consider nstenergy also? */
379 if (bUseLossyCompression)
381 gcd = ir->nstxout_compressed;
385 gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
386 gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
393 tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
396 /*! \libinternal \brief Set the data-writing intervals, and number of
397 * frames per frame set */
398 static void set_writing_intervals(tng_trajectory_t tng,
399 const gmx_bool bUseLossyCompression,
400 const t_inputrec *ir)
402 /* Define pointers to specific writing functions depending on if we
403 * write float or double data */
404 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
412 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
414 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
416 int xout, vout, fout;
417 int gcd = -1, lowest = -1;
420 tng_set_frames_per_frame_set(tng, bUseLossyCompression, ir);
422 if (bUseLossyCompression)
424 xout = ir->nstxout_compressed;
426 /* If there is no uncompressed coordinate output write forces
427 and velocities to the compressed tng file. */
438 compression = TNG_TNG_COMPRESSION;
445 compression = TNG_GZIP_COMPRESSION;
449 set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS,
450 "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
452 /* TODO: if/when we write energies to TNG also, reconsider how
453 * and when box information is written, because GROMACS
454 * behaviour pre-5.0 was to write the box with every
455 * trajectory frame and every energy frame, and probably
456 * people depend on this. */
458 gcd = greatest_common_divisor_if_positive(gcd, xout);
459 if (lowest < 0 || xout < lowest)
466 set_writing_interval(tng, vout, 3, TNG_TRAJ_VELOCITIES,
467 "VELOCITIES", TNG_PARTICLE_BLOCK_DATA,
470 gcd = greatest_common_divisor_if_positive(gcd, vout);
471 if (lowest < 0 || vout < lowest)
478 set_writing_interval(tng, fout, 3, TNG_TRAJ_FORCES,
479 "FORCES", TNG_PARTICLE_BLOCK_DATA,
480 TNG_GZIP_COMPRESSION);
482 gcd = greatest_common_divisor_if_positive(gcd, fout);
483 if (lowest < 0 || fout < lowest)
490 /* Lambdas and box shape written at an interval of the lowest common
491 denominator of other output */
492 set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA,
493 "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
494 TNG_GZIP_COMPRESSION);
496 set_writing_interval(tng, gcd, 9, TNG_TRAJ_BOX_SHAPE,
497 "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA,
498 TNG_GZIP_COMPRESSION);
499 if (gcd < lowest / 10)
501 gmx_warning("The lowest common denominator of trajectory output is "
502 "every %d step(s), whereas the shortest output interval "
503 "is every %d steps.", gcd, lowest);
509 void gmx_tng_prepare_md_writing(tng_trajectory_t tng,
510 const gmx_mtop_t *mtop,
511 const t_inputrec *ir)
514 gmx_tng_add_mtop(tng, mtop);
515 set_writing_intervals(tng, FALSE, ir);
516 tng_time_per_frame_set(tng, ir->delta_t * PICO);
518 GMX_UNUSED_VALUE(tng);
519 GMX_UNUSED_VALUE(mtop);
520 GMX_UNUSED_VALUE(ir);
525 /* Check if all atoms in the molecule system are selected
526 * by a selection group of type specified by gtype. */
527 static gmx_bool all_atoms_selected(const gmx_mtop_t *mtop,
530 const gmx_moltype_t *molType;
531 const t_atoms *atoms;
533 /* Iterate through all atoms in the molecule system and
534 * check if they belong to a selection group of the
536 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
538 molType = &mtop->moltype[mtop->molblock[molIndex].type];
540 atoms = &molType->atoms;
542 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
544 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
546 if (ggrpnr(&mtop->groups, gtype, i) != 0)
556 /* Create TNG molecules which will represent each of the selection
557 * groups for writing. But do that only if there is actually a
558 * specified selection group and it is not the whole system.
559 * TODO: Currently the only selection that is taken into account
560 * is egcCompressedX, but other selections should be added when
561 * e.g. writing energies is implemented.
563 static void add_selection_groups(tng_trajectory_t tng,
564 const gmx_mtop_t *mtop)
566 const gmx_moltype_t *molType;
567 const t_atoms *atoms;
569 const t_resinfo *resInfo;
570 const t_ilist *ilist;
573 tng_molecule_t mol, iterMol;
581 /* TODO: When the TNG molecules block is more flexible TNG selection
582 * groups should not need all atoms specified. It should be possible
583 * just to specify what molecules are selected (and/or which atoms in
584 * the molecule) and how many (if applicable). */
586 /* If no atoms are selected we do not need to create a
587 * TNG selection group molecule. */
588 if (mtop->groups.ngrpnr[egcCompressedX] == 0)
593 /* If all atoms are selected we do not have to create a selection
594 * group molecule in the TNG molecule system. */
595 if (all_atoms_selected(mtop, egcCompressedX))
600 /* The name of the TNG molecule containing the selection group is the
601 * same as the name of the selection group. */
602 nameIndex = *mtop->groups.grps[egcCompressedX].nm_ind;
603 groupName = *mtop->groups.grpname[nameIndex];
605 tng_molecule_alloc(tng, &mol);
606 tng_molecule_name_set(tng, mol, groupName);
607 tng_molecule_chain_add(tng, mol, "", &chain);
608 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
610 molType = &mtop->moltype[mtop->molblock[molIndex].type];
612 atoms = &molType->atoms;
614 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
616 bool bAtomsAdded = FALSE;
617 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
622 if (ggrpnr(&mtop->groups, egcCompressedX, i) != 0)
626 at = &atoms->atom[atomIndex];
629 resInfo = &atoms->resinfo[at->resind];
630 /* FIXME: When TNG supports both residue index and residue
631 * number the latter should be used. */
632 res_name = *resInfo->name;
633 res_id = at->resind + 1;
637 res_name = (char *)"";
640 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res)
643 /* Since there is ONE chain for selection groups do not keep the
644 * original residue IDs - otherwise there might be conflicts. */
645 tng_chain_residue_add(tng, chain, res_name, &res);
647 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
648 *(atoms->atomtype[atomIndex]),
649 atom_offset + atomIndex, &atom);
655 for (int k = 0; k < F_NRE; k++)
659 ilist = &molType->ilist[k];
663 while (l < ilist->nr)
666 atom1 = ilist->iatoms[l] + atom_offset;
667 atom2 = ilist->iatoms[l+1] + atom_offset;
668 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0 &&
669 ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
671 tng_molecule_bond_add(tng, mol, ilist->iatoms[l],
672 ilist->iatoms[l+1], &tngBond);
679 /* Settle is described using three atoms */
680 ilist = &molType->ilist[F_SETTLE];
684 while (l < ilist->nr)
686 int atom1, atom2, atom3;
687 atom1 = ilist->iatoms[l] + atom_offset;
688 atom2 = ilist->iatoms[l+1] + atom_offset;
689 atom3 = ilist->iatoms[l+2] + atom_offset;
690 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0)
692 if (ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
694 tng_molecule_bond_add(tng, mol, atom1,
697 if (ggrpnr(&mtop->groups, egcCompressedX, atom3) == 0)
699 tng_molecule_bond_add(tng, mol, atom1,
707 atom_offset += atoms->nr;
710 tng_molecule_existing_add(tng, &mol);
711 tng_molecule_cnt_set(tng, mol, 1);
712 tng_num_molecule_types_get(tng, &nMols);
713 for (gmx_int64_t k = 0; k < nMols; k++)
715 tng_molecule_of_index_get(tng, k, &iterMol);
720 tng_molecule_cnt_set(tng, iterMol, 0);
725 void gmx_tng_set_compression_precision(tng_trajectory_t tng,
729 tng_compression_precision_set(tng, prec);
731 GMX_UNUSED_VALUE(tng);
732 GMX_UNUSED_VALUE(prec);
736 void gmx_tng_prepare_low_prec_writing(tng_trajectory_t tng,
737 const gmx_mtop_t *mtop,
738 const t_inputrec *ir)
741 gmx_tng_add_mtop(tng, mtop);
742 add_selection_groups(tng, mtop);
743 set_writing_intervals(tng, TRUE, ir);
744 tng_time_per_frame_set(tng, ir->delta_t * PICO);
745 gmx_tng_set_compression_precision(tng, ir->x_compression_precision);
747 GMX_UNUSED_VALUE(tng);
748 GMX_UNUSED_VALUE(mtop);
749 GMX_UNUSED_VALUE(ir);
753 void gmx_fwrite_tng(tng_trajectory_t tng,
754 const gmx_bool bUseLossyCompression,
756 real elapsedPicoSeconds,
765 typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
775 static write_data_func_pointer write_data = tng_util_generic_with_time_double_write;
777 static write_data_func_pointer write_data = tng_util_generic_with_time_write;
779 double elapsedSeconds = elapsedPicoSeconds * PICO;
780 gmx_int64_t nParticles;
786 /* This function might get called when the type of the
787 compressed trajectory is actually XTC. So we exit and move
792 tng_num_particles_get(tng, &nParticles);
793 if (nAtoms != (int)nParticles)
795 tng_implicit_num_particles_set(tng, nAtoms);
798 if (bUseLossyCompression)
800 compression = TNG_TNG_COMPRESSION;
804 compression = TNG_GZIP_COMPRESSION;
807 /* The writing is done using write_data, which writes float or double
808 * depending on the GROMACS compilation. */
811 GMX_ASSERT(box, "Need a non-NULL box if positions are written");
813 if (write_data(tng, step, elapsedSeconds,
814 reinterpret_cast<const real *>(x),
815 3, TNG_TRAJ_POSITIONS, "POSITIONS",
816 TNG_PARTICLE_BLOCK_DATA,
817 compression) != TNG_SUCCESS)
819 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
825 if (write_data(tng, step, elapsedSeconds,
826 reinterpret_cast<const real *>(v),
827 3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
828 TNG_PARTICLE_BLOCK_DATA,
829 compression) != TNG_SUCCESS)
831 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
837 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
838 * compression for forces regardless of output mode */
839 if (write_data(tng, step, elapsedSeconds,
840 reinterpret_cast<const real *>(f),
841 3, TNG_TRAJ_FORCES, "FORCES",
842 TNG_PARTICLE_BLOCK_DATA,
843 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
845 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
849 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
850 * compression for lambdas and box shape regardless of output mode */
851 if (write_data(tng, step, elapsedSeconds,
852 reinterpret_cast<const real *>(box),
853 9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
854 TNG_NON_PARTICLE_BLOCK_DATA,
855 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
857 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
860 if (write_data(tng, step, elapsedSeconds,
861 reinterpret_cast<const real *>(&lambda),
862 1, TNG_GMX_LAMBDA, "LAMBDAS",
863 TNG_NON_PARTICLE_BLOCK_DATA,
864 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
866 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
869 GMX_UNUSED_VALUE(tng);
870 GMX_UNUSED_VALUE(bUseLossyCompression);
871 GMX_UNUSED_VALUE(step);
872 GMX_UNUSED_VALUE(elapsedPicoSeconds);
873 GMX_UNUSED_VALUE(lambda);
874 GMX_UNUSED_VALUE(box);
875 GMX_UNUSED_VALUE(nAtoms);
882 void fflush_tng(tng_trajectory_t tng)
889 tng_frame_set_premature_write(tng, TNG_USE_HASH);
891 GMX_UNUSED_VALUE(tng);
895 float gmx_tng_get_time_of_final_frame(tng_trajectory_t tng)
902 tng_num_frames_get(tng, &nFrames);
903 tng_util_time_of_frame_get(tng, nFrames - 1, &time);
908 GMX_UNUSED_VALUE(tng);
913 void gmx_prepare_tng_writing(const char *filename,
915 tng_trajectory_t *input,
916 tng_trajectory_t *output,
918 const gmx_mtop_t *mtop,
920 const char *indexGroupName)
923 /* FIXME after 5.0: Currently only standard block types are read */
924 const int defaultNumIds = 5;
925 static gmx_int64_t fallbackIds[defaultNumIds] =
927 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
928 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
931 static char fallbackNames[defaultNumIds][32] =
933 "BOX SHAPE", "POSITIONS", "VELOCITIES",
937 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
945 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
947 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
950 gmx_tng_open(filename, mode, output);
952 /* Do we have an input file in TNG format? If so, then there's
953 more data we can copy over, rather than having to improvise. */
956 /* Set parameters (compression, time per frame, molecule
957 * information, number of frames per frame set and writing
958 * intervals of positions, box shape and lambdas) of the
959 * output tng container based on their respective values int
960 * the input tng container */
961 double time, compression_precision;
962 gmx_int64_t n_frames_per_frame_set, interval = -1;
964 tng_compression_precision_get(*input, &compression_precision);
965 tng_compression_precision_set(*output, compression_precision);
966 // TODO make this configurable in a future version
967 char compression_type = TNG_TNG_COMPRESSION;
969 tng_molecule_system_copy(*input, *output);
971 tng_time_per_frame_get(*input, &time);
972 tng_time_per_frame_set(*output, time);
974 tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
975 tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
977 for (int i = 0; i < defaultNumIds; i++)
979 if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval)
982 switch (fallbackIds[i])
984 case TNG_TRAJ_POSITIONS:
985 case TNG_TRAJ_VELOCITIES:
986 set_writing_interval(*output, interval, 3, fallbackIds[i],
987 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
990 case TNG_TRAJ_FORCES:
991 set_writing_interval(*output, interval, 3, fallbackIds[i],
992 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
993 TNG_GZIP_COMPRESSION);
995 case TNG_TRAJ_BOX_SHAPE:
996 set_writing_interval(*output, interval, 9, fallbackIds[i],
997 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
998 TNG_GZIP_COMPRESSION);
1000 case TNG_GMX_LAMBDA:
1001 set_writing_interval(*output, interval, 1, fallbackIds[i],
1002 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
1003 TNG_GZIP_COMPRESSION);
1014 /* TODO after trjconv is modularized: fix this so the user can
1015 change precision when they are doing an operation where
1016 this makes sense, and not otherwise.
1018 char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
1019 gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
1021 gmx_tng_add_mtop(*output, mtop);
1022 tng_num_frames_per_frame_set_set(*output, 1);
1025 if (index && nAtoms > 0)
1027 gmx_tng_setup_atom_subgroup(*output, nAtoms, index, indexGroupName);
1030 /* If for some reason there are more requested atoms than there are atoms in the
1031 * molecular system create a number of implicit atoms (without atom data) to
1032 * compensate for that. */
1035 tng_implicit_num_particles_set(*output, nAtoms);
1038 GMX_UNUSED_VALUE(filename);
1039 GMX_UNUSED_VALUE(mode);
1040 GMX_UNUSED_VALUE(input);
1041 GMX_UNUSED_VALUE(output);
1042 GMX_UNUSED_VALUE(nAtoms);
1043 GMX_UNUSED_VALUE(mtop);
1044 GMX_UNUSED_VALUE(index);
1045 GMX_UNUSED_VALUE(indexGroupName);
1049 void gmx_write_tng_from_trxframe(tng_trajectory_t output,
1050 const t_trxframe *frame,
1054 if (frame->step > 0)
1056 double timePerFrame = frame->time * PICO / frame->step;
1057 tng_time_per_frame_set(output, timePerFrame);
1061 natoms = frame->natoms;
1063 gmx_fwrite_tng(output,
1074 GMX_UNUSED_VALUE(output);
1075 GMX_UNUSED_VALUE(frame);
1076 GMX_UNUSED_VALUE(natoms);
1085 convert_array_to_real_array(void *from,
1090 const char datatype)
1094 const bool useDouble = GMX_DOUBLE;
1097 case TNG_FLOAT_DATA:
1102 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1106 for (i = 0; i < nAtoms; i++)
1108 for (j = 0; j < nValues; j++)
1110 to[i*nValues+j] = reinterpret_cast<float *>(from)[i*nValues+j] * fact;
1117 for (i = 0; i < nAtoms; i++)
1119 for (j = 0; j < nValues; j++)
1121 to[i*nValues+j] = reinterpret_cast<float *>(from)[i*nValues+j] * fact;
1127 for (i = 0; i < nAtoms; i++)
1129 for (j = 0; j < nValues; j++)
1131 to[i*nValues+j] = reinterpret_cast<gmx_int64_t *>(from)[i*nValues+j] * fact;
1135 case TNG_DOUBLE_DATA:
1136 if (sizeof(real) == sizeof(double))
1140 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1144 for (i = 0; i < nAtoms; i++)
1146 for (j = 0; j < nValues; j++)
1148 to[i*nValues+j] = reinterpret_cast<double *>(from)[i*nValues+j] * fact;
1155 for (i = 0; i < nAtoms; i++)
1157 for (j = 0; j < nValues; j++)
1159 to[i*nValues+j] = reinterpret_cast<double *>(from)[i*nValues+j] * fact;
1165 gmx_incons("Illegal datatype when converting values to a real array!");
1171 real getDistanceScaleFactor(tng_trajectory_t in)
1173 gmx_int64_t exp = -1;
1174 real distanceScaleFactor;
1176 // TODO Hopefully, TNG 2.0 will do this kind of thing for us
1177 tng_distance_unit_exponential_get(in, &exp);
1179 // GROMACS expects distances in nm
1183 distanceScaleFactor = NANO/NANO;
1186 distanceScaleFactor = NANO/ANGSTROM;
1189 distanceScaleFactor = pow(10.0, exp + 9.0);
1192 return distanceScaleFactor;
1198 void gmx_tng_setup_atom_subgroup(tng_trajectory_t tng,
1204 gmx_int64_t nAtoms, cnt, nMols;
1205 tng_molecule_t mol, iterMol;
1209 tng_function_status stat;
1211 tng_num_particles_get(tng, &nAtoms);
1218 stat = tng_molecule_find(tng, name, -1, &mol);
1219 if (stat == TNG_SUCCESS)
1221 tng_molecule_num_atoms_get(tng, mol, &nAtoms);
1222 tng_molecule_cnt_get(tng, mol, &cnt);
1232 if (stat == TNG_FAILURE)
1234 /* The indexed atoms are added to one separate molecule. */
1235 tng_molecule_alloc(tng, &mol);
1236 tng_molecule_name_set(tng, mol, name);
1237 tng_molecule_chain_add(tng, mol, "", &chain);
1239 for (int i = 0; i < nind; i++)
1241 char temp_name[256], temp_type[256];
1243 /* Try to retrieve the residue name of the atom */
1244 stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1245 if (stat != TNG_SUCCESS)
1247 temp_name[0] = '\0';
1249 /* Check if the molecule of the selection already contains this residue */
1250 if (tng_chain_residue_find(tng, chain, temp_name, -1, &res)
1253 tng_chain_residue_add(tng, chain, temp_name, &res);
1255 /* Try to find the original name and type of the atom */
1256 stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1257 if (stat != TNG_SUCCESS)
1259 temp_name[0] = '\0';
1261 stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
1262 if (stat != TNG_SUCCESS)
1264 temp_type[0] = '\0';
1266 tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
1268 tng_molecule_existing_add(tng, &mol);
1270 /* Set the count of the molecule containing the selected atoms to 1 and all
1271 * other molecules to 0 */
1272 tng_molecule_cnt_set(tng, mol, 1);
1273 tng_num_molecule_types_get(tng, &nMols);
1274 for (gmx_int64_t k = 0; k < nMols; k++)
1276 tng_molecule_of_index_get(tng, k, &iterMol);
1281 tng_molecule_cnt_set(tng, iterMol, 0);
1284 GMX_UNUSED_VALUE(tng);
1285 GMX_UNUSED_VALUE(nind);
1286 GMX_UNUSED_VALUE(ind);
1287 GMX_UNUSED_VALUE(name);
1291 /* TODO: If/when TNG acquires the ability to copy data blocks without
1292 * uncompressing them, then this implemenation should be reconsidered.
1293 * Ideally, gmx trjconv -f a.tng -o b.tng -b 10 -e 20 would be fast
1294 * and lose no information. */
1295 gmx_bool gmx_read_next_tng_frame(tng_trajectory_t input,
1297 gmx_int64_t *requestedIds,
1298 int numRequestedIds)
1301 gmx_bool bOK = TRUE;
1302 tng_function_status stat;
1303 gmx_int64_t numberOfAtoms = -1, frameNumber = -1;
1304 gmx_int64_t nBlocks, blockId, *blockIds = nullptr, codecId;
1306 void *values = nullptr;
1307 double frameTime = -1.0;
1308 int size, blockDependency;
1310 const int defaultNumIds = 5;
1311 static gmx_int64_t fallbackRequestedIds[defaultNumIds] =
1313 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
1314 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
1321 fr->bLambda = FALSE;
1329 /* If no specific IDs were requested read all block types that can
1330 * currently be interpreted */
1331 if (!requestedIds || numRequestedIds == 0)
1333 numRequestedIds = defaultNumIds;
1334 requestedIds = fallbackRequestedIds;
1337 stat = tng_num_particles_get(input, &numberOfAtoms);
1338 if (stat != TNG_SUCCESS)
1340 gmx_file("Cannot determine number of atoms from TNG file.");
1342 fr->natoms = numberOfAtoms;
1344 bool nextFrameExists = gmx_get_tng_data_block_types_of_next_frame(input,
1351 gmx::unique_cptr<gmx_int64_t, gmx::free_wrapper> blockIdsGuard(blockIds);
1352 if (!nextFrameExists)
1362 for (gmx_int64_t i = 0; i < nBlocks; i++)
1364 blockId = blockIds[i];
1365 tng_data_block_dependency_get(input, blockId, &blockDependency);
1366 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1368 stat = tng_util_particle_data_next_frame_read(input,
1377 stat = tng_util_non_particle_data_next_frame_read(input,
1384 if (stat == TNG_CRITICAL)
1386 gmx_file("Cannot read positions from TNG file.");
1389 else if (stat == TNG_FAILURE)
1395 case TNG_TRAJ_BOX_SHAPE:
1399 size = sizeof(gmx_int64_t);
1401 case TNG_FLOAT_DATA:
1402 size = sizeof(float);
1404 case TNG_DOUBLE_DATA:
1405 size = sizeof(double);
1408 gmx_incons("Illegal datatype of box shape values!");
1410 for (int i = 0; i < DIM; i++)
1412 convert_array_to_real_array(reinterpret_cast<char *>(values) + size * i * DIM,
1413 reinterpret_cast<real *>(fr->box[i]),
1414 getDistanceScaleFactor(input),
1421 case TNG_TRAJ_POSITIONS:
1422 srenew(fr->x, fr->natoms);
1423 convert_array_to_real_array(values,
1424 reinterpret_cast<real *>(fr->x),
1425 getDistanceScaleFactor(input),
1430 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1431 /* This must be updated if/when more lossy compression methods are added */
1432 if (codecId == TNG_TNG_COMPRESSION)
1438 case TNG_TRAJ_VELOCITIES:
1439 srenew(fr->v, fr->natoms);
1440 convert_array_to_real_array(values,
1442 getDistanceScaleFactor(input),
1447 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1448 /* This must be updated if/when more lossy compression methods are added */
1449 if (codecId == TNG_TNG_COMPRESSION)
1455 case TNG_TRAJ_FORCES:
1456 srenew(fr->f, fr->natoms);
1457 convert_array_to_real_array(values,
1458 reinterpret_cast<real *>(fr->f),
1459 getDistanceScaleFactor(input),
1465 case TNG_GMX_LAMBDA:
1468 case TNG_FLOAT_DATA:
1469 fr->lambda = *(reinterpret_cast<float *>(values));
1471 case TNG_DOUBLE_DATA:
1472 fr->lambda = *(reinterpret_cast<double *>(values));
1475 gmx_incons("Illegal datatype lambda value!");
1480 gmx_warning("Illegal block type! Currently GROMACS tools can only handle certain data types. Skipping block.");
1482 /* values does not have to be freed before reading next frame. It will
1483 * be reallocated if it is not NULL. */
1486 fr->step = frameNumber;
1488 // Convert the time to ps
1489 fr->time = frameTime / PICO;
1492 // TODO This does not leak, but is not exception safe.
1493 /* values must be freed before leaving this function */
1498 GMX_UNUSED_VALUE(input);
1499 GMX_UNUSED_VALUE(fr);
1500 GMX_UNUSED_VALUE(requestedIds);
1501 GMX_UNUSED_VALUE(numRequestedIds);
1506 void gmx_print_tng_molecule_system(tng_trajectory_t input,
1510 gmx_int64_t nMolecules, nChains, nResidues, nAtoms, nFramesRead;
1511 gmx_int64_t strideLength, nParticlesRead, nValuesPerFrameRead, *molCntList;
1512 tng_molecule_t molecule;
1514 tng_residue_t residue;
1516 tng_function_status stat;
1520 void *data = nullptr;
1521 std::vector<real> atomCharges;
1522 std::vector<real> atomMasses;
1524 tng_num_molecule_types_get(input, &nMolecules);
1525 tng_molecule_cnt_list_get(input, &molCntList);
1526 /* Can the number of particles change in the trajectory or is it constant? */
1527 tng_num_particles_variable_get(input, &varNAtoms);
1529 for (gmx_int64_t i = 0; i < nMolecules; i++)
1531 tng_molecule_of_index_get(input, i, &molecule);
1532 tng_molecule_name_get(input, molecule, str, 256);
1533 if (varNAtoms == TNG_CONSTANT_N_ATOMS)
1535 if ((int)molCntList[i] == 0)
1539 fprintf(stream, "Molecule: %s, count: %d\n", str, (int)molCntList[i]);
1543 fprintf(stream, "Molecule: %s\n", str);
1545 tng_molecule_num_chains_get(input, molecule, &nChains);
1548 for (gmx_int64_t j = 0; j < nChains; j++)
1550 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
1551 tng_chain_name_get(input, chain, str, 256);
1552 fprintf(stream, "\tChain: %s\n", str);
1553 tng_chain_num_residues_get(input, chain, &nResidues);
1554 for (gmx_int64_t k = 0; k < nResidues; k++)
1556 tng_chain_residue_of_index_get(input, chain, k, &residue);
1557 tng_residue_name_get(input, residue, str, 256);
1558 fprintf(stream, "\t\tResidue: %s\n", str);
1559 tng_residue_num_atoms_get(input, residue, &nAtoms);
1560 for (gmx_int64_t l = 0; l < nAtoms; l++)
1562 tng_residue_atom_of_index_get(input, residue, l, &atom);
1563 tng_atom_name_get(input, atom, str, 256);
1564 fprintf(stream, "\t\t\tAtom: %s", str);
1565 tng_atom_type_get(input, atom, str, 256);
1566 fprintf(stream, " (%s)\n", str);
1571 /* It is possible to have a molecule without chains, in which case
1572 * residues in the molecule can be iterated through without going
1573 * through chains. */
1576 tng_molecule_num_residues_get(input, molecule, &nResidues);
1579 for (gmx_int64_t k = 0; k < nResidues; k++)
1581 tng_molecule_residue_of_index_get(input, molecule, k, &residue);
1582 tng_residue_name_get(input, residue, str, 256);
1583 fprintf(stream, "\t\tResidue: %s\n", str);
1584 tng_residue_num_atoms_get(input, residue, &nAtoms);
1585 for (gmx_int64_t l = 0; l < nAtoms; l++)
1587 tng_residue_atom_of_index_get(input, residue, l, &atom);
1588 tng_atom_name_get(input, atom, str, 256);
1589 fprintf(stream, "\t\t\tAtom: %s", str);
1590 tng_atom_type_get(input, atom, str, 256);
1591 fprintf(stream, " (%s)\n", str);
1597 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
1598 for (gmx_int64_t l = 0; l < nAtoms; l++)
1600 tng_molecule_atom_of_index_get(input, molecule, l, &atom);
1601 tng_atom_name_get(input, atom, str, 256);
1602 fprintf(stream, "\t\t\tAtom: %s", str);
1603 tng_atom_type_get(input, atom, str, 256);
1604 fprintf(stream, " (%s)\n", str);
1610 tng_num_particles_get(input, &nAtoms);
1611 stat = tng_particle_data_vector_get(input, TNG_TRAJ_PARTIAL_CHARGES, &data, &nFramesRead,
1612 &strideLength, &nParticlesRead,
1613 &nValuesPerFrameRead, &datatype);
1614 if (stat == TNG_SUCCESS)
1616 atomCharges.resize(nAtoms);
1617 convert_array_to_real_array(data,
1624 fprintf(stream, "Atom Charges (%d):\n", int(nAtoms));
1625 for (gmx_int64_t i = 0; i < nAtoms; i += 10)
1627 fprintf(stream, "Atom Charges [%8d-]=[", int(i));
1628 for (gmx_int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1630 fprintf(stream, " %12.5e", atomCharges[i + j]);
1632 fprintf(stream, "]\n");
1636 stat = tng_particle_data_vector_get(input, TNG_TRAJ_MASSES, &data, &nFramesRead,
1637 &strideLength, &nParticlesRead,
1638 &nValuesPerFrameRead, &datatype);
1639 if (stat == TNG_SUCCESS)
1641 atomMasses.resize(nAtoms);
1642 convert_array_to_real_array(data,
1649 fprintf(stream, "Atom Masses (%d):\n", int(nAtoms));
1650 for (gmx_int64_t i = 0; i < nAtoms; i += 10)
1652 fprintf(stream, "Atom Masses [%8d-]=[", int(i));
1653 for (gmx_int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1655 fprintf(stream, " %12.5e", atomMasses[i + j]);
1657 fprintf(stream, "]\n");
1663 GMX_UNUSED_VALUE(input);
1664 GMX_UNUSED_VALUE(stream);
1668 gmx_bool gmx_get_tng_data_block_types_of_next_frame(tng_trajectory_t input,
1671 gmx_int64_t *requestedIds,
1672 gmx_int64_t *nextFrame,
1673 gmx_int64_t *nBlocks,
1674 gmx_int64_t **blockIds)
1677 tng_function_status stat;
1679 stat = tng_util_trajectory_next_frame_present_data_blocks_find(input, frame,
1680 nRequestedIds, requestedIds,
1684 if (stat == TNG_CRITICAL)
1686 gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
1688 else if (stat == TNG_FAILURE)
1694 GMX_UNUSED_VALUE(input);
1695 GMX_UNUSED_VALUE(frame);
1696 GMX_UNUSED_VALUE(nRequestedIds);
1697 GMX_UNUSED_VALUE(requestedIds);
1698 GMX_UNUSED_VALUE(nextFrame);
1699 GMX_UNUSED_VALUE(nBlocks);
1700 GMX_UNUSED_VALUE(blockIds);
1705 gmx_bool gmx_get_tng_data_next_frame_of_block_type(tng_trajectory_t input,
1706 gmx_int64_t blockId,
1708 gmx_int64_t *frameNumber,
1710 gmx_int64_t *nValuesPerFrame,
1711 gmx_int64_t *nAtoms,
1718 tng_function_status stat;
1720 gmx_int64_t codecId;
1721 int blockDependency;
1722 void *data = nullptr;
1725 stat = tng_data_block_name_get(input, blockId, name, maxLen);
1726 if (stat != TNG_SUCCESS)
1728 gmx_file("Cannot read next frame of TNG file");
1730 stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
1731 if (stat != TNG_SUCCESS)
1733 gmx_file("Cannot read next frame of TNG file");
1735 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1737 tng_num_particles_get(input, nAtoms);
1738 stat = tng_util_particle_data_next_frame_read(input,
1747 *nAtoms = 1; /* There are not actually any atoms, but it is used for
1748 allocating memory */
1749 stat = tng_util_non_particle_data_next_frame_read(input,
1756 if (stat == TNG_CRITICAL)
1758 gmx_file("Cannot read next frame of TNG file");
1760 if (stat == TNG_FAILURE)
1766 stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
1767 if (stat != TNG_SUCCESS)
1769 gmx_file("Cannot read next frame of TNG file");
1771 srenew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
1772 convert_array_to_real_array(data,
1774 getDistanceScaleFactor(input),
1779 tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
1781 /* This must be updated if/when more lossy compression methods are added */
1782 if (codecId != TNG_TNG_COMPRESSION)
1796 GMX_UNUSED_VALUE(input);
1797 GMX_UNUSED_VALUE(blockId);
1798 GMX_UNUSED_VALUE(values);
1799 GMX_UNUSED_VALUE(frameNumber);
1800 GMX_UNUSED_VALUE(frameTime);
1801 GMX_UNUSED_VALUE(nValuesPerFrame);
1802 GMX_UNUSED_VALUE(nAtoms);
1803 GMX_UNUSED_VALUE(prec);
1804 GMX_UNUSED_VALUE(name);
1805 GMX_UNUSED_VALUE(maxLen);
1806 GMX_UNUSED_VALUE(bOK);