2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2018,2019,2020, 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(),
162 gmx_version(), precisionString);
165 tng_first_program_name_set(*tng, programInfo);
169 tng_last_program_name_set(*tng, programInfo);
173 if (!gmx_getusername(username, 256))
177 tng_first_user_name_set(*tng, username);
181 tng_last_user_name_set(*tng, username);
182 tng_file_headers_write(*tng, TNG_USE_HASH);
187 gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
188 GMX_UNUSED_VALUE(filename);
189 GMX_UNUSED_VALUE(mode);
190 GMX_UNUSED_VALUE(gmx_tng);
194 void gmx_tng_close(gmx_tng_trajectory_t* gmx_tng)
196 /* We have to check that tng is set because
197 * tng_util_trajectory_close wants to return a NULL in it, and
198 * gives a fatal error if it is NULL. */
200 if (gmx_tng == nullptr || *gmx_tng == nullptr)
204 tng_trajectory_t* tng = &(*gmx_tng)->tng;
208 tng_util_trajectory_close(tng);
214 GMX_UNUSED_VALUE(gmx_tng);
219 static void addTngMoleculeFromTopology(gmx_tng_trajectory_t gmx_tng,
220 const char* moleculeName,
221 const t_atoms* atoms,
222 int64_t numMolecules,
223 tng_molecule_t* tngMol)
225 tng_trajectory_t tng = gmx_tng->tng;
226 tng_chain_t tngChain = nullptr;
227 tng_residue_t tngRes = nullptr;
229 if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
231 gmx_file("Cannot add molecule to TNG molecular system.");
234 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++)
236 const t_atom* at = &atoms->atom[atomIndex];
237 /* FIXME: Currently the TNG API can only add atoms belonging to a
238 * residue and chain. Wait for TNG 2.0*/
241 const t_resinfo* resInfo = &atoms->resinfo[at->resind];
242 char chainName[2] = { resInfo->chainid, 0 };
243 tng_atom_t tngAtom = nullptr;
248 prevAtom = &atoms->atom[atomIndex - 1];
255 /* If this is the first atom or if the residue changed add the
256 * residue to the TNG molecular system. */
257 if (!prevAtom || resInfo != &atoms->resinfo[prevAtom->resind])
259 /* If this is the first atom or if the chain changed add
260 * the chain to the TNG molecular system. */
261 if (!prevAtom || resInfo->chainid != atoms->resinfo[prevAtom->resind].chainid)
263 tng_molecule_chain_add(tng, *tngMol, chainName, &tngChain);
265 /* FIXME: When TNG supports both residue index and residue
266 * number the latter should be used. Wait for TNG 2.0*/
267 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
269 tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIndex]),
270 *(atoms->atomtype[atomIndex]), &tngAtom);
273 tng_molecule_cnt_set(tng, *tngMol, numMolecules);
276 void gmx_tng_add_mtop(gmx_tng_trajectory_t gmx_tng, const gmx_mtop_t* mtop)
280 std::vector<real> atomCharges;
281 std::vector<real> atomMasses;
285 tng_trajectory_t tng = gmx_tng->tng;
289 /* No topology information available to add. */
294 datatype = TNG_DOUBLE_DATA;
296 datatype = TNG_FLOAT_DATA;
299 atomCharges.reserve(mtop->natoms);
300 atomMasses.reserve(mtop->natoms);
302 for (const gmx_molblock_t& molBlock : mtop->molblock)
304 tng_molecule_t tngMol = nullptr;
305 const gmx_moltype_t* molType = &mtop->moltype[molBlock.type];
307 /* Add a molecule to the TNG trajectory with the same name as the
308 * current molecule. */
309 addTngMoleculeFromTopology(gmx_tng, *(molType->name), &molType->atoms, molBlock.nmol, &tngMol);
311 /* Bonds have to be deduced from interactions (constraints etc). Different
312 * interactions have different sets of parameters. */
313 /* Constraints are specified using two atoms */
314 for (i = 0; i < F_NRE; i++)
318 const InteractionList& ilist = molType->ilist[i];
320 while (j < ilist.size())
322 tng_molecule_bond_add(tng, tngMol, ilist.iatoms[j], ilist.iatoms[j + 1], &tngBond);
327 /* Settle is described using three atoms */
328 const InteractionList& ilist = molType->ilist[F_SETTLE];
330 while (j < ilist.size())
332 tng_molecule_bond_add(tng, tngMol, ilist.iatoms[j], ilist.iatoms[j + 1], &tngBond);
333 tng_molecule_bond_add(tng, tngMol, ilist.iatoms[j], ilist.iatoms[j + 2], &tngBond);
336 /* First copy atom charges and masses, first atom by atom and then copy the memory for the molecule instances.
337 * FIXME: Atom B state data should also be written to TNG (v 2.0?) */
338 for (int atomCounter = 0; atomCounter < molType->atoms.nr; atomCounter++)
340 atomCharges.push_back(molType->atoms.atom[atomCounter].q);
341 atomMasses.push_back(molType->atoms.atom[atomCounter].m);
343 for (int molCounter = 1; molCounter < molBlock.nmol; molCounter++)
345 std::copy_n(atomCharges.end() - molType->atoms.nr, molType->atoms.nr,
346 std::back_inserter(atomCharges));
347 std::copy_n(atomMasses.end() - molType->atoms.nr, molType->atoms.nr,
348 std::back_inserter(atomMasses));
351 /* Write the TNG data blocks. */
352 tng_particle_data_block_add(tng, TNG_TRAJ_PARTIAL_CHARGES, "PARTIAL CHARGES", datatype,
353 TNG_NON_TRAJECTORY_BLOCK, 1, 1, 1, 0, mtop->natoms,
354 TNG_GZIP_COMPRESSION, atomCharges.data());
355 tng_particle_data_block_add(tng, TNG_TRAJ_MASSES, "ATOM MASSES", datatype, TNG_NON_TRAJECTORY_BLOCK,
356 1, 1, 1, 0, mtop->natoms, TNG_GZIP_COMPRESSION, atomMasses.data());
359 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
360 * if they are positive.
362 * If only one of n1 and n2 is positive, then return it.
363 * If neither n1 or n2 is positive, then return -1. */
364 static int greatest_common_divisor_if_positive(int n1, int n2)
368 return (0 >= n2) ? -1 : n2;
375 /* We have a non-trivial greatest common divisor to compute. */
376 return std::gcd(n1, n2);
379 /* By default try to write 100 frames (of actual output) in each frame set.
380 * This number is the number of outputs of the most frequently written data
381 * type per frame set.
382 * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
383 * setups regarding compression efficiency and compression time. Make this
384 * a hidden command-line option? */
385 const int defaultFramesPerFrameSet = 100;
387 /*! \libinternal \brief Set the number of frames per frame
388 * set according to output intervals.
389 * The default is that 100 frames are written of the data
390 * that is written most often. */
391 static void tng_set_frames_per_frame_set(gmx_tng_trajectory_t gmx_tng,
392 const gmx_bool bUseLossyCompression,
393 const t_inputrec* ir)
396 tng_trajectory_t tng = gmx_tng->tng;
398 /* Set the number of frames per frame set to contain at least
399 * defaultFramesPerFrameSet of the lowest common denominator of
400 * the writing interval of positions and velocities. */
401 /* FIXME after 5.0: consider nstenergy also? */
402 if (bUseLossyCompression)
404 gcd = ir->nstxout_compressed;
408 gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
409 gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
416 tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
419 /*! \libinternal \brief Set the data-writing intervals, and number of
420 * frames per frame set */
421 static void set_writing_intervals(gmx_tng_trajectory_t gmx_tng,
422 const gmx_bool bUseLossyCompression,
423 const t_inputrec* ir)
425 tng_trajectory_t tng = gmx_tng->tng;
427 /* Define pointers to specific writing functions depending on if we
428 * write float or double data */
429 typedef tng_function_status (*set_writing_interval_func_pointer)(
430 tng_trajectory_t, const int64_t, const int64_t, const int64_t, const char*, const char,
433 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
435 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
437 int xout, vout, fout;
438 int gcd = -1, lowest = -1;
441 tng_set_frames_per_frame_set(gmx_tng, bUseLossyCompression, ir);
443 if (bUseLossyCompression)
445 xout = ir->nstxout_compressed;
447 /* If there is no uncompressed coordinate output write forces
448 and velocities to the compressed tng file. */
459 compression = TNG_TNG_COMPRESSION;
466 compression = TNG_GZIP_COMPRESSION;
470 set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS, "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
472 /* TODO: if/when we write energies to TNG also, reconsider how
473 * and when box information is written, because GROMACS
474 * behaviour pre-5.0 was to write the box with every
475 * trajectory frame and every energy frame, and probably
476 * people depend on this. */
478 gcd = greatest_common_divisor_if_positive(gcd, xout);
479 if (lowest < 0 || xout < lowest)
486 set_writing_interval(tng, vout, 3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
487 TNG_PARTICLE_BLOCK_DATA, compression);
489 gcd = greatest_common_divisor_if_positive(gcd, vout);
490 if (lowest < 0 || vout < lowest)
497 set_writing_interval(tng, fout, 3, TNG_TRAJ_FORCES, "FORCES", TNG_PARTICLE_BLOCK_DATA,
498 TNG_GZIP_COMPRESSION);
500 gcd = greatest_common_divisor_if_positive(gcd, fout);
501 if (lowest < 0 || fout < lowest)
508 /* Lambdas and box shape written at an interval of the lowest common
509 denominator of other output */
510 set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA, "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
511 TNG_GZIP_COMPRESSION);
513 set_writing_interval(tng, gcd, 9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
514 TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION);
515 gmx_tng->lambdaOutputInterval = gcd;
516 gmx_tng->boxOutputInterval = gcd;
517 if (gcd < lowest / 10)
520 "The lowest common denominator of trajectory output is "
521 "every %d step(s), whereas the shortest output interval "
522 "is every %d steps.",
529 void gmx_tng_prepare_md_writing(gmx_tng_trajectory_t gmx_tng, const gmx_mtop_t* mtop, const t_inputrec* ir)
532 gmx_tng_add_mtop(gmx_tng, mtop);
533 set_writing_intervals(gmx_tng, FALSE, ir);
534 tng_time_per_frame_set(gmx_tng->tng, ir->delta_t * PICO);
535 gmx_tng->timePerFrameIsSet = true;
537 GMX_UNUSED_VALUE(gmx_tng);
538 GMX_UNUSED_VALUE(mtop);
539 GMX_UNUSED_VALUE(ir);
544 /* Check if all atoms in the molecule system are selected
545 * by a selection group of type specified by gtype. */
546 static gmx_bool all_atoms_selected(const gmx_mtop_t* mtop, const SimulationAtomGroupType gtype)
548 /* Iterate through all atoms in the molecule system and
549 * check if they belong to a selection group of the
552 for (const gmx_molblock_t& molBlock : mtop->molblock)
554 const gmx_moltype_t& molType = mtop->moltype[molBlock.type];
555 const t_atoms& atoms = molType.atoms;
557 for (int j = 0; j < molBlock.nmol; j++)
559 for (int atomIndex = 0; atomIndex < atoms.nr; atomIndex++, i++)
561 if (getGroupType(mtop->groups, gtype, i) != 0)
571 /* Create TNG molecules which will represent each of the selection
572 * groups for writing. But do that only if there is actually a
573 * specified selection group and it is not the whole system.
574 * TODO: Currently the only selection that is taken into account
575 * is egcCompressedX, but other selections should be added when
576 * e.g. writing energies is implemented.
578 static void add_selection_groups(gmx_tng_trajectory_t gmx_tng, const gmx_mtop_t* mtop)
580 const t_atoms* atoms;
582 const t_resinfo* resInfo;
585 tng_molecule_t mol, iterMol;
592 tng_trajectory_t tng = gmx_tng->tng;
594 /* TODO: When the TNG molecules block is more flexible TNG selection
595 * groups should not need all atoms specified. It should be possible
596 * just to specify what molecules are selected (and/or which atoms in
597 * the molecule) and how many (if applicable). */
599 /* If no atoms are selected we do not need to create a
600 * TNG selection group molecule. */
601 if (mtop->groups.numberOfGroupNumbers(SimulationAtomGroupType::CompressedPositionOutput) == 0)
606 /* If all atoms are selected we do not have to create a selection
607 * group molecule in the TNG molecule system. */
608 if (all_atoms_selected(mtop, SimulationAtomGroupType::CompressedPositionOutput))
613 /* The name of the TNG molecule containing the selection group is the
614 * same as the name of the selection group. */
615 nameIndex = mtop->groups.groups[SimulationAtomGroupType::CompressedPositionOutput][0];
616 groupName = *mtop->groups.groupNames[nameIndex];
618 tng_molecule_alloc(tng, &mol);
619 tng_molecule_name_set(tng, mol, groupName);
620 tng_molecule_chain_add(tng, mol, "", &chain);
622 for (const gmx_molblock_t& molBlock : mtop->molblock)
624 const gmx_moltype_t& molType = mtop->moltype[molBlock.type];
626 atoms = &molType.atoms;
628 for (int j = 0; j < molBlock.nmol; j++)
630 bool bAtomsAdded = FALSE;
631 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
636 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, i) != 0)
640 at = &atoms->atom[atomIndex];
643 resInfo = &atoms->resinfo[at->resind];
644 /* FIXME: When TNG supports both residue index and residue
645 * number the latter should be used. */
646 res_name = *resInfo->name;
647 res_id = at->resind + 1;
651 res_name = const_cast<char*>("");
654 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res) != TNG_SUCCESS)
656 /* Since there is ONE chain for selection groups do not keep the
657 * original residue IDs - otherwise there might be conflicts. */
658 tng_chain_residue_add(tng, chain, res_name, &res);
660 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
661 *(atoms->atomtype[atomIndex]), atom_offset + atomIndex, &atom);
667 for (int k = 0; k < F_NRE; k++)
671 const InteractionList& ilist = molType.ilist[k];
672 for (int l = 1; l < ilist.size(); l += 3)
675 atom1 = ilist.iatoms[l] + atom_offset;
676 atom2 = ilist.iatoms[l + 1] + atom_offset;
677 if (getGroupType(mtop->groups,
678 SimulationAtomGroupType::CompressedPositionOutput, atom1)
680 && getGroupType(mtop->groups,
681 SimulationAtomGroupType::CompressedPositionOutput, atom2)
684 tng_molecule_bond_add(tng, mol, ilist.iatoms[l],
685 ilist.iatoms[l + 1], &tngBond);
690 /* Settle is described using three atoms */
691 const InteractionList& ilist = molType.ilist[F_SETTLE];
692 for (int l = 1; l < ilist.size(); l += 4)
694 int atom1, atom2, atom3;
695 atom1 = ilist.iatoms[l] + atom_offset;
696 atom2 = ilist.iatoms[l + 1] + atom_offset;
697 atom3 = ilist.iatoms[l + 2] + atom_offset;
698 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, atom1)
701 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, atom2)
704 tng_molecule_bond_add(tng, mol, atom1, atom2, &tngBond);
706 if (getGroupType(mtop->groups, SimulationAtomGroupType::CompressedPositionOutput, atom3)
709 tng_molecule_bond_add(tng, mol, atom1, atom3, &tngBond);
714 atom_offset += atoms->nr;
717 tng_molecule_existing_add(tng, &mol);
718 tng_molecule_cnt_set(tng, mol, 1);
719 tng_num_molecule_types_get(tng, &nMols);
720 for (int64_t k = 0; k < nMols; k++)
722 tng_molecule_of_index_get(tng, k, &iterMol);
727 tng_molecule_cnt_set(tng, iterMol, 0);
732 void gmx_tng_set_compression_precision(gmx_tng_trajectory_t gmx_tng, real prec)
735 tng_compression_precision_set(gmx_tng->tng, prec);
737 GMX_UNUSED_VALUE(gmx_tng);
738 GMX_UNUSED_VALUE(prec);
742 void gmx_tng_prepare_low_prec_writing(gmx_tng_trajectory_t gmx_tng, const gmx_mtop_t* mtop, const t_inputrec* ir)
745 gmx_tng_add_mtop(gmx_tng, mtop);
746 add_selection_groups(gmx_tng, mtop);
747 set_writing_intervals(gmx_tng, TRUE, ir);
748 tng_time_per_frame_set(gmx_tng->tng, ir->delta_t * PICO);
749 gmx_tng->timePerFrameIsSet = true;
750 gmx_tng_set_compression_precision(gmx_tng, ir->x_compression_precision);
752 GMX_UNUSED_VALUE(gmx_tng);
753 GMX_UNUSED_VALUE(mtop);
754 GMX_UNUSED_VALUE(ir);
758 void gmx_fwrite_tng(gmx_tng_trajectory_t gmx_tng,
759 const gmx_bool bUseLossyCompression,
761 real elapsedPicoSeconds,
770 typedef tng_function_status (*write_data_func_pointer)(
771 tng_trajectory_t, const int64_t, const double, const real*, const int64_t,
772 const int64_t, const char*, const char, const char);
774 static write_data_func_pointer write_data = tng_util_generic_with_time_double_write;
776 static write_data_func_pointer write_data = tng_util_generic_with_time_write;
778 double elapsedSeconds = elapsedPicoSeconds * PICO;
785 /* This function might get called when the type of the
786 compressed trajectory is actually XTC. So we exit and move
790 tng_trajectory_t tng = gmx_tng->tng;
792 // While the GROMACS interface to this routine specifies 'step', TNG itself
793 // only uses 'frame index' internally, although it suggests that it's a good
794 // idea to let this represent the step rather than just frame numbers.
796 // However, the way the GROMACS interface works, there are cases where
797 // the step is simply not set, so all frames rather have step=0.
798 // If we call the lower-level TNG routines with the same frame index
799 // (which is set from the step) multiple times, things will go very
800 // wrong (overwritten data).
802 // To avoid this, we require the step value to always be larger than the
803 // previous value, and if this is not true we simply set it to a value
804 // one unit larger than the previous step.
806 // Note: We do allow the user to specify any crazy value the want for the
807 // first step. If it is "not set", GROMACS will have used the value 0.
808 if (gmx_tng->lastStepDataIsValid && step <= gmx_tng->lastStep)
810 step = gmx_tng->lastStep + 1;
813 // Now that we have fixed the potentially incorrect step, we can also set
814 // the time per frame if it was not already done, and we have
815 // step/time information from the last frame so it is possible to calculate it.
816 // Note that we do allow the time to be the same for all steps, or even
817 // decreasing. In that case we will simply say the time per step is 0
818 // or negative. A bit strange, but a correct representation of the data :-)
819 if (!gmx_tng->timePerFrameIsSet && gmx_tng->lastTimeDataIsValid && gmx_tng->lastStepDataIsValid)
821 double deltaTime = elapsedSeconds - gmx_tng->lastTime;
822 std::int64_t deltaStep = step - gmx_tng->lastStep;
823 tng_time_per_frame_set(tng, deltaTime / deltaStep);
824 gmx_tng->timePerFrameIsSet = true;
827 tng_num_particles_get(tng, &nParticles);
828 if (nAtoms != static_cast<int>(nParticles))
830 tng_implicit_num_particles_set(tng, nAtoms);
833 if (bUseLossyCompression)
835 compression = TNG_TNG_COMPRESSION;
839 compression = TNG_GZIP_COMPRESSION;
842 /* The writing is done using write_data, which writes float or double
843 * depending on the GROMACS compilation. */
846 GMX_ASSERT(box, "Need a non-NULL box if positions are written");
848 if (write_data(tng, step, elapsedSeconds, reinterpret_cast<const real*>(x), 3,
849 TNG_TRAJ_POSITIONS, "POSITIONS", TNG_PARTICLE_BLOCK_DATA, compression)
852 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
858 if (write_data(tng, step, elapsedSeconds, reinterpret_cast<const real*>(v), 3,
859 TNG_TRAJ_VELOCITIES, "VELOCITIES", TNG_PARTICLE_BLOCK_DATA, compression)
862 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
868 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
869 * compression for forces regardless of output mode */
870 if (write_data(tng, step, elapsedSeconds, reinterpret_cast<const real*>(f), 3,
871 TNG_TRAJ_FORCES, "FORCES", TNG_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION)
874 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
880 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
881 * compression for box shape regardless of output mode */
882 if (write_data(tng, step, elapsedSeconds, reinterpret_cast<const real*>(box), 9,
883 TNG_TRAJ_BOX_SHAPE, "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION)
886 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
892 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
893 * compression for lambda regardless of output mode */
894 if (write_data(tng, step, elapsedSeconds, reinterpret_cast<const real*>(&lambda), 1,
895 TNG_GMX_LAMBDA, "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION)
898 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
902 // Update the data in the wrapper
904 gmx_tng->lastStepDataIsValid = true;
905 gmx_tng->lastStep = step;
906 gmx_tng->lastTimeDataIsValid = true;
907 gmx_tng->lastTime = elapsedSeconds;
909 GMX_UNUSED_VALUE(gmx_tng);
910 GMX_UNUSED_VALUE(bUseLossyCompression);
911 GMX_UNUSED_VALUE(step);
912 GMX_UNUSED_VALUE(elapsedPicoSeconds);
913 GMX_UNUSED_VALUE(lambda);
914 GMX_UNUSED_VALUE(box);
915 GMX_UNUSED_VALUE(nAtoms);
922 void fflush_tng(gmx_tng_trajectory_t gmx_tng)
929 tng_frame_set_premature_write(gmx_tng->tng, TNG_USE_HASH);
931 GMX_UNUSED_VALUE(gmx_tng);
935 float gmx_tng_get_time_of_final_frame(gmx_tng_trajectory_t gmx_tng)
941 tng_trajectory_t tng = gmx_tng->tng;
943 tng_num_frames_get(tng, &nFrames);
944 tng_util_time_of_frame_get(tng, nFrames - 1, &time);
949 GMX_UNUSED_VALUE(gmx_tng);
954 void gmx_prepare_tng_writing(const char* filename,
956 gmx_tng_trajectory_t* gmx_tng_input,
957 gmx_tng_trajectory_t* gmx_tng_output,
959 const gmx_mtop_t* mtop,
960 gmx::ArrayRef<const int> index,
961 const char* indexGroupName)
964 tng_trajectory_t* input = (gmx_tng_input && *gmx_tng_input) ? &(*gmx_tng_input)->tng : nullptr;
965 /* FIXME after 5.0: Currently only standard block types are read */
966 const int defaultNumIds = 5;
967 static int64_t fallbackIds[defaultNumIds] = { TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
968 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES, TNG_GMX_LAMBDA };
969 static char fallbackNames[defaultNumIds][32] = { "BOX SHAPE", "POSITIONS", "VELOCITIES",
970 "FORCES", "LAMBDAS" };
972 typedef tng_function_status (*set_writing_interval_func_pointer)(
973 tng_trajectory_t, const int64_t, const int64_t, const int64_t, const char*, const char,
976 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
978 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
981 gmx_tng_open(filename, mode, gmx_tng_output);
982 tng_trajectory_t* output = &(*gmx_tng_output)->tng;
984 /* Do we have an input file in TNG format? If so, then there's
985 more data we can copy over, rather than having to improvise. */
986 if (gmx_tng_input && *gmx_tng_input)
988 /* Set parameters (compression, time per frame, molecule
989 * information, number of frames per frame set and writing
990 * intervals of positions, box shape and lambdas) of the
991 * output tng container based on their respective values int
992 * the input tng container */
993 double time, compression_precision;
994 int64_t n_frames_per_frame_set, interval = -1;
996 tng_compression_precision_get(*input, &compression_precision);
997 tng_compression_precision_set(*output, compression_precision);
998 // TODO make this configurable in a future version
999 char compression_type = TNG_TNG_COMPRESSION;
1001 tng_molecule_system_copy(*input, *output);
1003 tng_time_per_frame_get(*input, &time);
1004 tng_time_per_frame_set(*output, time);
1005 // Since we have copied the value from the input TNG we should not change it again
1006 (*gmx_tng_output)->timePerFrameIsSet = true;
1008 tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
1009 tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
1011 for (int i = 0; i < defaultNumIds; i++)
1013 if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval) == TNG_SUCCESS)
1015 switch (fallbackIds[i])
1017 case TNG_TRAJ_POSITIONS:
1018 case TNG_TRAJ_VELOCITIES:
1019 set_writing_interval(*output, interval, 3, fallbackIds[i], fallbackNames[i],
1020 TNG_PARTICLE_BLOCK_DATA, compression_type);
1022 case TNG_TRAJ_FORCES:
1023 set_writing_interval(*output, interval, 3, fallbackIds[i], fallbackNames[i],
1024 TNG_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION);
1026 case TNG_TRAJ_BOX_SHAPE:
1027 set_writing_interval(*output, interval, 9, fallbackIds[i], fallbackNames[i],
1028 TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION);
1029 (*gmx_tng_output)->boxOutputInterval = interval;
1031 case TNG_GMX_LAMBDA:
1032 set_writing_interval(*output, interval, 1, fallbackIds[i], fallbackNames[i],
1033 TNG_NON_PARTICLE_BLOCK_DATA, TNG_GZIP_COMPRESSION);
1034 (*gmx_tng_output)->lambdaOutputInterval = interval;
1043 /* TODO after trjconv is modularized: fix this so the user can
1044 change precision when they are doing an operation where
1045 this makes sense, and not otherwise.
1047 char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
1048 gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
1050 gmx_tng_add_mtop(*gmx_tng_output, mtop);
1051 tng_num_frames_per_frame_set_set(*output, 1);
1054 if ((!index.empty()) && nAtoms > 0)
1056 gmx_tng_setup_atom_subgroup(*gmx_tng_output, index, indexGroupName);
1059 /* If for some reason there are more requested atoms than there are atoms in the
1060 * molecular system create a number of implicit atoms (without atom data) to
1061 * compensate for that. */
1064 tng_implicit_num_particles_set(*output, nAtoms);
1067 GMX_UNUSED_VALUE(filename);
1068 GMX_UNUSED_VALUE(mode);
1069 GMX_UNUSED_VALUE(gmx_tng_input);
1070 GMX_UNUSED_VALUE(gmx_tng_output);
1071 GMX_UNUSED_VALUE(nAtoms);
1072 GMX_UNUSED_VALUE(mtop);
1073 GMX_UNUSED_VALUE(index);
1074 GMX_UNUSED_VALUE(indexGroupName);
1078 void gmx_write_tng_from_trxframe(gmx_tng_trajectory_t gmx_tng_output, const t_trxframe* frame, int natoms)
1083 natoms = frame->natoms;
1085 gmx_fwrite_tng(gmx_tng_output, TRUE, frame->step, frame->time, 0, frame->box, natoms, frame->x,
1086 frame->v, frame->f);
1088 GMX_UNUSED_VALUE(gmx_tng_output);
1089 GMX_UNUSED_VALUE(frame);
1090 GMX_UNUSED_VALUE(natoms);
1098 void convert_array_to_real_array(void* from,
1103 const char datatype)
1107 const bool useDouble = GMX_DOUBLE;
1110 case TNG_FLOAT_DATA:
1115 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1119 for (i = 0; i < nAtoms; i++)
1121 for (j = 0; j < nValues; j++)
1123 to[i * nValues + j] = reinterpret_cast<float*>(from)[i * nValues + j] * fact;
1130 for (i = 0; i < nAtoms; i++)
1132 for (j = 0; j < nValues; j++)
1134 to[i * nValues + j] = reinterpret_cast<float*>(from)[i * nValues + j] * fact;
1140 for (i = 0; i < nAtoms; i++)
1142 for (j = 0; j < nValues; j++)
1144 to[i * nValues + j] = reinterpret_cast<int64_t*>(from)[i * nValues + j] * fact;
1148 case TNG_DOUBLE_DATA:
1149 if (sizeof(real) == sizeof(double))
1153 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1157 for (i = 0; i < nAtoms; i++)
1159 for (j = 0; j < nValues; j++)
1161 to[i * nValues + j] = reinterpret_cast<double*>(from)[i * nValues + j] * fact;
1168 for (i = 0; i < nAtoms; i++)
1170 for (j = 0; j < nValues; j++)
1172 to[i * nValues + j] = reinterpret_cast<double*>(from)[i * nValues + j] * fact;
1177 default: gmx_incons("Illegal datatype when converting values to a real array!");
1181 real getDistanceScaleFactor(gmx_tng_trajectory_t in)
1184 real distanceScaleFactor;
1186 // TODO Hopefully, TNG 2.0 will do this kind of thing for us
1187 tng_distance_unit_exponential_get(in->tng, &exp);
1189 // GROMACS expects distances in nm
1192 case 9: distanceScaleFactor = NANO / NANO; break;
1193 case 10: distanceScaleFactor = NANO / ANGSTROM; break;
1194 default: distanceScaleFactor = pow(10.0, exp + 9.0);
1197 return distanceScaleFactor;
1203 void gmx_tng_setup_atom_subgroup(gmx_tng_trajectory_t gmx_tng, gmx::ArrayRef<const int> ind, const char* name)
1206 int64_t nAtoms, cnt, nMols;
1207 tng_molecule_t mol, iterMol;
1211 tng_function_status stat;
1212 tng_trajectory_t tng = gmx_tng->tng;
1214 tng_num_particles_get(tng, &nAtoms);
1216 if (nAtoms == ind.ssize())
1221 stat = tng_molecule_find(tng, name, -1, &mol);
1222 if (stat == TNG_SUCCESS)
1224 tng_molecule_num_atoms_get(tng, mol, &nAtoms);
1225 tng_molecule_cnt_get(tng, mol, &cnt);
1226 if (nAtoms == ind.ssize())
1235 if (stat == TNG_FAILURE)
1237 /* The indexed atoms are added to one separate molecule. */
1238 tng_molecule_alloc(tng, &mol);
1239 tng_molecule_name_set(tng, mol, name);
1240 tng_molecule_chain_add(tng, mol, "", &chain);
1242 for (gmx::index i = 0; i < ind.ssize(); i++)
1244 char temp_name[256], temp_type[256];
1246 /* Try to retrieve the residue name of the atom */
1247 stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1248 if (stat != TNG_SUCCESS)
1250 temp_name[0] = '\0';
1252 /* Check if the molecule of the selection already contains this residue */
1253 if (tng_chain_residue_find(tng, chain, temp_name, -1, &res) != TNG_SUCCESS)
1255 tng_chain_residue_add(tng, chain, temp_name, &res);
1257 /* Try to find the original name and type of the atom */
1258 stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1259 if (stat != TNG_SUCCESS)
1261 temp_name[0] = '\0';
1263 stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
1264 if (stat != TNG_SUCCESS)
1266 temp_type[0] = '\0';
1268 tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
1270 tng_molecule_existing_add(tng, &mol);
1272 /* Set the count of the molecule containing the selected atoms to 1 and all
1273 * other molecules to 0 */
1274 tng_molecule_cnt_set(tng, mol, 1);
1275 tng_num_molecule_types_get(tng, &nMols);
1276 for (int64_t k = 0; k < nMols; k++)
1278 tng_molecule_of_index_get(tng, k, &iterMol);
1283 tng_molecule_cnt_set(tng, iterMol, 0);
1286 GMX_UNUSED_VALUE(gmx_tng);
1287 GMX_UNUSED_VALUE(ind);
1288 GMX_UNUSED_VALUE(name);
1292 /* TODO: If/when TNG acquires the ability to copy data blocks without
1293 * uncompressing them, then this implemenation should be reconsidered.
1294 * Ideally, gmx trjconv -f a.tng -o b.tng -b 10 -e 20 would be fast
1295 * and lose no information. */
1296 gmx_bool gmx_read_next_tng_frame(gmx_tng_trajectory_t gmx_tng_input,
1298 int64_t* requestedIds,
1299 int numRequestedIds)
1302 tng_trajectory_t input = gmx_tng_input->tng;
1303 gmx_bool bOK = TRUE;
1304 tng_function_status stat;
1305 int64_t numberOfAtoms = -1, frameNumber = -1;
1306 int64_t nBlocks, blockId, *blockIds = nullptr, codecId;
1308 void* values = nullptr;
1309 double frameTime = -1.0;
1310 int size, blockDependency;
1312 const int defaultNumIds = 5;
1313 static int64_t fallbackRequestedIds[defaultNumIds] = { TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
1314 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
1320 fr->bLambda = FALSE;
1328 /* If no specific IDs were requested read all block types that can
1329 * currently be interpreted */
1330 if (!requestedIds || numRequestedIds == 0)
1332 numRequestedIds = defaultNumIds;
1333 requestedIds = fallbackRequestedIds;
1336 stat = tng_num_particles_get(input, &numberOfAtoms);
1337 if (stat != TNG_SUCCESS)
1339 gmx_file("Cannot determine number of atoms from TNG file.");
1341 fr->natoms = numberOfAtoms;
1343 bool nextFrameExists = gmx_get_tng_data_block_types_of_next_frame(
1344 gmx_tng_input, fr->step, numRequestedIds, requestedIds, &frameNumber, &nBlocks, &blockIds);
1345 gmx::unique_cptr<int64_t, gmx::free_wrapper> blockIdsGuard(blockIds);
1346 if (!nextFrameExists)
1356 for (int64_t i = 0; i < nBlocks; i++)
1358 blockId = blockIds[i];
1359 tng_data_block_dependency_get(input, blockId, &blockDependency);
1360 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1362 stat = tng_util_particle_data_next_frame_read(input, blockId, &values, &datatype,
1363 &frameNumber, &frameTime);
1367 stat = tng_util_non_particle_data_next_frame_read(input, blockId, &values, &datatype,
1368 &frameNumber, &frameTime);
1370 if (stat == TNG_CRITICAL)
1372 gmx_file("Cannot read positions from TNG file.");
1375 else if (stat == TNG_FAILURE)
1381 case TNG_TRAJ_BOX_SHAPE:
1384 case TNG_INT_DATA: size = sizeof(int64_t); break;
1385 case TNG_FLOAT_DATA: size = sizeof(float); break;
1386 case TNG_DOUBLE_DATA: size = sizeof(double); break;
1387 default: gmx_incons("Illegal datatype of box shape values!");
1389 for (int i = 0; i < DIM; i++)
1391 convert_array_to_real_array(reinterpret_cast<char*>(values) + size * i * DIM,
1392 reinterpret_cast<real*>(fr->box[i]),
1393 getDistanceScaleFactor(gmx_tng_input), 1, DIM, datatype);
1397 case TNG_TRAJ_POSITIONS:
1398 srenew(fr->x, fr->natoms);
1399 convert_array_to_real_array(values, reinterpret_cast<real*>(fr->x),
1400 getDistanceScaleFactor(gmx_tng_input), fr->natoms, DIM,
1403 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1404 /* This must be updated if/when more lossy compression methods are added */
1405 if (codecId == TNG_TNG_COMPRESSION)
1411 case TNG_TRAJ_VELOCITIES:
1412 srenew(fr->v, fr->natoms);
1413 convert_array_to_real_array(values, reinterpret_cast<real*>(fr->v),
1414 getDistanceScaleFactor(gmx_tng_input), fr->natoms, DIM,
1417 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1418 /* This must be updated if/when more lossy compression methods are added */
1419 if (codecId == TNG_TNG_COMPRESSION)
1425 case TNG_TRAJ_FORCES:
1426 srenew(fr->f, fr->natoms);
1427 convert_array_to_real_array(values, reinterpret_cast<real*>(fr->f),
1428 getDistanceScaleFactor(gmx_tng_input), fr->natoms, DIM,
1432 case TNG_GMX_LAMBDA:
1435 case TNG_FLOAT_DATA: fr->lambda = *(reinterpret_cast<float*>(values)); break;
1436 case TNG_DOUBLE_DATA: fr->lambda = *(reinterpret_cast<double*>(values)); break;
1437 default: gmx_incons("Illegal datatype lambda value!");
1443 "Illegal block type! Currently GROMACS tools can only handle certain data "
1444 "types. Skipping block.");
1446 /* values does not have to be freed before reading next frame. It will
1447 * be reallocated if it is not NULL. */
1450 fr->step = frameNumber;
1453 // Convert the time to ps
1454 fr->time = frameTime / PICO;
1455 fr->bTime = (frameTime > 0);
1457 // TODO This does not leak, but is not exception safe.
1458 /* values must be freed before leaving this function */
1463 GMX_UNUSED_VALUE(gmx_tng_input);
1464 GMX_UNUSED_VALUE(fr);
1465 GMX_UNUSED_VALUE(requestedIds);
1466 GMX_UNUSED_VALUE(numRequestedIds);
1471 void gmx_print_tng_molecule_system(gmx_tng_trajectory_t gmx_tng_input, FILE* stream)
1474 int64_t nMolecules, nChains, nResidues, nAtoms, nFramesRead;
1475 int64_t strideLength, nParticlesRead, nValuesPerFrameRead, *molCntList;
1476 tng_molecule_t molecule;
1478 tng_residue_t residue;
1480 tng_function_status stat;
1484 void* data = nullptr;
1485 std::vector<real> atomCharges;
1486 std::vector<real> atomMasses;
1487 tng_trajectory_t input = gmx_tng_input->tng;
1489 tng_num_molecule_types_get(input, &nMolecules);
1490 tng_molecule_cnt_list_get(input, &molCntList);
1491 /* Can the number of particles change in the trajectory or is it constant? */
1492 tng_num_particles_variable_get(input, &varNAtoms);
1494 for (int64_t i = 0; i < nMolecules; i++)
1496 tng_molecule_of_index_get(input, i, &molecule);
1497 tng_molecule_name_get(input, molecule, str, 256);
1498 if (varNAtoms == TNG_CONSTANT_N_ATOMS)
1500 if (static_cast<int>(molCntList[i]) == 0)
1504 fprintf(stream, "Molecule: %s, count: %d\n", str, static_cast<int>(molCntList[i]));
1508 fprintf(stream, "Molecule: %s\n", str);
1510 tng_molecule_num_chains_get(input, molecule, &nChains);
1513 for (int64_t j = 0; j < nChains; j++)
1515 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
1516 tng_chain_name_get(input, chain, str, 256);
1517 fprintf(stream, "\tChain: %s\n", str);
1518 tng_chain_num_residues_get(input, chain, &nResidues);
1519 for (int64_t k = 0; k < nResidues; k++)
1521 tng_chain_residue_of_index_get(input, chain, k, &residue);
1522 tng_residue_name_get(input, residue, str, 256);
1523 fprintf(stream, "\t\tResidue: %s\n", str);
1524 tng_residue_num_atoms_get(input, residue, &nAtoms);
1525 for (int64_t l = 0; l < nAtoms; l++)
1527 tng_residue_atom_of_index_get(input, residue, l, &atom);
1528 tng_atom_name_get(input, atom, str, 256);
1529 fprintf(stream, "\t\t\tAtom: %s", str);
1530 tng_atom_type_get(input, atom, str, 256);
1531 fprintf(stream, " (%s)\n", str);
1536 /* It is possible to have a molecule without chains, in which case
1537 * residues in the molecule can be iterated through without going
1538 * through chains. */
1541 tng_molecule_num_residues_get(input, molecule, &nResidues);
1544 for (int64_t k = 0; k < nResidues; k++)
1546 tng_molecule_residue_of_index_get(input, molecule, k, &residue);
1547 tng_residue_name_get(input, residue, str, 256);
1548 fprintf(stream, "\t\tResidue: %s\n", str);
1549 tng_residue_num_atoms_get(input, residue, &nAtoms);
1550 for (int64_t l = 0; l < nAtoms; l++)
1552 tng_residue_atom_of_index_get(input, residue, l, &atom);
1553 tng_atom_name_get(input, atom, str, 256);
1554 fprintf(stream, "\t\t\tAtom: %s", str);
1555 tng_atom_type_get(input, atom, str, 256);
1556 fprintf(stream, " (%s)\n", str);
1562 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
1563 for (int64_t l = 0; l < nAtoms; l++)
1565 tng_molecule_atom_of_index_get(input, molecule, l, &atom);
1566 tng_atom_name_get(input, atom, str, 256);
1567 fprintf(stream, "\t\t\tAtom: %s", str);
1568 tng_atom_type_get(input, atom, str, 256);
1569 fprintf(stream, " (%s)\n", str);
1575 tng_num_particles_get(input, &nAtoms);
1576 stat = tng_particle_data_vector_get(input, TNG_TRAJ_PARTIAL_CHARGES, &data, &nFramesRead,
1577 &strideLength, &nParticlesRead, &nValuesPerFrameRead, &datatype);
1578 if (stat == TNG_SUCCESS)
1580 atomCharges.resize(nAtoms);
1581 convert_array_to_real_array(data, atomCharges.data(), 1, nAtoms, 1, datatype);
1583 fprintf(stream, "Atom Charges (%d):\n", int(nAtoms));
1584 for (int64_t i = 0; i < nAtoms; i += 10)
1586 fprintf(stream, "Atom Charges [%8d-]=[", int(i));
1587 for (int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1589 fprintf(stream, " %12.5e", atomCharges[i + j]);
1591 fprintf(stream, "]\n");
1595 stat = tng_particle_data_vector_get(input, TNG_TRAJ_MASSES, &data, &nFramesRead, &strideLength,
1596 &nParticlesRead, &nValuesPerFrameRead, &datatype);
1597 if (stat == TNG_SUCCESS)
1599 atomMasses.resize(nAtoms);
1600 convert_array_to_real_array(data, atomMasses.data(), 1, nAtoms, 1, datatype);
1602 fprintf(stream, "Atom Masses (%d):\n", int(nAtoms));
1603 for (int64_t i = 0; i < nAtoms; i += 10)
1605 fprintf(stream, "Atom Masses [%8d-]=[", int(i));
1606 for (int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1608 fprintf(stream, " %12.5e", atomMasses[i + j]);
1610 fprintf(stream, "]\n");
1616 GMX_UNUSED_VALUE(gmx_tng_input);
1617 GMX_UNUSED_VALUE(stream);
1621 gmx_bool gmx_get_tng_data_block_types_of_next_frame(gmx_tng_trajectory_t gmx_tng_input,
1624 int64_t* requestedIds,
1630 tng_function_status stat;
1631 tng_trajectory_t input = gmx_tng_input->tng;
1633 stat = tng_util_trajectory_next_frame_present_data_blocks_find(
1634 input, frame, nRequestedIds, requestedIds, nextFrame, nBlocks, blockIds);
1636 if (stat == TNG_CRITICAL)
1638 gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
1640 else if (stat == TNG_FAILURE)
1646 GMX_UNUSED_VALUE(gmx_tng_input);
1647 GMX_UNUSED_VALUE(frame);
1648 GMX_UNUSED_VALUE(nRequestedIds);
1649 GMX_UNUSED_VALUE(requestedIds);
1650 GMX_UNUSED_VALUE(nextFrame);
1651 GMX_UNUSED_VALUE(nBlocks);
1652 GMX_UNUSED_VALUE(blockIds);
1657 gmx_bool gmx_get_tng_data_next_frame_of_block_type(gmx_tng_trajectory_t gmx_tng_input,
1660 int64_t* frameNumber,
1662 int64_t* nValuesPerFrame,
1670 tng_function_status stat;
1673 int blockDependency;
1674 void* data = nullptr;
1676 tng_trajectory_t input = gmx_tng_input->tng;
1678 stat = tng_data_block_name_get(input, blockId, name, maxLen);
1679 if (stat != TNG_SUCCESS)
1681 gmx_file("Cannot read next frame of TNG file");
1683 stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
1684 if (stat != TNG_SUCCESS)
1686 gmx_file("Cannot read next frame of TNG file");
1688 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1690 tng_num_particles_get(input, nAtoms);
1691 stat = tng_util_particle_data_next_frame_read(input, blockId, &data, &datatype, frameNumber,
1696 *nAtoms = 1; /* There are not actually any atoms, but it is used for
1697 allocating memory */
1698 stat = tng_util_non_particle_data_next_frame_read(input, blockId, &data, &datatype,
1699 frameNumber, frameTime);
1701 if (stat == TNG_CRITICAL)
1703 gmx_file("Cannot read next frame of TNG file");
1705 if (stat == TNG_FAILURE)
1711 stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
1712 if (stat != TNG_SUCCESS)
1714 gmx_file("Cannot read next frame of TNG file");
1716 srenew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
1717 convert_array_to_real_array(data, *values, getDistanceScaleFactor(gmx_tng_input), *nAtoms,
1718 *nValuesPerFrame, datatype);
1720 tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
1722 /* This must be updated if/when more lossy compression methods are added */
1723 if (codecId != TNG_TNG_COMPRESSION)
1737 GMX_UNUSED_VALUE(gmx_tng_input);
1738 GMX_UNUSED_VALUE(blockId);
1739 GMX_UNUSED_VALUE(values);
1740 GMX_UNUSED_VALUE(frameNumber);
1741 GMX_UNUSED_VALUE(frameTime);
1742 GMX_UNUSED_VALUE(nValuesPerFrame);
1743 GMX_UNUSED_VALUE(nAtoms);
1744 GMX_UNUSED_VALUE(prec);
1745 GMX_UNUSED_VALUE(name);
1746 GMX_UNUSED_VALUE(maxLen);
1747 GMX_UNUSED_VALUE(bOK);
1752 int gmx_tng_get_box_output_interval(gmx_tng_trajectory_t gmx_tng)
1755 return gmx_tng->boxOutputInterval;
1757 GMX_UNUSED_VALUE(gmx_tng);
1762 int gmx_tng_get_lambda_output_interval(gmx_tng_trajectory_t gmx_tng)
1765 return gmx_tng->lambdaOutputInterval;
1767 GMX_UNUSED_VALUE(gmx_tng);