2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2018,2019, 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"
69 using tng_trajectory_t = void*;
72 /*! \brief Gromacs Wrapper around tng datatype
74 * This could in principle hold any GROMACS-specific requirements not yet
75 * implemented in or not relevant to the TNG library itself. However, for now
76 * we only use it to handle some shortcomings we have discovered, where the TNG
77 * API itself is a bit fragile and can end up overwriting data if called several
78 * times with the same frame number.
79 * The logic to determine the time per step was also a bit fragile. This is not
80 * critical, but since we anyway need a wrapper for ensuring unique frame
81 * numbers, we can also use it to store the time of the first step and use that
82 * to derive a slightly better/safer estimate of the time per step.
84 * At some future point where we have a second-generation TNG API we should
85 * consider removing this again.
87 struct gmx_tng_trajectory
89 tng_trajectory_t tng; //!< Actual TNG handle (pointer)
90 bool lastStepDataIsValid; //!< True if lastStep has been set
91 std::int64_t lastStep; //!< Index/step used for last frame
92 bool lastTimeDataIsValid; //!< True if lastTime has been set
93 double lastTime; //!< Time of last frame (TNG unit is seconds)
94 bool timePerFrameIsSet; //!< True if we have set the time per frame
95 int boxOutputInterval; //!< Number of steps between the output of box size
96 int lambdaOutputInterval; //!< Number of steps between the output of lambdas
100 static const char* modeToVerb(char mode)
105 case 'r': p = "reading"; break;
106 case 'w': p = "writing"; break;
107 case 'a': p = "appending"; break;
108 default: gmx_fatal(FARGS, "Invalid file opening mode %c", mode);
114 void gmx_tng_open(const char* filename, char mode, gmx_tng_trajectory_t* gmx_tng)
117 /* First check whether we have to make a backup,
118 * only for writing, not for read or append.
122 make_backup(filename);
125 *gmx_tng = new gmx_tng_trajectory;
126 (*gmx_tng)->lastStepDataIsValid = false;
127 (*gmx_tng)->lastTimeDataIsValid = false;
128 (*gmx_tng)->timePerFrameIsSet = false;
129 tng_trajectory_t* tng = &(*gmx_tng)->tng;
131 /* tng must not be pointing at already allocated memory.
132 * Memory will be allocated by tng_util_trajectory_open() and must
133 * later on be freed by tng_util_trajectory_close(). */
134 if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
136 /* TNG does return more than one degree of error, but there is
137 no use case for GROMACS handling the non-fatal errors
139 gmx_fatal(FARGS, "File I/O error while opening %s for %s", filename, modeToVerb(mode));
142 if (mode == 'w' || mode == 'a')
145 gmx_gethostname(hostname, 256);
148 tng_first_computer_name_set(*tng, hostname);
152 tng_last_computer_name_set(*tng, hostname);
155 char programInfo[256];
156 const char* precisionString = "";
158 precisionString = " (double precision)";
160 sprintf(programInfo, "%.100s %.128s%.24s", gmx::getProgramContext().displayName(),
161 gmx_version(), precisionString);
164 tng_first_program_name_set(*tng, programInfo);
168 tng_last_program_name_set(*tng, programInfo);
172 if (!gmx_getusername(username, 256))
176 tng_first_user_name_set(*tng, username);
180 tng_last_user_name_set(*tng, username);
181 tng_file_headers_write(*tng, TNG_USE_HASH);
186 gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
187 GMX_UNUSED_VALUE(filename);
188 GMX_UNUSED_VALUE(mode);
189 GMX_UNUSED_VALUE(gmx_tng);
193 void gmx_tng_close(gmx_tng_trajectory_t* gmx_tng)
195 /* We have to check that tng is set because
196 * tng_util_trajectory_close wants to return a NULL in it, and
197 * gives a fatal error if it is NULL. */
199 if (gmx_tng == nullptr || *gmx_tng == nullptr)
203 tng_trajectory_t* tng = &(*gmx_tng)->tng;
207 tng_util_trajectory_close(tng);
213 GMX_UNUSED_VALUE(gmx_tng);
218 static void addTngMoleculeFromTopology(gmx_tng_trajectory_t gmx_tng,
219 const char* moleculeName,
220 const t_atoms* atoms,
221 int64_t numMolecules,
222 tng_molecule_t* tngMol)
224 tng_trajectory_t tng = gmx_tng->tng;
225 tng_chain_t tngChain = nullptr;
226 tng_residue_t tngRes = nullptr;
228 if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
230 gmx_file("Cannot add molecule to TNG molecular system.");
233 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++)
235 const t_atom* at = &atoms->atom[atomIndex];
236 /* FIXME: Currently the TNG API can only add atoms belonging to a
237 * residue and chain. Wait for TNG 2.0*/
240 const t_resinfo* resInfo = &atoms->resinfo[at->resind];
241 char chainName[2] = { resInfo->chainid, 0 };
242 tng_atom_t tngAtom = nullptr;
247 prevAtom = &atoms->atom[atomIndex - 1];
254 /* If this is the first atom or if the residue changed add the
255 * residue to the TNG molecular system. */
256 if (!prevAtom || resInfo != &atoms->resinfo[prevAtom->resind])
258 /* If this is the first atom or if the chain changed add
259 * the chain to the TNG molecular system. */
260 if (!prevAtom || resInfo->chainid != atoms->resinfo[prevAtom->resind].chainid)
262 tng_molecule_chain_add(tng, *tngMol, chainName, &tngChain);
264 /* FIXME: When TNG supports both residue index and residue
265 * number the latter should be used. Wait for TNG 2.0*/
266 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
268 tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIndex]),
269 *(atoms->atomtype[atomIndex]), &tngAtom);
272 tng_molecule_cnt_set(tng, *tngMol, numMolecules);
275 void gmx_tng_add_mtop(gmx_tng_trajectory_t gmx_tng, const gmx_mtop_t* mtop)
279 std::vector<real> atomCharges;
280 std::vector<real> atomMasses;
284 tng_trajectory_t tng = gmx_tng->tng;
288 /* No topology information available to add. */
293 datatype = TNG_DOUBLE_DATA;
295 datatype = TNG_FLOAT_DATA;
298 atomCharges.reserve(mtop->natoms);
299 atomMasses.reserve(mtop->natoms);
301 for (const gmx_molblock_t& molBlock : mtop->molblock)
303 tng_molecule_t tngMol = nullptr;
304 const gmx_moltype_t* molType = &mtop->moltype[molBlock.type];
306 /* Add a molecule to the TNG trajectory with the same name as the
307 * current molecule. */
308 addTngMoleculeFromTopology(gmx_tng, *(molType->name), &molType->atoms, molBlock.nmol, &tngMol);
310 /* Bonds have to be deduced from interactions (constraints etc). Different
311 * interactions have different sets of parameters. */
312 /* Constraints are specified using two atoms */
313 for (i = 0; i < F_NRE; i++)
317 const InteractionList& ilist = molType->ilist[i];
319 while (j < ilist.size())
321 tng_molecule_bond_add(tng, tngMol, ilist.iatoms[j], ilist.iatoms[j + 1], &tngBond);
326 /* Settle is described using three atoms */
327 const InteractionList& ilist = molType->ilist[F_SETTLE];
329 while (j < ilist.size())
331 tng_molecule_bond_add(tng, tngMol, ilist.iatoms[j], ilist.iatoms[j + 1], &tngBond);
332 tng_molecule_bond_add(tng, tngMol, ilist.iatoms[j], ilist.iatoms[j + 2], &tngBond);
335 /* First copy atom charges and masses, first atom by atom and then copy the memory for the molecule instances.
336 * FIXME: Atom B state data should also be written to TNG (v 2.0?) */
337 for (int atomCounter = 0; atomCounter < molType->atoms.nr; atomCounter++)
339 atomCharges.push_back(molType->atoms.atom[atomCounter].q);
340 atomMasses.push_back(molType->atoms.atom[atomCounter].m);
342 for (int molCounter = 1; molCounter < molBlock.nmol; molCounter++)
344 std::copy_n(atomCharges.end() - molType->atoms.nr, molType->atoms.nr,
345 std::back_inserter(atomCharges));
346 std::copy_n(atomMasses.end() - molType->atoms.nr, molType->atoms.nr,
347 std::back_inserter(atomMasses));
350 /* Write the TNG data blocks. */
351 tng_particle_data_block_add(tng, TNG_TRAJ_PARTIAL_CHARGES, "PARTIAL CHARGES", datatype,
352 TNG_NON_TRAJECTORY_BLOCK, 1, 1, 1, 0, mtop->natoms,
353 TNG_GZIP_COMPRESSION, atomCharges.data());
354 tng_particle_data_block_add(tng, TNG_TRAJ_MASSES, "ATOM MASSES", datatype, TNG_NON_TRAJECTORY_BLOCK,
355 1, 1, 1, 0, mtop->natoms, TNG_GZIP_COMPRESSION, atomMasses.data());
358 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
359 * if they are positive.
361 * If only one of n1 and n2 is positive, then return it.
362 * If neither n1 or n2 is positive, then return -1. */
363 static int greatest_common_divisor_if_positive(int n1, int n2)
367 return (0 >= n2) ? -1 : n2;
374 /* We have a non-trivial greatest common divisor to compute. */
375 return gmx_greatest_common_divisor(n1, n2);
378 /* By default try to write 100 frames (of actual output) in each frame set.
379 * This number is the number of outputs of the most frequently written data
380 * type per frame set.
381 * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
382 * setups regarding compression efficiency and compression time. Make this
383 * a hidden command-line option? */
384 const int defaultFramesPerFrameSet = 100;
386 /*! \libinternal \brief Set the number of frames per frame
387 * set according to output intervals.
388 * The default is that 100 frames are written of the data
389 * that is written most often. */
390 static void tng_set_frames_per_frame_set(gmx_tng_trajectory_t gmx_tng,
391 const gmx_bool bUseLossyCompression,
392 const t_inputrec* ir)
395 tng_trajectory_t tng = gmx_tng->tng;
397 /* Set the number of frames per frame set to contain at least
398 * defaultFramesPerFrameSet of the lowest common denominator of
399 * the writing interval of positions and velocities. */
400 /* FIXME after 5.0: consider nstenergy also? */
401 if (bUseLossyCompression)
403 gcd = ir->nstxout_compressed;
407 gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
408 gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
415 tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
418 /*! \libinternal \brief Set the data-writing intervals, and number of
419 * frames per frame set */
420 static void set_writing_intervals(gmx_tng_trajectory_t gmx_tng,
421 const gmx_bool bUseLossyCompression,
422 const t_inputrec* ir)
424 tng_trajectory_t tng = gmx_tng->tng;
426 /* Define pointers to specific writing functions depending on if we
427 * write float or double data */
428 typedef tng_function_status (*set_writing_interval_func_pointer)(
429 tng_trajectory_t, const int64_t, const int64_t, const int64_t, const char*, const char,
432 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
434 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
436 int xout, vout, fout;
437 int gcd = -1, lowest = -1;
440 tng_set_frames_per_frame_set(gmx_tng, bUseLossyCompression, ir);
442 if (bUseLossyCompression)
444 xout = ir->nstxout_compressed;
446 /* If there is no uncompressed coordinate output write forces
447 and velocities to the compressed tng file. */
458 compression = TNG_TNG_COMPRESSION;
465 compression = TNG_GZIP_COMPRESSION;
469 set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS, "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
471 /* TODO: if/when we write energies to TNG also, reconsider how
472 * and when box information is written, because GROMACS
473 * behaviour pre-5.0 was to write the box with every
474 * trajectory frame and every energy frame, and probably
475 * people depend on this. */
477 gcd = greatest_common_divisor_if_positive(gcd, xout);
478 if (lowest < 0 || xout < lowest)
485 set_writing_interval(tng, vout, 3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
486 TNG_PARTICLE_BLOCK_DATA, compression);
488 gcd = greatest_common_divisor_if_positive(gcd, vout);
489 if (lowest < 0 || vout < lowest)
496 set_writing_interval(tng, fout, 3, TNG_TRAJ_FORCES, "FORCES", TNG_PARTICLE_BLOCK_DATA,
497 TNG_GZIP_COMPRESSION);
499 gcd = greatest_common_divisor_if_positive(gcd, fout);
500 if (lowest < 0 || fout < lowest)
507 /* Lambdas and box shape written at an interval of the lowest common
508 denominator of other output */
509 set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA, "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
510 TNG_GZIP_COMPRESSION);
512 set_writing_interval(tng, gcd, 9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
513 TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION);
514 gmx_tng->lambdaOutputInterval = gcd;
515 gmx_tng->boxOutputInterval = gcd;
516 if (gcd < lowest / 10)
519 "The lowest common denominator of trajectory output is "
520 "every %d step(s), whereas the shortest output interval "
521 "is every %d steps.",
528 void gmx_tng_prepare_md_writing(gmx_tng_trajectory_t gmx_tng, const gmx_mtop_t* mtop, const t_inputrec* ir)
531 gmx_tng_add_mtop(gmx_tng, mtop);
532 set_writing_intervals(gmx_tng, FALSE, ir);
533 tng_time_per_frame_set(gmx_tng->tng, ir->delta_t * PICO);
534 gmx_tng->timePerFrameIsSet = true;
536 GMX_UNUSED_VALUE(gmx_tng);
537 GMX_UNUSED_VALUE(mtop);
538 GMX_UNUSED_VALUE(ir);
543 /* Check if all atoms in the molecule system are selected
544 * by a selection group of type specified by gtype. */
545 static gmx_bool all_atoms_selected(const gmx_mtop_t* mtop, const SimulationAtomGroupType gtype)
547 /* Iterate through all atoms in the molecule system and
548 * check if they belong to a selection group of the
551 for (const gmx_molblock_t& molBlock : mtop->molblock)
553 const gmx_moltype_t& molType = mtop->moltype[molBlock.type];
554 const t_atoms& atoms = molType.atoms;
556 for (int j = 0; j < molBlock.nmol; j++)
558 for (int atomIndex = 0; atomIndex < atoms.nr; atomIndex++, i++)
560 if (getGroupType(mtop->groups, gtype, i) != 0)
570 /* Create TNG molecules which will represent each of the selection
571 * groups for writing. But do that only if there is actually a
572 * specified selection group and it is not the whole system.
573 * TODO: Currently the only selection that is taken into account
574 * is egcCompressedX, but other selections should be added when
575 * e.g. writing energies is implemented.
577 static void add_selection_groups(gmx_tng_trajectory_t gmx_tng, const gmx_mtop_t* mtop)
579 const t_atoms* atoms;
581 const t_resinfo* resInfo;
584 tng_molecule_t mol, iterMol;
591 tng_trajectory_t tng = gmx_tng->tng;
593 /* TODO: When the TNG molecules block is more flexible TNG selection
594 * groups should not need all atoms specified. It should be possible
595 * just to specify what molecules are selected (and/or which atoms in
596 * the molecule) and how many (if applicable). */
598 /* If no atoms are selected we do not need to create a
599 * TNG selection group molecule. */
600 if (mtop->groups.numberOfGroupNumbers(SimulationAtomGroupType::CompressedPositionOutput) == 0)
605 /* If all atoms are selected we do not have to create a selection
606 * group molecule in the TNG molecule system. */
607 if (all_atoms_selected(mtop, SimulationAtomGroupType::CompressedPositionOutput))
612 /* The name of the TNG molecule containing the selection group is the
613 * same as the name of the selection group. */
614 nameIndex = mtop->groups.groups[SimulationAtomGroupType::CompressedPositionOutput][0];
615 groupName = *mtop->groups.groupNames[nameIndex];
617 tng_molecule_alloc(tng, &mol);
618 tng_molecule_name_set(tng, mol, groupName);
619 tng_molecule_chain_add(tng, mol, "", &chain);
621 for (const gmx_molblock_t& molBlock : mtop->molblock)
623 const gmx_moltype_t& molType = mtop->moltype[molBlock.type];
625 atoms = &molType.atoms;
627 for (int j = 0; j < molBlock.nmol; j++)
629 bool bAtomsAdded = FALSE;
630 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
635 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, i) != 0)
639 at = &atoms->atom[atomIndex];
642 resInfo = &atoms->resinfo[at->resind];
643 /* FIXME: When TNG supports both residue index and residue
644 * number the latter should be used. */
645 res_name = *resInfo->name;
646 res_id = at->resind + 1;
650 res_name = const_cast<char*>("");
653 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res) != TNG_SUCCESS)
655 /* Since there is ONE chain for selection groups do not keep the
656 * original residue IDs - otherwise there might be conflicts. */
657 tng_chain_residue_add(tng, chain, res_name, &res);
659 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
660 *(atoms->atomtype[atomIndex]), atom_offset + atomIndex, &atom);
666 for (int k = 0; k < F_NRE; k++)
670 const InteractionList& ilist = molType.ilist[k];
671 for (int l = 1; l < ilist.size(); l += 3)
674 atom1 = ilist.iatoms[l] + atom_offset;
675 atom2 = ilist.iatoms[l + 1] + atom_offset;
676 if (getGroupType(mtop->groups,
677 SimulationAtomGroupType::CompressedPositionOutput, atom1)
679 && getGroupType(mtop->groups,
680 SimulationAtomGroupType::CompressedPositionOutput, atom2)
683 tng_molecule_bond_add(tng, mol, ilist.iatoms[l],
684 ilist.iatoms[l + 1], &tngBond);
689 /* Settle is described using three atoms */
690 const InteractionList& ilist = molType.ilist[F_SETTLE];
691 for (int l = 1; l < ilist.size(); l += 4)
693 int atom1, atom2, atom3;
694 atom1 = ilist.iatoms[l] + atom_offset;
695 atom2 = ilist.iatoms[l + 1] + atom_offset;
696 atom3 = ilist.iatoms[l + 2] + atom_offset;
697 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, atom1)
700 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, atom2)
703 tng_molecule_bond_add(tng, mol, atom1, atom2, &tngBond);
705 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, atom3)
708 tng_molecule_bond_add(tng, mol, atom1, atom3, &tngBond);
713 atom_offset += atoms->nr;
716 tng_molecule_existing_add(tng, &mol);
717 tng_molecule_cnt_set(tng, mol, 1);
718 tng_num_molecule_types_get(tng, &nMols);
719 for (int64_t k = 0; k < nMols; k++)
721 tng_molecule_of_index_get(tng, k, &iterMol);
726 tng_molecule_cnt_set(tng, iterMol, 0);
731 void gmx_tng_set_compression_precision(gmx_tng_trajectory_t gmx_tng, real prec)
734 tng_compression_precision_set(gmx_tng->tng, prec);
736 GMX_UNUSED_VALUE(gmx_tng);
737 GMX_UNUSED_VALUE(prec);
741 void gmx_tng_prepare_low_prec_writing(gmx_tng_trajectory_t gmx_tng, const gmx_mtop_t* mtop, const t_inputrec* ir)
744 gmx_tng_add_mtop(gmx_tng, mtop);
745 add_selection_groups(gmx_tng, mtop);
746 set_writing_intervals(gmx_tng, TRUE, ir);
747 tng_time_per_frame_set(gmx_tng->tng, ir->delta_t * PICO);
748 gmx_tng->timePerFrameIsSet = true;
749 gmx_tng_set_compression_precision(gmx_tng, ir->x_compression_precision);
751 GMX_UNUSED_VALUE(gmx_tng);
752 GMX_UNUSED_VALUE(mtop);
753 GMX_UNUSED_VALUE(ir);
757 void gmx_fwrite_tng(gmx_tng_trajectory_t gmx_tng,
758 const gmx_bool bUseLossyCompression,
760 real elapsedPicoSeconds,
769 typedef tng_function_status (*write_data_func_pointer)(
770 tng_trajectory_t, const int64_t, const double, const real*, const int64_t,
771 const int64_t, const char*, const char, const char);
773 static write_data_func_pointer write_data = tng_util_generic_with_time_double_write;
775 static write_data_func_pointer write_data = tng_util_generic_with_time_write;
777 double elapsedSeconds = elapsedPicoSeconds * PICO;
784 /* This function might get called when the type of the
785 compressed trajectory is actually XTC. So we exit and move
789 tng_trajectory_t tng = gmx_tng->tng;
791 // While the GROMACS interface to this routine specifies 'step', TNG itself
792 // only uses 'frame index' internally, although it suggests that it's a good
793 // idea to let this represent the step rather than just frame numbers.
795 // However, the way the GROMACS interface works, there are cases where
796 // the step is simply not set, so all frames rather have step=0.
797 // If we call the lower-level TNG routines with the same frame index
798 // (which is set from the step) multiple times, things will go very
799 // wrong (overwritten data).
801 // To avoid this, we require the step value to always be larger than the
802 // previous value, and if this is not true we simply set it to a value
803 // one unit larger than the previous step.
805 // Note: We do allow the user to specify any crazy value the want for the
806 // first step. If it is "not set", GROMACS will have used the value 0.
807 if (gmx_tng->lastStepDataIsValid && step <= gmx_tng->lastStep)
809 step = gmx_tng->lastStep + 1;
812 // Now that we have fixed the potentially incorrect step, we can also set
813 // the time per frame if it was not already done, and we have
814 // step/time information from the last frame so it is possible to calculate it.
815 // Note that we do allow the time to be the same for all steps, or even
816 // decreasing. In that case we will simply say the time per step is 0
817 // or negative. A bit strange, but a correct representation of the data :-)
818 if (!gmx_tng->timePerFrameIsSet && gmx_tng->lastTimeDataIsValid && gmx_tng->lastStepDataIsValid)
820 double deltaTime = elapsedSeconds - gmx_tng->lastTime;
821 std::int64_t deltaStep = step - gmx_tng->lastStep;
822 tng_time_per_frame_set(tng, deltaTime / deltaStep);
823 gmx_tng->timePerFrameIsSet = true;
826 tng_num_particles_get(tng, &nParticles);
827 if (nAtoms != static_cast<int>(nParticles))
829 tng_implicit_num_particles_set(tng, nAtoms);
832 if (bUseLossyCompression)
834 compression = TNG_TNG_COMPRESSION;
838 compression = TNG_GZIP_COMPRESSION;
841 /* The writing is done using write_data, which writes float or double
842 * depending on the GROMACS compilation. */
845 GMX_ASSERT(box, "Need a non-NULL box if positions are written");
847 if (write_data(tng, step, elapsedSeconds, reinterpret_cast<const real*>(x), 3,
848 TNG_TRAJ_POSITIONS, "POSITIONS", TNG_PARTICLE_BLOCK_DATA, compression)
851 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
857 if (write_data(tng, step, elapsedSeconds, reinterpret_cast<const real*>(v), 3,
858 TNG_TRAJ_VELOCITIES, "VELOCITIES", TNG_PARTICLE_BLOCK_DATA, compression)
861 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
867 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
868 * compression for forces regardless of output mode */
869 if (write_data(tng, step, elapsedSeconds, reinterpret_cast<const real*>(f), 3,
870 TNG_TRAJ_FORCES, "FORCES", TNG_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION)
873 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
879 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
880 * compression for box shape regardless of output mode */
881 if (write_data(tng, step, elapsedSeconds, reinterpret_cast<const real*>(box), 9,
882 TNG_TRAJ_BOX_SHAPE, "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION)
885 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
891 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
892 * compression for lambda regardless of output mode */
893 if (write_data(tng, step, elapsedSeconds, reinterpret_cast<const real*>(&lambda), 1,
894 TNG_GMX_LAMBDA, "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION)
897 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
901 // Update the data in the wrapper
903 gmx_tng->lastStepDataIsValid = true;
904 gmx_tng->lastStep = step;
905 gmx_tng->lastTimeDataIsValid = true;
906 gmx_tng->lastTime = elapsedSeconds;
908 GMX_UNUSED_VALUE(gmx_tng);
909 GMX_UNUSED_VALUE(bUseLossyCompression);
910 GMX_UNUSED_VALUE(step);
911 GMX_UNUSED_VALUE(elapsedPicoSeconds);
912 GMX_UNUSED_VALUE(lambda);
913 GMX_UNUSED_VALUE(box);
914 GMX_UNUSED_VALUE(nAtoms);
921 void fflush_tng(gmx_tng_trajectory_t gmx_tng)
928 tng_frame_set_premature_write(gmx_tng->tng, TNG_USE_HASH);
930 GMX_UNUSED_VALUE(gmx_tng);
934 float gmx_tng_get_time_of_final_frame(gmx_tng_trajectory_t gmx_tng)
940 tng_trajectory_t tng = gmx_tng->tng;
942 tng_num_frames_get(tng, &nFrames);
943 tng_util_time_of_frame_get(tng, nFrames - 1, &time);
948 GMX_UNUSED_VALUE(gmx_tng);
953 void gmx_prepare_tng_writing(const char* filename,
955 gmx_tng_trajectory_t* gmx_tng_input,
956 gmx_tng_trajectory_t* gmx_tng_output,
958 const gmx_mtop_t* mtop,
959 gmx::ArrayRef<const int> index,
960 const char* indexGroupName)
963 tng_trajectory_t* input = (gmx_tng_input && *gmx_tng_input) ? &(*gmx_tng_input)->tng : nullptr;
964 /* FIXME after 5.0: Currently only standard block types are read */
965 const int defaultNumIds = 5;
966 static int64_t fallbackIds[defaultNumIds] = { TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
967 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES, TNG_GMX_LAMBDA };
968 static char fallbackNames[defaultNumIds][32] = { "BOX SHAPE", "POSITIONS", "VELOCITIES",
969 "FORCES", "LAMBDAS" };
971 typedef tng_function_status (*set_writing_interval_func_pointer)(
972 tng_trajectory_t, const int64_t, const int64_t, const int64_t, const char*, const char,
975 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
977 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
980 gmx_tng_open(filename, mode, gmx_tng_output);
981 tng_trajectory_t* output = &(*gmx_tng_output)->tng;
983 /* Do we have an input file in TNG format? If so, then there's
984 more data we can copy over, rather than having to improvise. */
985 if (gmx_tng_input && *gmx_tng_input)
987 /* Set parameters (compression, time per frame, molecule
988 * information, number of frames per frame set and writing
989 * intervals of positions, box shape and lambdas) of the
990 * output tng container based on their respective values int
991 * the input tng container */
992 double time, compression_precision;
993 int64_t n_frames_per_frame_set, interval = -1;
995 tng_compression_precision_get(*input, &compression_precision);
996 tng_compression_precision_set(*output, compression_precision);
997 // TODO make this configurable in a future version
998 char compression_type = TNG_TNG_COMPRESSION;
1000 tng_molecule_system_copy(*input, *output);
1002 tng_time_per_frame_get(*input, &time);
1003 tng_time_per_frame_set(*output, time);
1004 // Since we have copied the value from the input TNG we should not change it again
1005 (*gmx_tng_output)->timePerFrameIsSet = true;
1007 tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
1008 tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
1010 for (int i = 0; i < defaultNumIds; i++)
1012 if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval) == TNG_SUCCESS)
1014 switch (fallbackIds[i])
1016 case TNG_TRAJ_POSITIONS:
1017 case TNG_TRAJ_VELOCITIES:
1018 set_writing_interval(*output, interval, 3, fallbackIds[i], fallbackNames[i],
1019 TNG_PARTICLE_BLOCK_DATA, compression_type);
1021 case TNG_TRAJ_FORCES:
1022 set_writing_interval(*output, interval, 3, fallbackIds[i], fallbackNames[i],
1023 TNG_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION);
1025 case TNG_TRAJ_BOX_SHAPE:
1026 set_writing_interval(*output, interval, 9, fallbackIds[i], fallbackNames[i],
1027 TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION);
1028 (*gmx_tng_output)->boxOutputInterval = interval;
1030 case TNG_GMX_LAMBDA:
1031 set_writing_interval(*output, interval, 1, fallbackIds[i], fallbackNames[i],
1032 TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION);
1033 (*gmx_tng_output)->lambdaOutputInterval = interval;
1042 /* TODO after trjconv is modularized: fix this so the user can
1043 change precision when they are doing an operation where
1044 this makes sense, and not otherwise.
1046 char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
1047 gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
1049 gmx_tng_add_mtop(*gmx_tng_output, mtop);
1050 tng_num_frames_per_frame_set_set(*output, 1);
1053 if ((!index.empty()) && nAtoms > 0)
1055 gmx_tng_setup_atom_subgroup(*gmx_tng_output, index, indexGroupName);
1058 /* If for some reason there are more requested atoms than there are atoms in the
1059 * molecular system create a number of implicit atoms (without atom data) to
1060 * compensate for that. */
1063 tng_implicit_num_particles_set(*output, nAtoms);
1066 GMX_UNUSED_VALUE(filename);
1067 GMX_UNUSED_VALUE(mode);
1068 GMX_UNUSED_VALUE(gmx_tng_input);
1069 GMX_UNUSED_VALUE(gmx_tng_output);
1070 GMX_UNUSED_VALUE(nAtoms);
1071 GMX_UNUSED_VALUE(mtop);
1072 GMX_UNUSED_VALUE(index);
1073 GMX_UNUSED_VALUE(indexGroupName);
1077 void gmx_write_tng_from_trxframe(gmx_tng_trajectory_t gmx_tng_output, const t_trxframe* frame, int natoms)
1082 natoms = frame->natoms;
1084 gmx_fwrite_tng(gmx_tng_output, TRUE, frame->step, frame->time, 0, frame->box, natoms, frame->x,
1085 frame->v, frame->f);
1087 GMX_UNUSED_VALUE(gmx_tng_output);
1088 GMX_UNUSED_VALUE(frame);
1089 GMX_UNUSED_VALUE(natoms);
1097 void convert_array_to_real_array(void* from,
1102 const char datatype)
1106 const bool useDouble = GMX_DOUBLE;
1109 case TNG_FLOAT_DATA:
1114 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1118 for (i = 0; i < nAtoms; i++)
1120 for (j = 0; j < nValues; j++)
1122 to[i * nValues + j] = reinterpret_cast<float*>(from)[i * nValues + j] * fact;
1129 for (i = 0; i < nAtoms; i++)
1131 for (j = 0; j < nValues; j++)
1133 to[i * nValues + j] = reinterpret_cast<float*>(from)[i * nValues + j] * fact;
1139 for (i = 0; i < nAtoms; i++)
1141 for (j = 0; j < nValues; j++)
1143 to[i * nValues + j] = reinterpret_cast<int64_t*>(from)[i * nValues + j] * fact;
1147 case TNG_DOUBLE_DATA:
1148 if (sizeof(real) == sizeof(double))
1152 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1156 for (i = 0; i < nAtoms; i++)
1158 for (j = 0; j < nValues; j++)
1160 to[i * nValues + j] = reinterpret_cast<double*>(from)[i * nValues + j] * fact;
1167 for (i = 0; i < nAtoms; i++)
1169 for (j = 0; j < nValues; j++)
1171 to[i * nValues + j] = reinterpret_cast<double*>(from)[i * nValues + j] * fact;
1176 default: gmx_incons("Illegal datatype when converting values to a real array!");
1180 real getDistanceScaleFactor(gmx_tng_trajectory_t in)
1183 real distanceScaleFactor;
1185 // TODO Hopefully, TNG 2.0 will do this kind of thing for us
1186 tng_distance_unit_exponential_get(in->tng, &exp);
1188 // GROMACS expects distances in nm
1191 case 9: distanceScaleFactor = NANO / NANO; break;
1192 case 10: distanceScaleFactor = NANO / ANGSTROM; break;
1193 default: distanceScaleFactor = pow(10.0, exp + 9.0);
1196 return distanceScaleFactor;
1202 void gmx_tng_setup_atom_subgroup(gmx_tng_trajectory_t gmx_tng, gmx::ArrayRef<const int> ind, const char* name)
1205 int64_t nAtoms, cnt, nMols;
1206 tng_molecule_t mol, iterMol;
1210 tng_function_status stat;
1211 tng_trajectory_t tng = gmx_tng->tng;
1213 tng_num_particles_get(tng, &nAtoms);
1215 if (nAtoms == ind.ssize())
1220 stat = tng_molecule_find(tng, name, -1, &mol);
1221 if (stat == TNG_SUCCESS)
1223 tng_molecule_num_atoms_get(tng, mol, &nAtoms);
1224 tng_molecule_cnt_get(tng, mol, &cnt);
1225 if (nAtoms == ind.ssize())
1234 if (stat == TNG_FAILURE)
1236 /* The indexed atoms are added to one separate molecule. */
1237 tng_molecule_alloc(tng, &mol);
1238 tng_molecule_name_set(tng, mol, name);
1239 tng_molecule_chain_add(tng, mol, "", &chain);
1241 for (gmx::index i = 0; i < ind.ssize(); i++)
1243 char temp_name[256], temp_type[256];
1245 /* Try to retrieve the residue name of the atom */
1246 stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1247 if (stat != TNG_SUCCESS)
1249 temp_name[0] = '\0';
1251 /* Check if the molecule of the selection already contains this residue */
1252 if (tng_chain_residue_find(tng, chain, temp_name, -1, &res) != TNG_SUCCESS)
1254 tng_chain_residue_add(tng, chain, temp_name, &res);
1256 /* Try to find the original name and type of the atom */
1257 stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1258 if (stat != TNG_SUCCESS)
1260 temp_name[0] = '\0';
1262 stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
1263 if (stat != TNG_SUCCESS)
1265 temp_type[0] = '\0';
1267 tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
1269 tng_molecule_existing_add(tng, &mol);
1271 /* Set the count of the molecule containing the selected atoms to 1 and all
1272 * other molecules to 0 */
1273 tng_molecule_cnt_set(tng, mol, 1);
1274 tng_num_molecule_types_get(tng, &nMols);
1275 for (int64_t k = 0; k < nMols; k++)
1277 tng_molecule_of_index_get(tng, k, &iterMol);
1282 tng_molecule_cnt_set(tng, iterMol, 0);
1285 GMX_UNUSED_VALUE(gmx_tng);
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(gmx_tng_trajectory_t gmx_tng_input,
1297 int64_t* requestedIds,
1298 int numRequestedIds)
1301 tng_trajectory_t input = gmx_tng_input->tng;
1302 gmx_bool bOK = TRUE;
1303 tng_function_status stat;
1304 int64_t numberOfAtoms = -1, frameNumber = -1;
1305 int64_t nBlocks, blockId, *blockIds = nullptr, codecId;
1307 void* values = nullptr;
1308 double frameTime = -1.0;
1309 int size, blockDependency;
1311 const int defaultNumIds = 5;
1312 static int64_t fallbackRequestedIds[defaultNumIds] = { TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
1313 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
1319 fr->bLambda = FALSE;
1327 /* If no specific IDs were requested read all block types that can
1328 * currently be interpreted */
1329 if (!requestedIds || numRequestedIds == 0)
1331 numRequestedIds = defaultNumIds;
1332 requestedIds = fallbackRequestedIds;
1335 stat = tng_num_particles_get(input, &numberOfAtoms);
1336 if (stat != TNG_SUCCESS)
1338 gmx_file("Cannot determine number of atoms from TNG file.");
1340 fr->natoms = numberOfAtoms;
1342 bool nextFrameExists = gmx_get_tng_data_block_types_of_next_frame(
1343 gmx_tng_input, fr->step, numRequestedIds, requestedIds, &frameNumber, &nBlocks, &blockIds);
1344 gmx::unique_cptr<int64_t, gmx::free_wrapper> blockIdsGuard(blockIds);
1345 if (!nextFrameExists)
1355 for (int64_t i = 0; i < nBlocks; i++)
1357 blockId = blockIds[i];
1358 tng_data_block_dependency_get(input, blockId, &blockDependency);
1359 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1361 stat = tng_util_particle_data_next_frame_read(input, blockId, &values, &datatype,
1362 &frameNumber, &frameTime);
1366 stat = tng_util_non_particle_data_next_frame_read(input, blockId, &values, &datatype,
1367 &frameNumber, &frameTime);
1369 if (stat == TNG_CRITICAL)
1371 gmx_file("Cannot read positions from TNG file.");
1374 else if (stat == TNG_FAILURE)
1380 case TNG_TRAJ_BOX_SHAPE:
1383 case TNG_INT_DATA: size = sizeof(int64_t); break;
1384 case TNG_FLOAT_DATA: size = sizeof(float); break;
1385 case TNG_DOUBLE_DATA: size = sizeof(double); break;
1386 default: gmx_incons("Illegal datatype of box shape values!");
1388 for (int i = 0; i < DIM; i++)
1390 convert_array_to_real_array(reinterpret_cast<char*>(values) + size * i * DIM,
1391 reinterpret_cast<real*>(fr->box[i]),
1392 getDistanceScaleFactor(gmx_tng_input), 1, DIM, datatype);
1396 case TNG_TRAJ_POSITIONS:
1397 srenew(fr->x, fr->natoms);
1398 convert_array_to_real_array(values, reinterpret_cast<real*>(fr->x),
1399 getDistanceScaleFactor(gmx_tng_input), fr->natoms, DIM,
1402 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1403 /* This must be updated if/when more lossy compression methods are added */
1404 if (codecId == TNG_TNG_COMPRESSION)
1410 case TNG_TRAJ_VELOCITIES:
1411 srenew(fr->v, fr->natoms);
1412 convert_array_to_real_array(values, reinterpret_cast<real*>(fr->v),
1413 getDistanceScaleFactor(gmx_tng_input), fr->natoms, DIM,
1416 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1417 /* This must be updated if/when more lossy compression methods are added */
1418 if (codecId == TNG_TNG_COMPRESSION)
1424 case TNG_TRAJ_FORCES:
1425 srenew(fr->f, fr->natoms);
1426 convert_array_to_real_array(values, reinterpret_cast<real*>(fr->f),
1427 getDistanceScaleFactor(gmx_tng_input), fr->natoms, DIM,
1431 case TNG_GMX_LAMBDA:
1434 case TNG_FLOAT_DATA: fr->lambda = *(reinterpret_cast<float*>(values)); break;
1435 case TNG_DOUBLE_DATA: fr->lambda = *(reinterpret_cast<double*>(values)); break;
1436 default: gmx_incons("Illegal datatype lambda value!");
1442 "Illegal block type! Currently GROMACS tools can only handle certain data "
1443 "types. Skipping block.");
1445 /* values does not have to be freed before reading next frame. It will
1446 * be reallocated if it is not NULL. */
1449 fr->step = frameNumber;
1452 // Convert the time to ps
1453 fr->time = frameTime / PICO;
1454 fr->bTime = (frameTime > 0);
1456 // TODO This does not leak, but is not exception safe.
1457 /* values must be freed before leaving this function */
1462 GMX_UNUSED_VALUE(gmx_tng_input);
1463 GMX_UNUSED_VALUE(fr);
1464 GMX_UNUSED_VALUE(requestedIds);
1465 GMX_UNUSED_VALUE(numRequestedIds);
1470 void gmx_print_tng_molecule_system(gmx_tng_trajectory_t gmx_tng_input, FILE* stream)
1473 int64_t nMolecules, nChains, nResidues, nAtoms, nFramesRead;
1474 int64_t strideLength, nParticlesRead, nValuesPerFrameRead, *molCntList;
1475 tng_molecule_t molecule;
1477 tng_residue_t residue;
1479 tng_function_status stat;
1483 void* data = nullptr;
1484 std::vector<real> atomCharges;
1485 std::vector<real> atomMasses;
1486 tng_trajectory_t input = gmx_tng_input->tng;
1488 tng_num_molecule_types_get(input, &nMolecules);
1489 tng_molecule_cnt_list_get(input, &molCntList);
1490 /* Can the number of particles change in the trajectory or is it constant? */
1491 tng_num_particles_variable_get(input, &varNAtoms);
1493 for (int64_t i = 0; i < nMolecules; i++)
1495 tng_molecule_of_index_get(input, i, &molecule);
1496 tng_molecule_name_get(input, molecule, str, 256);
1497 if (varNAtoms == TNG_CONSTANT_N_ATOMS)
1499 if (static_cast<int>(molCntList[i]) == 0)
1503 fprintf(stream, "Molecule: %s, count: %d\n", str, static_cast<int>(molCntList[i]));
1507 fprintf(stream, "Molecule: %s\n", str);
1509 tng_molecule_num_chains_get(input, molecule, &nChains);
1512 for (int64_t j = 0; j < nChains; j++)
1514 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
1515 tng_chain_name_get(input, chain, str, 256);
1516 fprintf(stream, "\tChain: %s\n", str);
1517 tng_chain_num_residues_get(input, chain, &nResidues);
1518 for (int64_t k = 0; k < nResidues; k++)
1520 tng_chain_residue_of_index_get(input, chain, k, &residue);
1521 tng_residue_name_get(input, residue, str, 256);
1522 fprintf(stream, "\t\tResidue: %s\n", str);
1523 tng_residue_num_atoms_get(input, residue, &nAtoms);
1524 for (int64_t l = 0; l < nAtoms; l++)
1526 tng_residue_atom_of_index_get(input, residue, l, &atom);
1527 tng_atom_name_get(input, atom, str, 256);
1528 fprintf(stream, "\t\t\tAtom: %s", str);
1529 tng_atom_type_get(input, atom, str, 256);
1530 fprintf(stream, " (%s)\n", str);
1535 /* It is possible to have a molecule without chains, in which case
1536 * residues in the molecule can be iterated through without going
1537 * through chains. */
1540 tng_molecule_num_residues_get(input, molecule, &nResidues);
1543 for (int64_t k = 0; k < nResidues; k++)
1545 tng_molecule_residue_of_index_get(input, molecule, k, &residue);
1546 tng_residue_name_get(input, residue, str, 256);
1547 fprintf(stream, "\t\tResidue: %s\n", str);
1548 tng_residue_num_atoms_get(input, residue, &nAtoms);
1549 for (int64_t l = 0; l < nAtoms; l++)
1551 tng_residue_atom_of_index_get(input, residue, l, &atom);
1552 tng_atom_name_get(input, atom, str, 256);
1553 fprintf(stream, "\t\t\tAtom: %s", str);
1554 tng_atom_type_get(input, atom, str, 256);
1555 fprintf(stream, " (%s)\n", str);
1561 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
1562 for (int64_t l = 0; l < nAtoms; l++)
1564 tng_molecule_atom_of_index_get(input, molecule, l, &atom);
1565 tng_atom_name_get(input, atom, str, 256);
1566 fprintf(stream, "\t\t\tAtom: %s", str);
1567 tng_atom_type_get(input, atom, str, 256);
1568 fprintf(stream, " (%s)\n", str);
1574 tng_num_particles_get(input, &nAtoms);
1575 stat = tng_particle_data_vector_get(input, TNG_TRAJ_PARTIAL_CHARGES, &data, &nFramesRead,
1576 &strideLength, &nParticlesRead, &nValuesPerFrameRead, &datatype);
1577 if (stat == TNG_SUCCESS)
1579 atomCharges.resize(nAtoms);
1580 convert_array_to_real_array(data, atomCharges.data(), 1, nAtoms, 1, datatype);
1582 fprintf(stream, "Atom Charges (%d):\n", int(nAtoms));
1583 for (int64_t i = 0; i < nAtoms; i += 10)
1585 fprintf(stream, "Atom Charges [%8d-]=[", int(i));
1586 for (int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1588 fprintf(stream, " %12.5e", atomCharges[i + j]);
1590 fprintf(stream, "]\n");
1594 stat = tng_particle_data_vector_get(input, TNG_TRAJ_MASSES, &data, &nFramesRead, &strideLength,
1595 &nParticlesRead, &nValuesPerFrameRead, &datatype);
1596 if (stat == TNG_SUCCESS)
1598 atomMasses.resize(nAtoms);
1599 convert_array_to_real_array(data, atomMasses.data(), 1, nAtoms, 1, datatype);
1601 fprintf(stream, "Atom Masses (%d):\n", int(nAtoms));
1602 for (int64_t i = 0; i < nAtoms; i += 10)
1604 fprintf(stream, "Atom Masses [%8d-]=[", int(i));
1605 for (int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1607 fprintf(stream, " %12.5e", atomMasses[i + j]);
1609 fprintf(stream, "]\n");
1615 GMX_UNUSED_VALUE(gmx_tng_input);
1616 GMX_UNUSED_VALUE(stream);
1620 gmx_bool gmx_get_tng_data_block_types_of_next_frame(gmx_tng_trajectory_t gmx_tng_input,
1623 int64_t* requestedIds,
1629 tng_function_status stat;
1630 tng_trajectory_t input = gmx_tng_input->tng;
1632 stat = tng_util_trajectory_next_frame_present_data_blocks_find(
1633 input, frame, nRequestedIds, requestedIds, nextFrame, nBlocks, blockIds);
1635 if (stat == TNG_CRITICAL)
1637 gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
1639 else if (stat == TNG_FAILURE)
1645 GMX_UNUSED_VALUE(gmx_tng_input);
1646 GMX_UNUSED_VALUE(frame);
1647 GMX_UNUSED_VALUE(nRequestedIds);
1648 GMX_UNUSED_VALUE(requestedIds);
1649 GMX_UNUSED_VALUE(nextFrame);
1650 GMX_UNUSED_VALUE(nBlocks);
1651 GMX_UNUSED_VALUE(blockIds);
1656 gmx_bool gmx_get_tng_data_next_frame_of_block_type(gmx_tng_trajectory_t gmx_tng_input,
1659 int64_t* frameNumber,
1661 int64_t* nValuesPerFrame,
1669 tng_function_status stat;
1672 int blockDependency;
1673 void* data = nullptr;
1675 tng_trajectory_t input = gmx_tng_input->tng;
1677 stat = tng_data_block_name_get(input, blockId, name, maxLen);
1678 if (stat != TNG_SUCCESS)
1680 gmx_file("Cannot read next frame of TNG file");
1682 stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
1683 if (stat != TNG_SUCCESS)
1685 gmx_file("Cannot read next frame of TNG file");
1687 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1689 tng_num_particles_get(input, nAtoms);
1690 stat = tng_util_particle_data_next_frame_read(input, blockId, &data, &datatype, frameNumber,
1695 *nAtoms = 1; /* There are not actually any atoms, but it is used for
1696 allocating memory */
1697 stat = tng_util_non_particle_data_next_frame_read(input, blockId, &data, &datatype,
1698 frameNumber, frameTime);
1700 if (stat == TNG_CRITICAL)
1702 gmx_file("Cannot read next frame of TNG file");
1704 if (stat == TNG_FAILURE)
1710 stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
1711 if (stat != TNG_SUCCESS)
1713 gmx_file("Cannot read next frame of TNG file");
1715 srenew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
1716 convert_array_to_real_array(data, *values, getDistanceScaleFactor(gmx_tng_input), *nAtoms,
1717 *nValuesPerFrame, datatype);
1719 tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
1721 /* This must be updated if/when more lossy compression methods are added */
1722 if (codecId != TNG_TNG_COMPRESSION)
1736 GMX_UNUSED_VALUE(gmx_tng_input);
1737 GMX_UNUSED_VALUE(blockId);
1738 GMX_UNUSED_VALUE(values);
1739 GMX_UNUSED_VALUE(frameNumber);
1740 GMX_UNUSED_VALUE(frameTime);
1741 GMX_UNUSED_VALUE(nValuesPerFrame);
1742 GMX_UNUSED_VALUE(nAtoms);
1743 GMX_UNUSED_VALUE(prec);
1744 GMX_UNUSED_VALUE(name);
1745 GMX_UNUSED_VALUE(maxLen);
1746 GMX_UNUSED_VALUE(bOK);
1751 int gmx_tng_get_box_output_interval(gmx_tng_trajectory_t gmx_tng)
1754 return gmx_tng->boxOutputInterval;
1756 GMX_UNUSED_VALUE(gmx_tng);
1761 int gmx_tng_get_lambda_output_interval(gmx_tng_trajectory_t gmx_tng)
1764 return gmx_tng->lambdaOutputInterval;
1766 GMX_UNUSED_VALUE(gmx_tng);