2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2018,2019,2020,2021, 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.
50 # include "tng/tng_io.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/mdtypes/inputrec.h"
56 #include "gromacs/topology/ifunc.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/trajectory/trajectoryframe.h"
59 #include "gromacs/utility/basedefinitions.h"
60 #include "gromacs/utility/baseversion.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/futil.h"
63 #include "gromacs/utility/gmxassert.h"
64 #include "gromacs/utility/programcontext.h"
65 #include "gromacs/utility/smalloc.h"
66 #include "gromacs/utility/sysinfo.h"
67 #include "gromacs/utility/unique_cptr.h"
70 using tng_trajectory_t = void*;
73 /*! \brief Gromacs Wrapper around tng datatype
75 * This could in principle hold any GROMACS-specific requirements not yet
76 * implemented in or not relevant to the TNG library itself. However, for now
77 * we only use it to handle some shortcomings we have discovered, where the TNG
78 * API itself is a bit fragile and can end up overwriting data if called several
79 * times with the same frame number.
80 * The logic to determine the time per step was also a bit fragile. This is not
81 * critical, but since we anyway need a wrapper for ensuring unique frame
82 * numbers, we can also use it to store the time of the first step and use that
83 * to derive a slightly better/safer estimate of the time per step.
85 * At some future point where we have a second-generation TNG API we should
86 * consider removing this again.
88 struct gmx_tng_trajectory
90 tng_trajectory_t tng; //!< Actual TNG handle (pointer)
91 bool lastStepDataIsValid; //!< True if lastStep has been set
92 std::int64_t lastStep; //!< Index/step used for last frame
93 bool lastTimeDataIsValid; //!< True if lastTime has been set
94 double lastTime; //!< Time of last frame (TNG unit is seconds)
95 bool timePerFrameIsSet; //!< True if we have set the time per frame
96 int boxOutputInterval; //!< Number of steps between the output of box size
97 int lambdaOutputInterval; //!< Number of steps between the output of lambdas
101 static const char* modeToVerb(char mode)
106 case 'r': p = "reading"; break;
107 case 'w': p = "writing"; break;
108 case 'a': p = "appending"; break;
109 default: gmx_fatal(FARGS, "Invalid file opening mode %c", mode);
115 void gmx_tng_open(const char* filename, char mode, gmx_tng_trajectory_t* gmx_tng)
118 /* First check whether we have to make a backup,
119 * only for writing, not for read or append.
123 make_backup(filename);
126 *gmx_tng = new gmx_tng_trajectory;
127 (*gmx_tng)->lastStepDataIsValid = false;
128 (*gmx_tng)->lastTimeDataIsValid = false;
129 (*gmx_tng)->timePerFrameIsSet = false;
130 tng_trajectory_t* tng = &(*gmx_tng)->tng;
132 /* tng must not be pointing at already allocated memory.
133 * Memory will be allocated by tng_util_trajectory_open() and must
134 * later on be freed by tng_util_trajectory_close(). */
135 if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
137 /* TNG does return more than one degree of error, but there is
138 no use case for GROMACS handling the non-fatal errors
140 gmx_fatal(FARGS, "File I/O error while opening %s for %s", filename, modeToVerb(mode));
143 if (mode == 'w' || mode == 'a')
146 gmx_gethostname(hostname, 256);
149 tng_first_computer_name_set(*tng, hostname);
153 tng_last_computer_name_set(*tng, hostname);
156 char programInfo[256];
157 const char* precisionString = "";
159 precisionString = " (double precision)";
161 sprintf(programInfo, "%.100s %.128s%.24s", gmx::getProgramContext().displayName(), 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(
269 tng, tngRes, *(atoms->atomname[atomIndex]), *(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,
346 std::back_inserter(atomCharges));
347 std::copy_n(atomMasses.end() - molType->atoms.nr, molType->atoms.nr, std::back_inserter(atomMasses));
350 /* Write the TNG data blocks. */
351 tng_particle_data_block_add(tng,
352 TNG_TRAJ_PARTIAL_CHARGES,
355 TNG_NON_TRAJECTORY_BLOCK,
361 TNG_GZIP_COMPRESSION,
363 tng_particle_data_block_add(tng,
367 TNG_NON_TRAJECTORY_BLOCK,
373 TNG_GZIP_COMPRESSION,
377 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
378 * if they are positive.
380 * If only one of n1 and n2 is positive, then return it.
381 * If neither n1 or n2 is positive, then return -1. */
382 static int greatest_common_divisor_if_positive(int n1, int n2)
386 return (0 >= n2) ? -1 : n2;
393 /* We have a non-trivial greatest common divisor to compute. */
394 return std::gcd(n1, n2);
397 /* By default try to write 100 frames (of actual output) in each frame set.
398 * This number is the number of outputs of the most frequently written data
399 * type per frame set.
400 * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
401 * setups regarding compression efficiency and compression time. Make this
402 * a hidden command-line option? */
403 const int defaultFramesPerFrameSet = 100;
405 /*! \libinternal \brief Set the number of frames per frame
406 * set according to output intervals.
407 * The default is that 100 frames are written of the data
408 * that is written most often. */
409 static void tng_set_frames_per_frame_set(gmx_tng_trajectory_t gmx_tng,
410 const gmx_bool bUseLossyCompression,
411 const t_inputrec* ir)
414 tng_trajectory_t tng = gmx_tng->tng;
416 /* Set the number of frames per frame set to contain at least
417 * defaultFramesPerFrameSet of the lowest common denominator of
418 * the writing interval of positions and velocities. */
419 /* FIXME after 5.0: consider nstenergy also? */
420 if (bUseLossyCompression)
422 gcd = ir->nstxout_compressed;
426 gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
427 gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
434 tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
437 /*! \libinternal \brief Set the data-writing intervals, and number of
438 * frames per frame set */
439 static void set_writing_intervals(gmx_tng_trajectory_t gmx_tng,
440 const gmx_bool bUseLossyCompression,
441 const t_inputrec* ir)
443 tng_trajectory_t tng = gmx_tng->tng;
445 /* Define pointers to specific writing functions depending on if we
446 * write float or double data */
447 typedef tng_function_status (*set_writing_interval_func_pointer)(
448 tng_trajectory_t, const int64_t, const int64_t, const int64_t, const char*, const char, const char);
450 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
452 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
454 int xout, vout, fout;
455 int gcd = -1, lowest = -1;
458 tng_set_frames_per_frame_set(gmx_tng, bUseLossyCompression, ir);
460 if (bUseLossyCompression)
462 xout = ir->nstxout_compressed;
464 /* If there is no uncompressed coordinate output write forces
465 and velocities to the compressed tng file. */
476 compression = TNG_TNG_COMPRESSION;
483 compression = TNG_GZIP_COMPRESSION;
487 set_writing_interval(
488 tng, xout, 3, TNG_TRAJ_POSITIONS, "POSITIONS", TNG_PARTICLE_BLOCK_DATA, compression);
489 /* TODO: if/when we write energies to TNG also, reconsider how
490 * and when box information is written, because GROMACS
491 * behaviour pre-5.0 was to write the box with every
492 * trajectory frame and every energy frame, and probably
493 * people depend on this. */
495 gcd = greatest_common_divisor_if_positive(gcd, xout);
496 if (lowest < 0 || xout < lowest)
503 set_writing_interval(
504 tng, vout, 3, TNG_TRAJ_VELOCITIES, "VELOCITIES", TNG_PARTICLE_BLOCK_DATA, compression);
506 gcd = greatest_common_divisor_if_positive(gcd, vout);
507 if (lowest < 0 || vout < lowest)
514 set_writing_interval(
515 tng, fout, 3, TNG_TRAJ_FORCES, "FORCES", TNG_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION);
517 gcd = greatest_common_divisor_if_positive(gcd, fout);
518 if (lowest < 0 || fout < lowest)
525 /* Lambdas and box shape written at an interval of the lowest common
526 denominator of other output */
527 set_writing_interval(
528 tng, gcd, 1, TNG_GMX_LAMBDA, "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION);
530 set_writing_interval(
531 tng, gcd, 9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION);
532 gmx_tng->lambdaOutputInterval = gcd;
533 gmx_tng->boxOutputInterval = gcd;
534 if (gcd < lowest / 10)
537 "The lowest common denominator of trajectory output is "
538 "every %d step(s), whereas the shortest output interval "
539 "is every %d steps.",
547 void gmx_tng_prepare_md_writing(gmx_tng_trajectory_t gmx_tng, const gmx_mtop_t* mtop, const t_inputrec* ir)
550 gmx_tng_add_mtop(gmx_tng, mtop);
551 set_writing_intervals(gmx_tng, FALSE, ir);
552 tng_time_per_frame_set(gmx_tng->tng, ir->delta_t * gmx::c_pico);
553 gmx_tng->timePerFrameIsSet = true;
555 GMX_UNUSED_VALUE(gmx_tng);
556 GMX_UNUSED_VALUE(mtop);
557 GMX_UNUSED_VALUE(ir);
562 /* Check if all atoms in the molecule system are selected
563 * by a selection group of type specified by gtype. */
564 static gmx_bool all_atoms_selected(const gmx_mtop_t* mtop, const SimulationAtomGroupType gtype)
566 /* Iterate through all atoms in the molecule system and
567 * check if they belong to a selection group of the
570 for (const gmx_molblock_t& molBlock : mtop->molblock)
572 const gmx_moltype_t& molType = mtop->moltype[molBlock.type];
573 const t_atoms& atoms = molType.atoms;
575 for (int j = 0; j < molBlock.nmol; j++)
577 for (int atomIndex = 0; atomIndex < atoms.nr; atomIndex++, i++)
579 if (getGroupType(mtop->groups, gtype, i) != 0)
589 /* Create TNG molecules which will represent each of the selection
590 * groups for writing. But do that only if there is actually a
591 * specified selection group and it is not the whole system.
592 * TODO: Currently the only selection that is taken into account
593 * is egcCompressedX, but other selections should be added when
594 * e.g. writing energies is implemented.
596 static void add_selection_groups(gmx_tng_trajectory_t gmx_tng, const gmx_mtop_t* mtop)
598 const t_atoms* atoms;
600 const t_resinfo* resInfo;
603 tng_molecule_t mol, iterMol;
610 tng_trajectory_t tng = gmx_tng->tng;
612 /* TODO: When the TNG molecules block is more flexible TNG selection
613 * groups should not need all atoms specified. It should be possible
614 * just to specify what molecules are selected (and/or which atoms in
615 * the molecule) and how many (if applicable). */
617 /* If no atoms are selected we do not need to create a
618 * TNG selection group molecule. */
619 if (mtop->groups.numberOfGroupNumbers(SimulationAtomGroupType::CompressedPositionOutput) == 0)
624 /* If all atoms are selected we do not have to create a selection
625 * group molecule in the TNG molecule system. */
626 if (all_atoms_selected(mtop, SimulationAtomGroupType::CompressedPositionOutput))
631 /* The name of the TNG molecule containing the selection group is the
632 * same as the name of the selection group. */
633 nameIndex = mtop->groups.groups[SimulationAtomGroupType::CompressedPositionOutput][0];
634 groupName = *mtop->groups.groupNames[nameIndex];
636 tng_molecule_alloc(tng, &mol);
637 tng_molecule_name_set(tng, mol, groupName);
638 tng_molecule_chain_add(tng, mol, "", &chain);
640 for (const gmx_molblock_t& molBlock : mtop->molblock)
642 const gmx_moltype_t& molType = mtop->moltype[molBlock.type];
644 atoms = &molType.atoms;
646 for (int j = 0; j < molBlock.nmol; j++)
648 bool bAtomsAdded = FALSE;
649 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
654 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, i) != 0)
658 at = &atoms->atom[atomIndex];
661 resInfo = &atoms->resinfo[at->resind];
662 /* FIXME: When TNG supports both residue index and residue
663 * number the latter should be used. */
664 res_name = *resInfo->name;
665 res_id = at->resind + 1;
669 res_name = const_cast<char*>("");
672 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res) != TNG_SUCCESS)
674 /* Since there is ONE chain for selection groups do not keep the
675 * original residue IDs - otherwise there might be conflicts. */
676 tng_chain_residue_add(tng, chain, res_name, &res);
678 tng_residue_atom_w_id_add(tng,
680 *(atoms->atomname[atomIndex]),
681 *(atoms->atomtype[atomIndex]),
682 atom_offset + atomIndex,
689 for (int k = 0; k < F_NRE; k++)
693 const InteractionList& ilist = molType.ilist[k];
694 for (int l = 1; l < ilist.size(); l += 3)
697 atom1 = ilist.iatoms[l] + atom_offset;
698 atom2 = ilist.iatoms[l + 1] + atom_offset;
699 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, atom1)
701 && getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, atom2)
704 tng_molecule_bond_add(
705 tng, mol, ilist.iatoms[l], ilist.iatoms[l + 1], &tngBond);
710 /* Settle is described using three atoms */
711 const InteractionList& ilist = molType.ilist[F_SETTLE];
712 for (int l = 1; l < ilist.size(); l += 4)
714 int atom1, atom2, atom3;
715 atom1 = ilist.iatoms[l] + atom_offset;
716 atom2 = ilist.iatoms[l + 1] + atom_offset;
717 atom3 = ilist.iatoms[l + 2] + atom_offset;
718 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, atom1)
721 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, atom2)
724 tng_molecule_bond_add(tng, mol, atom1, atom2, &tngBond);
726 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, atom3)
729 tng_molecule_bond_add(tng, mol, atom1, atom3, &tngBond);
734 atom_offset += atoms->nr;
737 tng_molecule_existing_add(tng, &mol);
738 tng_molecule_cnt_set(tng, mol, 1);
739 tng_num_molecule_types_get(tng, &nMols);
740 for (int64_t k = 0; k < nMols; k++)
742 tng_molecule_of_index_get(tng, k, &iterMol);
747 tng_molecule_cnt_set(tng, iterMol, 0);
752 void gmx_tng_set_compression_precision(gmx_tng_trajectory_t gmx_tng, real prec)
755 tng_compression_precision_set(gmx_tng->tng, prec);
757 GMX_UNUSED_VALUE(gmx_tng);
758 GMX_UNUSED_VALUE(prec);
762 void gmx_tng_prepare_low_prec_writing(gmx_tng_trajectory_t gmx_tng, const gmx_mtop_t* mtop, const t_inputrec* ir)
765 gmx_tng_add_mtop(gmx_tng, mtop);
766 add_selection_groups(gmx_tng, mtop);
767 set_writing_intervals(gmx_tng, TRUE, ir);
768 tng_time_per_frame_set(gmx_tng->tng, ir->delta_t * gmx::c_pico);
769 gmx_tng->timePerFrameIsSet = true;
770 gmx_tng_set_compression_precision(gmx_tng, ir->x_compression_precision);
772 GMX_UNUSED_VALUE(gmx_tng);
773 GMX_UNUSED_VALUE(mtop);
774 GMX_UNUSED_VALUE(ir);
778 void gmx_fwrite_tng(gmx_tng_trajectory_t gmx_tng,
779 const gmx_bool bUseLossyCompression,
781 real elapsedPicoSeconds,
790 typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
800 static write_data_func_pointer write_data = tng_util_generic_with_time_double_write;
802 static write_data_func_pointer write_data = tng_util_generic_with_time_write;
804 double elapsedSeconds = elapsedPicoSeconds * gmx::c_pico;
811 /* This function might get called when the type of the
812 compressed trajectory is actually XTC. So we exit and move
816 tng_trajectory_t tng = gmx_tng->tng;
818 // While the GROMACS interface to this routine specifies 'step', TNG itself
819 // only uses 'frame index' internally, although it suggests that it's a good
820 // idea to let this represent the step rather than just frame numbers.
822 // However, the way the GROMACS interface works, there are cases where
823 // the step is simply not set, so all frames rather have step=0.
824 // If we call the lower-level TNG routines with the same frame index
825 // (which is set from the step) multiple times, things will go very
826 // wrong (overwritten data).
828 // To avoid this, we require the step value to always be larger than the
829 // previous value, and if this is not true we simply set it to a value
830 // one unit larger than the previous step.
832 // Note: We do allow the user to specify any crazy value the want for the
833 // first step. If it is "not set", GROMACS will have used the value 0.
834 if (gmx_tng->lastStepDataIsValid && step <= gmx_tng->lastStep)
836 step = gmx_tng->lastStep + 1;
839 // Now that we have fixed the potentially incorrect step, we can also set
840 // the time per frame if it was not already done, and we have
841 // step/time information from the last frame so it is possible to calculate it.
842 // Note that we do allow the time to be the same for all steps, or even
843 // decreasing. In that case we will simply say the time per step is 0
844 // or negative. A bit strange, but a correct representation of the data :-)
845 if (!gmx_tng->timePerFrameIsSet && gmx_tng->lastTimeDataIsValid && gmx_tng->lastStepDataIsValid)
847 double deltaTime = elapsedSeconds - gmx_tng->lastTime;
848 std::int64_t deltaStep = step - gmx_tng->lastStep;
849 tng_time_per_frame_set(tng, deltaTime / deltaStep);
850 gmx_tng->timePerFrameIsSet = true;
853 tng_num_particles_get(tng, &nParticles);
854 if (nAtoms != static_cast<int>(nParticles))
856 tng_implicit_num_particles_set(tng, nAtoms);
859 if (bUseLossyCompression)
861 compression = TNG_TNG_COMPRESSION;
865 compression = TNG_GZIP_COMPRESSION;
868 /* The writing is done using write_data, which writes float or double
869 * depending on the GROMACS compilation. */
872 GMX_ASSERT(box, "Need a non-NULL box if positions are written");
877 reinterpret_cast<const real*>(x),
881 TNG_PARTICLE_BLOCK_DATA,
885 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
894 reinterpret_cast<const real*>(v),
898 TNG_PARTICLE_BLOCK_DATA,
902 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
908 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
909 * compression for forces regardless of output mode */
913 reinterpret_cast<const real*>(f),
917 TNG_PARTICLE_BLOCK_DATA,
918 TNG_GZIP_COMPRESSION)
921 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
927 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
928 * compression for box shape regardless of output mode */
932 reinterpret_cast<const real*>(box),
936 TNG_NON_PARTICLE_BLOCK_DATA,
937 TNG_GZIP_COMPRESSION)
940 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
946 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
947 * compression for lambda regardless of output mode */
951 reinterpret_cast<const real*>(&lambda),
955 TNG_NON_PARTICLE_BLOCK_DATA,
956 TNG_GZIP_COMPRESSION)
959 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
963 // Update the data in the wrapper
965 gmx_tng->lastStepDataIsValid = true;
966 gmx_tng->lastStep = step;
967 gmx_tng->lastTimeDataIsValid = true;
968 gmx_tng->lastTime = elapsedSeconds;
970 GMX_UNUSED_VALUE(gmx_tng);
971 GMX_UNUSED_VALUE(bUseLossyCompression);
972 GMX_UNUSED_VALUE(step);
973 GMX_UNUSED_VALUE(elapsedPicoSeconds);
974 GMX_UNUSED_VALUE(lambda);
975 GMX_UNUSED_VALUE(box);
976 GMX_UNUSED_VALUE(nAtoms);
983 void fflush_tng(gmx_tng_trajectory_t gmx_tng)
990 tng_frame_set_premature_write(gmx_tng->tng, TNG_USE_HASH);
992 GMX_UNUSED_VALUE(gmx_tng);
996 float gmx_tng_get_time_of_final_frame(gmx_tng_trajectory_t gmx_tng)
1002 tng_trajectory_t tng = gmx_tng->tng;
1004 tng_num_frames_get(tng, &nFrames);
1005 tng_util_time_of_frame_get(tng, nFrames - 1, &time);
1007 fTime = time / gmx::c_pico;
1010 GMX_UNUSED_VALUE(gmx_tng);
1015 void gmx_prepare_tng_writing(const char* filename,
1017 gmx_tng_trajectory_t* gmx_tng_input,
1018 gmx_tng_trajectory_t* gmx_tng_output,
1020 const gmx_mtop_t* mtop,
1021 gmx::ArrayRef<const int> index,
1022 const char* indexGroupName)
1025 tng_trajectory_t* input = (gmx_tng_input && *gmx_tng_input) ? &(*gmx_tng_input)->tng : nullptr;
1026 /* FIXME after 5.0: Currently only standard block types are read */
1027 const int defaultNumIds = 5;
1028 static int64_t fallbackIds[defaultNumIds] = {
1029 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS, TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES, TNG_GMX_LAMBDA
1031 static char fallbackNames[defaultNumIds][32] = {
1032 "BOX SHAPE", "POSITIONS", "VELOCITIES", "FORCES", "LAMBDAS"
1035 typedef tng_function_status (*set_writing_interval_func_pointer)(
1036 tng_trajectory_t, const int64_t, const int64_t, const int64_t, const char*, const char, const char);
1038 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
1040 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
1043 gmx_tng_open(filename, mode, gmx_tng_output);
1044 tng_trajectory_t* output = &(*gmx_tng_output)->tng;
1046 /* Do we have an input file in TNG format? If so, then there's
1047 more data we can copy over, rather than having to improvise. */
1048 if (gmx_tng_input && *gmx_tng_input)
1050 /* Set parameters (compression, time per frame, molecule
1051 * information, number of frames per frame set and writing
1052 * intervals of positions, box shape and lambdas) of the
1053 * output tng container based on their respective values int
1054 * the input tng container */
1055 double time, compression_precision;
1056 int64_t n_frames_per_frame_set, interval = -1;
1058 tng_compression_precision_get(*input, &compression_precision);
1059 tng_compression_precision_set(*output, compression_precision);
1060 // TODO make this configurable in a future version
1061 char compression_type = TNG_TNG_COMPRESSION;
1063 tng_molecule_system_copy(*input, *output);
1065 tng_time_per_frame_get(*input, &time);
1066 tng_time_per_frame_set(*output, time);
1067 // Since we have copied the value from the input TNG we should not change it again
1068 (*gmx_tng_output)->timePerFrameIsSet = true;
1070 tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
1071 tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
1073 for (int i = 0; i < defaultNumIds; i++)
1075 if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval) == TNG_SUCCESS)
1077 switch (fallbackIds[i])
1079 case TNG_TRAJ_POSITIONS:
1080 case TNG_TRAJ_VELOCITIES:
1081 set_writing_interval(*output,
1086 TNG_PARTICLE_BLOCK_DATA,
1089 case TNG_TRAJ_FORCES:
1090 set_writing_interval(*output,
1095 TNG_PARTICLE_BLOCK_DATA,
1096 TNG_GZIP_COMPRESSION);
1098 case TNG_TRAJ_BOX_SHAPE:
1099 set_writing_interval(*output,
1104 TNG_NON_PARTICLE_BLOCK_DATA,
1105 TNG_GZIP_COMPRESSION);
1106 (*gmx_tng_output)->boxOutputInterval = interval;
1108 case TNG_GMX_LAMBDA:
1109 set_writing_interval(*output,
1114 TNG_NON_PARTICLE_BLOCK_DATA,
1115 TNG_GZIP_COMPRESSION);
1116 (*gmx_tng_output)->lambdaOutputInterval = interval;
1125 /* TODO after trjconv is modularized: fix this so the user can
1126 change precision when they are doing an operation where
1127 this makes sense, and not otherwise.
1129 char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
1130 gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
1132 gmx_tng_add_mtop(*gmx_tng_output, mtop);
1133 tng_num_frames_per_frame_set_set(*output, 1);
1136 if ((!index.empty()) && nAtoms > 0)
1138 gmx_tng_setup_atom_subgroup(*gmx_tng_output, index, indexGroupName);
1141 /* If for some reason there are more requested atoms than there are atoms in the
1142 * molecular system create a number of implicit atoms (without atom data) to
1143 * compensate for that. */
1146 tng_implicit_num_particles_set(*output, nAtoms);
1149 GMX_UNUSED_VALUE(filename);
1150 GMX_UNUSED_VALUE(mode);
1151 GMX_UNUSED_VALUE(gmx_tng_input);
1152 GMX_UNUSED_VALUE(gmx_tng_output);
1153 GMX_UNUSED_VALUE(nAtoms);
1154 GMX_UNUSED_VALUE(mtop);
1155 GMX_UNUSED_VALUE(index);
1156 GMX_UNUSED_VALUE(indexGroupName);
1160 void gmx_write_tng_from_trxframe(gmx_tng_trajectory_t gmx_tng_output, const t_trxframe* frame, int natoms)
1165 natoms = frame->natoms;
1168 gmx_tng_output, TRUE, frame->step, frame->time, 0, frame->box, natoms, frame->x, frame->v, frame->f);
1170 GMX_UNUSED_VALUE(gmx_tng_output);
1171 GMX_UNUSED_VALUE(frame);
1172 GMX_UNUSED_VALUE(natoms);
1180 void convert_array_to_real_array(void* from,
1185 const char datatype)
1189 const bool useDouble = GMX_DOUBLE;
1192 case TNG_FLOAT_DATA:
1197 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1201 for (i = 0; i < nAtoms; i++)
1203 for (j = 0; j < nValues; j++)
1205 to[i * nValues + j] = reinterpret_cast<float*>(from)[i * nValues + j] * fact;
1212 for (i = 0; i < nAtoms; i++)
1214 for (j = 0; j < nValues; j++)
1216 to[i * nValues + j] = reinterpret_cast<float*>(from)[i * nValues + j] * fact;
1222 for (i = 0; i < nAtoms; i++)
1224 for (j = 0; j < nValues; j++)
1226 to[i * nValues + j] = reinterpret_cast<int64_t*>(from)[i * nValues + j] * fact;
1230 case TNG_DOUBLE_DATA:
1231 if (sizeof(real) == sizeof(double))
1235 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1239 for (i = 0; i < nAtoms; i++)
1241 for (j = 0; j < nValues; j++)
1243 to[i * nValues + j] = reinterpret_cast<double*>(from)[i * nValues + j] * fact;
1250 for (i = 0; i < nAtoms; i++)
1252 for (j = 0; j < nValues; j++)
1254 to[i * nValues + j] = reinterpret_cast<double*>(from)[i * nValues + j] * fact;
1259 default: gmx_incons("Illegal datatype when converting values to a real array!");
1263 real getDistanceScaleFactor(gmx_tng_trajectory_t in)
1266 real distanceScaleFactor;
1268 // TODO Hopefully, TNG 2.0 will do this kind of thing for us
1269 tng_distance_unit_exponential_get(in->tng, &exp);
1271 // GROMACS expects distances in nm
1274 case 9: distanceScaleFactor = gmx::c_nano / gmx::c_nano; break;
1275 case 10: distanceScaleFactor = gmx::c_nano / gmx::c_angstrom; break;
1276 default: distanceScaleFactor = pow(10.0, exp + 9.0);
1279 return distanceScaleFactor;
1285 void gmx_tng_setup_atom_subgroup(gmx_tng_trajectory_t gmx_tng, gmx::ArrayRef<const int> ind, const char* name)
1288 int64_t nAtoms, cnt, nMols;
1289 tng_molecule_t mol, iterMol;
1293 tng_function_status stat;
1294 tng_trajectory_t tng = gmx_tng->tng;
1296 tng_num_particles_get(tng, &nAtoms);
1298 if (nAtoms == ind.ssize())
1303 stat = tng_molecule_find(tng, name, -1, &mol);
1304 if (stat == TNG_SUCCESS)
1306 tng_molecule_num_atoms_get(tng, mol, &nAtoms);
1307 tng_molecule_cnt_get(tng, mol, &cnt);
1308 if (nAtoms == ind.ssize())
1317 if (stat == TNG_FAILURE)
1319 /* The indexed atoms are added to one separate molecule. */
1320 tng_molecule_alloc(tng, &mol);
1321 tng_molecule_name_set(tng, mol, name);
1322 tng_molecule_chain_add(tng, mol, "", &chain);
1324 for (gmx::index i = 0; i < ind.ssize(); i++)
1326 char temp_name[256], temp_type[256];
1328 /* Try to retrieve the residue name of the atom */
1329 stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1330 if (stat != TNG_SUCCESS)
1332 temp_name[0] = '\0';
1334 /* Check if the molecule of the selection already contains this residue */
1335 if (tng_chain_residue_find(tng, chain, temp_name, -1, &res) != TNG_SUCCESS)
1337 tng_chain_residue_add(tng, chain, temp_name, &res);
1339 /* Try to find the original name and type of the atom */
1340 stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1341 if (stat != TNG_SUCCESS)
1343 temp_name[0] = '\0';
1345 stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
1346 if (stat != TNG_SUCCESS)
1348 temp_type[0] = '\0';
1350 tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
1352 tng_molecule_existing_add(tng, &mol);
1354 /* Set the count of the molecule containing the selected atoms to 1 and all
1355 * other molecules to 0 */
1356 tng_molecule_cnt_set(tng, mol, 1);
1357 tng_num_molecule_types_get(tng, &nMols);
1358 for (int64_t k = 0; k < nMols; k++)
1360 tng_molecule_of_index_get(tng, k, &iterMol);
1365 tng_molecule_cnt_set(tng, iterMol, 0);
1368 GMX_UNUSED_VALUE(gmx_tng);
1369 GMX_UNUSED_VALUE(ind);
1370 GMX_UNUSED_VALUE(name);
1374 /* TODO: If/when TNG acquires the ability to copy data blocks without
1375 * uncompressing them, then this implemenation should be reconsidered.
1376 * Ideally, gmx trjconv -f a.tng -o b.tng -b 10 -e 20 would be fast
1377 * and lose no information. */
1378 gmx_bool gmx_read_next_tng_frame(gmx_tng_trajectory_t gmx_tng_input,
1380 int64_t* requestedIds,
1381 int numRequestedIds)
1384 tng_trajectory_t input = gmx_tng_input->tng;
1385 gmx_bool bOK = TRUE;
1386 tng_function_status stat;
1387 int64_t numberOfAtoms = -1, frameNumber = -1;
1388 int64_t nBlocks, blockId, *blockIds = nullptr, codecId;
1390 void* values = nullptr;
1391 double frameTime = -1.0;
1392 int size, blockDependency;
1394 const int defaultNumIds = 5;
1395 static int64_t fallbackRequestedIds[defaultNumIds] = {
1396 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS, TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES, TNG_GMX_LAMBDA
1402 fr->bLambda = FALSE;
1410 /* If no specific IDs were requested read all block types that can
1411 * currently be interpreted */
1412 if (!requestedIds || numRequestedIds == 0)
1414 numRequestedIds = defaultNumIds;
1415 requestedIds = fallbackRequestedIds;
1418 stat = tng_num_particles_get(input, &numberOfAtoms);
1419 if (stat != TNG_SUCCESS)
1421 gmx_file("Cannot determine number of atoms from TNG file.");
1423 fr->natoms = numberOfAtoms;
1425 bool nextFrameExists = gmx_get_tng_data_block_types_of_next_frame(
1426 gmx_tng_input, fr->step, numRequestedIds, requestedIds, &frameNumber, &nBlocks, &blockIds);
1427 gmx::unique_cptr<int64_t, gmx::free_wrapper> blockIdsGuard(blockIds);
1428 if (!nextFrameExists)
1438 for (int64_t i = 0; i < nBlocks; i++)
1440 blockId = blockIds[i];
1441 tng_data_block_dependency_get(input, blockId, &blockDependency);
1442 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1444 stat = tng_util_particle_data_next_frame_read(
1445 input, blockId, &values, &datatype, &frameNumber, &frameTime);
1449 stat = tng_util_non_particle_data_next_frame_read(
1450 input, blockId, &values, &datatype, &frameNumber, &frameTime);
1452 if (stat == TNG_CRITICAL)
1454 gmx_file("Cannot read positions from TNG file.");
1457 else if (stat == TNG_FAILURE)
1463 case TNG_TRAJ_BOX_SHAPE:
1466 case TNG_INT_DATA: size = sizeof(int64_t); break;
1467 case TNG_FLOAT_DATA: size = sizeof(float); break;
1468 case TNG_DOUBLE_DATA: size = sizeof(double); break;
1469 default: gmx_incons("Illegal datatype of box shape values!");
1471 for (int i = 0; i < DIM; i++)
1473 convert_array_to_real_array(reinterpret_cast<char*>(values) + size * i * DIM,
1474 reinterpret_cast<real*>(fr->box[i]),
1475 getDistanceScaleFactor(gmx_tng_input),
1482 case TNG_TRAJ_POSITIONS:
1483 srenew(fr->x, fr->natoms);
1484 convert_array_to_real_array(values,
1485 reinterpret_cast<real*>(fr->x),
1486 getDistanceScaleFactor(gmx_tng_input),
1491 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1492 /* This must be updated if/when more lossy compression methods are added */
1493 if (codecId == TNG_TNG_COMPRESSION)
1499 case TNG_TRAJ_VELOCITIES:
1500 srenew(fr->v, fr->natoms);
1501 convert_array_to_real_array(values,
1502 reinterpret_cast<real*>(fr->v),
1503 getDistanceScaleFactor(gmx_tng_input),
1508 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1509 /* This must be updated if/when more lossy compression methods are added */
1510 if (codecId == TNG_TNG_COMPRESSION)
1516 case TNG_TRAJ_FORCES:
1517 srenew(fr->f, fr->natoms);
1518 convert_array_to_real_array(values,
1519 reinterpret_cast<real*>(fr->f),
1520 getDistanceScaleFactor(gmx_tng_input),
1526 case TNG_GMX_LAMBDA:
1529 case TNG_FLOAT_DATA: fr->lambda = *(reinterpret_cast<float*>(values)); break;
1530 case TNG_DOUBLE_DATA: fr->lambda = *(reinterpret_cast<double*>(values)); break;
1531 default: gmx_incons("Illegal datatype lambda value!");
1537 "Illegal block type! Currently GROMACS tools can only handle certain data "
1538 "types. Skipping block.");
1540 /* values does not have to be freed before reading next frame. It will
1541 * be reallocated if it is not NULL. */
1544 fr->step = frameNumber;
1547 // Convert the time to ps
1548 fr->time = frameTime / gmx::c_pico;
1549 fr->bTime = (frameTime > 0);
1551 // TODO This does not leak, but is not exception safe.
1552 /* values must be freed before leaving this function */
1557 GMX_UNUSED_VALUE(gmx_tng_input);
1558 GMX_UNUSED_VALUE(fr);
1559 GMX_UNUSED_VALUE(requestedIds);
1560 GMX_UNUSED_VALUE(numRequestedIds);
1565 void gmx_print_tng_molecule_system(gmx_tng_trajectory_t gmx_tng_input, FILE* stream)
1568 int64_t nMolecules, nChains, nResidues, nAtoms, nFramesRead;
1569 int64_t strideLength, nParticlesRead, nValuesPerFrameRead, *molCntList;
1570 tng_molecule_t molecule;
1572 tng_residue_t residue;
1574 tng_function_status stat;
1578 void* data = nullptr;
1579 std::vector<real> atomCharges;
1580 std::vector<real> atomMasses;
1581 tng_trajectory_t input = gmx_tng_input->tng;
1583 tng_num_molecule_types_get(input, &nMolecules);
1584 tng_molecule_cnt_list_get(input, &molCntList);
1585 /* Can the number of particles change in the trajectory or is it constant? */
1586 tng_num_particles_variable_get(input, &varNAtoms);
1588 for (int64_t i = 0; i < nMolecules; i++)
1590 tng_molecule_of_index_get(input, i, &molecule);
1591 tng_molecule_name_get(input, molecule, str, 256);
1592 if (varNAtoms == TNG_CONSTANT_N_ATOMS)
1594 if (static_cast<int>(molCntList[i]) == 0)
1598 fprintf(stream, "Molecule: %s, count: %d\n", str, static_cast<int>(molCntList[i]));
1602 fprintf(stream, "Molecule: %s\n", str);
1604 tng_molecule_num_chains_get(input, molecule, &nChains);
1607 for (int64_t j = 0; j < nChains; j++)
1609 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
1610 tng_chain_name_get(input, chain, str, 256);
1611 fprintf(stream, "\tChain: %s\n", str);
1612 tng_chain_num_residues_get(input, chain, &nResidues);
1613 for (int64_t k = 0; k < nResidues; k++)
1615 tng_chain_residue_of_index_get(input, chain, k, &residue);
1616 tng_residue_name_get(input, residue, str, 256);
1617 fprintf(stream, "\t\tResidue: %s\n", str);
1618 tng_residue_num_atoms_get(input, residue, &nAtoms);
1619 for (int64_t l = 0; l < nAtoms; l++)
1621 tng_residue_atom_of_index_get(input, residue, l, &atom);
1622 tng_atom_name_get(input, atom, str, 256);
1623 fprintf(stream, "\t\t\tAtom: %s", str);
1624 tng_atom_type_get(input, atom, str, 256);
1625 fprintf(stream, " (%s)\n", str);
1630 /* It is possible to have a molecule without chains, in which case
1631 * residues in the molecule can be iterated through without going
1632 * through chains. */
1635 tng_molecule_num_residues_get(input, molecule, &nResidues);
1638 for (int64_t k = 0; k < nResidues; k++)
1640 tng_molecule_residue_of_index_get(input, molecule, k, &residue);
1641 tng_residue_name_get(input, residue, str, 256);
1642 fprintf(stream, "\t\tResidue: %s\n", str);
1643 tng_residue_num_atoms_get(input, residue, &nAtoms);
1644 for (int64_t l = 0; l < nAtoms; l++)
1646 tng_residue_atom_of_index_get(input, residue, l, &atom);
1647 tng_atom_name_get(input, atom, str, 256);
1648 fprintf(stream, "\t\t\tAtom: %s", str);
1649 tng_atom_type_get(input, atom, str, 256);
1650 fprintf(stream, " (%s)\n", str);
1656 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
1657 for (int64_t l = 0; l < nAtoms; l++)
1659 tng_molecule_atom_of_index_get(input, molecule, l, &atom);
1660 tng_atom_name_get(input, atom, str, 256);
1661 fprintf(stream, "\t\t\tAtom: %s", str);
1662 tng_atom_type_get(input, atom, str, 256);
1663 fprintf(stream, " (%s)\n", str);
1669 tng_num_particles_get(input, &nAtoms);
1670 stat = tng_particle_data_vector_get(input,
1671 TNG_TRAJ_PARTIAL_CHARGES,
1676 &nValuesPerFrameRead,
1678 if (stat == TNG_SUCCESS)
1680 atomCharges.resize(nAtoms);
1681 convert_array_to_real_array(data, atomCharges.data(), 1, nAtoms, 1, datatype);
1683 fprintf(stream, "Atom Charges (%d):\n", int(nAtoms));
1684 for (int64_t i = 0; i < nAtoms; i += 10)
1686 fprintf(stream, "Atom Charges [%8d-]=[", int(i));
1687 for (int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1689 fprintf(stream, " %12.5e", atomCharges[i + j]);
1691 fprintf(stream, "]\n");
1695 stat = tng_particle_data_vector_get(
1696 input, TNG_TRAJ_MASSES, &data, &nFramesRead, &strideLength, &nParticlesRead, &nValuesPerFrameRead, &datatype);
1697 if (stat == TNG_SUCCESS)
1699 atomMasses.resize(nAtoms);
1700 convert_array_to_real_array(data, atomMasses.data(), 1, nAtoms, 1, datatype);
1702 fprintf(stream, "Atom Masses (%d):\n", int(nAtoms));
1703 for (int64_t i = 0; i < nAtoms; i += 10)
1705 fprintf(stream, "Atom Masses [%8d-]=[", int(i));
1706 for (int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1708 fprintf(stream, " %12.5e", atomMasses[i + j]);
1710 fprintf(stream, "]\n");
1716 GMX_UNUSED_VALUE(gmx_tng_input);
1717 GMX_UNUSED_VALUE(stream);
1721 gmx_bool gmx_get_tng_data_block_types_of_next_frame(gmx_tng_trajectory_t gmx_tng_input,
1724 int64_t* requestedIds,
1730 tng_function_status stat;
1731 tng_trajectory_t input = gmx_tng_input->tng;
1733 stat = tng_util_trajectory_next_frame_present_data_blocks_find(
1734 input, frame, nRequestedIds, requestedIds, nextFrame, nBlocks, blockIds);
1736 if (stat == TNG_CRITICAL)
1738 gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
1740 else if (stat == TNG_FAILURE)
1746 GMX_UNUSED_VALUE(gmx_tng_input);
1747 GMX_UNUSED_VALUE(frame);
1748 GMX_UNUSED_VALUE(nRequestedIds);
1749 GMX_UNUSED_VALUE(requestedIds);
1750 GMX_UNUSED_VALUE(nextFrame);
1751 GMX_UNUSED_VALUE(nBlocks);
1752 GMX_UNUSED_VALUE(blockIds);
1757 gmx_bool gmx_get_tng_data_next_frame_of_block_type(gmx_tng_trajectory_t gmx_tng_input,
1760 int64_t* frameNumber,
1762 int64_t* nValuesPerFrame,
1770 tng_function_status stat;
1773 int blockDependency;
1774 void* data = nullptr;
1776 tng_trajectory_t input = gmx_tng_input->tng;
1778 stat = tng_data_block_name_get(input, blockId, name, maxLen);
1779 if (stat != TNG_SUCCESS)
1781 gmx_file("Cannot read next frame of TNG file");
1783 stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
1784 if (stat != TNG_SUCCESS)
1786 gmx_file("Cannot read next frame of TNG file");
1788 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1790 tng_num_particles_get(input, nAtoms);
1791 stat = tng_util_particle_data_next_frame_read(
1792 input, blockId, &data, &datatype, frameNumber, frameTime);
1796 *nAtoms = 1; /* There are not actually any atoms, but it is used for
1797 allocating memory */
1798 stat = tng_util_non_particle_data_next_frame_read(
1799 input, blockId, &data, &datatype, frameNumber, frameTime);
1801 if (stat == TNG_CRITICAL)
1803 gmx_file("Cannot read next frame of TNG file");
1805 if (stat == TNG_FAILURE)
1811 stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
1812 if (stat != TNG_SUCCESS)
1814 gmx_file("Cannot read next frame of TNG file");
1816 srenew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
1817 convert_array_to_real_array(
1818 data, *values, getDistanceScaleFactor(gmx_tng_input), *nAtoms, *nValuesPerFrame, datatype);
1820 tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
1822 /* This must be updated if/when more lossy compression methods are added */
1823 if (codecId != TNG_TNG_COMPRESSION)
1837 GMX_UNUSED_VALUE(gmx_tng_input);
1838 GMX_UNUSED_VALUE(blockId);
1839 GMX_UNUSED_VALUE(values);
1840 GMX_UNUSED_VALUE(frameNumber);
1841 GMX_UNUSED_VALUE(frameTime);
1842 GMX_UNUSED_VALUE(nValuesPerFrame);
1843 GMX_UNUSED_VALUE(nAtoms);
1844 GMX_UNUSED_VALUE(prec);
1845 GMX_UNUSED_VALUE(name);
1846 GMX_UNUSED_VALUE(maxLen);
1847 GMX_UNUSED_VALUE(bOK);
1852 int gmx_tng_get_box_output_interval(gmx_tng_trajectory_t gmx_tng)
1855 return gmx_tng->boxOutputInterval;
1857 GMX_UNUSED_VALUE(gmx_tng);
1862 int gmx_tng_get_lambda_output_interval(gmx_tng_trajectory_t gmx_tng)
1865 return gmx_tng->lambdaOutputInterval;
1867 GMX_UNUSED_VALUE(gmx_tng);