clang-tidy: readibility redundant
[alexxy/gromacs.git] / src / gromacs / fileio / tngio.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "tngio.h"
38
39 #include "config.h"
40
41 #include <cmath>
42
43 #include <algorithm>
44 #include <iterator>
45 #include <memory>
46 #include <vector>
47
48 #if GMX_USE_TNG
49 #include "tng/tng_io.h"
50 #endif
51
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"
67
68 /*! \brief Gromacs Wrapper around tng datatype
69  *
70  * This could in principle hold any GROMACS-specific requirements not yet
71  * implemented in or not relevant to the TNG library itself. However, for now
72  * we only use it to handle some shortcomings we have discovered, where the TNG
73  * API itself is a bit fragile and can end up overwriting data if called several
74  * times with the same frame number.
75  * The logic to determine the time per step was also a bit fragile. This is not
76  * critical, but since we anyway need a wrapper for ensuring unique frame
77  * numbers, we can also use it to store the time of the first step and use that
78  * to derive a slightly better/safer estimate of the time per step.
79  *
80  * At some future point where we have a second-generation TNG API we should
81  * consider removing this again.
82  */
83 struct gmx_tng_trajectory
84 {
85     tng_trajectory_t tng;                  //!< Actual TNG handle (pointer)
86     bool             lastStepDataIsValid;  //!< True if lastStep has been set
87     std::int64_t     lastStep;             //!< Index/step used for last frame
88     bool             lastTimeDataIsValid;  //!< True if lastTime has been set
89     double           lastTime;             //!< Time of last frame (TNG unit is seconds)
90     bool             timePerFrameIsSet;    //!< True if we have set the time per frame
91     int              boxOutputInterval;    //!< Number of steps between the output of box size
92     int              lambdaOutputInterval; //!< Number of steps between the output of lambdas
93 };
94
95 static const char *modeToVerb(char mode)
96 {
97     const char *p;
98     switch (mode)
99     {
100         case 'r':
101             p = "reading";
102             break;
103         case 'w':
104             p = "writing";
105             break;
106         case 'a':
107             p = "appending";
108             break;
109         default:
110             gmx_fatal(FARGS, "Invalid file opening mode %c", mode);
111             p = "";
112             break;
113     }
114     return p;
115 }
116
117 void gmx_tng_open(const char           *filename,
118                   char                  mode,
119                   gmx_tng_trajectory_t *gmx_tng)
120 {
121 #if GMX_USE_TNG
122     /* First check whether we have to make a backup,
123      * only for writing, not for read or append.
124      */
125     if (mode == 'w')
126     {
127         make_backup(filename);
128     }
129
130     *gmx_tng                        = new gmx_tng_trajectory;
131     (*gmx_tng)->lastStepDataIsValid = false;
132     (*gmx_tng)->lastTimeDataIsValid = false;
133     (*gmx_tng)->timePerFrameIsSet   = false;
134     tng_trajectory_t * tng = &(*gmx_tng)->tng;
135
136     /* tng must not be pointing at already allocated memory.
137      * Memory will be allocated by tng_util_trajectory_open() and must
138      * later on be freed by tng_util_trajectory_close(). */
139     if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
140     {
141         /* TNG does return more than one degree of error, but there is
142            no use case for GROMACS handling the non-fatal errors
143            gracefully. */
144         gmx_fatal(FARGS,
145                   "File I/O error while opening %s for %s",
146                   filename,
147                   modeToVerb(mode));
148     }
149
150     if (mode == 'w' || mode == 'a')
151     {
152         char hostname[256];
153         gmx_gethostname(hostname, 256);
154         if (mode == 'w')
155         {
156             tng_first_computer_name_set(*tng, hostname);
157         }
158         else
159         {
160             tng_last_computer_name_set(*tng, hostname);
161         }
162
163         char        programInfo[256];
164         const char *precisionString = "";
165 #if GMX_DOUBLE
166         precisionString = " (double precision)";
167 #endif
168         sprintf(programInfo, "%.100s %.128s%.24s",
169                 gmx::getProgramContext().displayName(),
170                 gmx_version(), precisionString);
171         if (mode == 'w')
172         {
173             tng_first_program_name_set(*tng, programInfo);
174         }
175         else
176         {
177             tng_last_program_name_set(*tng, programInfo);
178         }
179
180         char username[256];
181         if (!gmx_getusername(username, 256))
182         {
183             if (mode == 'w')
184             {
185                 tng_first_user_name_set(*tng, username);
186             }
187             else
188             {
189                 tng_last_user_name_set(*tng, username);
190                 tng_file_headers_write(*tng, TNG_USE_HASH);
191             }
192         }
193     }
194 #else
195     gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
196     GMX_UNUSED_VALUE(filename);
197     GMX_UNUSED_VALUE(mode);
198     GMX_UNUSED_VALUE(gmx_tng);
199 #endif
200 }
201
202 void gmx_tng_close(gmx_tng_trajectory_t *gmx_tng)
203 {
204     /* We have to check that tng is set because
205      * tng_util_trajectory_close wants to return a NULL in it, and
206      * gives a fatal error if it is NULL. */
207 #if GMX_USE_TNG
208     if (gmx_tng == nullptr || *gmx_tng == nullptr)
209     {
210         return;
211     }
212     tng_trajectory_t * tng = &(*gmx_tng)->tng;
213
214     if (tng)
215     {
216         tng_util_trajectory_close(tng);
217     }
218     delete *gmx_tng;
219     *gmx_tng = nullptr;
220
221 #else
222     GMX_UNUSED_VALUE(gmx_tng);
223 #endif
224 }
225
226 #if GMX_USE_TNG
227 static void addTngMoleculeFromTopology(gmx_tng_trajectory_t gmx_tng,
228                                        const char          *moleculeName,
229                                        const t_atoms       *atoms,
230                                        gmx_int64_t          numMolecules,
231                                        tng_molecule_t      *tngMol)
232 {
233     tng_trajectory_t tng      = gmx_tng->tng;
234     tng_chain_t      tngChain = nullptr;
235     tng_residue_t    tngRes   = nullptr;
236
237     if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
238     {
239         gmx_file("Cannot add molecule to TNG molecular system.");
240     }
241
242     for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++)
243     {
244         const t_atom *at = &atoms->atom[atomIndex];
245         /* FIXME: Currently the TNG API can only add atoms belonging to a
246          * residue and chain. Wait for TNG 2.0*/
247         if (atoms->nres > 0)
248         {
249             const t_resinfo *resInfo        = &atoms->resinfo[at->resind];
250             char             chainName[2]   = {resInfo->chainid, 0};
251             tng_atom_t       tngAtom        = nullptr;
252             t_atom          *prevAtom;
253
254             if (atomIndex > 0)
255             {
256                 prevAtom = &atoms->atom[atomIndex - 1];
257             }
258             else
259             {
260                 prevAtom = nullptr;
261             }
262
263             /* If this is the first atom or if the residue changed add the
264              * residue to the TNG molecular system. */
265             if (!prevAtom || resInfo != &atoms->resinfo[prevAtom->resind])
266             {
267                 /* If this is the first atom or if the chain changed add
268                  * the chain to the TNG molecular system. */
269                 if (!prevAtom || resInfo->chainid !=
270                     atoms->resinfo[prevAtom->resind].chainid)
271                 {
272                     tng_molecule_chain_add(tng, *tngMol, chainName,
273                                            &tngChain);
274                 }
275                 /* FIXME: When TNG supports both residue index and residue
276                  * number the latter should be used. Wait for TNG 2.0*/
277                 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
278             }
279             tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIndex]), *(atoms->atomtype[atomIndex]), &tngAtom);
280         }
281     }
282     tng_molecule_cnt_set(tng, *tngMol, numMolecules);
283 }
284
285 void gmx_tng_add_mtop(gmx_tng_trajectory_t  gmx_tng,
286                       const gmx_mtop_t     *mtop)
287 {
288     int                i;
289     int                j;
290     std::vector<real>  atomCharges;
291     std::vector<real>  atomMasses;
292     const t_ilist     *ilist;
293     tng_bond_t         tngBond;
294     char               datatype;
295
296     tng_trajectory_t   tng = gmx_tng->tng;
297
298     if (!mtop)
299     {
300         /* No topology information available to add. */
301         return;
302     }
303
304 #if GMX_DOUBLE
305     datatype = TNG_DOUBLE_DATA;
306 #else
307     datatype = TNG_FLOAT_DATA;
308 #endif
309
310     atomCharges.reserve(mtop->natoms);
311     atomMasses.reserve(mtop->natoms);
312
313     for (const gmx_molblock_t &molBlock :  mtop->molblock)
314     {
315         tng_molecule_t       tngMol  = nullptr;
316         const gmx_moltype_t *molType = &mtop->moltype[molBlock.type];
317
318         /* Add a molecule to the TNG trajectory with the same name as the
319          * current molecule. */
320         addTngMoleculeFromTopology(gmx_tng,
321                                    *(molType->name),
322                                    &molType->atoms,
323                                    molBlock.nmol,
324                                    &tngMol);
325
326         /* Bonds have to be deduced from interactions (constraints etc). Different
327          * interactions have different sets of parameters. */
328         /* Constraints are specified using two atoms */
329         for (i = 0; i < F_NRE; i++)
330         {
331             if (IS_CHEMBOND(i))
332             {
333                 ilist = &molType->ilist[i];
334                 if (ilist)
335                 {
336                     j = 1;
337                     while (j < ilist->nr)
338                     {
339                         tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
340                         j += 3;
341                     }
342                 }
343             }
344         }
345         /* Settle is described using three atoms */
346         ilist = &molType->ilist[F_SETTLE];
347         if (ilist)
348         {
349             j = 1;
350             while (j < ilist->nr)
351             {
352                 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
353                 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+2], &tngBond);
354                 j += 4;
355             }
356         }
357         /* First copy atom charges and masses, first atom by atom and then copy the memory for the molecule instances.
358          * FIXME: Atom B state data should also be written to TNG (v 2.0?) */
359         for (int atomCounter = 0; atomCounter < molType->atoms.nr; atomCounter++)
360         {
361             atomCharges.push_back(molType->atoms.atom[atomCounter].q);
362             atomMasses.push_back(molType->atoms.atom[atomCounter].m);
363         }
364         for (int molCounter = 1; molCounter < molBlock.nmol; molCounter++)
365         {
366             std::copy_n(atomCharges.end() - molType->atoms.nr, molType->atoms.nr, std::back_inserter(atomCharges));
367             std::copy_n(atomMasses.end()  - molType->atoms.nr, molType->atoms.nr, std::back_inserter(atomMasses));
368         }
369     }
370     /* Write the TNG data blocks. */
371     tng_particle_data_block_add(tng, TNG_TRAJ_PARTIAL_CHARGES, "PARTIAL CHARGES",
372                                 datatype, TNG_NON_TRAJECTORY_BLOCK,
373                                 1, 1, 1, 0, mtop->natoms,
374                                 TNG_GZIP_COMPRESSION, atomCharges.data());
375     tng_particle_data_block_add(tng, TNG_TRAJ_MASSES, "ATOM MASSES",
376                                 datatype, TNG_NON_TRAJECTORY_BLOCK,
377                                 1, 1, 1, 0, mtop->natoms,
378                                 TNG_GZIP_COMPRESSION, atomMasses.data());
379 }
380
381 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
382  * if they are positive.
383  *
384  * If only one of n1 and n2 is positive, then return it.
385  * If neither n1 or n2 is positive, then return -1. */
386 static int
387 greatest_common_divisor_if_positive(int n1, int n2)
388 {
389     if (0 >= n1)
390     {
391         return (0 >= n2) ? -1 : n2;
392     }
393     if (0 >= n2)
394     {
395         return n1;
396     }
397
398     /* We have a non-trivial greatest common divisor to compute. */
399     return gmx_greatest_common_divisor(n1, n2);
400 }
401
402 /* By default try to write 100 frames (of actual output) in each frame set.
403  * This number is the number of outputs of the most frequently written data
404  * type per frame set.
405  * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
406  * setups regarding compression efficiency and compression time. Make this
407  * a hidden command-line option? */
408 const int defaultFramesPerFrameSet = 100;
409
410 /*! \libinternal \brief  Set the number of frames per frame
411  * set according to output intervals.
412  * The default is that 100 frames are written of the data
413  * that is written most often. */
414 static void tng_set_frames_per_frame_set(gmx_tng_trajectory_t  gmx_tng,
415                                          const gmx_bool        bUseLossyCompression,
416                                          const t_inputrec     *ir)
417 {
418     int              gcd = -1;
419     tng_trajectory_t tng = gmx_tng->tng;
420
421     /* Set the number of frames per frame set to contain at least
422      * defaultFramesPerFrameSet of the lowest common denominator of
423      * the writing interval of positions and velocities. */
424     /* FIXME after 5.0: consider nstenergy also? */
425     if (bUseLossyCompression)
426     {
427         gcd = ir->nstxout_compressed;
428     }
429     else
430     {
431         gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
432         gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
433     }
434     if (0 >= gcd)
435     {
436         return;
437     }
438
439     tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
440 }
441
442 /*! \libinternal \brief Set the data-writing intervals, and number of
443  * frames per frame set */
444 static void set_writing_intervals(gmx_tng_trajectory_t  gmx_tng,
445                                   const gmx_bool        bUseLossyCompression,
446                                   const t_inputrec     *ir)
447 {
448     tng_trajectory_t tng = gmx_tng->tng;
449
450     /* Define pointers to specific writing functions depending on if we
451      * write float or double data */
452     typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
453                                                                      const gmx_int64_t,
454                                                                      const gmx_int64_t,
455                                                                      const gmx_int64_t,
456                                                                      const char*,
457                                                                      const char,
458                                                                      const char);
459 #if GMX_DOUBLE
460     set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
461 #else
462     set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
463 #endif
464     int  xout, vout, fout;
465     int  gcd = -1, lowest = -1;
466     char compression;
467
468     tng_set_frames_per_frame_set(gmx_tng, bUseLossyCompression, ir);
469
470     if (bUseLossyCompression)
471     {
472         xout        = ir->nstxout_compressed;
473
474         /* If there is no uncompressed coordinate output write forces
475            and velocities to the compressed tng file. */
476         if (ir->nstxout)
477         {
478             vout        = 0;
479             fout        = 0;
480         }
481         else
482         {
483             vout        = ir->nstvout;
484             fout        = ir->nstfout;
485         }
486         compression = TNG_TNG_COMPRESSION;
487     }
488     else
489     {
490         xout        = ir->nstxout;
491         vout        = ir->nstvout;
492         fout        = ir->nstfout;
493         compression = TNG_GZIP_COMPRESSION;
494     }
495     if (xout)
496     {
497         set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS,
498                              "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
499                              compression);
500         /* TODO: if/when we write energies to TNG also, reconsider how
501          * and when box information is written, because GROMACS
502          * behaviour pre-5.0 was to write the box with every
503          * trajectory frame and every energy frame, and probably
504          * people depend on this. */
505
506         gcd = greatest_common_divisor_if_positive(gcd, xout);
507         if (lowest < 0 || xout < lowest)
508         {
509             lowest = xout;
510         }
511     }
512     if (vout)
513     {
514         set_writing_interval(tng, vout, 3, TNG_TRAJ_VELOCITIES,
515                              "VELOCITIES", TNG_PARTICLE_BLOCK_DATA,
516                              compression);
517
518         gcd = greatest_common_divisor_if_positive(gcd, vout);
519         if (lowest < 0 || vout < lowest)
520         {
521             lowest = vout;
522         }
523     }
524     if (fout)
525     {
526         set_writing_interval(tng, fout, 3, TNG_TRAJ_FORCES,
527                              "FORCES", TNG_PARTICLE_BLOCK_DATA,
528                              TNG_GZIP_COMPRESSION);
529
530         gcd = greatest_common_divisor_if_positive(gcd, fout);
531         if (lowest < 0 || fout < lowest)
532         {
533             lowest = fout;
534         }
535     }
536     if (gcd > 0)
537     {
538         /* Lambdas and box shape written at an interval of the lowest common
539            denominator of other output */
540         set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA,
541                              "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
542                              TNG_GZIP_COMPRESSION);
543
544         set_writing_interval(tng, gcd, 9, TNG_TRAJ_BOX_SHAPE,
545                              "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA,
546                              TNG_GZIP_COMPRESSION);
547         gmx_tng->lambdaOutputInterval = gcd;
548         gmx_tng->boxOutputInterval    = gcd;
549         if (gcd < lowest / 10)
550         {
551             gmx_warning("The lowest common denominator of trajectory output is "
552                         "every %d step(s), whereas the shortest output interval "
553                         "is every %d steps.", gcd, lowest);
554         }
555     }
556 }
557 #endif
558
559 void gmx_tng_prepare_md_writing(gmx_tng_trajectory_t  gmx_tng,
560                                 const gmx_mtop_t     *mtop,
561                                 const t_inputrec     *ir)
562 {
563 #if GMX_USE_TNG
564     gmx_tng_add_mtop(gmx_tng, mtop);
565     set_writing_intervals(gmx_tng, FALSE, ir);
566     tng_time_per_frame_set(gmx_tng->tng, ir->delta_t * PICO);
567     gmx_tng->timePerFrameIsSet = true;
568 #else
569     GMX_UNUSED_VALUE(gmx_tng);
570     GMX_UNUSED_VALUE(mtop);
571     GMX_UNUSED_VALUE(ir);
572 #endif
573 }
574
575 #if GMX_USE_TNG
576 /* Check if all atoms in the molecule system are selected
577  * by a selection group of type specified by gtype. */
578 static gmx_bool all_atoms_selected(const gmx_mtop_t *mtop,
579                                    const int         gtype)
580 {
581     /* Iterate through all atoms in the molecule system and
582      * check if they belong to a selection group of the
583      * requested type. */
584     int i = 0;
585     for (const gmx_molblock_t &molBlock : mtop->molblock)
586     {
587         const gmx_moltype_t &molType = mtop->moltype[molBlock.type];
588         const t_atoms       &atoms   = molType.atoms;
589
590         for (int j = 0; j < molBlock.nmol; j++)
591         {
592             for (int atomIndex = 0; atomIndex < atoms.nr; atomIndex++, i++)
593             {
594                 if (ggrpnr(&mtop->groups, gtype, i) != 0)
595                 {
596                     return FALSE;
597                 }
598             }
599         }
600     }
601     return TRUE;
602 }
603
604 /* Create TNG molecules which will represent each of the selection
605  * groups for writing. But do that only if there is actually a
606  * specified selection group and it is not the whole system.
607  * TODO: Currently the only selection that is taken into account
608  * is egcCompressedX, but other selections should be added when
609  * e.g. writing energies is implemented.
610  */
611 static void add_selection_groups(gmx_tng_trajectory_t  gmx_tng,
612                                  const gmx_mtop_t     *mtop)
613 {
614     const t_atoms           *atoms;
615     const t_atom            *at;
616     const t_resinfo         *resInfo;
617     const t_ilist           *ilist;
618     int                      nameIndex;
619     int                      atom_offset = 0;
620     tng_molecule_t           mol, iterMol;
621     tng_chain_t              chain;
622     tng_residue_t            res;
623     tng_atom_t               atom;
624     tng_bond_t               tngBond;
625     gmx_int64_t              nMols;
626     char                    *groupName;
627     tng_trajectory_t         tng = gmx_tng->tng;
628
629     /* TODO: When the TNG molecules block is more flexible TNG selection
630      * groups should not need all atoms specified. It should be possible
631      * just to specify what molecules are selected (and/or which atoms in
632      * the molecule) and how many (if applicable). */
633
634     /* If no atoms are selected we do not need to create a
635      * TNG selection group molecule. */
636     if (mtop->groups.ngrpnr[egcCompressedX] == 0)
637     {
638         return;
639     }
640
641     /* If all atoms are selected we do not have to create a selection
642      * group molecule in the TNG molecule system. */
643     if (all_atoms_selected(mtop, egcCompressedX))
644     {
645         return;
646     }
647
648     /* The name of the TNG molecule containing the selection group is the
649      * same as the name of the selection group. */
650     nameIndex = *mtop->groups.grps[egcCompressedX].nm_ind;
651     groupName = *mtop->groups.grpname[nameIndex];
652
653     tng_molecule_alloc(tng, &mol);
654     tng_molecule_name_set(tng, mol, groupName);
655     tng_molecule_chain_add(tng, mol, "", &chain);
656     int i = 0;
657     for (const gmx_molblock_t &molBlock :  mtop->molblock)
658     {
659         const gmx_moltype_t &molType = mtop->moltype[molBlock.type];
660
661         atoms = &molType.atoms;
662
663         for (int j = 0; j < molBlock.nmol; j++)
664         {
665             bool bAtomsAdded = FALSE;
666             for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
667             {
668                 char *res_name;
669                 int   res_id;
670
671                 if (ggrpnr(&mtop->groups, egcCompressedX, i) != 0)
672                 {
673                     continue;
674                 }
675                 at = &atoms->atom[atomIndex];
676                 if (atoms->nres > 0)
677                 {
678                     resInfo = &atoms->resinfo[at->resind];
679                     /* FIXME: When TNG supports both residue index and residue
680                      * number the latter should be used. */
681                     res_name = *resInfo->name;
682                     res_id   = at->resind + 1;
683                 }
684                 else
685                 {
686                     res_name = (char *)"";
687                     res_id   = 0;
688                 }
689                 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res)
690                     != TNG_SUCCESS)
691                 {
692                     /* Since there is ONE chain for selection groups do not keep the
693                      * original residue IDs - otherwise there might be conflicts. */
694                     tng_chain_residue_add(tng, chain, res_name, &res);
695                 }
696                 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
697                                           *(atoms->atomtype[atomIndex]),
698                                           atom_offset + atomIndex, &atom);
699                 bAtomsAdded = TRUE;
700             }
701             /* Add bonds. */
702             if (bAtomsAdded)
703             {
704                 for (int k = 0; k < F_NRE; k++)
705                 {
706                     if (IS_CHEMBOND(k))
707                     {
708                         ilist = &molType.ilist[k];
709                         if (ilist)
710                         {
711                             int l = 1;
712                             while (l < ilist->nr)
713                             {
714                                 int atom1, atom2;
715                                 atom1 = ilist->iatoms[l] + atom_offset;
716                                 atom2 = ilist->iatoms[l+1] + atom_offset;
717                                 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0 &&
718                                     ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
719                                 {
720                                     tng_molecule_bond_add(tng, mol, ilist->iatoms[l],
721                                                           ilist->iatoms[l+1], &tngBond);
722                                 }
723                                 l += 3;
724                             }
725                         }
726                     }
727                 }
728                 /* Settle is described using three atoms */
729                 ilist = &molType.ilist[F_SETTLE];
730                 if (ilist)
731                 {
732                     int l = 1;
733                     while (l < ilist->nr)
734                     {
735                         int atom1, atom2, atom3;
736                         atom1 = ilist->iatoms[l] + atom_offset;
737                         atom2 = ilist->iatoms[l+1] + atom_offset;
738                         atom3 = ilist->iatoms[l+2] + atom_offset;
739                         if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0)
740                         {
741                             if (ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
742                             {
743                                 tng_molecule_bond_add(tng, mol, atom1,
744                                                       atom2, &tngBond);
745                             }
746                             if (ggrpnr(&mtop->groups, egcCompressedX, atom3) == 0)
747                             {
748                                 tng_molecule_bond_add(tng, mol, atom1,
749                                                       atom3, &tngBond);
750                             }
751                         }
752                         l += 4;
753                     }
754                 }
755             }
756             atom_offset += atoms->nr;
757         }
758     }
759     tng_molecule_existing_add(tng, &mol);
760     tng_molecule_cnt_set(tng, mol, 1);
761     tng_num_molecule_types_get(tng, &nMols);
762     for (gmx_int64_t k = 0; k < nMols; k++)
763     {
764         tng_molecule_of_index_get(tng, k, &iterMol);
765         if (iterMol == mol)
766         {
767             continue;
768         }
769         tng_molecule_cnt_set(tng, iterMol, 0);
770     }
771 }
772 #endif
773
774 void gmx_tng_set_compression_precision(gmx_tng_trajectory_t gmx_tng,
775                                        real                 prec)
776 {
777 #if GMX_USE_TNG
778     tng_compression_precision_set(gmx_tng->tng, prec);
779 #else
780     GMX_UNUSED_VALUE(gmx_tng);
781     GMX_UNUSED_VALUE(prec);
782 #endif
783 }
784
785 void gmx_tng_prepare_low_prec_writing(gmx_tng_trajectory_t  gmx_tng,
786                                       const gmx_mtop_t     *mtop,
787                                       const t_inputrec     *ir)
788 {
789 #if GMX_USE_TNG
790     gmx_tng_add_mtop(gmx_tng, mtop);
791     add_selection_groups(gmx_tng, mtop);
792     set_writing_intervals(gmx_tng, TRUE, ir);
793     tng_time_per_frame_set(gmx_tng->tng, ir->delta_t * PICO);
794     gmx_tng->timePerFrameIsSet = true;
795     gmx_tng_set_compression_precision(gmx_tng, ir->x_compression_precision);
796 #else
797     GMX_UNUSED_VALUE(gmx_tng);
798     GMX_UNUSED_VALUE(mtop);
799     GMX_UNUSED_VALUE(ir);
800 #endif
801 }
802
803 void gmx_fwrite_tng(gmx_tng_trajectory_t gmx_tng,
804                     const gmx_bool       bUseLossyCompression,
805                     gmx_int64_t          step,
806                     real                 elapsedPicoSeconds,
807                     real                 lambda,
808                     const rvec          *box,
809                     int                  nAtoms,
810                     const rvec          *x,
811                     const rvec          *v,
812                     const rvec          *f)
813 {
814 #if GMX_USE_TNG
815     typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
816                                                            const gmx_int64_t,
817                                                            const double,
818                                                            const real*,
819                                                            const gmx_int64_t,
820                                                            const gmx_int64_t,
821                                                            const char*,
822                                                            const char,
823                                                            const char);
824 #if GMX_DOUBLE
825     static write_data_func_pointer           write_data           = tng_util_generic_with_time_double_write;
826 #else
827     static write_data_func_pointer           write_data           = tng_util_generic_with_time_write;
828 #endif
829     double                                   elapsedSeconds = elapsedPicoSeconds * PICO;
830     gmx_int64_t                              nParticles;
831     char                                     compression;
832
833
834     if (!gmx_tng)
835     {
836         /* This function might get called when the type of the
837            compressed trajectory is actually XTC. So we exit and move
838            on. */
839         return;
840     }
841     tng_trajectory_t tng = gmx_tng->tng;
842
843     // While the GROMACS interface to this routine specifies 'step', TNG itself
844     // only uses 'frame index' internally, although it suggests that it's a good
845     // idea to let this represent the step rather than just frame numbers.
846     //
847     // However, the way the GROMACS interface works, there are cases where
848     // the step is simply not set, so all frames rather have step=0.
849     // If we call the lower-level TNG routines with the same frame index
850     // (which is set from the step) multiple times, things will go very
851     // wrong (overwritten data).
852     //
853     // To avoid this, we require the step value to always be larger than the
854     // previous value, and if this is not true we simply set it to a value
855     // one unit larger than the previous step.
856     //
857     // Note: We do allow the user to specify any crazy value the want for the
858     // first step. If it is "not set", GROMACS will have used the value 0.
859     if (gmx_tng->lastStepDataIsValid && step <= gmx_tng->lastStep)
860     {
861         step = gmx_tng->lastStep + 1;
862     }
863
864     // Now that we have fixed the potentially incorrect step, we can also set
865     // the time per frame if it was not already done, and we have
866     // step/time information from the last frame so it is possible to calculate it.
867     // Note that we do allow the time to be the same for all steps, or even
868     // decreasing. In that case we will simply say the time per step is 0
869     // or negative. A bit strange, but a correct representation of the data :-)
870     if (!gmx_tng->timePerFrameIsSet && gmx_tng->lastTimeDataIsValid && gmx_tng->lastStepDataIsValid)
871     {
872         double       deltaTime = elapsedSeconds - gmx_tng->lastTime;
873         std::int64_t deltaStep = step - gmx_tng->lastStep;
874         tng_time_per_frame_set(tng, deltaTime / deltaStep );
875         gmx_tng->timePerFrameIsSet = true;
876     }
877
878     tng_num_particles_get(tng, &nParticles);
879     if (nAtoms != (int)nParticles)
880     {
881         tng_implicit_num_particles_set(tng, nAtoms);
882     }
883
884     if (bUseLossyCompression)
885     {
886         compression = TNG_TNG_COMPRESSION;
887     }
888     else
889     {
890         compression = TNG_GZIP_COMPRESSION;
891     }
892
893     /* The writing is done using write_data, which writes float or double
894      * depending on the GROMACS compilation. */
895     if (x)
896     {
897         GMX_ASSERT(box, "Need a non-NULL box if positions are written");
898
899         if (write_data(tng, step, elapsedSeconds,
900                        reinterpret_cast<const real *>(x),
901                        3, TNG_TRAJ_POSITIONS, "POSITIONS",
902                        TNG_PARTICLE_BLOCK_DATA,
903                        compression) != TNG_SUCCESS)
904         {
905             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
906         }
907     }
908
909     if (v)
910     {
911         if (write_data(tng, step, elapsedSeconds,
912                        reinterpret_cast<const real *>(v),
913                        3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
914                        TNG_PARTICLE_BLOCK_DATA,
915                        compression) != TNG_SUCCESS)
916         {
917             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
918         }
919     }
920
921     if (f)
922     {
923         /* TNG-MF1 compression only compresses positions and velocities. Use lossless
924          * compression for forces regardless of output mode */
925         if (write_data(tng, step, elapsedSeconds,
926                        reinterpret_cast<const real *>(f),
927                        3, TNG_TRAJ_FORCES, "FORCES",
928                        TNG_PARTICLE_BLOCK_DATA,
929                        TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
930         {
931             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
932         }
933     }
934
935     if (box)
936     {
937         /* TNG-MF1 compression only compresses positions and velocities. Use lossless
938          * compression for box shape regardless of output mode */
939         if (write_data(tng, step, elapsedSeconds,
940                        reinterpret_cast<const real *>(box),
941                        9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
942                        TNG_NON_PARTICLE_BLOCK_DATA,
943                        TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
944         {
945             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
946         }
947     }
948
949     if (lambda >= 0)
950     {
951         /* TNG-MF1 compression only compresses positions and velocities. Use lossless
952          * compression for lambda regardless of output mode */
953         if (write_data(tng, step, elapsedSeconds,
954                        reinterpret_cast<const real *>(&lambda),
955                        1, TNG_GMX_LAMBDA, "LAMBDAS",
956                        TNG_NON_PARTICLE_BLOCK_DATA,
957                        TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
958         {
959             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
960         }
961     }
962
963     // Update the data in the wrapper
964
965     gmx_tng->lastStepDataIsValid = true;
966     gmx_tng->lastStep            = step;
967     gmx_tng->lastTimeDataIsValid = true;
968     gmx_tng->lastTime            = elapsedSeconds;
969 #else
970     GMX_UNUSED_VALUE(gmx_tng);
971     GMX_UNUSED_VALUE(bUseLossyCompression);
972     GMX_UNUSED_VALUE(step);
973     GMX_UNUSED_VALUE(elapsedPicoSeconds);
974     GMX_UNUSED_VALUE(lambda);
975     GMX_UNUSED_VALUE(box);
976     GMX_UNUSED_VALUE(nAtoms);
977     GMX_UNUSED_VALUE(x);
978     GMX_UNUSED_VALUE(v);
979     GMX_UNUSED_VALUE(f);
980 #endif
981 }
982
983 void fflush_tng(gmx_tng_trajectory_t gmx_tng)
984 {
985 #if GMX_USE_TNG
986     if (!gmx_tng)
987     {
988         return;
989     }
990     tng_frame_set_premature_write(gmx_tng->tng, TNG_USE_HASH);
991 #else
992     GMX_UNUSED_VALUE(gmx_tng);
993 #endif
994 }
995
996 float gmx_tng_get_time_of_final_frame(gmx_tng_trajectory_t gmx_tng)
997 {
998 #if GMX_USE_TNG
999     gmx_int64_t      nFrames;
1000     double           time;
1001     float            fTime;
1002     tng_trajectory_t tng = gmx_tng->tng;
1003
1004     tng_num_frames_get(tng, &nFrames);
1005     tng_util_time_of_frame_get(tng, nFrames - 1, &time);
1006
1007     fTime = time / PICO;
1008     return fTime;
1009 #else
1010     GMX_UNUSED_VALUE(gmx_tng);
1011     return -1.0;
1012 #endif
1013 }
1014
1015 void gmx_prepare_tng_writing(const char              *filename,
1016                              char                     mode,
1017                              gmx_tng_trajectory_t    *gmx_tng_input,
1018                              gmx_tng_trajectory_t    *gmx_tng_output,
1019                              int                      nAtoms,
1020                              const gmx_mtop_t        *mtop,
1021                              const int               *index,
1022                              const char              *indexGroupName)
1023 {
1024 #if GMX_USE_TNG
1025     tng_trajectory_t   *input  = (gmx_tng_input && *gmx_tng_input) ? &(*gmx_tng_input)->tng : nullptr;
1026     /* FIXME after 5.0: Currently only standard block types are read */
1027     const int           defaultNumIds              = 5;
1028     static gmx_int64_t  fallbackIds[defaultNumIds] =
1029     {
1030         TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
1031         TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
1032         TNG_GMX_LAMBDA
1033     };
1034     static char         fallbackNames[defaultNumIds][32] =
1035     {
1036         "BOX SHAPE", "POSITIONS", "VELOCITIES",
1037         "FORCES", "LAMBDAS"
1038     };
1039
1040     typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
1041                                                                      const gmx_int64_t,
1042                                                                      const gmx_int64_t,
1043                                                                      const gmx_int64_t,
1044                                                                      const char*,
1045                                                                      const char,
1046                                                                      const char);
1047 #if GMX_DOUBLE
1048     set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
1049 #else
1050     set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
1051 #endif
1052
1053     gmx_tng_open(filename, mode, gmx_tng_output);
1054     tng_trajectory_t *output = &(*gmx_tng_output)->tng;
1055
1056     /* Do we have an input file in TNG format? If so, then there's
1057        more data we can copy over, rather than having to improvise. */
1058     if (gmx_tng_input && *gmx_tng_input)
1059     {
1060         /* Set parameters (compression, time per frame, molecule
1061          * information, number of frames per frame set and writing
1062          * intervals of positions, box shape and lambdas) of the
1063          * output tng container based on their respective values int
1064          * the input tng container */
1065         double      time, compression_precision;
1066         gmx_int64_t n_frames_per_frame_set, interval = -1;
1067
1068         tng_compression_precision_get(*input, &compression_precision);
1069         tng_compression_precision_set(*output, compression_precision);
1070         // TODO make this configurable in a future version
1071         char compression_type = TNG_TNG_COMPRESSION;
1072
1073         tng_molecule_system_copy(*input, *output);
1074
1075         tng_time_per_frame_get(*input, &time);
1076         tng_time_per_frame_set(*output, time);
1077         // Since we have copied the value from the input TNG we should not change it again
1078         (*gmx_tng_output)->timePerFrameIsSet = true;
1079
1080         tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
1081         tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
1082
1083         for (int i = 0; i < defaultNumIds; i++)
1084         {
1085             if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval)
1086                 == TNG_SUCCESS)
1087             {
1088                 switch (fallbackIds[i])
1089                 {
1090                     case TNG_TRAJ_POSITIONS:
1091                     case TNG_TRAJ_VELOCITIES:
1092                         set_writing_interval(*output, interval, 3, fallbackIds[i],
1093                                              fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
1094                                              compression_type);
1095                         break;
1096                     case TNG_TRAJ_FORCES:
1097                         set_writing_interval(*output, interval, 3, fallbackIds[i],
1098                                              fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
1099                                              TNG_GZIP_COMPRESSION);
1100                         break;
1101                     case TNG_TRAJ_BOX_SHAPE:
1102                         set_writing_interval(*output, interval, 9, fallbackIds[i],
1103                                              fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
1104                                              TNG_GZIP_COMPRESSION);
1105                         (*gmx_tng_output)->boxOutputInterval = interval;
1106                         break;
1107                     case TNG_GMX_LAMBDA:
1108                         set_writing_interval(*output, interval, 1, fallbackIds[i],
1109                                              fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
1110                                              TNG_GZIP_COMPRESSION);
1111                         (*gmx_tng_output)->lambdaOutputInterval = interval;
1112                         break;
1113                     default:
1114                         continue;
1115                 }
1116             }
1117         }
1118
1119     }
1120     else
1121     {
1122         /* TODO after trjconv is modularized: fix this so the user can
1123            change precision when they are doing an operation where
1124            this makes sense, and not otherwise.
1125
1126            char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
1127            gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
1128          */
1129         gmx_tng_add_mtop(*gmx_tng_output, mtop);
1130         tng_num_frames_per_frame_set_set(*output, 1);
1131     }
1132
1133     if (index && nAtoms > 0)
1134     {
1135         gmx_tng_setup_atom_subgroup(*gmx_tng_output, nAtoms, index, indexGroupName);
1136     }
1137
1138     /* If for some reason there are more requested atoms than there are atoms in the
1139      * molecular system create a number of implicit atoms (without atom data) to
1140      * compensate for that. */
1141     if (nAtoms >= 0)
1142     {
1143         tng_implicit_num_particles_set(*output, nAtoms);
1144     }
1145 #else
1146     GMX_UNUSED_VALUE(filename);
1147     GMX_UNUSED_VALUE(mode);
1148     GMX_UNUSED_VALUE(gmx_tng_input);
1149     GMX_UNUSED_VALUE(gmx_tng_output);
1150     GMX_UNUSED_VALUE(nAtoms);
1151     GMX_UNUSED_VALUE(mtop);
1152     GMX_UNUSED_VALUE(index);
1153     GMX_UNUSED_VALUE(indexGroupName);
1154 #endif
1155 }
1156
1157 void gmx_write_tng_from_trxframe(gmx_tng_trajectory_t    gmx_tng_output,
1158                                  const t_trxframe       *frame,
1159                                  int                     natoms)
1160 {
1161 #if GMX_USE_TNG
1162     if (natoms < 0)
1163     {
1164         natoms = frame->natoms;
1165     }
1166     gmx_fwrite_tng(gmx_tng_output,
1167                    TRUE,
1168                    frame->step,
1169                    frame->time,
1170                    0,
1171                    frame->box,
1172                    natoms,
1173                    frame->x,
1174                    frame->v,
1175                    frame->f);
1176 #else
1177     GMX_UNUSED_VALUE(gmx_tng_output);
1178     GMX_UNUSED_VALUE(frame);
1179     GMX_UNUSED_VALUE(natoms);
1180 #endif
1181 }
1182
1183 namespace
1184 {
1185
1186 #if GMX_USE_TNG
1187 void
1188 convert_array_to_real_array(void       *from,
1189                             real       *to,
1190                             const float fact,
1191                             const int   nAtoms,
1192                             const int   nValues,
1193                             const char  datatype)
1194 {
1195     int        i, j;
1196
1197     const bool useDouble = GMX_DOUBLE;
1198     switch (datatype)
1199     {
1200         case TNG_FLOAT_DATA:
1201             if (!useDouble)
1202             {
1203                 if (fact == 1)
1204                 {
1205                     memcpy(to, from, nValues * sizeof(real) * nAtoms);
1206                 }
1207                 else
1208                 {
1209                     for (i = 0; i < nAtoms; i++)
1210                     {
1211                         for (j = 0; j < nValues; j++)
1212                         {
1213                             to[i*nValues+j] = reinterpret_cast<float *>(from)[i*nValues+j] * fact;
1214                         }
1215                     }
1216                 }
1217             }
1218             else
1219             {
1220                 for (i = 0; i < nAtoms; i++)
1221                 {
1222                     for (j = 0; j < nValues; j++)
1223                     {
1224                         to[i*nValues+j] = reinterpret_cast<float *>(from)[i*nValues+j] * fact;
1225                     }
1226                 }
1227             }
1228             break;
1229         case TNG_INT_DATA:
1230             for (i = 0; i < nAtoms; i++)
1231             {
1232                 for (j = 0; j < nValues; j++)
1233                 {
1234                     to[i*nValues+j] = reinterpret_cast<gmx_int64_t *>(from)[i*nValues+j] * fact;
1235                 }
1236             }
1237             break;
1238         case TNG_DOUBLE_DATA:
1239             if (sizeof(real) == sizeof(double))
1240             {
1241                 if (fact == 1)
1242                 {
1243                     memcpy(to, from, nValues * sizeof(real) * nAtoms);
1244                 }
1245                 else
1246                 {
1247                     for (i = 0; i < nAtoms; i++)
1248                     {
1249                         for (j = 0; j < nValues; j++)
1250                         {
1251                             to[i*nValues+j] = reinterpret_cast<double *>(from)[i*nValues+j] * fact;
1252                         }
1253                     }
1254                 }
1255             }
1256             else
1257             {
1258                 for (i = 0; i < nAtoms; i++)
1259                 {
1260                     for (j = 0; j < nValues; j++)
1261                     {
1262                         to[i*nValues+j] = reinterpret_cast<double *>(from)[i*nValues+j] * fact;
1263                     }
1264                 }
1265             }
1266             break;
1267         default:
1268             gmx_incons("Illegal datatype when converting values to a real array!");
1269             return;
1270     }
1271 }
1272
1273 real getDistanceScaleFactor(gmx_tng_trajectory_t in)
1274 {
1275     gmx_int64_t exp = -1;
1276     real        distanceScaleFactor;
1277
1278     // TODO Hopefully, TNG 2.0 will do this kind of thing for us
1279     tng_distance_unit_exponential_get(in->tng, &exp);
1280
1281     // GROMACS expects distances in nm
1282     switch (exp)
1283     {
1284         case 9:
1285             distanceScaleFactor = NANO/NANO;
1286             break;
1287         case 10:
1288             distanceScaleFactor = NANO/ANGSTROM;
1289             break;
1290         default:
1291             distanceScaleFactor = pow(10.0, exp + 9.0);
1292     }
1293
1294     return distanceScaleFactor;
1295 }
1296 #endif
1297
1298 } // namespace
1299
1300 void gmx_tng_setup_atom_subgroup(gmx_tng_trajectory_t gmx_tng,
1301                                  const int            nind,
1302                                  const int           *ind,
1303                                  const char          *name)
1304 {
1305 #if GMX_USE_TNG
1306     gmx_int64_t              nAtoms, cnt, nMols;
1307     tng_molecule_t           mol, iterMol;
1308     tng_chain_t              chain;
1309     tng_residue_t            res;
1310     tng_atom_t               atom;
1311     tng_function_status      stat;
1312     tng_trajectory_t         tng = gmx_tng->tng;
1313
1314     tng_num_particles_get(tng, &nAtoms);
1315
1316     if (nAtoms == nind)
1317     {
1318         return;
1319     }
1320
1321     stat = tng_molecule_find(tng, name, -1, &mol);
1322     if (stat == TNG_SUCCESS)
1323     {
1324         tng_molecule_num_atoms_get(tng, mol, &nAtoms);
1325         tng_molecule_cnt_get(tng, mol, &cnt);
1326         if (nAtoms == nind)
1327         {
1328             stat = TNG_SUCCESS;
1329         }
1330         else
1331         {
1332             stat = TNG_FAILURE;
1333         }
1334     }
1335     if (stat == TNG_FAILURE)
1336     {
1337         /* The indexed atoms are added to one separate molecule. */
1338         tng_molecule_alloc(tng, &mol);
1339         tng_molecule_name_set(tng, mol, name);
1340         tng_molecule_chain_add(tng, mol, "", &chain);
1341
1342         for (int i = 0; i < nind; i++)
1343         {
1344             char        temp_name[256], temp_type[256];
1345
1346             /* Try to retrieve the residue name of the atom */
1347             stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1348             if (stat != TNG_SUCCESS)
1349             {
1350                 temp_name[0] = '\0';
1351             }
1352             /* Check if the molecule of the selection already contains this residue */
1353             if (tng_chain_residue_find(tng, chain, temp_name, -1, &res)
1354                 != TNG_SUCCESS)
1355             {
1356                 tng_chain_residue_add(tng, chain, temp_name, &res);
1357             }
1358             /* Try to find the original name and type of the atom */
1359             stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1360             if (stat != TNG_SUCCESS)
1361             {
1362                 temp_name[0] = '\0';
1363             }
1364             stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
1365             if (stat != TNG_SUCCESS)
1366             {
1367                 temp_type[0] = '\0';
1368             }
1369             tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
1370         }
1371         tng_molecule_existing_add(tng, &mol);
1372     }
1373     /* Set the count of the molecule containing the selected atoms to 1 and all
1374      * other molecules to 0 */
1375     tng_molecule_cnt_set(tng, mol, 1);
1376     tng_num_molecule_types_get(tng, &nMols);
1377     for (gmx_int64_t k = 0; k < nMols; k++)
1378     {
1379         tng_molecule_of_index_get(tng, k, &iterMol);
1380         if (iterMol == mol)
1381         {
1382             continue;
1383         }
1384         tng_molecule_cnt_set(tng, iterMol, 0);
1385     }
1386 #else
1387     GMX_UNUSED_VALUE(gmx_tng);
1388     GMX_UNUSED_VALUE(nind);
1389     GMX_UNUSED_VALUE(ind);
1390     GMX_UNUSED_VALUE(name);
1391 #endif
1392 }
1393
1394 /* TODO: If/when TNG acquires the ability to copy data blocks without
1395  * uncompressing them, then this implemenation should be reconsidered.
1396  * Ideally, gmx trjconv -f a.tng -o b.tng -b 10 -e 20 would be fast
1397  * and lose no information. */
1398 gmx_bool gmx_read_next_tng_frame(gmx_tng_trajectory_t        gmx_tng_input,
1399                                  t_trxframe                 *fr,
1400                                  gmx_int64_t                *requestedIds,
1401                                  int                         numRequestedIds)
1402 {
1403 #if GMX_USE_TNG
1404     tng_trajectory_t        input = gmx_tng_input->tng;
1405     gmx_bool                bOK   = TRUE;
1406     tng_function_status     stat;
1407     gmx_int64_t             numberOfAtoms = -1, frameNumber = -1;
1408     gmx_int64_t             nBlocks, blockId, *blockIds = nullptr, codecId;
1409     char                    datatype      = -1;
1410     void                   *values        = nullptr;
1411     double                  frameTime     = -1.0;
1412     int                     size, blockDependency;
1413     double                  prec;
1414     const int               defaultNumIds = 5;
1415     static gmx_int64_t      fallbackRequestedIds[defaultNumIds] =
1416     {
1417         TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
1418         TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
1419         TNG_GMX_LAMBDA
1420     };
1421
1422
1423     fr->bStep     = FALSE;
1424     fr->bTime     = FALSE;
1425     fr->bLambda   = FALSE;
1426     fr->bAtoms    = FALSE;
1427     fr->bPrec     = FALSE;
1428     fr->bX        = FALSE;
1429     fr->bV        = FALSE;
1430     fr->bF        = FALSE;
1431     fr->bBox      = FALSE;
1432
1433     /* If no specific IDs were requested read all block types that can
1434      * currently be interpreted */
1435     if (!requestedIds || numRequestedIds == 0)
1436     {
1437         numRequestedIds = defaultNumIds;
1438         requestedIds    = fallbackRequestedIds;
1439     }
1440
1441     stat = tng_num_particles_get(input, &numberOfAtoms);
1442     if (stat != TNG_SUCCESS)
1443     {
1444         gmx_file("Cannot determine number of atoms from TNG file.");
1445     }
1446     fr->natoms = numberOfAtoms;
1447
1448     bool nextFrameExists = gmx_get_tng_data_block_types_of_next_frame(gmx_tng_input,
1449                                                                       fr->step,
1450                                                                       numRequestedIds,
1451                                                                       requestedIds,
1452                                                                       &frameNumber,
1453                                                                       &nBlocks,
1454                                                                       &blockIds);
1455     gmx::unique_cptr<gmx_int64_t, gmx::free_wrapper> blockIdsGuard(blockIds);
1456     if (!nextFrameExists)
1457     {
1458         return FALSE;
1459     }
1460
1461     if (nBlocks == 0)
1462     {
1463         return FALSE;
1464     }
1465
1466     for (gmx_int64_t i = 0; i < nBlocks; i++)
1467     {
1468         blockId = blockIds[i];
1469         tng_data_block_dependency_get(input, blockId, &blockDependency);
1470         if (blockDependency & TNG_PARTICLE_DEPENDENT)
1471         {
1472             stat = tng_util_particle_data_next_frame_read(input,
1473                                                           blockId,
1474                                                           &values,
1475                                                           &datatype,
1476                                                           &frameNumber,
1477                                                           &frameTime);
1478         }
1479         else
1480         {
1481             stat = tng_util_non_particle_data_next_frame_read(input,
1482                                                               blockId,
1483                                                               &values,
1484                                                               &datatype,
1485                                                               &frameNumber,
1486                                                               &frameTime);
1487         }
1488         if (stat == TNG_CRITICAL)
1489         {
1490             gmx_file("Cannot read positions from TNG file.");
1491             return FALSE;
1492         }
1493         else if (stat == TNG_FAILURE)
1494         {
1495             continue;
1496         }
1497         switch (blockId)
1498         {
1499             case TNG_TRAJ_BOX_SHAPE:
1500                 switch (datatype)
1501                 {
1502                     case TNG_INT_DATA:
1503                         size = sizeof(gmx_int64_t);
1504                         break;
1505                     case TNG_FLOAT_DATA:
1506                         size = sizeof(float);
1507                         break;
1508                     case TNG_DOUBLE_DATA:
1509                         size = sizeof(double);
1510                         break;
1511                     default:
1512                         gmx_incons("Illegal datatype of box shape values!");
1513                 }
1514                 for (int i = 0; i < DIM; i++)
1515                 {
1516                     convert_array_to_real_array(reinterpret_cast<char *>(values) + size * i * DIM,
1517                                                 reinterpret_cast<real *>(fr->box[i]),
1518                                                 getDistanceScaleFactor(gmx_tng_input),
1519                                                 1,
1520                                                 DIM,
1521                                                 datatype);
1522                 }
1523                 fr->bBox = TRUE;
1524                 break;
1525             case TNG_TRAJ_POSITIONS:
1526                 srenew(fr->x, fr->natoms);
1527                 convert_array_to_real_array(values,
1528                                             reinterpret_cast<real *>(fr->x),
1529                                             getDistanceScaleFactor(gmx_tng_input),
1530                                             fr->natoms,
1531                                             DIM,
1532                                             datatype);
1533                 fr->bX = TRUE;
1534                 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1535                 /* This must be updated if/when more lossy compression methods are added */
1536                 if (codecId == TNG_TNG_COMPRESSION)
1537                 {
1538                     fr->prec  = prec;
1539                     fr->bPrec = TRUE;
1540                 }
1541                 break;
1542             case TNG_TRAJ_VELOCITIES:
1543                 srenew(fr->v, fr->natoms);
1544                 convert_array_to_real_array(values,
1545                                             (real *) fr->v,
1546                                             getDistanceScaleFactor(gmx_tng_input),
1547                                             fr->natoms,
1548                                             DIM,
1549                                             datatype);
1550                 fr->bV = TRUE;
1551                 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1552                 /* This must be updated if/when more lossy compression methods are added */
1553                 if (codecId == TNG_TNG_COMPRESSION)
1554                 {
1555                     fr->prec  = prec;
1556                     fr->bPrec = TRUE;
1557                 }
1558                 break;
1559             case TNG_TRAJ_FORCES:
1560                 srenew(fr->f, fr->natoms);
1561                 convert_array_to_real_array(values,
1562                                             reinterpret_cast<real *>(fr->f),
1563                                             getDistanceScaleFactor(gmx_tng_input),
1564                                             fr->natoms,
1565                                             DIM,
1566                                             datatype);
1567                 fr->bF = TRUE;
1568                 break;
1569             case TNG_GMX_LAMBDA:
1570                 switch (datatype)
1571                 {
1572                     case TNG_FLOAT_DATA:
1573                         fr->lambda = *(reinterpret_cast<float *>(values));
1574                         break;
1575                     case TNG_DOUBLE_DATA:
1576                         fr->lambda = *(reinterpret_cast<double *>(values));
1577                         break;
1578                     default:
1579                         gmx_incons("Illegal datatype lambda value!");
1580                 }
1581                 fr->bLambda = TRUE;
1582                 break;
1583             default:
1584                 gmx_warning("Illegal block type! Currently GROMACS tools can only handle certain data types. Skipping block.");
1585         }
1586         /* values does not have to be freed before reading next frame. It will
1587          * be reallocated if it is not NULL. */
1588     }
1589
1590     fr->step  = frameNumber;
1591     fr->bStep = TRUE;
1592
1593     // Convert the time to ps
1594     fr->time  = frameTime / PICO;
1595     fr->bTime = (frameTime > 0);
1596
1597     // TODO This does not leak, but is not exception safe.
1598     /* values must be freed before leaving this function */
1599     sfree(values);
1600
1601     return bOK;
1602 #else
1603     GMX_UNUSED_VALUE(gmx_tng_input);
1604     GMX_UNUSED_VALUE(fr);
1605     GMX_UNUSED_VALUE(requestedIds);
1606     GMX_UNUSED_VALUE(numRequestedIds);
1607     return FALSE;
1608 #endif
1609 }
1610
1611 void gmx_print_tng_molecule_system(gmx_tng_trajectory_t gmx_tng_input,
1612                                    FILE                *stream)
1613 {
1614 #if GMX_USE_TNG
1615     gmx_int64_t         nMolecules, nChains, nResidues, nAtoms, nFramesRead;
1616     gmx_int64_t         strideLength, nParticlesRead, nValuesPerFrameRead, *molCntList;
1617     tng_molecule_t      molecule;
1618     tng_chain_t         chain;
1619     tng_residue_t       residue;
1620     tng_atom_t          atom;
1621     tng_function_status stat;
1622     char                str[256];
1623     char                varNAtoms;
1624     char                datatype;
1625     void               *data = nullptr;
1626     std::vector<real>   atomCharges;
1627     std::vector<real>   atomMasses;
1628     tng_trajectory_t    input = gmx_tng_input->tng;
1629
1630     tng_num_molecule_types_get(input, &nMolecules);
1631     tng_molecule_cnt_list_get(input, &molCntList);
1632     /* Can the number of particles change in the trajectory or is it constant? */
1633     tng_num_particles_variable_get(input, &varNAtoms);
1634
1635     for (gmx_int64_t i = 0; i < nMolecules; i++)
1636     {
1637         tng_molecule_of_index_get(input, i, &molecule);
1638         tng_molecule_name_get(input, molecule, str, 256);
1639         if (varNAtoms == TNG_CONSTANT_N_ATOMS)
1640         {
1641             if ((int)molCntList[i] == 0)
1642             {
1643                 continue;
1644             }
1645             fprintf(stream, "Molecule: %s, count: %d\n", str, (int)molCntList[i]);
1646         }
1647         else
1648         {
1649             fprintf(stream, "Molecule: %s\n", str);
1650         }
1651         tng_molecule_num_chains_get(input, molecule, &nChains);
1652         if (nChains > 0)
1653         {
1654             for (gmx_int64_t j = 0; j < nChains; j++)
1655             {
1656                 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
1657                 tng_chain_name_get(input, chain, str, 256);
1658                 fprintf(stream, "\tChain: %s\n", str);
1659                 tng_chain_num_residues_get(input, chain, &nResidues);
1660                 for (gmx_int64_t k = 0; k < nResidues; k++)
1661                 {
1662                     tng_chain_residue_of_index_get(input, chain, k, &residue);
1663                     tng_residue_name_get(input, residue, str, 256);
1664                     fprintf(stream, "\t\tResidue: %s\n", str);
1665                     tng_residue_num_atoms_get(input, residue, &nAtoms);
1666                     for (gmx_int64_t l = 0; l < nAtoms; l++)
1667                     {
1668                         tng_residue_atom_of_index_get(input, residue, l, &atom);
1669                         tng_atom_name_get(input, atom, str, 256);
1670                         fprintf(stream, "\t\t\tAtom: %s", str);
1671                         tng_atom_type_get(input, atom, str, 256);
1672                         fprintf(stream, " (%s)\n", str);
1673                     }
1674                 }
1675             }
1676         }
1677         /* It is possible to have a molecule without chains, in which case
1678          * residues in the molecule can be iterated through without going
1679          * through chains. */
1680         else
1681         {
1682             tng_molecule_num_residues_get(input, molecule, &nResidues);
1683             if (nResidues > 0)
1684             {
1685                 for (gmx_int64_t k = 0; k < nResidues; k++)
1686                 {
1687                     tng_molecule_residue_of_index_get(input, molecule, k, &residue);
1688                     tng_residue_name_get(input, residue, str, 256);
1689                     fprintf(stream, "\t\tResidue: %s\n", str);
1690                     tng_residue_num_atoms_get(input, residue, &nAtoms);
1691                     for (gmx_int64_t l = 0; l < nAtoms; l++)
1692                     {
1693                         tng_residue_atom_of_index_get(input, residue, l, &atom);
1694                         tng_atom_name_get(input, atom, str, 256);
1695                         fprintf(stream, "\t\t\tAtom: %s", str);
1696                         tng_atom_type_get(input, atom, str, 256);
1697                         fprintf(stream, " (%s)\n", str);
1698                     }
1699                 }
1700             }
1701             else
1702             {
1703                 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
1704                 for (gmx_int64_t l = 0; l < nAtoms; l++)
1705                 {
1706                     tng_molecule_atom_of_index_get(input, molecule, l, &atom);
1707                     tng_atom_name_get(input, atom, str, 256);
1708                     fprintf(stream, "\t\t\tAtom: %s", str);
1709                     tng_atom_type_get(input, atom, str, 256);
1710                     fprintf(stream, " (%s)\n", str);
1711                 }
1712             }
1713         }
1714     }
1715
1716     tng_num_particles_get(input, &nAtoms);
1717     stat = tng_particle_data_vector_get(input, TNG_TRAJ_PARTIAL_CHARGES, &data, &nFramesRead,
1718                                         &strideLength, &nParticlesRead,
1719                                         &nValuesPerFrameRead, &datatype);
1720     if (stat == TNG_SUCCESS)
1721     {
1722         atomCharges.resize(nAtoms);
1723         convert_array_to_real_array(data,
1724                                     atomCharges.data(),
1725                                     1,
1726                                     nAtoms,
1727                                     1,
1728                                     datatype);
1729
1730         fprintf(stream, "Atom Charges (%d):\n", int(nAtoms));
1731         for (gmx_int64_t i = 0; i < nAtoms; i += 10)
1732         {
1733             fprintf(stream, "Atom Charges [%8d-]=[", int(i));
1734             for (gmx_int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1735             {
1736                 fprintf(stream, " %12.5e", atomCharges[i + j]);
1737             }
1738             fprintf(stream, "]\n");
1739         }
1740     }
1741
1742     stat = tng_particle_data_vector_get(input, TNG_TRAJ_MASSES, &data, &nFramesRead,
1743                                         &strideLength, &nParticlesRead,
1744                                         &nValuesPerFrameRead, &datatype);
1745     if (stat == TNG_SUCCESS)
1746     {
1747         atomMasses.resize(nAtoms);
1748         convert_array_to_real_array(data,
1749                                     atomMasses.data(),
1750                                     1,
1751                                     nAtoms,
1752                                     1,
1753                                     datatype);
1754
1755         fprintf(stream, "Atom Masses (%d):\n", int(nAtoms));
1756         for (gmx_int64_t i = 0; i < nAtoms; i += 10)
1757         {
1758             fprintf(stream, "Atom Masses [%8d-]=[", int(i));
1759             for (gmx_int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1760             {
1761                 fprintf(stream, " %12.5e", atomMasses[i + j]);
1762             }
1763             fprintf(stream, "]\n");
1764         }
1765     }
1766
1767     sfree(data);
1768 #else
1769     GMX_UNUSED_VALUE(gmx_tng_input);
1770     GMX_UNUSED_VALUE(stream);
1771 #endif
1772 }
1773
1774 gmx_bool gmx_get_tng_data_block_types_of_next_frame(gmx_tng_trajectory_t gmx_tng_input,
1775                                                     int                  frame,
1776                                                     int                  nRequestedIds,
1777                                                     gmx_int64_t         *requestedIds,
1778                                                     gmx_int64_t         *nextFrame,
1779                                                     gmx_int64_t         *nBlocks,
1780                                                     gmx_int64_t        **blockIds)
1781 {
1782 #if GMX_USE_TNG
1783     tng_function_status stat;
1784     tng_trajectory_t    input = gmx_tng_input->tng;
1785
1786     stat = tng_util_trajectory_next_frame_present_data_blocks_find(input, frame,
1787                                                                    nRequestedIds, requestedIds,
1788                                                                    nextFrame,
1789                                                                    nBlocks, blockIds);
1790
1791     if (stat == TNG_CRITICAL)
1792     {
1793         gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
1794     }
1795     else if (stat == TNG_FAILURE)
1796     {
1797         return FALSE;
1798     }
1799     return TRUE;
1800 #else
1801     GMX_UNUSED_VALUE(gmx_tng_input);
1802     GMX_UNUSED_VALUE(frame);
1803     GMX_UNUSED_VALUE(nRequestedIds);
1804     GMX_UNUSED_VALUE(requestedIds);
1805     GMX_UNUSED_VALUE(nextFrame);
1806     GMX_UNUSED_VALUE(nBlocks);
1807     GMX_UNUSED_VALUE(blockIds);
1808     return FALSE;
1809 #endif
1810 }
1811
1812 gmx_bool gmx_get_tng_data_next_frame_of_block_type(gmx_tng_trajectory_t gmx_tng_input,
1813                                                    gmx_int64_t          blockId,
1814                                                    real               **values,
1815                                                    gmx_int64_t         *frameNumber,
1816                                                    double              *frameTime,
1817                                                    gmx_int64_t         *nValuesPerFrame,
1818                                                    gmx_int64_t         *nAtoms,
1819                                                    real                *prec,
1820                                                    char                *name,
1821                                                    int                  maxLen,
1822                                                    gmx_bool            *bOK)
1823 {
1824 #if GMX_USE_TNG
1825     tng_function_status stat;
1826     char                datatype = -1;
1827     gmx_int64_t         codecId;
1828     int                 blockDependency;
1829     void               *data = nullptr;
1830     double              localPrec;
1831     tng_trajectory_t    input = gmx_tng_input->tng;
1832
1833     stat = tng_data_block_name_get(input, blockId, name, maxLen);
1834     if (stat != TNG_SUCCESS)
1835     {
1836         gmx_file("Cannot read next frame of TNG file");
1837     }
1838     stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
1839     if (stat != TNG_SUCCESS)
1840     {
1841         gmx_file("Cannot read next frame of TNG file");
1842     }
1843     if (blockDependency & TNG_PARTICLE_DEPENDENT)
1844     {
1845         tng_num_particles_get(input, nAtoms);
1846         stat = tng_util_particle_data_next_frame_read(input,
1847                                                       blockId,
1848                                                       &data,
1849                                                       &datatype,
1850                                                       frameNumber,
1851                                                       frameTime);
1852     }
1853     else
1854     {
1855         *nAtoms = 1; /* There are not actually any atoms, but it is used for
1856                         allocating memory */
1857         stat    = tng_util_non_particle_data_next_frame_read(input,
1858                                                              blockId,
1859                                                              &data,
1860                                                              &datatype,
1861                                                              frameNumber,
1862                                                              frameTime);
1863     }
1864     if (stat == TNG_CRITICAL)
1865     {
1866         gmx_file("Cannot read next frame of TNG file");
1867     }
1868     if (stat == TNG_FAILURE)
1869     {
1870         *bOK = TRUE;
1871         return FALSE;
1872     }
1873
1874     stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
1875     if (stat != TNG_SUCCESS)
1876     {
1877         gmx_file("Cannot read next frame of TNG file");
1878     }
1879     srenew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
1880     convert_array_to_real_array(data,
1881                                 *values,
1882                                 getDistanceScaleFactor(gmx_tng_input),
1883                                 *nAtoms,
1884                                 *nValuesPerFrame,
1885                                 datatype);
1886
1887     tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
1888
1889     /* This must be updated if/when more lossy compression methods are added */
1890     if (codecId != TNG_TNG_COMPRESSION)
1891     {
1892         *prec = -1.0;
1893     }
1894     else
1895     {
1896         *prec = localPrec;
1897     }
1898
1899     sfree(data);
1900
1901     *bOK = TRUE;
1902     return TRUE;
1903 #else
1904     GMX_UNUSED_VALUE(gmx_tng_input);
1905     GMX_UNUSED_VALUE(blockId);
1906     GMX_UNUSED_VALUE(values);
1907     GMX_UNUSED_VALUE(frameNumber);
1908     GMX_UNUSED_VALUE(frameTime);
1909     GMX_UNUSED_VALUE(nValuesPerFrame);
1910     GMX_UNUSED_VALUE(nAtoms);
1911     GMX_UNUSED_VALUE(prec);
1912     GMX_UNUSED_VALUE(name);
1913     GMX_UNUSED_VALUE(maxLen);
1914     GMX_UNUSED_VALUE(bOK);
1915     return FALSE;
1916 #endif
1917 }
1918
1919 int gmx_tng_get_box_output_interval(gmx_tng_trajectory_t gmx_tng)
1920 {
1921 #if GMX_USE_TNG
1922     return gmx_tng->boxOutputInterval;
1923 #else
1924     GMX_UNUSED_VALUE(gmx_tng);
1925 #endif
1926 }
1927
1928 int gmx_tng_get_lambda_output_interval(gmx_tng_trajectory_t gmx_tng)
1929 {
1930 #if GMX_USE_TNG
1931     return gmx_tng->lambdaOutputInterval;
1932 #else
1933     GMX_UNUSED_VALUE(gmx_tng);
1934 #endif
1935 }