2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014, 2015 by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
46 #include "tng/tng_io.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/topology/ifunc.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/trajectory/trajectoryframe.h"
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/baseversion.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/programcontext.h"
61 #include "gromacs/utility/smalloc.h"
62 #include "gromacs/utility/sysinfo.h"
64 static const char *modeToVerb(char mode)
79 gmx_fatal(FARGS, "Invalid file opening mode %c", mode);
86 void gmx_tng_open(const char *filename,
88 tng_trajectory_t *tng)
91 /* First check whether we have to make a backup,
92 * only for writing, not for read or append.
96 make_backup(filename);
99 /* tng must not be pointing at already allocated memory.
100 * Memory will be allocated by tng_util_trajectory_open() and must
101 * later on be freed by tng_util_trajectory_close(). */
102 if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
104 /* TNG does return more than one degree of error, but there is
105 no use case for GROMACS handling the non-fatal errors
108 "File I/O error while opening %s for %s",
113 if (mode == 'w' || mode == 'a')
116 gmx_gethostname(hostname, 256);
119 tng_first_computer_name_set(*tng, hostname);
123 tng_last_computer_name_set(*tng, hostname);
126 char programInfo[256];
127 const char *precisionString = "";
129 precisionString = " (double precision)";
131 sprintf(programInfo, "%.100s %.128s%.24s",
132 gmx::getProgramContext().displayName(),
133 gmx_version(), precisionString);
136 tng_first_program_name_set(*tng, programInfo);
140 tng_last_program_name_set(*tng, programInfo);
144 if (!gmx_getusername(username, 256))
148 tng_first_user_name_set(*tng, username);
152 tng_last_user_name_set(*tng, username);
153 tng_file_headers_write(*tng, TNG_USE_HASH);
158 gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
159 GMX_UNUSED_VALUE(filename);
160 GMX_UNUSED_VALUE(mode);
161 GMX_UNUSED_VALUE(tng);
165 void gmx_tng_close(tng_trajectory_t *tng)
167 /* We have to check that tng is set because
168 * tng_util_trajectory_close wants to return a NULL in it, and
169 * gives a fatal error if it is NULL. */
173 tng_util_trajectory_close(tng);
176 GMX_UNUSED_VALUE(tng);
181 static void addTngMoleculeFromTopology(tng_trajectory_t tng,
182 const char *moleculeName,
183 const t_atoms *atoms,
184 gmx_int64_t numMolecules,
185 tng_molecule_t *tngMol)
187 tng_chain_t tngChain = nullptr;
188 tng_residue_t tngRes = nullptr;
190 if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
192 gmx_file("Cannot add molecule to TNG molecular system.");
195 /* FIXME: The TNG atoms should contain mass and atomB info (for free
196 * energy calculations), i.e. in when it's available in TNG (2.0). */
197 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++)
199 const t_atom *at = &atoms->atom[atomIndex];
200 /* FIXME: Currently the TNG API can only add atoms belonging to a
201 * residue and chain. Wait for TNG 2.0*/
204 const t_resinfo *resInfo = &atoms->resinfo[at->resind];
205 char chainName[2] = {resInfo->chainid, 0};
206 tng_atom_t tngAtom = nullptr;
211 prevAtom = &atoms->atom[atomIndex - 1];
218 /* If this is the first atom or if the residue changed add the
219 * residue to the TNG molecular system. */
220 if (!prevAtom || resInfo != &atoms->resinfo[prevAtom->resind])
222 /* If this is the first atom or if the chain changed add
223 * the chain to the TNG molecular system. */
224 if (!prevAtom || resInfo->chainid !=
225 atoms->resinfo[prevAtom->resind].chainid)
227 tng_molecule_chain_add(tng, *tngMol, chainName,
230 /* FIXME: When TNG supports both residue index and residue
231 * number the latter should be used. Wait for TNG 2.0*/
232 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
234 tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIndex]), *(atoms->atomtype[atomIndex]), &tngAtom);
237 tng_molecule_cnt_set(tng, *tngMol, numMolecules);
240 void gmx_tng_add_mtop(tng_trajectory_t tng,
241 const gmx_mtop_t *mtop)
244 const t_ilist *ilist;
249 /* No topology information available to add. */
253 for (int molIndex = 0; molIndex < mtop->nmolblock; molIndex++)
255 tng_molecule_t tngMol = nullptr;
256 const gmx_moltype_t *molType =
257 &mtop->moltype[mtop->molblock[molIndex].type];
259 /* Add a molecule to the TNG trajectory with the same name as the
260 * current molecule. */
261 addTngMoleculeFromTopology(tng,
264 mtop->molblock[molIndex].nmol,
267 /* Bonds have to be deduced from interactions (constraints etc). Different
268 * interactions have different sets of parameters. */
269 /* Constraints are specified using two atoms */
270 for (i = 0; i < F_NRE; i++)
274 ilist = &molType->ilist[i];
278 while (j < ilist->nr)
280 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
286 /* Settle is described using three atoms */
287 ilist = &molType->ilist[F_SETTLE];
291 while (j < ilist->nr)
293 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
294 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+2], &tngBond);
301 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
302 * if they are positive.
304 * If only one of n1 and n2 is positive, then return it.
305 * If neither n1 or n2 is positive, then return -1. */
307 greatest_common_divisor_if_positive(int n1, int n2)
311 return (0 >= n2) ? -1 : n2;
318 /* We have a non-trivial greatest common divisor to compute. */
319 return gmx_greatest_common_divisor(n1, n2);
322 /* By default try to write 100 frames (of actual output) in each frame set.
323 * This number is the number of outputs of the most frequently written data
324 * type per frame set.
325 * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
326 * setups regarding compression efficiency and compression time. Make this
327 * a hidden command-line option? */
328 const int defaultFramesPerFrameSet = 100;
330 /*! \libinternal \brief Set the number of frames per frame
331 * set according to output intervals.
332 * The default is that 100 frames are written of the data
333 * that is written most often. */
334 static void tng_set_frames_per_frame_set(tng_trajectory_t tng,
335 const gmx_bool bUseLossyCompression,
336 const t_inputrec *ir)
340 /* Set the number of frames per frame set to contain at least
341 * defaultFramesPerFrameSet of the lowest common denominator of
342 * the writing interval of positions and velocities. */
343 /* FIXME after 5.0: consider nstenergy also? */
344 if (bUseLossyCompression)
346 gcd = ir->nstxout_compressed;
350 gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
351 gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
358 tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
361 /*! \libinternal \brief Set the data-writing intervals, and number of
362 * frames per frame set */
363 static void set_writing_intervals(tng_trajectory_t tng,
364 const gmx_bool bUseLossyCompression,
365 const t_inputrec *ir)
367 /* Define pointers to specific writing functions depending on if we
368 * write float or double data */
369 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
377 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
379 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
381 int xout, vout, fout;
382 int gcd = -1, lowest = -1;
385 tng_set_frames_per_frame_set(tng, bUseLossyCompression, ir);
387 if (bUseLossyCompression)
389 xout = ir->nstxout_compressed;
391 /* If there is no uncompressed coordinate output write forces
392 and velocities to the compressed tng file. */
403 compression = TNG_TNG_COMPRESSION;
410 compression = TNG_GZIP_COMPRESSION;
414 set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS,
415 "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
417 /* TODO: if/when we write energies to TNG also, reconsider how
418 * and when box information is written, because GROMACS
419 * behaviour pre-5.0 was to write the box with every
420 * trajectory frame and every energy frame, and probably
421 * people depend on this. */
423 gcd = greatest_common_divisor_if_positive(gcd, xout);
424 if (lowest < 0 || xout < lowest)
431 set_writing_interval(tng, vout, 3, TNG_TRAJ_VELOCITIES,
432 "VELOCITIES", TNG_PARTICLE_BLOCK_DATA,
435 gcd = greatest_common_divisor_if_positive(gcd, vout);
436 if (lowest < 0 || vout < lowest)
443 set_writing_interval(tng, fout, 3, TNG_TRAJ_FORCES,
444 "FORCES", TNG_PARTICLE_BLOCK_DATA,
445 TNG_GZIP_COMPRESSION);
447 gcd = greatest_common_divisor_if_positive(gcd, fout);
448 if (lowest < 0 || fout < lowest)
455 /* Lambdas and box shape written at an interval of the lowest common
456 denominator of other output */
457 set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA,
458 "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
459 TNG_GZIP_COMPRESSION);
461 set_writing_interval(tng, gcd, 9, TNG_TRAJ_BOX_SHAPE,
462 "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA,
463 TNG_GZIP_COMPRESSION);
464 if (gcd < lowest / 10)
466 gmx_warning("The lowest common denominator of trajectory output is "
467 "every %d step(s), whereas the shortest output interval "
468 "is every %d steps.", gcd, lowest);
474 void gmx_tng_prepare_md_writing(tng_trajectory_t tng,
475 const gmx_mtop_t *mtop,
476 const t_inputrec *ir)
479 gmx_tng_add_mtop(tng, mtop);
480 set_writing_intervals(tng, FALSE, ir);
481 tng_time_per_frame_set(tng, ir->delta_t * PICO);
483 GMX_UNUSED_VALUE(tng);
484 GMX_UNUSED_VALUE(mtop);
485 GMX_UNUSED_VALUE(ir);
490 /* Check if all atoms in the molecule system are selected
491 * by a selection group of type specified by gtype. */
492 static gmx_bool all_atoms_selected(const gmx_mtop_t *mtop,
495 const gmx_moltype_t *molType;
496 const t_atoms *atoms;
498 /* Iterate through all atoms in the molecule system and
499 * check if they belong to a selection group of the
501 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
503 molType = &mtop->moltype[mtop->molblock[molIndex].type];
505 atoms = &molType->atoms;
507 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
509 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
511 if (ggrpnr(&mtop->groups, gtype, i) != 0)
521 /* Create TNG molecules which will represent each of the selection
522 * groups for writing. But do that only if there is actually a
523 * specified selection group and it is not the whole system.
524 * TODO: Currently the only selection that is taken into account
525 * is egcCompressedX, but other selections should be added when
526 * e.g. writing energies is implemented.
528 static void add_selection_groups(tng_trajectory_t tng,
529 const gmx_mtop_t *mtop)
531 const gmx_moltype_t *molType;
532 const t_atoms *atoms;
534 const t_resinfo *resInfo;
535 const t_ilist *ilist;
538 tng_molecule_t mol, iterMol;
546 /* TODO: When the TNG molecules block is more flexible TNG selection
547 * groups should not need all atoms specified. It should be possible
548 * just to specify what molecules are selected (and/or which atoms in
549 * the molecule) and how many (if applicable). */
551 /* If no atoms are selected we do not need to create a
552 * TNG selection group molecule. */
553 if (mtop->groups.ngrpnr[egcCompressedX] == 0)
558 /* If all atoms are selected we do not have to create a selection
559 * group molecule in the TNG molecule system. */
560 if (all_atoms_selected(mtop, egcCompressedX))
565 /* The name of the TNG molecule containing the selection group is the
566 * same as the name of the selection group. */
567 nameIndex = *mtop->groups.grps[egcCompressedX].nm_ind;
568 groupName = *mtop->groups.grpname[nameIndex];
570 tng_molecule_alloc(tng, &mol);
571 tng_molecule_name_set(tng, mol, groupName);
572 tng_molecule_chain_add(tng, mol, "", &chain);
573 for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
575 molType = &mtop->moltype[mtop->molblock[molIndex].type];
577 atoms = &molType->atoms;
579 for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
581 bool bAtomsAdded = FALSE;
582 for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
587 if (ggrpnr(&mtop->groups, egcCompressedX, i) != 0)
591 at = &atoms->atom[atomIndex];
594 resInfo = &atoms->resinfo[at->resind];
595 /* FIXME: When TNG supports both residue index and residue
596 * number the latter should be used. */
597 res_name = *resInfo->name;
598 res_id = at->resind + 1;
602 res_name = (char *)"";
605 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res)
608 /* Since there is ONE chain for selection groups do not keep the
609 * original residue IDs - otherwise there might be conflicts. */
610 tng_chain_residue_add(tng, chain, res_name, &res);
612 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
613 *(atoms->atomtype[atomIndex]),
614 atom_offset + atomIndex, &atom);
620 for (int k = 0; k < F_NRE; k++)
624 ilist = &molType->ilist[k];
628 while (l < ilist->nr)
631 atom1 = ilist->iatoms[l] + atom_offset;
632 atom2 = ilist->iatoms[l+1] + atom_offset;
633 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0 &&
634 ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
636 tng_molecule_bond_add(tng, mol, ilist->iatoms[l],
637 ilist->iatoms[l+1], &tngBond);
644 /* Settle is described using three atoms */
645 ilist = &molType->ilist[F_SETTLE];
649 while (l < ilist->nr)
651 int atom1, atom2, atom3;
652 atom1 = ilist->iatoms[l] + atom_offset;
653 atom2 = ilist->iatoms[l+1] + atom_offset;
654 atom3 = ilist->iatoms[l+2] + atom_offset;
655 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0)
657 if (ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
659 tng_molecule_bond_add(tng, mol, atom1,
662 if (ggrpnr(&mtop->groups, egcCompressedX, atom3) == 0)
664 tng_molecule_bond_add(tng, mol, atom1,
672 atom_offset += atoms->nr;
675 tng_molecule_existing_add(tng, &mol);
676 tng_molecule_cnt_set(tng, mol, 1);
677 tng_num_molecule_types_get(tng, &nMols);
678 for (gmx_int64_t k = 0; k < nMols; k++)
680 tng_molecule_of_index_get(tng, k, &iterMol);
685 tng_molecule_cnt_set(tng, iterMol, 0);
690 void gmx_tng_set_compression_precision(tng_trajectory_t tng,
694 tng_compression_precision_set(tng, prec);
696 GMX_UNUSED_VALUE(tng);
697 GMX_UNUSED_VALUE(prec);
701 void gmx_tng_prepare_low_prec_writing(tng_trajectory_t tng,
702 const gmx_mtop_t *mtop,
703 const t_inputrec *ir)
706 gmx_tng_add_mtop(tng, mtop);
707 add_selection_groups(tng, mtop);
708 set_writing_intervals(tng, TRUE, ir);
709 tng_time_per_frame_set(tng, ir->delta_t * PICO);
710 gmx_tng_set_compression_precision(tng, ir->x_compression_precision);
712 GMX_UNUSED_VALUE(tng);
713 GMX_UNUSED_VALUE(mtop);
714 GMX_UNUSED_VALUE(ir);
718 void gmx_fwrite_tng(tng_trajectory_t tng,
719 const gmx_bool bUseLossyCompression,
721 real elapsedPicoSeconds,
730 typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
740 static write_data_func_pointer write_data = tng_util_generic_with_time_double_write;
742 static write_data_func_pointer write_data = tng_util_generic_with_time_write;
744 double elapsedSeconds = elapsedPicoSeconds * PICO;
745 gmx_int64_t nParticles;
751 /* This function might get called when the type of the
752 compressed trajectory is actually XTC. So we exit and move
757 tng_num_particles_get(tng, &nParticles);
758 if (nAtoms != (int)nParticles)
760 tng_implicit_num_particles_set(tng, nAtoms);
763 if (bUseLossyCompression)
765 compression = TNG_TNG_COMPRESSION;
769 compression = TNG_GZIP_COMPRESSION;
772 /* The writing is done using write_data, which writes float or double
773 * depending on the GROMACS compilation. */
776 GMX_ASSERT(box, "Need a non-NULL box if positions are written");
778 if (write_data(tng, step, elapsedSeconds,
779 reinterpret_cast<const real *>(x),
780 3, TNG_TRAJ_POSITIONS, "POSITIONS",
781 TNG_PARTICLE_BLOCK_DATA,
782 compression) != TNG_SUCCESS)
784 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
790 if (write_data(tng, step, elapsedSeconds,
791 reinterpret_cast<const real *>(v),
792 3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
793 TNG_PARTICLE_BLOCK_DATA,
794 compression) != TNG_SUCCESS)
796 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
802 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
803 * compression for forces regardless of output mode */
804 if (write_data(tng, step, elapsedSeconds,
805 reinterpret_cast<const real *>(f),
806 3, TNG_TRAJ_FORCES, "FORCES",
807 TNG_PARTICLE_BLOCK_DATA,
808 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
810 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
814 /* TNG-MF1 compression only compresses positions and velocities. Use lossless
815 * compression for lambdas and box shape regardless of output mode */
816 if (write_data(tng, step, elapsedSeconds,
817 reinterpret_cast<const real *>(box),
818 9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
819 TNG_NON_PARTICLE_BLOCK_DATA,
820 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
822 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
825 if (write_data(tng, step, elapsedSeconds,
826 reinterpret_cast<const real *>(&lambda),
827 1, TNG_GMX_LAMBDA, "LAMBDAS",
828 TNG_NON_PARTICLE_BLOCK_DATA,
829 TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
831 gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
834 GMX_UNUSED_VALUE(tng);
835 GMX_UNUSED_VALUE(bUseLossyCompression);
836 GMX_UNUSED_VALUE(step);
837 GMX_UNUSED_VALUE(elapsedPicoSeconds);
838 GMX_UNUSED_VALUE(lambda);
839 GMX_UNUSED_VALUE(box);
840 GMX_UNUSED_VALUE(nAtoms);
847 void fflush_tng(tng_trajectory_t tng)
854 tng_frame_set_premature_write(tng, TNG_USE_HASH);
856 GMX_UNUSED_VALUE(tng);
860 float gmx_tng_get_time_of_final_frame(tng_trajectory_t tng)
867 tng_num_frames_get(tng, &nFrames);
868 tng_util_time_of_frame_get(tng, nFrames - 1, &time);
873 GMX_UNUSED_VALUE(tng);
878 void gmx_prepare_tng_writing(const char *filename,
880 tng_trajectory_t *input,
881 tng_trajectory_t *output,
883 const gmx_mtop_t *mtop,
885 const char *indexGroupName)
888 /* FIXME after 5.0: Currently only standard block types are read */
889 const int defaultNumIds = 5;
890 static gmx_int64_t fallbackIds[defaultNumIds] =
892 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
893 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
896 static char fallbackNames[defaultNumIds][32] =
898 "BOX SHAPE", "POSITIONS", "VELOCITIES",
902 typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
910 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
912 set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
915 gmx_tng_open(filename, mode, output);
917 /* Do we have an input file in TNG format? If so, then there's
918 more data we can copy over, rather than having to improvise. */
921 /* Set parameters (compression, time per frame, molecule
922 * information, number of frames per frame set and writing
923 * intervals of positions, box shape and lambdas) of the
924 * output tng container based on their respective values int
925 * the input tng container */
926 double time, compression_precision;
927 gmx_int64_t n_frames_per_frame_set, interval = -1;
929 tng_compression_precision_get(*input, &compression_precision);
930 tng_compression_precision_set(*output, compression_precision);
931 // TODO make this configurable in a future version
932 char compression_type = TNG_TNG_COMPRESSION;
934 tng_molecule_system_copy(*input, *output);
936 tng_time_per_frame_get(*input, &time);
937 tng_time_per_frame_set(*output, time);
939 tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
940 tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
942 for (int i = 0; i < defaultNumIds; i++)
944 if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval)
947 switch (fallbackIds[i])
949 case TNG_TRAJ_POSITIONS:
950 case TNG_TRAJ_VELOCITIES:
951 set_writing_interval(*output, interval, 3, fallbackIds[i],
952 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
955 case TNG_TRAJ_FORCES:
956 set_writing_interval(*output, interval, 3, fallbackIds[i],
957 fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
958 TNG_GZIP_COMPRESSION);
960 case TNG_TRAJ_BOX_SHAPE:
961 set_writing_interval(*output, interval, 9, fallbackIds[i],
962 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
963 TNG_GZIP_COMPRESSION);
966 set_writing_interval(*output, interval, 1, fallbackIds[i],
967 fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
968 TNG_GZIP_COMPRESSION);
979 /* TODO after trjconv is modularized: fix this so the user can
980 change precision when they are doing an operation where
981 this makes sense, and not otherwise.
983 char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
984 gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
986 gmx_tng_add_mtop(*output, mtop);
987 tng_num_frames_per_frame_set_set(*output, 1);
990 if (index && nAtoms > 0)
992 gmx_tng_setup_atom_subgroup(*output, nAtoms, index, indexGroupName);
995 /* If for some reason there are more requested atoms than there are atoms in the
996 * molecular system create a number of implicit atoms (without atom data) to
997 * compensate for that. */
1000 tng_implicit_num_particles_set(*output, nAtoms);
1003 GMX_UNUSED_VALUE(filename);
1004 GMX_UNUSED_VALUE(mode);
1005 GMX_UNUSED_VALUE(input);
1006 GMX_UNUSED_VALUE(output);
1007 GMX_UNUSED_VALUE(nAtoms);
1008 GMX_UNUSED_VALUE(mtop);
1009 GMX_UNUSED_VALUE(index);
1010 GMX_UNUSED_VALUE(indexGroupName);
1014 void gmx_write_tng_from_trxframe(tng_trajectory_t output,
1015 const t_trxframe *frame,
1019 if (frame->step > 0)
1021 double timePerFrame = frame->time * PICO / frame->step;
1022 tng_time_per_frame_set(output, timePerFrame);
1026 natoms = frame->natoms;
1028 gmx_fwrite_tng(output,
1039 GMX_UNUSED_VALUE(output);
1040 GMX_UNUSED_VALUE(frame);
1041 GMX_UNUSED_VALUE(natoms);
1050 convert_array_to_real_array(void *from,
1055 const char datatype)
1059 const bool useDouble = GMX_DOUBLE;
1062 case TNG_FLOAT_DATA:
1067 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1071 for (i = 0; i < nAtoms; i++)
1073 for (j = 0; j < nValues; j++)
1075 to[i*nValues+j] = reinterpret_cast<float *>(from)[i*nValues+j] * fact;
1082 for (i = 0; i < nAtoms; i++)
1084 for (j = 0; j < nValues; j++)
1086 to[i*nValues+j] = reinterpret_cast<float *>(from)[i*nValues+j] * fact;
1092 for (i = 0; i < nAtoms; i++)
1094 for (j = 0; j < nValues; j++)
1096 to[i*nValues+j] = reinterpret_cast<gmx_int64_t *>(from)[i*nValues+j] * fact;
1100 case TNG_DOUBLE_DATA:
1101 if (sizeof(real) == sizeof(double))
1105 memcpy(to, from, nValues * sizeof(real) * nAtoms);
1109 for (i = 0; i < nAtoms; i++)
1111 for (j = 0; j < nValues; j++)
1113 to[i*nValues+j] = reinterpret_cast<double *>(from)[i*nValues+j] * fact;
1120 for (i = 0; i < nAtoms; i++)
1122 for (j = 0; j < nValues; j++)
1124 to[i*nValues+j] = reinterpret_cast<double *>(from)[i*nValues+j] * fact;
1130 gmx_incons("Illegal datatype when converting values to a real array!");
1136 real getDistanceScaleFactor(tng_trajectory_t in)
1138 gmx_int64_t exp = -1;
1139 real distanceScaleFactor;
1141 // TODO Hopefully, TNG 2.0 will do this kind of thing for us
1142 tng_distance_unit_exponential_get(in, &exp);
1144 // GROMACS expects distances in nm
1148 distanceScaleFactor = NANO/NANO;
1151 distanceScaleFactor = NANO/ANGSTROM;
1154 distanceScaleFactor = pow(10.0, exp + 9.0);
1157 return distanceScaleFactor;
1163 void gmx_tng_setup_atom_subgroup(tng_trajectory_t tng,
1169 gmx_int64_t nAtoms, cnt, nMols;
1170 tng_molecule_t mol, iterMol;
1174 tng_function_status stat;
1176 tng_num_particles_get(tng, &nAtoms);
1183 stat = tng_molecule_find(tng, name, -1, &mol);
1184 if (stat == TNG_SUCCESS)
1186 tng_molecule_num_atoms_get(tng, mol, &nAtoms);
1187 tng_molecule_cnt_get(tng, mol, &cnt);
1197 if (stat == TNG_FAILURE)
1199 /* The indexed atoms are added to one separate molecule. */
1200 tng_molecule_alloc(tng, &mol);
1201 tng_molecule_name_set(tng, mol, name);
1202 tng_molecule_chain_add(tng, mol, "", &chain);
1204 for (int i = 0; i < nind; i++)
1206 char temp_name[256], temp_type[256];
1208 /* Try to retrieve the residue name of the atom */
1209 stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1210 if (stat != TNG_SUCCESS)
1212 temp_name[0] = '\0';
1214 /* Check if the molecule of the selection already contains this residue */
1215 if (tng_chain_residue_find(tng, chain, temp_name, -1, &res)
1218 tng_chain_residue_add(tng, chain, temp_name, &res);
1220 /* Try to find the original name and type of the atom */
1221 stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1222 if (stat != TNG_SUCCESS)
1224 temp_name[0] = '\0';
1226 stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
1227 if (stat != TNG_SUCCESS)
1229 temp_type[0] = '\0';
1231 tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
1233 tng_molecule_existing_add(tng, &mol);
1235 /* Set the count of the molecule containing the selected atoms to 1 and all
1236 * other molecules to 0 */
1237 tng_molecule_cnt_set(tng, mol, 1);
1238 tng_num_molecule_types_get(tng, &nMols);
1239 for (gmx_int64_t k = 0; k < nMols; k++)
1241 tng_molecule_of_index_get(tng, k, &iterMol);
1246 tng_molecule_cnt_set(tng, iterMol, 0);
1249 GMX_UNUSED_VALUE(tng);
1250 GMX_UNUSED_VALUE(nind);
1251 GMX_UNUSED_VALUE(ind);
1252 GMX_UNUSED_VALUE(name);
1256 /* TODO: If/when TNG acquires the ability to copy data blocks without
1257 * uncompressing them, then this implemenation should be reconsidered.
1258 * Ideally, gmx trjconv -f a.tng -o b.tng -b 10 -e 20 would be fast
1259 * and lose no information. */
1260 gmx_bool gmx_read_next_tng_frame(tng_trajectory_t input,
1262 gmx_int64_t *requestedIds,
1263 int numRequestedIds)
1266 gmx_bool bOK = TRUE;
1267 tng_function_status stat;
1268 gmx_int64_t numberOfAtoms = -1, frameNumber = -1;
1269 gmx_int64_t nBlocks, blockId, *blockIds = nullptr, codecId;
1271 void *values = nullptr;
1272 double frameTime = -1.0;
1273 int size, blockDependency;
1275 const int defaultNumIds = 5;
1276 static gmx_int64_t fallbackRequestedIds[defaultNumIds] =
1278 TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
1279 TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
1286 fr->bLambda = FALSE;
1294 /* If no specific IDs were requested read all block types that can
1295 * currently be interpreted */
1296 if (!requestedIds || numRequestedIds == 0)
1298 numRequestedIds = defaultNumIds;
1299 requestedIds = fallbackRequestedIds;
1302 stat = tng_num_particles_get(input, &numberOfAtoms);
1303 if (stat != TNG_SUCCESS)
1305 gmx_file("Cannot determine number of atoms from TNG file.");
1307 fr->natoms = numberOfAtoms;
1309 if (!gmx_get_tng_data_block_types_of_next_frame(input,
1325 for (gmx_int64_t i = 0; i < nBlocks; i++)
1327 blockId = blockIds[i];
1328 tng_data_block_dependency_get(input, blockId, &blockDependency);
1329 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1331 stat = tng_util_particle_data_next_frame_read(input,
1340 stat = tng_util_non_particle_data_next_frame_read(input,
1347 if (stat == TNG_CRITICAL)
1349 gmx_file("Cannot read positions from TNG file.");
1352 else if (stat == TNG_FAILURE)
1358 case TNG_TRAJ_BOX_SHAPE:
1362 size = sizeof(gmx_int64_t);
1364 case TNG_FLOAT_DATA:
1365 size = sizeof(float);
1367 case TNG_DOUBLE_DATA:
1368 size = sizeof(double);
1371 gmx_incons("Illegal datatype of box shape values!");
1373 for (int i = 0; i < DIM; i++)
1375 convert_array_to_real_array(reinterpret_cast<char *>(values) + size * i * DIM,
1376 reinterpret_cast<real *>(fr->box[i]),
1377 getDistanceScaleFactor(input),
1384 case TNG_TRAJ_POSITIONS:
1385 srenew(fr->x, fr->natoms);
1386 convert_array_to_real_array(values,
1387 reinterpret_cast<real *>(fr->x),
1388 getDistanceScaleFactor(input),
1393 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1394 /* This must be updated if/when more lossy compression methods are added */
1395 if (codecId == TNG_TNG_COMPRESSION)
1401 case TNG_TRAJ_VELOCITIES:
1402 srenew(fr->v, fr->natoms);
1403 convert_array_to_real_array(values,
1405 getDistanceScaleFactor(input),
1410 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1411 /* This must be updated if/when more lossy compression methods are added */
1412 if (codecId == TNG_TNG_COMPRESSION)
1418 case TNG_TRAJ_FORCES:
1419 srenew(fr->f, fr->natoms);
1420 convert_array_to_real_array(values,
1421 reinterpret_cast<real *>(fr->f),
1422 getDistanceScaleFactor(input),
1428 case TNG_GMX_LAMBDA:
1431 case TNG_FLOAT_DATA:
1432 fr->lambda = *(reinterpret_cast<float *>(values));
1434 case TNG_DOUBLE_DATA:
1435 fr->lambda = *(reinterpret_cast<double *>(values));
1438 gmx_incons("Illegal datatype lambda value!");
1443 gmx_warning("Illegal block type! Currently GROMACS tools can only handle certain data types. Skipping block.");
1445 /* values does not have to be freed before reading next frame. It will
1446 * be reallocated if it is not NULL. */
1449 fr->step = frameNumber;
1451 // Convert the time to ps
1452 fr->time = frameTime / PICO;
1455 /* values must be freed before leaving this function */
1460 GMX_UNUSED_VALUE(input);
1461 GMX_UNUSED_VALUE(fr);
1462 GMX_UNUSED_VALUE(requestedIds);
1463 GMX_UNUSED_VALUE(numRequestedIds);
1468 void gmx_print_tng_molecule_system(tng_trajectory_t input,
1472 gmx_int64_t nMolecules, nChains, nResidues, nAtoms, *molCntList;
1473 tng_molecule_t molecule;
1475 tng_residue_t residue;
1477 char str[256], varNAtoms;
1479 tng_num_molecule_types_get(input, &nMolecules);
1480 tng_molecule_cnt_list_get(input, &molCntList);
1481 /* Can the number of particles change in the trajectory or is it constant? */
1482 tng_num_particles_variable_get(input, &varNAtoms);
1484 for (gmx_int64_t i = 0; i < nMolecules; i++)
1486 tng_molecule_of_index_get(input, i, &molecule);
1487 tng_molecule_name_get(input, molecule, str, 256);
1488 if (varNAtoms == TNG_CONSTANT_N_ATOMS)
1490 if ((int)molCntList[i] == 0)
1494 fprintf(stream, "Molecule: %s, count: %d\n", str, (int)molCntList[i]);
1498 fprintf(stream, "Molecule: %s\n", str);
1500 tng_molecule_num_chains_get(input, molecule, &nChains);
1503 for (gmx_int64_t j = 0; j < nChains; j++)
1505 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
1506 tng_chain_name_get(input, chain, str, 256);
1507 fprintf(stream, "\tChain: %s\n", str);
1508 tng_chain_num_residues_get(input, chain, &nResidues);
1509 for (gmx_int64_t k = 0; k < nResidues; k++)
1511 tng_chain_residue_of_index_get(input, chain, k, &residue);
1512 tng_residue_name_get(input, residue, str, 256);
1513 fprintf(stream, "\t\tResidue: %s\n", str);
1514 tng_residue_num_atoms_get(input, residue, &nAtoms);
1515 for (gmx_int64_t l = 0; l < nAtoms; l++)
1517 tng_residue_atom_of_index_get(input, residue, l, &atom);
1518 tng_atom_name_get(input, atom, str, 256);
1519 fprintf(stream, "\t\t\tAtom: %s", str);
1520 tng_atom_type_get(input, atom, str, 256);
1521 fprintf(stream, " (%s)\n", str);
1526 /* It is possible to have a molecule without chains, in which case
1527 * residues in the molecule can be iterated through without going
1528 * through chains. */
1531 tng_molecule_num_residues_get(input, molecule, &nResidues);
1534 for (gmx_int64_t k = 0; k < nResidues; k++)
1536 tng_molecule_residue_of_index_get(input, molecule, k, &residue);
1537 tng_residue_name_get(input, residue, str, 256);
1538 fprintf(stream, "\t\tResidue: %s\n", str);
1539 tng_residue_num_atoms_get(input, residue, &nAtoms);
1540 for (gmx_int64_t l = 0; l < nAtoms; l++)
1542 tng_residue_atom_of_index_get(input, residue, l, &atom);
1543 tng_atom_name_get(input, atom, str, 256);
1544 fprintf(stream, "\t\t\tAtom: %s", str);
1545 tng_atom_type_get(input, atom, str, 256);
1546 fprintf(stream, " (%s)\n", str);
1552 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
1553 for (gmx_int64_t l = 0; l < nAtoms; l++)
1555 tng_molecule_atom_of_index_get(input, molecule, l, &atom);
1556 tng_atom_name_get(input, atom, str, 256);
1557 fprintf(stream, "\t\t\tAtom: %s", str);
1558 tng_atom_type_get(input, atom, str, 256);
1559 fprintf(stream, " (%s)\n", str);
1565 GMX_UNUSED_VALUE(input);
1566 GMX_UNUSED_VALUE(stream);
1570 gmx_bool gmx_get_tng_data_block_types_of_next_frame(tng_trajectory_t input,
1573 gmx_int64_t *requestedIds,
1574 gmx_int64_t *nextFrame,
1575 gmx_int64_t *nBlocks,
1576 gmx_int64_t **blockIds)
1579 tng_function_status stat;
1581 stat = tng_util_trajectory_next_frame_present_data_blocks_find(input, frame,
1582 nRequestedIds, requestedIds,
1586 if (stat == TNG_CRITICAL)
1588 gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
1590 else if (stat == TNG_FAILURE)
1596 GMX_UNUSED_VALUE(input);
1597 GMX_UNUSED_VALUE(frame);
1598 GMX_UNUSED_VALUE(nRequestedIds);
1599 GMX_UNUSED_VALUE(requestedIds);
1600 GMX_UNUSED_VALUE(nextFrame);
1601 GMX_UNUSED_VALUE(nBlocks);
1602 GMX_UNUSED_VALUE(blockIds);
1607 gmx_bool gmx_get_tng_data_next_frame_of_block_type(tng_trajectory_t input,
1608 gmx_int64_t blockId,
1610 gmx_int64_t *frameNumber,
1612 gmx_int64_t *nValuesPerFrame,
1613 gmx_int64_t *nAtoms,
1620 tng_function_status stat;
1622 gmx_int64_t codecId;
1623 int blockDependency;
1624 void *data = nullptr;
1627 stat = tng_data_block_name_get(input, blockId, name, maxLen);
1628 if (stat != TNG_SUCCESS)
1630 gmx_file("Cannot read next frame of TNG file");
1632 stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
1633 if (stat != TNG_SUCCESS)
1635 gmx_file("Cannot read next frame of TNG file");
1637 if (blockDependency & TNG_PARTICLE_DEPENDENT)
1639 tng_num_particles_get(input, nAtoms);
1640 stat = tng_util_particle_data_next_frame_read(input,
1649 *nAtoms = 1; /* There are not actually any atoms, but it is used for
1650 allocating memory */
1651 stat = tng_util_non_particle_data_next_frame_read(input,
1658 if (stat == TNG_CRITICAL)
1660 gmx_file("Cannot read next frame of TNG file");
1662 if (stat == TNG_FAILURE)
1668 stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
1669 if (stat != TNG_SUCCESS)
1671 gmx_file("Cannot read next frame of TNG file");
1673 srenew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
1674 convert_array_to_real_array(data,
1676 getDistanceScaleFactor(input),
1681 tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
1683 /* This must be updated if/when more lossy compression methods are added */
1684 if (codecId != TNG_TNG_COMPRESSION)
1696 GMX_UNUSED_VALUE(input);
1697 GMX_UNUSED_VALUE(blockId);
1698 GMX_UNUSED_VALUE(values);
1699 GMX_UNUSED_VALUE(frameNumber);
1700 GMX_UNUSED_VALUE(frameTime);
1701 GMX_UNUSED_VALUE(nValuesPerFrame);
1702 GMX_UNUSED_VALUE(nAtoms);
1703 GMX_UNUSED_VALUE(prec);
1704 GMX_UNUSED_VALUE(name);
1705 GMX_UNUSED_VALUE(maxLen);
1706 GMX_UNUSED_VALUE(bOK);