2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
49 #include "tng/tng_io.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/utilities.h"
54 #include "gromacs/mdtypes/inputrec.h"
55 #include "gromacs/topology/ifunc.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/trajectory/trajectoryframe.h"
58 #include "gromacs/utility/basedefinitions.h"
59 #include "gromacs/utility/baseversion.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/gmxassert.h"
63 #include "gromacs/utility/programcontext.h"
64 #include "gromacs/utility/smalloc.h"
65 #include "gromacs/utility/sysinfo.h"
66 #include "gromacs/utility/unique_cptr.h"
69 using tng_trajectory_t = void *;
72 /*! \brief Gromacs Wrapper around tng datatype
74 * This could in principle hold any GROMACS-specific requirements not yet
75 * implemented in or not relevant to the TNG library itself. However, for now
76 * we only use it to handle some shortcomings we have discovered, where the TNG
77 * API itself is a bit fragile and can end up overwriting data if called several
78 * times with the same frame number.
79 * The logic to determine the time per step was also a bit fragile. This is not
80 * critical, but since we anyway need a wrapper for ensuring unique frame
81 * numbers, we can also use it to store the time of the first step and use that
82 * to derive a slightly better/safer estimate of the time per step.
84 * At some future point where we have a second-generation TNG API we should
85 * consider removing this again.
87 struct gmx_tng_trajectory
89 tng_trajectory_t tng; //!< Actual TNG handle (pointer)
90 bool lastStepDataIsValid; //!< True if lastStep has been set
91 std::int64_t lastStep; //!< Index/step used for last frame
92 bool lastTimeDataIsValid; //!< True if lastTime has been set
93 double lastTime; //!< Time of last frame (TNG unit is seconds)
94 bool timePerFrameIsSet; //!< True if we have set the time per frame
95 int boxOutputInterval; //!< Number of steps between the output of box size
96 int lambdaOutputInterval; //!< Number of steps between the output of lambdas
100 static const char *modeToVerb(char mode)
115 gmx_fatal(FARGS, "Invalid file opening mode %c", mode);
121 void gmx_tng_open(const char *filename,
123 gmx_tng_trajectory_t *gmx_tng)
126 /* First check whether we have to make a backup,
127 * only for writing, not for read or append.
131 make_backup(filename);
134 *gmx_tng = new gmx_tng_trajectory;
135 (*gmx_tng)->lastStepDataIsValid = false;
136 (*gmx_tng)->lastTimeDataIsValid = false;
137 (*gmx_tng)->timePerFrameIsSet = false;
138 tng_trajectory_t * tng = &(*gmx_tng)->tng;
140 /* tng must not be pointing at already allocated memory.
141 * Memory will be allocated by tng_util_trajectory_open() and must
142 * later on be freed by tng_util_trajectory_close(). */
143 if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
145 /* TNG does return more than one degree of error, but there is
146 no use case for GROMACS handling the non-fatal errors
149 "File I/O error while opening %s for %s",
154 if (mode == 'w' || mode == 'a')
157 gmx_gethostname(hostname, 256);
160 tng_first_computer_name_set(*tng, hostname);
164 tng_last_computer_name_set(*tng, hostname);
167 char programInfo[256];
168 const char *precisionString = "";
170 precisionString = " (double precision)";
172 sprintf(programInfo, "%.100s %.128s%.24s",
173 gmx::getProgramContext().displayName(),
174 gmx_version(), precisionString);
177 tng_first_program_name_set(*tng, programInfo);
181 tng_last_program_name_set(*tng, programInfo);
185 if (!gmx_getusername(username, 256))
189 tng_first_user_name_set(*tng, username);
193 tng_last_user_name_set(*tng, username);
194 tng_file_headers_write(*tng, TNG_USE_HASH);
199 gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
200 GMX_UNUSED_VALUE(filename);
201 GMX_UNUSED_VALUE(mode);
202 GMX_UNUSED_VALUE(gmx_tng);
206 void gmx_tng_close(gmx_tng_trajectory_t *gmx_tng)
208 /* We have to check that tng is set because
209 * tng_util_trajectory_close wants to return a NULL in it, and
210 * gives a fatal error if it is NULL. */
212 if (gmx_tng == nullptr || *gmx_tng == nullptr)
216 tng_trajectory_t * tng = &(*gmx_tng)->tng;
220 tng_util_trajectory_close(tng);
226 GMX_UNUSED_VALUE(gmx_tng);
231 static void addTngMoleculeFromTopology(gmx_tng_trajectory_t gmx_tng,
232 const char *moleculeName,
233 const t_atoms *atoms,
234 int64_t numMolecules,
235 tng_molecule_t *tngMol)
237 tng_trajectory_t tng = gmx_tng->tng;
238 tng_chain_t tngChain = nullptr;
239 tng_residue_t tngRes = nullptr;
241 if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
243 gmx_file("Cannot add molecule to TNG molecular system.");
246 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++)
248 const t_atom *at = &atoms->atom[atomIndex];
249 /* FIXME: Currently the TNG API can only add atoms belonging to a
250 * residue and chain. Wait for TNG 2.0*/
253 const t_resinfo *resInfo = &atoms->resinfo[at->resind];
254 char chainName[2] = {resInfo->chainid, 0};
255 tng_atom_t tngAtom = nullptr;
260 prevAtom = &atoms->atom[atomIndex - 1];
267 /* If this is the first atom or if the residue changed add the
268 * residue to the TNG molecular system. */
269 if (!prevAtom || resInfo != &atoms->resinfo[prevAtom->resind])
271 /* If this is the first atom or if the chain changed add
272 * the chain to the TNG molecular system. */
273 if (!prevAtom || resInfo->chainid !=
274 atoms->resinfo[prevAtom->resind].chainid)
276 tng_molecule_chain_add(tng, *tngMol, chainName,
279 /* FIXME: When TNG supports both residue index and residue
280 * number the latter should be used. Wait for TNG 2.0*/
281 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
283 tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIndex]), *(atoms->atomtype[atomIndex]), &tngAtom);
286 tng_molecule_cnt_set(tng, *tngMol, numMolecules);
289 void gmx_tng_add_mtop(gmx_tng_trajectory_t gmx_tng,
290 const gmx_mtop_t *mtop)
294 std::vector<real> atomCharges;
295 std::vector<real> atomMasses;
299 tng_trajectory_t tng = gmx_tng->tng;
303 /* No topology information available to add. */
308 datatype = TNG_DOUBLE_DATA;
310 datatype = TNG_FLOAT_DATA;
313 atomCharges.reserve(mtop->natoms);
314 atomMasses.reserve(mtop->natoms);
316 for (const gmx_molblock_t &molBlock : mtop->molblock)
318 tng_molecule_t tngMol = nullptr;
319 const gmx_moltype_t *molType = &mtop->moltype[molBlock.type];
321 /* Add a molecule to the TNG trajectory with the same name as the
322 * current molecule. */
323 addTngMoleculeFromTopology(gmx_tng,
329 /* Bonds have to be deduced from interactions (constraints etc). Different
330 * interactions have different sets of parameters. */
331 /* Constraints are specified using two atoms */
332 for (i = 0; i < F_NRE; i++)
336 const InteractionList &ilist = molType->ilist[i];
338 while (j < ilist.size())
340 tng_molecule_bond_add(tng, tngMol, ilist.iatoms[j], ilist.iatoms[j + 1], &tngBond);
345 /* Settle is described using three atoms */
346 const InteractionList &ilist = molType->ilist[F_SETTLE];
348 while (j < ilist.size())
350 tng_molecule_bond_add(tng, tngMol, ilist.iatoms[j], ilist.iatoms[j + 1], &tngBond);
351 tng_molecule_bond_add(tng, tngMol, ilist.iatoms[j], ilist.iatoms[j + 2], &tngBond);
354 /* First copy atom charges and masses, first atom by atom and then copy the memory for the molecule instances.
355 * FIXME: Atom B state data should also be written to TNG (v 2.0?) */
356 for (int atomCounter = 0; atomCounter < molType->atoms.nr; atomCounter++)
358 atomCharges.push_back(molType->atoms.atom[atomCounter].q);
359 atomMasses.push_back(molType->atoms.atom[atomCounter].m);
361 for (int molCounter = 1; molCounter < molBlock.nmol; molCounter++)
363 std::copy_n(atomCharges.end() - molType->atoms.nr, molType->atoms.nr, std::back_inserter(atomCharges));
364 std::copy_n(atomMasses.end() - molType->atoms.nr, molType->atoms.nr, std::back_inserter(atomMasses));
367 /* Write the TNG data blocks. */
368 tng_particle_data_block_add(tng, TNG_TRAJ_PARTIAL_CHARGES, "PARTIAL CHARGES",
369 datatype, TNG_NON_TRAJECTORY_BLOCK,
370 1, 1, 1, 0, mtop->natoms,
371 TNG_GZIP_COMPRESSION, atomCharges.data());
372 tng_particle_data_block_add(tng, TNG_TRAJ_MASSES, "ATOM MASSES",
373 datatype, TNG_NON_TRAJECTORY_BLOCK,
374 1, 1, 1, 0, mtop->natoms,
375 TNG_GZIP_COMPRESSION, atomMasses.data());
378 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
379 * if they are positive.
381 * If only one of n1 and n2 is positive, then return it.
382 * If neither n1 or n2 is positive, then return -1. */
384 greatest_common_divisor_if_positive(int n1, int n2)
388 return (0 >= n2) ? -1 : n2;
395 /* We have a non-trivial greatest common divisor to compute. */
396 return gmx_greatest_common_divisor(n1, n2);
399 /* By default try to write 100 frames (of actual output) in each frame set.
400 * This number is the number of outputs of the most frequently written data
401 * type per frame set.
402 * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
403 * setups regarding compression efficiency and compression time. Make this
404 * a hidden command-line option? */
405 const int defaultFramesPerFrameSet = 100;
407 /*! \libinternal \brief Set the number of frames per frame
408 * set according to output intervals.
409 * The default is that 100 frames are written of the data
410 * that is written most often. */
411 static void tng_set_frames_per_frame_set(gmx_tng_trajectory_t gmx_tng,
412 const gmx_bool bUseLossyCompression,
413 const t_inputrec *ir)
416 tng_trajectory_t tng = gmx_tng->tng;
418 /* Set the number of frames per frame set to contain at least
419 * defaultFramesPerFrameSet of the lowest common denominator of
420 * the writing interval of positions and velocities. */
421 /* FIXME after 5.0: consider nstenergy also? */
422 if (bUseLossyCompression)
424 gcd = ir->nstxout_compressed;
428 gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
429 gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
436 tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
439 /*! \libinternal \brief Set the data-writing intervals, and number of
440 * frames per frame set */
441 static void set_writing_intervals(gmx_tng_trajectory_t gmx_tng,
442 const gmx_bool bUseLossyCompression,
443 const t_inputrec *ir)
445 tng_trajectory_t tng = gmx_tng->tng;
447 /* Define pointers to specific writing functions depending on if we
448 * write float or double data */
449 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
457 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
459 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
461 int xout, vout, fout;
462 int gcd = -1, lowest = -1;
465 tng_set_frames_per_frame_set(gmx_tng, bUseLossyCompression, ir);
467 if (bUseLossyCompression)
469 xout = ir->nstxout_compressed;
471 /* If there is no uncompressed coordinate output write forces
472 and velocities to the compressed tng file. */
483 compression = TNG_TNG_COMPRESSION;
490 compression = TNG_GZIP_COMPRESSION;
494 set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS,
495 "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
497 /* TODO: if/when we write energies to TNG also, reconsider how
498 * and when box information is written, because GROMACS
499 * behaviour pre-5.0 was to write the box with every
500 * trajectory frame and every energy frame, and probably
501 * people depend on this. */
503 gcd = greatest_common_divisor_if_positive(gcd, xout);
504 if (lowest < 0 || xout < lowest)
511 set_writing_interval(tng, vout, 3, TNG_TRAJ_VELOCITIES,
512 "VELOCITIES", TNG_PARTICLE_BLOCK_DATA,
515 gcd = greatest_common_divisor_if_positive(gcd, vout);
516 if (lowest < 0 || vout < lowest)
523 set_writing_interval(tng, fout, 3, TNG_TRAJ_FORCES,
524 "FORCES", TNG_PARTICLE_BLOCK_DATA,
525 TNG_GZIP_COMPRESSION);
527 gcd = greatest_common_divisor_if_positive(gcd, fout);
528 if (lowest < 0 || fout < lowest)
535 /* Lambdas and box shape written at an interval of the lowest common
536 denominator of other output */
537 set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA,
538 "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
539 TNG_GZIP_COMPRESSION);
541 set_writing_interval(tng, gcd, 9, TNG_TRAJ_BOX_SHAPE,
542 "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA,
543 TNG_GZIP_COMPRESSION);
544 gmx_tng->lambdaOutputInterval = gcd;
545 gmx_tng->boxOutputInterval = gcd;
546 if (gcd < lowest / 10)
548 gmx_warning("The lowest common denominator of trajectory output is "
549 "every %d step(s), whereas the shortest output interval "
550 "is every %d steps.", gcd, lowest);
556 void gmx_tng_prepare_md_writing(gmx_tng_trajectory_t gmx_tng,
557 const gmx_mtop_t *mtop,
558 const t_inputrec *ir)
561 gmx_tng_add_mtop(gmx_tng, mtop);
562 set_writing_intervals(gmx_tng, FALSE, ir);
563 tng_time_per_frame_set(gmx_tng->tng, ir->delta_t * PICO);
564 gmx_tng->timePerFrameIsSet = true;
566 GMX_UNUSED_VALUE(gmx_tng);
567 GMX_UNUSED_VALUE(mtop);
568 GMX_UNUSED_VALUE(ir);
573 /* Check if all atoms in the molecule system are selected
574 * by a selection group of type specified by gtype. */
575 static gmx_bool all_atoms_selected(const gmx_mtop_t *mtop,
578 /* Iterate through all atoms in the molecule system and
579 * check if they belong to a selection group of the
582 for (const gmx_molblock_t &molBlock : mtop->molblock)
584 const gmx_moltype_t &molType = mtop->moltype[molBlock.type];
585 const t_atoms &atoms = molType.atoms;
587 for (int j = 0; j < molBlock.nmol; j++)
589 for (int atomIndex = 0; atomIndex < atoms.nr; atomIndex++, i++)
591 if (getGroupType(&mtop->groups, gtype, i) != 0)
601 /* Create TNG molecules which will represent each of the selection
602 * groups for writing. But do that only if there is actually a
603 * specified selection group and it is not the whole system.
604 * TODO: Currently the only selection that is taken into account
605 * is egcCompressedX, but other selections should be added when
606 * e.g. writing energies is implemented.
608 static void add_selection_groups(gmx_tng_trajectory_t gmx_tng,
609 const gmx_mtop_t *mtop)
611 const t_atoms *atoms;
613 const t_resinfo *resInfo;
616 tng_molecule_t mol, iterMol;
623 tng_trajectory_t tng = gmx_tng->tng;
625 /* TODO: When the TNG molecules block is more flexible TNG selection
626 * groups should not need all atoms specified. It should be possible
627 * just to specify what molecules are selected (and/or which atoms in
628 * the molecule) and how many (if applicable). */
630 /* If no atoms are selected we do not need to create a
631 * TNG selection group molecule. */
632 if (mtop->groups.ngrpnr[egcCompressedX] == 0)
637 /* If all atoms are selected we do not have to create a selection
638 * group molecule in the TNG molecule system. */
639 if (all_atoms_selected(mtop, egcCompressedX))
644 /* The name of the TNG molecule containing the selection group is the
645 * same as the name of the selection group. */
646 nameIndex = *mtop->groups.grps[egcCompressedX].nm_ind;
647 groupName = *mtop->groups.grpname[nameIndex];
649 tng_molecule_alloc(tng, &mol);
650 tng_molecule_name_set(tng, mol, groupName);
651 tng_molecule_chain_add(tng, mol, "", &chain);
653 for (const gmx_molblock_t &molBlock : mtop->molblock)
655 const gmx_moltype_t &molType = mtop->moltype[molBlock.type];
657 atoms = &molType.atoms;
659 for (int j = 0; j < molBlock.nmol; j++)
661 bool bAtomsAdded = FALSE;
662 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
667 if (getGroupType(&mtop->groups, egcCompressedX, i) != 0)
671 at = &atoms->atom[atomIndex];
674 resInfo = &atoms->resinfo[at->resind];
675 /* FIXME: When TNG supports both residue index and residue
676 * number the latter should be used. */
677 res_name = *resInfo->name;
678 res_id = at->resind + 1;
682 res_name = const_cast<char*>("");
685 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res)
688 /* Since there is ONE chain for selection groups do not keep the
689 * original residue IDs - otherwise there might be conflicts. */
690 tng_chain_residue_add(tng, chain, res_name, &res);
692 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
693 *(atoms->atomtype[atomIndex]),
694 atom_offset + atomIndex, &atom);
700 for (int k = 0; k < F_NRE; k++)
704 const InteractionList &ilist = molType.ilist[k];
705 for (int l = 1; l < ilist.size(); l += 3)
708 atom1 = ilist.iatoms[l] + atom_offset;
709 atom2 = ilist.iatoms[l + 1] + atom_offset;
710 if (getGroupType(&mtop->groups, egcCompressedX, atom1) == 0 &&
711 getGroupType(&mtop->groups, egcCompressedX, atom2) == 0)
713 tng_molecule_bond_add(tng, mol, ilist.iatoms[l],
714 ilist.iatoms[l + 1], &tngBond);
719 /* Settle is described using three atoms */
720 const InteractionList &ilist = molType.ilist[F_SETTLE];
721 for (int l = 1; l < ilist.size(); l += 4)
723 int atom1, atom2, atom3;
724 atom1 = ilist.iatoms[l] + atom_offset;
725 atom2 = ilist.iatoms[l + 1] + atom_offset;
726 atom3 = ilist.iatoms[l + 2] + atom_offset;
727 if (getGroupType(&mtop->groups, egcCompressedX, atom1) == 0)
729 if (getGroupType(&mtop->groups, egcCompressedX, atom2) == 0)
731 tng_molecule_bond_add(tng, mol, atom1,
734 if (getGroupType(&mtop->groups, egcCompressedX, atom3) == 0)
736 tng_molecule_bond_add(tng, mol, atom1,
742 atom_offset += atoms->nr;
745 tng_molecule_existing_add(tng, &mol);
746 tng_molecule_cnt_set(tng, mol, 1);
747 tng_num_molecule_types_get(tng, &nMols);
748 for (int64_t k = 0; k < nMols; k++)
750 tng_molecule_of_index_get(tng, k, &iterMol);
755 tng_molecule_cnt_set(tng, iterMol, 0);
760 void gmx_tng_set_compression_precision(gmx_tng_trajectory_t gmx_tng,
764 tng_compression_precision_set(gmx_tng->tng, prec);
766 GMX_UNUSED_VALUE(gmx_tng);
767 GMX_UNUSED_VALUE(prec);
771 void gmx_tng_prepare_low_prec_writing(gmx_tng_trajectory_t gmx_tng,
772 const gmx_mtop_t *mtop,
773 const t_inputrec *ir)
776 gmx_tng_add_mtop(gmx_tng, mtop);
777 add_selection_groups(gmx_tng, mtop);
778 set_writing_intervals(gmx_tng, TRUE, ir);
779 tng_time_per_frame_set(gmx_tng->tng, ir->delta_t * PICO);
780 gmx_tng->timePerFrameIsSet = true;
781 gmx_tng_set_compression_precision(gmx_tng, ir->x_compression_precision);
783 GMX_UNUSED_VALUE(gmx_tng);
784 GMX_UNUSED_VALUE(mtop);
785 GMX_UNUSED_VALUE(ir);
789 void gmx_fwrite_tng(gmx_tng_trajectory_t gmx_tng,
790 const gmx_bool bUseLossyCompression,
792 real elapsedPicoSeconds,
801 typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
811 static write_data_func_pointer write_data = tng_util_generic_with_time_double_write;
813 static write_data_func_pointer write_data = tng_util_generic_with_time_write;
815 double elapsedSeconds = elapsedPicoSeconds * PICO;
822 /* This function might get called when the type of the
823 compressed trajectory is actually XTC. So we exit and move
827 tng_trajectory_t tng = gmx_tng->tng;
829 // While the GROMACS interface to this routine specifies 'step', TNG itself
830 // only uses 'frame index' internally, although it suggests that it's a good
831 // idea to let this represent the step rather than just frame numbers.
833 // However, the way the GROMACS interface works, there are cases where
834 // the step is simply not set, so all frames rather have step=0.
835 // If we call the lower-level TNG routines with the same frame index
836 // (which is set from the step) multiple times, things will go very
837 // wrong (overwritten data).
839 // To avoid this, we require the step value to always be larger than the
840 // previous value, and if this is not true we simply set it to a value
841 // one unit larger than the previous step.
843 // Note: We do allow the user to specify any crazy value the want for the
844 // first step. If it is "not set", GROMACS will have used the value 0.
845 if (gmx_tng->lastStepDataIsValid && step <= gmx_tng->lastStep)
847 step = gmx_tng->lastStep + 1;
850 // Now that we have fixed the potentially incorrect step, we can also set
851 // the time per frame if it was not already done, and we have
852 // step/time information from the last frame so it is possible to calculate it.
853 // Note that we do allow the time to be the same for all steps, or even
854 // decreasing. In that case we will simply say the time per step is 0
855 // or negative. A bit strange, but a correct representation of the data :-)
856 if (!gmx_tng->timePerFrameIsSet && gmx_tng->lastTimeDataIsValid && gmx_tng->lastStepDataIsValid)
858 double deltaTime = elapsedSeconds - gmx_tng->lastTime;
859 std::int64_t deltaStep = step - gmx_tng->lastStep;
860 tng_time_per_frame_set(tng, deltaTime / deltaStep );
861 gmx_tng->timePerFrameIsSet = true;
864 tng_num_particles_get(tng, &nParticles);
865 if (nAtoms != static_cast<int>(nParticles))
867 tng_implicit_num_particles_set(tng, nAtoms);
870 if (bUseLossyCompression)
872 compression = TNG_TNG_COMPRESSION;
876 compression = TNG_GZIP_COMPRESSION;
879 /* The writing is done using write_data, which writes float or double
880 * depending on the GROMACS compilation. */
883 GMX_ASSERT(box, "Need a non-NULL box if positions are written");
885 if (write_data(tng, step, elapsedSeconds,
886 reinterpret_cast<const real *>(x),
887 3, TNG_TRAJ_POSITIONS, "POSITIONS",
888 TNG_PARTICLE_BLOCK_DATA,
889 compression) != TNG_SUCCESS)
891 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
897 if (write_data(tng, step, elapsedSeconds,
898 reinterpret_cast<const real *>(v),
899 3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
900 TNG_PARTICLE_BLOCK_DATA,
901 compression) != TNG_SUCCESS)
903 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
909 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
910 * compression for forces regardless of output mode */
911 if (write_data(tng, step, elapsedSeconds,
912 reinterpret_cast<const real *>(f),
913 3, TNG_TRAJ_FORCES, "FORCES",
914 TNG_PARTICLE_BLOCK_DATA,
915 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
917 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
923 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
924 * compression for box shape regardless of output mode */
925 if (write_data(tng, step, elapsedSeconds,
926 reinterpret_cast<const real *>(box),
927 9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
928 TNG_NON_PARTICLE_BLOCK_DATA,
929 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
931 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
937 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
938 * compression for lambda regardless of output mode */
939 if (write_data(tng, step, elapsedSeconds,
940 reinterpret_cast<const real *>(&lambda),
941 1, TNG_GMX_LAMBDA, "LAMBDAS",
942 TNG_NON_PARTICLE_BLOCK_DATA,
943 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
945 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
949 // Update the data in the wrapper
951 gmx_tng->lastStepDataIsValid = true;
952 gmx_tng->lastStep = step;
953 gmx_tng->lastTimeDataIsValid = true;
954 gmx_tng->lastTime = elapsedSeconds;
956 GMX_UNUSED_VALUE(gmx_tng);
957 GMX_UNUSED_VALUE(bUseLossyCompression);
958 GMX_UNUSED_VALUE(step);
959 GMX_UNUSED_VALUE(elapsedPicoSeconds);
960 GMX_UNUSED_VALUE(lambda);
961 GMX_UNUSED_VALUE(box);
962 GMX_UNUSED_VALUE(nAtoms);
969 void fflush_tng(gmx_tng_trajectory_t gmx_tng)
976 tng_frame_set_premature_write(gmx_tng->tng, TNG_USE_HASH);
978 GMX_UNUSED_VALUE(gmx_tng);
982 float gmx_tng_get_time_of_final_frame(gmx_tng_trajectory_t gmx_tng)
988 tng_trajectory_t tng = gmx_tng->tng;
990 tng_num_frames_get(tng, &nFrames);
991 tng_util_time_of_frame_get(tng, nFrames - 1, &time);
996 GMX_UNUSED_VALUE(gmx_tng);
1001 void gmx_prepare_tng_writing(const char *filename,
1003 gmx_tng_trajectory_t *gmx_tng_input,
1004 gmx_tng_trajectory_t *gmx_tng_output,
1006 const gmx_mtop_t *mtop,
1008 const char *indexGroupName)
1011 tng_trajectory_t *input = (gmx_tng_input && *gmx_tng_input) ? &(*gmx_tng_input)->tng : nullptr;
1012 /* FIXME after 5.0: Currently only standard block types are read */
1013 const int defaultNumIds = 5;
1014 static int64_t fallbackIds[defaultNumIds] =
1016 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
1017 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
1020 static char fallbackNames[defaultNumIds][32] =
1022 "BOX SHAPE", "POSITIONS", "VELOCITIES",
1026 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
1034 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
1036 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
1039 gmx_tng_open(filename, mode, gmx_tng_output);
1040 tng_trajectory_t *output = &(*gmx_tng_output)->tng;
1042 /* Do we have an input file in TNG format? If so, then there's
1043 more data we can copy over, rather than having to improvise. */
1044 if (gmx_tng_input && *gmx_tng_input)
1046 /* Set parameters (compression, time per frame, molecule
1047 * information, number of frames per frame set and writing
1048 * intervals of positions, box shape and lambdas) of the
1049 * output tng container based on their respective values int
1050 * the input tng container */
1051 double time, compression_precision;
1052 int64_t n_frames_per_frame_set, interval = -1;
1054 tng_compression_precision_get(*input, &compression_precision);
1055 tng_compression_precision_set(*output, compression_precision);
1056 // TODO make this configurable in a future version
1057 char compression_type = TNG_TNG_COMPRESSION;
1059 tng_molecule_system_copy(*input, *output);
1061 tng_time_per_frame_get(*input, &time);
1062 tng_time_per_frame_set(*output, time);
1063 // Since we have copied the value from the input TNG we should not change it again
1064 (*gmx_tng_output)->timePerFrameIsSet = true;
1066 tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
1067 tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
1069 for (int i = 0; i < defaultNumIds; i++)
1071 if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval)
1074 switch (fallbackIds[i])
1076 case TNG_TRAJ_POSITIONS:
1077 case TNG_TRAJ_VELOCITIES:
1078 set_writing_interval(*output, interval, 3, fallbackIds[i],
1079 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
1082 case TNG_TRAJ_FORCES:
1083 set_writing_interval(*output, interval, 3, fallbackIds[i],
1084 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
1085 TNG_GZIP_COMPRESSION);
1087 case TNG_TRAJ_BOX_SHAPE:
1088 set_writing_interval(*output, interval, 9, fallbackIds[i],
1089 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
1090 TNG_GZIP_COMPRESSION);
1091 (*gmx_tng_output)->boxOutputInterval = interval;
1093 case TNG_GMX_LAMBDA:
1094 set_writing_interval(*output, interval, 1, fallbackIds[i],
1095 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
1096 TNG_GZIP_COMPRESSION);
1097 (*gmx_tng_output)->lambdaOutputInterval = interval;
1108 /* TODO after trjconv is modularized: fix this so the user can
1109 change precision when they are doing an operation where
1110 this makes sense, and not otherwise.
1112 char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
1113 gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
1115 gmx_tng_add_mtop(*gmx_tng_output, mtop);
1116 tng_num_frames_per_frame_set_set(*output, 1);
1119 if (index && nAtoms > 0)
1121 gmx_tng_setup_atom_subgroup(*gmx_tng_output, nAtoms, index, indexGroupName);
1124 /* If for some reason there are more requested atoms than there are atoms in the
1125 * molecular system create a number of implicit atoms (without atom data) to
1126 * compensate for that. */
1129 tng_implicit_num_particles_set(*output, nAtoms);
1132 GMX_UNUSED_VALUE(filename);
1133 GMX_UNUSED_VALUE(mode);
1134 GMX_UNUSED_VALUE(gmx_tng_input);
1135 GMX_UNUSED_VALUE(gmx_tng_output);
1136 GMX_UNUSED_VALUE(nAtoms);
1137 GMX_UNUSED_VALUE(mtop);
1138 GMX_UNUSED_VALUE(index);
1139 GMX_UNUSED_VALUE(indexGroupName);
1143 void gmx_write_tng_from_trxframe(gmx_tng_trajectory_t gmx_tng_output,
1144 const t_trxframe *frame,
1150 natoms = frame->natoms;
1152 gmx_fwrite_tng(gmx_tng_output,
1163 GMX_UNUSED_VALUE(gmx_tng_output);
1164 GMX_UNUSED_VALUE(frame);
1165 GMX_UNUSED_VALUE(natoms);
1174 convert_array_to_real_array(void *from,
1179 const char datatype)
1183 const bool useDouble = GMX_DOUBLE;
1186 case TNG_FLOAT_DATA:
1191 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1195 for (i = 0; i < nAtoms; i++)
1197 for (j = 0; j < nValues; j++)
1199 to[i*nValues+j] = reinterpret_cast<float *>(from)[i*nValues+j] * fact;
1206 for (i = 0; i < nAtoms; i++)
1208 for (j = 0; j < nValues; j++)
1210 to[i*nValues+j] = reinterpret_cast<float *>(from)[i*nValues+j] * fact;
1216 for (i = 0; i < nAtoms; i++)
1218 for (j = 0; j < nValues; j++)
1220 to[i*nValues+j] = reinterpret_cast<int64_t *>(from)[i*nValues+j] * fact;
1224 case TNG_DOUBLE_DATA:
1225 if (sizeof(real) == sizeof(double))
1229 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1233 for (i = 0; i < nAtoms; i++)
1235 for (j = 0; j < nValues; j++)
1237 to[i*nValues+j] = reinterpret_cast<double *>(from)[i*nValues+j] * fact;
1244 for (i = 0; i < nAtoms; i++)
1246 for (j = 0; j < nValues; j++)
1248 to[i*nValues+j] = reinterpret_cast<double *>(from)[i*nValues+j] * fact;
1254 gmx_incons("Illegal datatype when converting values to a real array!");
1258 real getDistanceScaleFactor(gmx_tng_trajectory_t in)
1261 real distanceScaleFactor;
1263 // TODO Hopefully, TNG 2.0 will do this kind of thing for us
1264 tng_distance_unit_exponential_get(in->tng, &exp);
1266 // GROMACS expects distances in nm
1270 distanceScaleFactor = NANO/NANO;
1273 distanceScaleFactor = NANO/ANGSTROM;
1276 distanceScaleFactor = pow(10.0, exp + 9.0);
1279 return distanceScaleFactor;
1285 void gmx_tng_setup_atom_subgroup(gmx_tng_trajectory_t gmx_tng,
1291 int64_t nAtoms, cnt, nMols;
1292 tng_molecule_t mol, iterMol;
1296 tng_function_status stat;
1297 tng_trajectory_t tng = gmx_tng->tng;
1299 tng_num_particles_get(tng, &nAtoms);
1306 stat = tng_molecule_find(tng, name, -1, &mol);
1307 if (stat == TNG_SUCCESS)
1309 tng_molecule_num_atoms_get(tng, mol, &nAtoms);
1310 tng_molecule_cnt_get(tng, mol, &cnt);
1320 if (stat == TNG_FAILURE)
1322 /* The indexed atoms are added to one separate molecule. */
1323 tng_molecule_alloc(tng, &mol);
1324 tng_molecule_name_set(tng, mol, name);
1325 tng_molecule_chain_add(tng, mol, "", &chain);
1327 for (int i = 0; i < nind; i++)
1329 char temp_name[256], temp_type[256];
1331 /* Try to retrieve the residue name of the atom */
1332 stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1333 if (stat != TNG_SUCCESS)
1335 temp_name[0] = '\0';
1337 /* Check if the molecule of the selection already contains this residue */
1338 if (tng_chain_residue_find(tng, chain, temp_name, -1, &res)
1341 tng_chain_residue_add(tng, chain, temp_name, &res);
1343 /* Try to find the original name and type of the atom */
1344 stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1345 if (stat != TNG_SUCCESS)
1347 temp_name[0] = '\0';
1349 stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
1350 if (stat != TNG_SUCCESS)
1352 temp_type[0] = '\0';
1354 tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
1356 tng_molecule_existing_add(tng, &mol);
1358 /* Set the count of the molecule containing the selected atoms to 1 and all
1359 * other molecules to 0 */
1360 tng_molecule_cnt_set(tng, mol, 1);
1361 tng_num_molecule_types_get(tng, &nMols);
1362 for (int64_t k = 0; k < nMols; k++)
1364 tng_molecule_of_index_get(tng, k, &iterMol);
1369 tng_molecule_cnt_set(tng, iterMol, 0);
1372 GMX_UNUSED_VALUE(gmx_tng);
1373 GMX_UNUSED_VALUE(nind);
1374 GMX_UNUSED_VALUE(ind);
1375 GMX_UNUSED_VALUE(name);
1379 /* TODO: If/when TNG acquires the ability to copy data blocks without
1380 * uncompressing them, then this implemenation should be reconsidered.
1381 * Ideally, gmx trjconv -f a.tng -o b.tng -b 10 -e 20 would be fast
1382 * and lose no information. */
1383 gmx_bool gmx_read_next_tng_frame(gmx_tng_trajectory_t gmx_tng_input,
1385 int64_t *requestedIds,
1386 int numRequestedIds)
1389 tng_trajectory_t input = gmx_tng_input->tng;
1390 gmx_bool bOK = TRUE;
1391 tng_function_status stat;
1392 int64_t numberOfAtoms = -1, frameNumber = -1;
1393 int64_t nBlocks, blockId, *blockIds = nullptr, codecId;
1395 void *values = nullptr;
1396 double frameTime = -1.0;
1397 int size, blockDependency;
1399 const int defaultNumIds = 5;
1400 static int64_t fallbackRequestedIds[defaultNumIds] =
1402 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
1403 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
1410 fr->bLambda = FALSE;
1418 /* If no specific IDs were requested read all block types that can
1419 * currently be interpreted */
1420 if (!requestedIds || numRequestedIds == 0)
1422 numRequestedIds = defaultNumIds;
1423 requestedIds = fallbackRequestedIds;
1426 stat = tng_num_particles_get(input, &numberOfAtoms);
1427 if (stat != TNG_SUCCESS)
1429 gmx_file("Cannot determine number of atoms from TNG file.");
1431 fr->natoms = numberOfAtoms;
1433 bool nextFrameExists = gmx_get_tng_data_block_types_of_next_frame(gmx_tng_input,
1440 gmx::unique_cptr<int64_t, gmx::free_wrapper> blockIdsGuard(blockIds);
1441 if (!nextFrameExists)
1451 for (int64_t i = 0; i < nBlocks; i++)
1453 blockId = blockIds[i];
1454 tng_data_block_dependency_get(input, blockId, &blockDependency);
1455 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1457 stat = tng_util_particle_data_next_frame_read(input,
1466 stat = tng_util_non_particle_data_next_frame_read(input,
1473 if (stat == TNG_CRITICAL)
1475 gmx_file("Cannot read positions from TNG file.");
1478 else if (stat == TNG_FAILURE)
1484 case TNG_TRAJ_BOX_SHAPE:
1488 size = sizeof(int64_t);
1490 case TNG_FLOAT_DATA:
1491 size = sizeof(float);
1493 case TNG_DOUBLE_DATA:
1494 size = sizeof(double);
1497 gmx_incons("Illegal datatype of box shape values!");
1499 for (int i = 0; i < DIM; i++)
1501 convert_array_to_real_array(reinterpret_cast<char *>(values) + size * i * DIM,
1502 reinterpret_cast<real *>(fr->box[i]),
1503 getDistanceScaleFactor(gmx_tng_input),
1510 case TNG_TRAJ_POSITIONS:
1511 srenew(fr->x, fr->natoms);
1512 convert_array_to_real_array(values,
1513 reinterpret_cast<real *>(fr->x),
1514 getDistanceScaleFactor(gmx_tng_input),
1519 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1520 /* This must be updated if/when more lossy compression methods are added */
1521 if (codecId == TNG_TNG_COMPRESSION)
1527 case TNG_TRAJ_VELOCITIES:
1528 srenew(fr->v, fr->natoms);
1529 convert_array_to_real_array(values,
1530 reinterpret_cast<real *>(fr->v),
1531 getDistanceScaleFactor(gmx_tng_input),
1536 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1537 /* This must be updated if/when more lossy compression methods are added */
1538 if (codecId == TNG_TNG_COMPRESSION)
1544 case TNG_TRAJ_FORCES:
1545 srenew(fr->f, fr->natoms);
1546 convert_array_to_real_array(values,
1547 reinterpret_cast<real *>(fr->f),
1548 getDistanceScaleFactor(gmx_tng_input),
1554 case TNG_GMX_LAMBDA:
1557 case TNG_FLOAT_DATA:
1558 fr->lambda = *(reinterpret_cast<float *>(values));
1560 case TNG_DOUBLE_DATA:
1561 fr->lambda = *(reinterpret_cast<double *>(values));
1564 gmx_incons("Illegal datatype lambda value!");
1569 gmx_warning("Illegal block type! Currently GROMACS tools can only handle certain data types. Skipping block.");
1571 /* values does not have to be freed before reading next frame. It will
1572 * be reallocated if it is not NULL. */
1575 fr->step = frameNumber;
1578 // Convert the time to ps
1579 fr->time = frameTime / PICO;
1580 fr->bTime = (frameTime > 0);
1582 // TODO This does not leak, but is not exception safe.
1583 /* values must be freed before leaving this function */
1588 GMX_UNUSED_VALUE(gmx_tng_input);
1589 GMX_UNUSED_VALUE(fr);
1590 GMX_UNUSED_VALUE(requestedIds);
1591 GMX_UNUSED_VALUE(numRequestedIds);
1596 void gmx_print_tng_molecule_system(gmx_tng_trajectory_t gmx_tng_input,
1600 int64_t nMolecules, nChains, nResidues, nAtoms, nFramesRead;
1601 int64_t strideLength, nParticlesRead, nValuesPerFrameRead, *molCntList;
1602 tng_molecule_t molecule;
1604 tng_residue_t residue;
1606 tng_function_status stat;
1610 void *data = nullptr;
1611 std::vector<real> atomCharges;
1612 std::vector<real> atomMasses;
1613 tng_trajectory_t input = gmx_tng_input->tng;
1615 tng_num_molecule_types_get(input, &nMolecules);
1616 tng_molecule_cnt_list_get(input, &molCntList);
1617 /* Can the number of particles change in the trajectory or is it constant? */
1618 tng_num_particles_variable_get(input, &varNAtoms);
1620 for (int64_t i = 0; i < nMolecules; i++)
1622 tng_molecule_of_index_get(input, i, &molecule);
1623 tng_molecule_name_get(input, molecule, str, 256);
1624 if (varNAtoms == TNG_CONSTANT_N_ATOMS)
1626 if (static_cast<int>(molCntList[i]) == 0)
1630 fprintf(stream, "Molecule: %s, count: %d\n", str, static_cast<int>(molCntList[i]));
1634 fprintf(stream, "Molecule: %s\n", str);
1636 tng_molecule_num_chains_get(input, molecule, &nChains);
1639 for (int64_t j = 0; j < nChains; j++)
1641 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
1642 tng_chain_name_get(input, chain, str, 256);
1643 fprintf(stream, "\tChain: %s\n", str);
1644 tng_chain_num_residues_get(input, chain, &nResidues);
1645 for (int64_t k = 0; k < nResidues; k++)
1647 tng_chain_residue_of_index_get(input, chain, k, &residue);
1648 tng_residue_name_get(input, residue, str, 256);
1649 fprintf(stream, "\t\tResidue: %s\n", str);
1650 tng_residue_num_atoms_get(input, residue, &nAtoms);
1651 for (int64_t l = 0; l < nAtoms; l++)
1653 tng_residue_atom_of_index_get(input, residue, l, &atom);
1654 tng_atom_name_get(input, atom, str, 256);
1655 fprintf(stream, "\t\t\tAtom: %s", str);
1656 tng_atom_type_get(input, atom, str, 256);
1657 fprintf(stream, " (%s)\n", str);
1662 /* It is possible to have a molecule without chains, in which case
1663 * residues in the molecule can be iterated through without going
1664 * through chains. */
1667 tng_molecule_num_residues_get(input, molecule, &nResidues);
1670 for (int64_t k = 0; k < nResidues; k++)
1672 tng_molecule_residue_of_index_get(input, molecule, k, &residue);
1673 tng_residue_name_get(input, residue, str, 256);
1674 fprintf(stream, "\t\tResidue: %s\n", str);
1675 tng_residue_num_atoms_get(input, residue, &nAtoms);
1676 for (int64_t l = 0; l < nAtoms; l++)
1678 tng_residue_atom_of_index_get(input, residue, l, &atom);
1679 tng_atom_name_get(input, atom, str, 256);
1680 fprintf(stream, "\t\t\tAtom: %s", str);
1681 tng_atom_type_get(input, atom, str, 256);
1682 fprintf(stream, " (%s)\n", str);
1688 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
1689 for (int64_t l = 0; l < nAtoms; l++)
1691 tng_molecule_atom_of_index_get(input, molecule, l, &atom);
1692 tng_atom_name_get(input, atom, str, 256);
1693 fprintf(stream, "\t\t\tAtom: %s", str);
1694 tng_atom_type_get(input, atom, str, 256);
1695 fprintf(stream, " (%s)\n", str);
1701 tng_num_particles_get(input, &nAtoms);
1702 stat = tng_particle_data_vector_get(input, TNG_TRAJ_PARTIAL_CHARGES, &data, &nFramesRead,
1703 &strideLength, &nParticlesRead,
1704 &nValuesPerFrameRead, &datatype);
1705 if (stat == TNG_SUCCESS)
1707 atomCharges.resize(nAtoms);
1708 convert_array_to_real_array(data,
1715 fprintf(stream, "Atom Charges (%d):\n", int(nAtoms));
1716 for (int64_t i = 0; i < nAtoms; i += 10)
1718 fprintf(stream, "Atom Charges [%8d-]=[", int(i));
1719 for (int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1721 fprintf(stream, " %12.5e", atomCharges[i + j]);
1723 fprintf(stream, "]\n");
1727 stat = tng_particle_data_vector_get(input, TNG_TRAJ_MASSES, &data, &nFramesRead,
1728 &strideLength, &nParticlesRead,
1729 &nValuesPerFrameRead, &datatype);
1730 if (stat == TNG_SUCCESS)
1732 atomMasses.resize(nAtoms);
1733 convert_array_to_real_array(data,
1740 fprintf(stream, "Atom Masses (%d):\n", int(nAtoms));
1741 for (int64_t i = 0; i < nAtoms; i += 10)
1743 fprintf(stream, "Atom Masses [%8d-]=[", int(i));
1744 for (int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1746 fprintf(stream, " %12.5e", atomMasses[i + j]);
1748 fprintf(stream, "]\n");
1754 GMX_UNUSED_VALUE(gmx_tng_input);
1755 GMX_UNUSED_VALUE(stream);
1759 gmx_bool gmx_get_tng_data_block_types_of_next_frame(gmx_tng_trajectory_t gmx_tng_input,
1762 int64_t *requestedIds,
1768 tng_function_status stat;
1769 tng_trajectory_t input = gmx_tng_input->tng;
1771 stat = tng_util_trajectory_next_frame_present_data_blocks_find(input, frame,
1772 nRequestedIds, requestedIds,
1776 if (stat == TNG_CRITICAL)
1778 gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
1780 else if (stat == TNG_FAILURE)
1786 GMX_UNUSED_VALUE(gmx_tng_input);
1787 GMX_UNUSED_VALUE(frame);
1788 GMX_UNUSED_VALUE(nRequestedIds);
1789 GMX_UNUSED_VALUE(requestedIds);
1790 GMX_UNUSED_VALUE(nextFrame);
1791 GMX_UNUSED_VALUE(nBlocks);
1792 GMX_UNUSED_VALUE(blockIds);
1797 gmx_bool gmx_get_tng_data_next_frame_of_block_type(gmx_tng_trajectory_t gmx_tng_input,
1800 int64_t *frameNumber,
1802 int64_t *nValuesPerFrame,
1810 tng_function_status stat;
1813 int blockDependency;
1814 void *data = nullptr;
1816 tng_trajectory_t input = gmx_tng_input->tng;
1818 stat = tng_data_block_name_get(input, blockId, name, maxLen);
1819 if (stat != TNG_SUCCESS)
1821 gmx_file("Cannot read next frame of TNG file");
1823 stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
1824 if (stat != TNG_SUCCESS)
1826 gmx_file("Cannot read next frame of TNG file");
1828 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1830 tng_num_particles_get(input, nAtoms);
1831 stat = tng_util_particle_data_next_frame_read(input,
1840 *nAtoms = 1; /* There are not actually any atoms, but it is used for
1841 allocating memory */
1842 stat = tng_util_non_particle_data_next_frame_read(input,
1849 if (stat == TNG_CRITICAL)
1851 gmx_file("Cannot read next frame of TNG file");
1853 if (stat == TNG_FAILURE)
1859 stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
1860 if (stat != TNG_SUCCESS)
1862 gmx_file("Cannot read next frame of TNG file");
1864 srenew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
1865 convert_array_to_real_array(data,
1867 getDistanceScaleFactor(gmx_tng_input),
1872 tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
1874 /* This must be updated if/when more lossy compression methods are added */
1875 if (codecId != TNG_TNG_COMPRESSION)
1889 GMX_UNUSED_VALUE(gmx_tng_input);
1890 GMX_UNUSED_VALUE(blockId);
1891 GMX_UNUSED_VALUE(values);
1892 GMX_UNUSED_VALUE(frameNumber);
1893 GMX_UNUSED_VALUE(frameTime);
1894 GMX_UNUSED_VALUE(nValuesPerFrame);
1895 GMX_UNUSED_VALUE(nAtoms);
1896 GMX_UNUSED_VALUE(prec);
1897 GMX_UNUSED_VALUE(name);
1898 GMX_UNUSED_VALUE(maxLen);
1899 GMX_UNUSED_VALUE(bOK);
1904 int gmx_tng_get_box_output_interval(gmx_tng_trajectory_t gmx_tng)
1907 return gmx_tng->boxOutputInterval;
1909 GMX_UNUSED_VALUE(gmx_tng);
1914 int gmx_tng_get_lambda_output_interval(gmx_tng_trajectory_t gmx_tng)
1917 return gmx_tng->lambdaOutputInterval;
1919 GMX_UNUSED_VALUE(gmx_tng);