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