Fix range bug in add_selection_groups in tngio.cpp
[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 (int molIndex = 0; molIndex < mtop->nmolblock; molIndex++)
312     {
313         tng_molecule_t       tngMol  = nullptr;
314         const gmx_moltype_t *molType = &mtop->moltype[mtop->molblock[molIndex].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                                    mtop->molblock[molIndex].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 < mtop->molblock[molIndex].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     const gmx_moltype_t     *molType;
578     const t_atoms           *atoms;
579
580     /* Iterate through all atoms in the molecule system and
581      * check if they belong to a selection group of the
582      * requested type. */
583     for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
584     {
585         molType = &mtop->moltype[mtop->molblock[molIndex].type];
586
587         atoms = &molType->atoms;
588
589         for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
590         {
591             for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
592             {
593                 if (ggrpnr(&mtop->groups, gtype, i) != 0)
594                 {
595                     return FALSE;
596                 }
597             }
598         }
599     }
600     return TRUE;
601 }
602
603 /* Create TNG molecules which will represent each of the selection
604  * groups for writing. But do that only if there is actually a
605  * specified selection group and it is not the whole system.
606  * TODO: Currently the only selection that is taken into account
607  * is egcCompressedX, but other selections should be added when
608  * e.g. writing energies is implemented.
609  */
610 static void add_selection_groups(gmx_tng_trajectory_t  gmx_tng,
611                                  const gmx_mtop_t     *mtop)
612 {
613     const gmx_moltype_t     *molType;
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     for (int molIndex = 0, i = 0; molIndex < mtop->nmolblock; molIndex++)
657     {
658         molType = &mtop->moltype[mtop->molblock[molIndex].type];
659
660         atoms = &molType->atoms;
661
662         for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
663         {
664             bool bAtomsAdded = FALSE;
665             for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
666             {
667                 char *res_name;
668                 int   res_id;
669
670                 if (ggrpnr(&mtop->groups, egcCompressedX, i) != 0)
671                 {
672                     continue;
673                 }
674                 at = &atoms->atom[atomIndex];
675                 if (atoms->nres > 0)
676                 {
677                     resInfo = &atoms->resinfo[at->resind];
678                     /* FIXME: When TNG supports both residue index and residue
679                      * number the latter should be used. */
680                     res_name = *resInfo->name;
681                     res_id   = at->resind + 1;
682                 }
683                 else
684                 {
685                     res_name = (char *)"";
686                     res_id   = 0;
687                 }
688                 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res)
689                     != TNG_SUCCESS)
690                 {
691                     /* Since there is ONE chain for selection groups do not keep the
692                      * original residue IDs - otherwise there might be conflicts. */
693                     tng_chain_residue_add(tng, chain, res_name, &res);
694                 }
695                 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
696                                           *(atoms->atomtype[atomIndex]),
697                                           atom_offset + atomIndex, &atom);
698                 bAtomsAdded = TRUE;
699             }
700             /* Add bonds. */
701             if (bAtomsAdded)
702             {
703                 for (int k = 0; k < F_NRE; k++)
704                 {
705                     if (IS_CHEMBOND(k))
706                     {
707                         ilist = &molType->ilist[k];
708                         if (ilist)
709                         {
710                             int l = 1;
711                             while (l < ilist->nr)
712                             {
713                                 int atom1, atom2;
714                                 atom1 = ilist->iatoms[l] + atom_offset;
715                                 atom2 = ilist->iatoms[l+1] + atom_offset;
716                                 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0 &&
717                                     ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
718                                 {
719                                     tng_molecule_bond_add(tng, mol, ilist->iatoms[l],
720                                                           ilist->iatoms[l+1], &tngBond);
721                                 }
722                                 l += 3;
723                             }
724                         }
725                     }
726                 }
727                 /* Settle is described using three atoms */
728                 ilist = &molType->ilist[F_SETTLE];
729                 if (ilist)
730                 {
731                     int l = 1;
732                     while (l < ilist->nr)
733                     {
734                         int atom1, atom2, atom3;
735                         atom1 = ilist->iatoms[l] + atom_offset;
736                         atom2 = ilist->iatoms[l+1] + atom_offset;
737                         atom3 = ilist->iatoms[l+2] + atom_offset;
738                         if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0)
739                         {
740                             if (ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
741                             {
742                                 tng_molecule_bond_add(tng, mol, atom1,
743                                                       atom2, &tngBond);
744                             }
745                             if (ggrpnr(&mtop->groups, egcCompressedX, atom3) == 0)
746                             {
747                                 tng_molecule_bond_add(tng, mol, atom1,
748                                                       atom3, &tngBond);
749                             }
750                         }
751                         l += 4;
752                     }
753                 }
754             }
755             atom_offset += atoms->nr;
756         }
757     }
758     tng_molecule_existing_add(tng, &mol);
759     tng_molecule_cnt_set(tng, mol, 1);
760     tng_num_molecule_types_get(tng, &nMols);
761     for (gmx_int64_t k = 0; k < nMols; k++)
762     {
763         tng_molecule_of_index_get(tng, k, &iterMol);
764         if (iterMol == mol)
765         {
766             continue;
767         }
768         tng_molecule_cnt_set(tng, iterMol, 0);
769     }
770 }
771 #endif
772
773 void gmx_tng_set_compression_precision(gmx_tng_trajectory_t gmx_tng,
774                                        real                 prec)
775 {
776 #if GMX_USE_TNG
777     tng_compression_precision_set(gmx_tng->tng, prec);
778 #else
779     GMX_UNUSED_VALUE(gmx_tng);
780     GMX_UNUSED_VALUE(prec);
781 #endif
782 }
783
784 void gmx_tng_prepare_low_prec_writing(gmx_tng_trajectory_t  gmx_tng,
785                                       const gmx_mtop_t     *mtop,
786                                       const t_inputrec     *ir)
787 {
788 #if GMX_USE_TNG
789     gmx_tng_add_mtop(gmx_tng, mtop);
790     add_selection_groups(gmx_tng, mtop);
791     set_writing_intervals(gmx_tng, TRUE, ir);
792     tng_time_per_frame_set(gmx_tng->tng, ir->delta_t * PICO);
793     gmx_tng->timePerFrameIsSet = true;
794     gmx_tng_set_compression_precision(gmx_tng, ir->x_compression_precision);
795 #else
796     GMX_UNUSED_VALUE(gmx_tng);
797     GMX_UNUSED_VALUE(mtop);
798     GMX_UNUSED_VALUE(ir);
799 #endif
800 }
801
802 void gmx_fwrite_tng(gmx_tng_trajectory_t gmx_tng,
803                     const gmx_bool       bUseLossyCompression,
804                     gmx_int64_t          step,
805                     real                 elapsedPicoSeconds,
806                     real                 lambda,
807                     const rvec          *box,
808                     int                  nAtoms,
809                     const rvec          *x,
810                     const rvec          *v,
811                     const rvec          *f)
812 {
813 #if GMX_USE_TNG
814     typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
815                                                            const gmx_int64_t,
816                                                            const double,
817                                                            const real*,
818                                                            const gmx_int64_t,
819                                                            const gmx_int64_t,
820                                                            const char*,
821                                                            const char,
822                                                            const char);
823 #if GMX_DOUBLE
824     static write_data_func_pointer           write_data           = tng_util_generic_with_time_double_write;
825 #else
826     static write_data_func_pointer           write_data           = tng_util_generic_with_time_write;
827 #endif
828     double                                   elapsedSeconds = elapsedPicoSeconds * PICO;
829     gmx_int64_t                              nParticles;
830     char                                     compression;
831
832
833     if (!gmx_tng)
834     {
835         /* This function might get called when the type of the
836            compressed trajectory is actually XTC. So we exit and move
837            on. */
838         return;
839     }
840     tng_trajectory_t tng = gmx_tng->tng;
841
842     // While the GROMACS interface to this routine specifies 'step', TNG itself
843     // only uses 'frame index' internally, although it suggests that it's a good
844     // idea to let this represent the step rather than just frame numbers.
845     //
846     // However, the way the GROMACS interface works, there are cases where
847     // the step is simply not set, so all frames rather have step=0.
848     // If we call the lower-level TNG routines with the same frame index
849     // (which is set from the step) multiple times, things will go very
850     // wrong (overwritten data).
851     //
852     // To avoid this, we require the step value to always be larger than the
853     // previous value, and if this is not true we simply set it to a value
854     // one unit larger than the previous step.
855     //
856     // Note: We do allow the user to specify any crazy value the want for the
857     // first step. If it is "not set", GROMACS will have used the value 0.
858     if (gmx_tng->lastStepDataIsValid && step <= gmx_tng->lastStep)
859     {
860         step = gmx_tng->lastStep + 1;
861     }
862
863     // Now that we have fixed the potentially incorrect step, we can also set
864     // the time per frame if it was not already done, and we have
865     // step/time information from the last frame so it is possible to calculate it.
866     // Note that we do allow the time to be the same for all steps, or even
867     // decreasing. In that case we will simply say the time per step is 0
868     // or negative. A bit strange, but a correct representation of the data :-)
869     if (!gmx_tng->timePerFrameIsSet && gmx_tng->lastTimeDataIsValid && gmx_tng->lastStepDataIsValid)
870     {
871         double       deltaTime = elapsedSeconds - gmx_tng->lastTime;
872         std::int64_t deltaStep = step - gmx_tng->lastStep;
873         tng_time_per_frame_set(tng, deltaTime / deltaStep );
874         gmx_tng->timePerFrameIsSet = true;
875     }
876
877     tng_num_particles_get(tng, &nParticles);
878     if (nAtoms != (int)nParticles)
879     {
880         tng_implicit_num_particles_set(tng, nAtoms);
881     }
882
883     if (bUseLossyCompression)
884     {
885         compression = TNG_TNG_COMPRESSION;
886     }
887     else
888     {
889         compression = TNG_GZIP_COMPRESSION;
890     }
891
892     /* The writing is done using write_data, which writes float or double
893      * depending on the GROMACS compilation. */
894     if (x)
895     {
896         GMX_ASSERT(box, "Need a non-NULL box if positions are written");
897
898         if (write_data(tng, step, elapsedSeconds,
899                        reinterpret_cast<const real *>(x),
900                        3, TNG_TRAJ_POSITIONS, "POSITIONS",
901                        TNG_PARTICLE_BLOCK_DATA,
902                        compression) != TNG_SUCCESS)
903         {
904             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
905         }
906     }
907
908     if (v)
909     {
910         if (write_data(tng, step, elapsedSeconds,
911                        reinterpret_cast<const real *>(v),
912                        3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
913                        TNG_PARTICLE_BLOCK_DATA,
914                        compression) != TNG_SUCCESS)
915         {
916             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
917         }
918     }
919
920     if (f)
921     {
922         /* TNG-MF1 compression only compresses positions and velocities. Use lossless
923          * compression for forces regardless of output mode */
924         if (write_data(tng, step, elapsedSeconds,
925                        reinterpret_cast<const real *>(f),
926                        3, TNG_TRAJ_FORCES, "FORCES",
927                        TNG_PARTICLE_BLOCK_DATA,
928                        TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
929         {
930             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
931         }
932     }
933
934     /* TNG-MF1 compression only compresses positions and velocities. Use lossless
935      * compression for lambdas and box shape regardless of output mode */
936     if (write_data(tng, step, elapsedSeconds,
937                    reinterpret_cast<const real *>(box),
938                    9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
939                    TNG_NON_PARTICLE_BLOCK_DATA,
940                    TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
941     {
942         gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
943     }
944
945     if (write_data(tng, step, elapsedSeconds,
946                    reinterpret_cast<const real *>(&lambda),
947                    1, TNG_GMX_LAMBDA, "LAMBDAS",
948                    TNG_NON_PARTICLE_BLOCK_DATA,
949                    TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
950     {
951         gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
952     }
953
954     // Update the data in the wrapper
955
956     gmx_tng->lastStepDataIsValid = true;
957     gmx_tng->lastStep            = step;
958     gmx_tng->lastTimeDataIsValid = true;
959     gmx_tng->lastTime            = elapsedSeconds;
960 #else
961     GMX_UNUSED_VALUE(gmx_tng);
962     GMX_UNUSED_VALUE(bUseLossyCompression);
963     GMX_UNUSED_VALUE(step);
964     GMX_UNUSED_VALUE(elapsedPicoSeconds);
965     GMX_UNUSED_VALUE(lambda);
966     GMX_UNUSED_VALUE(box);
967     GMX_UNUSED_VALUE(nAtoms);
968     GMX_UNUSED_VALUE(x);
969     GMX_UNUSED_VALUE(v);
970     GMX_UNUSED_VALUE(f);
971 #endif
972 }
973
974 void fflush_tng(gmx_tng_trajectory_t gmx_tng)
975 {
976 #if GMX_USE_TNG
977     if (!gmx_tng)
978     {
979         return;
980     }
981     tng_frame_set_premature_write(gmx_tng->tng, TNG_USE_HASH);
982 #else
983     GMX_UNUSED_VALUE(gmx_tng);
984 #endif
985 }
986
987 float gmx_tng_get_time_of_final_frame(gmx_tng_trajectory_t gmx_tng)
988 {
989 #if GMX_USE_TNG
990     gmx_int64_t      nFrames;
991     double           time;
992     float            fTime;
993     tng_trajectory_t tng = gmx_tng->tng;
994
995     tng_num_frames_get(tng, &nFrames);
996     tng_util_time_of_frame_get(tng, nFrames - 1, &time);
997
998     fTime = time / PICO;
999     return fTime;
1000 #else
1001     GMX_UNUSED_VALUE(gmx_tng);
1002     return -1.0;
1003 #endif
1004 }
1005
1006 void gmx_prepare_tng_writing(const char              *filename,
1007                              char                     mode,
1008                              gmx_tng_trajectory_t    *gmx_tng_input,
1009                              gmx_tng_trajectory_t    *gmx_tng_output,
1010                              int                      nAtoms,
1011                              const gmx_mtop_t        *mtop,
1012                              const int               *index,
1013                              const char              *indexGroupName)
1014 {
1015 #if GMX_USE_TNG
1016     tng_trajectory_t   *input  = (gmx_tng_input && *gmx_tng_input) ? &(*gmx_tng_input)->tng : nullptr;
1017     /* FIXME after 5.0: Currently only standard block types are read */
1018     const int           defaultNumIds              = 5;
1019     static gmx_int64_t  fallbackIds[defaultNumIds] =
1020     {
1021         TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
1022         TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
1023         TNG_GMX_LAMBDA
1024     };
1025     static char         fallbackNames[defaultNumIds][32] =
1026     {
1027         "BOX SHAPE", "POSITIONS", "VELOCITIES",
1028         "FORCES", "LAMBDAS"
1029     };
1030
1031     typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
1032                                                                      const gmx_int64_t,
1033                                                                      const gmx_int64_t,
1034                                                                      const gmx_int64_t,
1035                                                                      const char*,
1036                                                                      const char,
1037                                                                      const char);
1038 #if GMX_DOUBLE
1039     set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
1040 #else
1041     set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
1042 #endif
1043
1044     gmx_tng_open(filename, mode, gmx_tng_output);
1045     tng_trajectory_t *output = &(*gmx_tng_output)->tng;
1046
1047     /* Do we have an input file in TNG format? If so, then there's
1048        more data we can copy over, rather than having to improvise. */
1049     if (gmx_tng_input && *gmx_tng_input)
1050     {
1051         /* Set parameters (compression, time per frame, molecule
1052          * information, number of frames per frame set and writing
1053          * intervals of positions, box shape and lambdas) of the
1054          * output tng container based on their respective values int
1055          * the input tng container */
1056         double      time, compression_precision;
1057         gmx_int64_t n_frames_per_frame_set, interval = -1;
1058
1059         tng_compression_precision_get(*input, &compression_precision);
1060         tng_compression_precision_set(*output, compression_precision);
1061         // TODO make this configurable in a future version
1062         char compression_type = TNG_TNG_COMPRESSION;
1063
1064         tng_molecule_system_copy(*input, *output);
1065
1066         tng_time_per_frame_get(*input, &time);
1067         tng_time_per_frame_set(*output, time);
1068         // Since we have copied the value from the input TNG we should not change it again
1069         (*gmx_tng_output)->timePerFrameIsSet = true;
1070
1071         tng_num_frames_per_frame_set_get(*input, &n_frames_per_frame_set);
1072         tng_num_frames_per_frame_set_set(*output, n_frames_per_frame_set);
1073
1074         for (int i = 0; i < defaultNumIds; i++)
1075         {
1076             if (tng_data_get_stride_length(*input, fallbackIds[i], -1, &interval)
1077                 == TNG_SUCCESS)
1078             {
1079                 switch (fallbackIds[i])
1080                 {
1081                     case TNG_TRAJ_POSITIONS:
1082                     case TNG_TRAJ_VELOCITIES:
1083                         set_writing_interval(*output, interval, 3, fallbackIds[i],
1084                                              fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
1085                                              compression_type);
1086                         break;
1087                     case TNG_TRAJ_FORCES:
1088                         set_writing_interval(*output, interval, 3, fallbackIds[i],
1089                                              fallbackNames[i], TNG_PARTICLE_BLOCK_DATA,
1090                                              TNG_GZIP_COMPRESSION);
1091                         break;
1092                     case TNG_TRAJ_BOX_SHAPE:
1093                         set_writing_interval(*output, interval, 9, fallbackIds[i],
1094                                              fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
1095                                              TNG_GZIP_COMPRESSION);
1096                         break;
1097                     case TNG_GMX_LAMBDA:
1098                         set_writing_interval(*output, interval, 1, fallbackIds[i],
1099                                              fallbackNames[i], TNG_NON_PARTICLE_BLOCK_DATA,
1100                                              TNG_GZIP_COMPRESSION);
1101                         break;
1102                     default:
1103                         continue;
1104                 }
1105             }
1106         }
1107
1108     }
1109     else
1110     {
1111         /* TODO after trjconv is modularized: fix this so the user can
1112            change precision when they are doing an operation where
1113            this makes sense, and not otherwise.
1114
1115            char compression = bUseLossyCompression ? TNG_TNG_COMPRESSION : TNG_GZIP_COMPRESSION;
1116            gmx_tng_set_compression_precision(*output, ndec2prec(nDecimalsOfPrecision));
1117          */
1118         gmx_tng_add_mtop(*gmx_tng_output, mtop);
1119         tng_num_frames_per_frame_set_set(*output, 1);
1120     }
1121
1122     if (index && nAtoms > 0)
1123     {
1124         gmx_tng_setup_atom_subgroup(*gmx_tng_output, nAtoms, index, indexGroupName);
1125     }
1126
1127     /* If for some reason there are more requested atoms than there are atoms in the
1128      * molecular system create a number of implicit atoms (without atom data) to
1129      * compensate for that. */
1130     if (nAtoms >= 0)
1131     {
1132         tng_implicit_num_particles_set(*output, nAtoms);
1133     }
1134 #else
1135     GMX_UNUSED_VALUE(filename);
1136     GMX_UNUSED_VALUE(mode);
1137     GMX_UNUSED_VALUE(gmx_tng_input);
1138     GMX_UNUSED_VALUE(gmx_tng_output);
1139     GMX_UNUSED_VALUE(nAtoms);
1140     GMX_UNUSED_VALUE(mtop);
1141     GMX_UNUSED_VALUE(index);
1142     GMX_UNUSED_VALUE(indexGroupName);
1143 #endif
1144 }
1145
1146 void gmx_write_tng_from_trxframe(gmx_tng_trajectory_t    gmx_tng_output,
1147                                  const t_trxframe       *frame,
1148                                  int                     natoms)
1149 {
1150 #if GMX_USE_TNG
1151     if (natoms < 0)
1152     {
1153         natoms = frame->natoms;
1154     }
1155     gmx_fwrite_tng(gmx_tng_output,
1156                    TRUE,
1157                    frame->step,
1158                    frame->time,
1159                    0,
1160                    frame->box,
1161                    natoms,
1162                    frame->x,
1163                    frame->v,
1164                    frame->f);
1165 #else
1166     GMX_UNUSED_VALUE(gmx_tng_output);
1167     GMX_UNUSED_VALUE(frame);
1168     GMX_UNUSED_VALUE(natoms);
1169 #endif
1170 }
1171
1172 namespace
1173 {
1174
1175 #if GMX_USE_TNG
1176 void
1177 convert_array_to_real_array(void       *from,
1178                             real       *to,
1179                             const float fact,
1180                             const int   nAtoms,
1181                             const int   nValues,
1182                             const char  datatype)
1183 {
1184     int        i, j;
1185
1186     const bool useDouble = GMX_DOUBLE;
1187     switch (datatype)
1188     {
1189         case TNG_FLOAT_DATA:
1190             if (!useDouble)
1191             {
1192                 if (fact == 1)
1193                 {
1194                     memcpy(to, from, nValues * sizeof(real) * nAtoms);
1195                 }
1196                 else
1197                 {
1198                     for (i = 0; i < nAtoms; i++)
1199                     {
1200                         for (j = 0; j < nValues; j++)
1201                         {
1202                             to[i*nValues+j] = reinterpret_cast<float *>(from)[i*nValues+j] * fact;
1203                         }
1204                     }
1205                 }
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             break;
1218         case TNG_INT_DATA:
1219             for (i = 0; i < nAtoms; i++)
1220             {
1221                 for (j = 0; j < nValues; j++)
1222                 {
1223                     to[i*nValues+j] = reinterpret_cast<gmx_int64_t *>(from)[i*nValues+j] * fact;
1224                 }
1225             }
1226             break;
1227         case TNG_DOUBLE_DATA:
1228             if (sizeof(real) == sizeof(double))
1229             {
1230                 if (fact == 1)
1231                 {
1232                     memcpy(to, from, nValues * sizeof(real) * nAtoms);
1233                 }
1234                 else
1235                 {
1236                     for (i = 0; i < nAtoms; i++)
1237                     {
1238                         for (j = 0; j < nValues; j++)
1239                         {
1240                             to[i*nValues+j] = reinterpret_cast<double *>(from)[i*nValues+j] * fact;
1241                         }
1242                     }
1243                 }
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             break;
1256         default:
1257             gmx_incons("Illegal datatype when converting values to a real array!");
1258             return;
1259     }
1260     return;
1261 }
1262
1263 real getDistanceScaleFactor(gmx_tng_trajectory_t in)
1264 {
1265     gmx_int64_t exp = -1;
1266     real        distanceScaleFactor;
1267
1268     // TODO Hopefully, TNG 2.0 will do this kind of thing for us
1269     tng_distance_unit_exponential_get(in->tng, &exp);
1270
1271     // GROMACS expects distances in nm
1272     switch (exp)
1273     {
1274         case 9:
1275             distanceScaleFactor = NANO/NANO;
1276             break;
1277         case 10:
1278             distanceScaleFactor = NANO/ANGSTROM;
1279             break;
1280         default:
1281             distanceScaleFactor = pow(10.0, exp + 9.0);
1282     }
1283
1284     return distanceScaleFactor;
1285 }
1286 #endif
1287
1288 } // namespace
1289
1290 void gmx_tng_setup_atom_subgroup(gmx_tng_trajectory_t gmx_tng,
1291                                  const int            nind,
1292                                  const int           *ind,
1293                                  const char          *name)
1294 {
1295 #if GMX_USE_TNG
1296     gmx_int64_t              nAtoms, cnt, nMols;
1297     tng_molecule_t           mol, iterMol;
1298     tng_chain_t              chain;
1299     tng_residue_t            res;
1300     tng_atom_t               atom;
1301     tng_function_status      stat;
1302     tng_trajectory_t         tng = gmx_tng->tng;
1303
1304     tng_num_particles_get(tng, &nAtoms);
1305
1306     if (nAtoms == nind)
1307     {
1308         return;
1309     }
1310
1311     stat = tng_molecule_find(tng, name, -1, &mol);
1312     if (stat == TNG_SUCCESS)
1313     {
1314         tng_molecule_num_atoms_get(tng, mol, &nAtoms);
1315         tng_molecule_cnt_get(tng, mol, &cnt);
1316         if (nAtoms == nind)
1317         {
1318             stat = TNG_SUCCESS;
1319         }
1320         else
1321         {
1322             stat = TNG_FAILURE;
1323         }
1324     }
1325     if (stat == TNG_FAILURE)
1326     {
1327         /* The indexed atoms are added to one separate molecule. */
1328         tng_molecule_alloc(tng, &mol);
1329         tng_molecule_name_set(tng, mol, name);
1330         tng_molecule_chain_add(tng, mol, "", &chain);
1331
1332         for (int i = 0; i < nind; i++)
1333         {
1334             char        temp_name[256], temp_type[256];
1335
1336             /* Try to retrieve the residue name of the atom */
1337             stat = tng_residue_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1338             if (stat != TNG_SUCCESS)
1339             {
1340                 temp_name[0] = '\0';
1341             }
1342             /* Check if the molecule of the selection already contains this residue */
1343             if (tng_chain_residue_find(tng, chain, temp_name, -1, &res)
1344                 != TNG_SUCCESS)
1345             {
1346                 tng_chain_residue_add(tng, chain, temp_name, &res);
1347             }
1348             /* Try to find the original name and type of the atom */
1349             stat = tng_atom_name_of_particle_nr_get(tng, ind[i], temp_name, 256);
1350             if (stat != TNG_SUCCESS)
1351             {
1352                 temp_name[0] = '\0';
1353             }
1354             stat = tng_atom_type_of_particle_nr_get(tng, ind[i], temp_type, 256);
1355             if (stat != TNG_SUCCESS)
1356             {
1357                 temp_type[0] = '\0';
1358             }
1359             tng_residue_atom_w_id_add(tng, res, temp_name, temp_type, ind[i], &atom);
1360         }
1361         tng_molecule_existing_add(tng, &mol);
1362     }
1363     /* Set the count of the molecule containing the selected atoms to 1 and all
1364      * other molecules to 0 */
1365     tng_molecule_cnt_set(tng, mol, 1);
1366     tng_num_molecule_types_get(tng, &nMols);
1367     for (gmx_int64_t k = 0; k < nMols; k++)
1368     {
1369         tng_molecule_of_index_get(tng, k, &iterMol);
1370         if (iterMol == mol)
1371         {
1372             continue;
1373         }
1374         tng_molecule_cnt_set(tng, iterMol, 0);
1375     }
1376 #else
1377     GMX_UNUSED_VALUE(gmx_tng);
1378     GMX_UNUSED_VALUE(nind);
1379     GMX_UNUSED_VALUE(ind);
1380     GMX_UNUSED_VALUE(name);
1381 #endif
1382 }
1383
1384 /* TODO: If/when TNG acquires the ability to copy data blocks without
1385  * uncompressing them, then this implemenation should be reconsidered.
1386  * Ideally, gmx trjconv -f a.tng -o b.tng -b 10 -e 20 would be fast
1387  * and lose no information. */
1388 gmx_bool gmx_read_next_tng_frame(gmx_tng_trajectory_t        gmx_tng_input,
1389                                  t_trxframe                 *fr,
1390                                  gmx_int64_t                *requestedIds,
1391                                  int                         numRequestedIds)
1392 {
1393 #if GMX_USE_TNG
1394     tng_trajectory_t        input = gmx_tng_input->tng;
1395     gmx_bool                bOK   = TRUE;
1396     tng_function_status     stat;
1397     gmx_int64_t             numberOfAtoms = -1, frameNumber = -1;
1398     gmx_int64_t             nBlocks, blockId, *blockIds = nullptr, codecId;
1399     char                    datatype      = -1;
1400     void                   *values        = nullptr;
1401     double                  frameTime     = -1.0;
1402     int                     size, blockDependency;
1403     double                  prec;
1404     const int               defaultNumIds = 5;
1405     static gmx_int64_t      fallbackRequestedIds[defaultNumIds] =
1406     {
1407         TNG_TRAJ_BOX_SHAPE, TNG_TRAJ_POSITIONS,
1408         TNG_TRAJ_VELOCITIES, TNG_TRAJ_FORCES,
1409         TNG_GMX_LAMBDA
1410     };
1411
1412
1413     fr->bStep     = FALSE;
1414     fr->bTime     = FALSE;
1415     fr->bLambda   = FALSE;
1416     fr->bAtoms    = FALSE;
1417     fr->bPrec     = FALSE;
1418     fr->bX        = FALSE;
1419     fr->bV        = FALSE;
1420     fr->bF        = FALSE;
1421     fr->bBox      = FALSE;
1422
1423     /* If no specific IDs were requested read all block types that can
1424      * currently be interpreted */
1425     if (!requestedIds || numRequestedIds == 0)
1426     {
1427         numRequestedIds = defaultNumIds;
1428         requestedIds    = fallbackRequestedIds;
1429     }
1430
1431     stat = tng_num_particles_get(input, &numberOfAtoms);
1432     if (stat != TNG_SUCCESS)
1433     {
1434         gmx_file("Cannot determine number of atoms from TNG file.");
1435     }
1436     fr->natoms = numberOfAtoms;
1437
1438     bool nextFrameExists = gmx_get_tng_data_block_types_of_next_frame(gmx_tng_input,
1439                                                                       fr->step,
1440                                                                       numRequestedIds,
1441                                                                       requestedIds,
1442                                                                       &frameNumber,
1443                                                                       &nBlocks,
1444                                                                       &blockIds);
1445     gmx::unique_cptr<gmx_int64_t, gmx::free_wrapper> blockIdsGuard(blockIds);
1446     if (!nextFrameExists)
1447     {
1448         return FALSE;
1449     }
1450
1451     if (nBlocks == 0)
1452     {
1453         return FALSE;
1454     }
1455
1456     for (gmx_int64_t i = 0; i < nBlocks; i++)
1457     {
1458         blockId = blockIds[i];
1459         tng_data_block_dependency_get(input, blockId, &blockDependency);
1460         if (blockDependency & TNG_PARTICLE_DEPENDENT)
1461         {
1462             stat = tng_util_particle_data_next_frame_read(input,
1463                                                           blockId,
1464                                                           &values,
1465                                                           &datatype,
1466                                                           &frameNumber,
1467                                                           &frameTime);
1468         }
1469         else
1470         {
1471             stat = tng_util_non_particle_data_next_frame_read(input,
1472                                                               blockId,
1473                                                               &values,
1474                                                               &datatype,
1475                                                               &frameNumber,
1476                                                               &frameTime);
1477         }
1478         if (stat == TNG_CRITICAL)
1479         {
1480             gmx_file("Cannot read positions from TNG file.");
1481             return FALSE;
1482         }
1483         else if (stat == TNG_FAILURE)
1484         {
1485             continue;
1486         }
1487         switch (blockId)
1488         {
1489             case TNG_TRAJ_BOX_SHAPE:
1490                 switch (datatype)
1491                 {
1492                     case TNG_INT_DATA:
1493                         size = sizeof(gmx_int64_t);
1494                         break;
1495                     case TNG_FLOAT_DATA:
1496                         size = sizeof(float);
1497                         break;
1498                     case TNG_DOUBLE_DATA:
1499                         size = sizeof(double);
1500                         break;
1501                     default:
1502                         gmx_incons("Illegal datatype of box shape values!");
1503                 }
1504                 for (int i = 0; i < DIM; i++)
1505                 {
1506                     convert_array_to_real_array(reinterpret_cast<char *>(values) + size * i * DIM,
1507                                                 reinterpret_cast<real *>(fr->box[i]),
1508                                                 getDistanceScaleFactor(gmx_tng_input),
1509                                                 1,
1510                                                 DIM,
1511                                                 datatype);
1512                 }
1513                 fr->bBox = TRUE;
1514                 break;
1515             case TNG_TRAJ_POSITIONS:
1516                 srenew(fr->x, fr->natoms);
1517                 convert_array_to_real_array(values,
1518                                             reinterpret_cast<real *>(fr->x),
1519                                             getDistanceScaleFactor(gmx_tng_input),
1520                                             fr->natoms,
1521                                             DIM,
1522                                             datatype);
1523                 fr->bX = TRUE;
1524                 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1525                 /* This must be updated if/when more lossy compression methods are added */
1526                 if (codecId == TNG_TNG_COMPRESSION)
1527                 {
1528                     fr->prec  = prec;
1529                     fr->bPrec = TRUE;
1530                 }
1531                 break;
1532             case TNG_TRAJ_VELOCITIES:
1533                 srenew(fr->v, fr->natoms);
1534                 convert_array_to_real_array(values,
1535                                             (real *) fr->v,
1536                                             getDistanceScaleFactor(gmx_tng_input),
1537                                             fr->natoms,
1538                                             DIM,
1539                                             datatype);
1540                 fr->bV = TRUE;
1541                 tng_util_frame_current_compression_get(input, blockId, &codecId, &prec);
1542                 /* This must be updated if/when more lossy compression methods are added */
1543                 if (codecId == TNG_TNG_COMPRESSION)
1544                 {
1545                     fr->prec  = prec;
1546                     fr->bPrec = TRUE;
1547                 }
1548                 break;
1549             case TNG_TRAJ_FORCES:
1550                 srenew(fr->f, fr->natoms);
1551                 convert_array_to_real_array(values,
1552                                             reinterpret_cast<real *>(fr->f),
1553                                             getDistanceScaleFactor(gmx_tng_input),
1554                                             fr->natoms,
1555                                             DIM,
1556                                             datatype);
1557                 fr->bF = TRUE;
1558                 break;
1559             case TNG_GMX_LAMBDA:
1560                 switch (datatype)
1561                 {
1562                     case TNG_FLOAT_DATA:
1563                         fr->lambda = *(reinterpret_cast<float *>(values));
1564                         break;
1565                     case TNG_DOUBLE_DATA:
1566                         fr->lambda = *(reinterpret_cast<double *>(values));
1567                         break;
1568                     default:
1569                         gmx_incons("Illegal datatype lambda value!");
1570                 }
1571                 fr->bLambda = TRUE;
1572                 break;
1573             default:
1574                 gmx_warning("Illegal block type! Currently GROMACS tools can only handle certain data types. Skipping block.");
1575         }
1576         /* values does not have to be freed before reading next frame. It will
1577          * be reallocated if it is not NULL. */
1578     }
1579
1580     fr->step  = frameNumber;
1581     fr->bStep = TRUE;
1582
1583     // Convert the time to ps
1584     fr->time  = frameTime / PICO;
1585     fr->bTime = (frameTime > 0);
1586
1587     // TODO This does not leak, but is not exception safe.
1588     /* values must be freed before leaving this function */
1589     sfree(values);
1590
1591     return bOK;
1592 #else
1593     GMX_UNUSED_VALUE(gmx_tng_input);
1594     GMX_UNUSED_VALUE(fr);
1595     GMX_UNUSED_VALUE(requestedIds);
1596     GMX_UNUSED_VALUE(numRequestedIds);
1597     return FALSE;
1598 #endif
1599 }
1600
1601 void gmx_print_tng_molecule_system(gmx_tng_trajectory_t gmx_tng_input,
1602                                    FILE                *stream)
1603 {
1604 #if GMX_USE_TNG
1605     gmx_int64_t         nMolecules, nChains, nResidues, nAtoms, nFramesRead;
1606     gmx_int64_t         strideLength, nParticlesRead, nValuesPerFrameRead, *molCntList;
1607     tng_molecule_t      molecule;
1608     tng_chain_t         chain;
1609     tng_residue_t       residue;
1610     tng_atom_t          atom;
1611     tng_function_status stat;
1612     char                str[256];
1613     char                varNAtoms;
1614     char                datatype;
1615     void               *data = nullptr;
1616     std::vector<real>   atomCharges;
1617     std::vector<real>   atomMasses;
1618     tng_trajectory_t    input = gmx_tng_input->tng;
1619
1620     tng_num_molecule_types_get(input, &nMolecules);
1621     tng_molecule_cnt_list_get(input, &molCntList);
1622     /* Can the number of particles change in the trajectory or is it constant? */
1623     tng_num_particles_variable_get(input, &varNAtoms);
1624
1625     for (gmx_int64_t i = 0; i < nMolecules; i++)
1626     {
1627         tng_molecule_of_index_get(input, i, &molecule);
1628         tng_molecule_name_get(input, molecule, str, 256);
1629         if (varNAtoms == TNG_CONSTANT_N_ATOMS)
1630         {
1631             if ((int)molCntList[i] == 0)
1632             {
1633                 continue;
1634             }
1635             fprintf(stream, "Molecule: %s, count: %d\n", str, (int)molCntList[i]);
1636         }
1637         else
1638         {
1639             fprintf(stream, "Molecule: %s\n", str);
1640         }
1641         tng_molecule_num_chains_get(input, molecule, &nChains);
1642         if (nChains > 0)
1643         {
1644             for (gmx_int64_t j = 0; j < nChains; j++)
1645             {
1646                 tng_molecule_chain_of_index_get(input, molecule, j, &chain);
1647                 tng_chain_name_get(input, chain, str, 256);
1648                 fprintf(stream, "\tChain: %s\n", str);
1649                 tng_chain_num_residues_get(input, chain, &nResidues);
1650                 for (gmx_int64_t k = 0; k < nResidues; k++)
1651                 {
1652                     tng_chain_residue_of_index_get(input, chain, k, &residue);
1653                     tng_residue_name_get(input, residue, str, 256);
1654                     fprintf(stream, "\t\tResidue: %s\n", str);
1655                     tng_residue_num_atoms_get(input, residue, &nAtoms);
1656                     for (gmx_int64_t l = 0; l < nAtoms; l++)
1657                     {
1658                         tng_residue_atom_of_index_get(input, residue, l, &atom);
1659                         tng_atom_name_get(input, atom, str, 256);
1660                         fprintf(stream, "\t\t\tAtom: %s", str);
1661                         tng_atom_type_get(input, atom, str, 256);
1662                         fprintf(stream, " (%s)\n", str);
1663                     }
1664                 }
1665             }
1666         }
1667         /* It is possible to have a molecule without chains, in which case
1668          * residues in the molecule can be iterated through without going
1669          * through chains. */
1670         else
1671         {
1672             tng_molecule_num_residues_get(input, molecule, &nResidues);
1673             if (nResidues > 0)
1674             {
1675                 for (gmx_int64_t k = 0; k < nResidues; k++)
1676                 {
1677                     tng_molecule_residue_of_index_get(input, molecule, k, &residue);
1678                     tng_residue_name_get(input, residue, str, 256);
1679                     fprintf(stream, "\t\tResidue: %s\n", str);
1680                     tng_residue_num_atoms_get(input, residue, &nAtoms);
1681                     for (gmx_int64_t l = 0; l < nAtoms; l++)
1682                     {
1683                         tng_residue_atom_of_index_get(input, residue, l, &atom);
1684                         tng_atom_name_get(input, atom, str, 256);
1685                         fprintf(stream, "\t\t\tAtom: %s", str);
1686                         tng_atom_type_get(input, atom, str, 256);
1687                         fprintf(stream, " (%s)\n", str);
1688                     }
1689                 }
1690             }
1691             else
1692             {
1693                 tng_molecule_num_atoms_get(input, molecule, &nAtoms);
1694                 for (gmx_int64_t l = 0; l < nAtoms; l++)
1695                 {
1696                     tng_molecule_atom_of_index_get(input, molecule, l, &atom);
1697                     tng_atom_name_get(input, atom, str, 256);
1698                     fprintf(stream, "\t\t\tAtom: %s", str);
1699                     tng_atom_type_get(input, atom, str, 256);
1700                     fprintf(stream, " (%s)\n", str);
1701                 }
1702             }
1703         }
1704     }
1705
1706     tng_num_particles_get(input, &nAtoms);
1707     stat = tng_particle_data_vector_get(input, TNG_TRAJ_PARTIAL_CHARGES, &data, &nFramesRead,
1708                                         &strideLength, &nParticlesRead,
1709                                         &nValuesPerFrameRead, &datatype);
1710     if (stat == TNG_SUCCESS)
1711     {
1712         atomCharges.resize(nAtoms);
1713         convert_array_to_real_array(data,
1714                                     atomCharges.data(),
1715                                     1,
1716                                     nAtoms,
1717                                     1,
1718                                     datatype);
1719
1720         fprintf(stream, "Atom Charges (%d):\n", int(nAtoms));
1721         for (gmx_int64_t i = 0; i < nAtoms; i += 10)
1722         {
1723             fprintf(stream, "Atom Charges [%8d-]=[", int(i));
1724             for (gmx_int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1725             {
1726                 fprintf(stream, " %12.5e", atomCharges[i + j]);
1727             }
1728             fprintf(stream, "]\n");
1729         }
1730     }
1731
1732     stat = tng_particle_data_vector_get(input, TNG_TRAJ_MASSES, &data, &nFramesRead,
1733                                         &strideLength, &nParticlesRead,
1734                                         &nValuesPerFrameRead, &datatype);
1735     if (stat == TNG_SUCCESS)
1736     {
1737         atomMasses.resize(nAtoms);
1738         convert_array_to_real_array(data,
1739                                     atomMasses.data(),
1740                                     1,
1741                                     nAtoms,
1742                                     1,
1743                                     datatype);
1744
1745         fprintf(stream, "Atom Masses (%d):\n", int(nAtoms));
1746         for (gmx_int64_t i = 0; i < nAtoms; i += 10)
1747         {
1748             fprintf(stream, "Atom Masses [%8d-]=[", int(i));
1749             for (gmx_int64_t j = 0; (j < 10 && i + j < nAtoms); j++)
1750             {
1751                 fprintf(stream, " %12.5e", atomMasses[i + j]);
1752             }
1753             fprintf(stream, "]\n");
1754         }
1755     }
1756
1757     sfree(data);
1758 #else
1759     GMX_UNUSED_VALUE(gmx_tng_input);
1760     GMX_UNUSED_VALUE(stream);
1761 #endif
1762 }
1763
1764 gmx_bool gmx_get_tng_data_block_types_of_next_frame(gmx_tng_trajectory_t gmx_tng_input,
1765                                                     int                  frame,
1766                                                     int                  nRequestedIds,
1767                                                     gmx_int64_t         *requestedIds,
1768                                                     gmx_int64_t         *nextFrame,
1769                                                     gmx_int64_t         *nBlocks,
1770                                                     gmx_int64_t        **blockIds)
1771 {
1772 #if GMX_USE_TNG
1773     tng_function_status stat;
1774     tng_trajectory_t    input = gmx_tng_input->tng;
1775
1776     stat = tng_util_trajectory_next_frame_present_data_blocks_find(input, frame,
1777                                                                    nRequestedIds, requestedIds,
1778                                                                    nextFrame,
1779                                                                    nBlocks, blockIds);
1780
1781     if (stat == TNG_CRITICAL)
1782     {
1783         gmx_file("Cannot read TNG file. Cannot find data blocks of next frame.");
1784     }
1785     else if (stat == TNG_FAILURE)
1786     {
1787         return FALSE;
1788     }
1789     return TRUE;
1790 #else
1791     GMX_UNUSED_VALUE(gmx_tng_input);
1792     GMX_UNUSED_VALUE(frame);
1793     GMX_UNUSED_VALUE(nRequestedIds);
1794     GMX_UNUSED_VALUE(requestedIds);
1795     GMX_UNUSED_VALUE(nextFrame);
1796     GMX_UNUSED_VALUE(nBlocks);
1797     GMX_UNUSED_VALUE(blockIds);
1798     return FALSE;
1799 #endif
1800 }
1801
1802 gmx_bool gmx_get_tng_data_next_frame_of_block_type(gmx_tng_trajectory_t gmx_tng_input,
1803                                                    gmx_int64_t          blockId,
1804                                                    real               **values,
1805                                                    gmx_int64_t         *frameNumber,
1806                                                    double              *frameTime,
1807                                                    gmx_int64_t         *nValuesPerFrame,
1808                                                    gmx_int64_t         *nAtoms,
1809                                                    real                *prec,
1810                                                    char                *name,
1811                                                    int                  maxLen,
1812                                                    gmx_bool            *bOK)
1813 {
1814 #if GMX_USE_TNG
1815     tng_function_status stat;
1816     char                datatype = -1;
1817     gmx_int64_t         codecId;
1818     int                 blockDependency;
1819     void               *data = nullptr;
1820     double              localPrec;
1821     tng_trajectory_t    input = gmx_tng_input->tng;
1822
1823     stat = tng_data_block_name_get(input, blockId, name, maxLen);
1824     if (stat != TNG_SUCCESS)
1825     {
1826         gmx_file("Cannot read next frame of TNG file");
1827     }
1828     stat = tng_data_block_dependency_get(input, blockId, &blockDependency);
1829     if (stat != TNG_SUCCESS)
1830     {
1831         gmx_file("Cannot read next frame of TNG file");
1832     }
1833     if (blockDependency & TNG_PARTICLE_DEPENDENT)
1834     {
1835         tng_num_particles_get(input, nAtoms);
1836         stat = tng_util_particle_data_next_frame_read(input,
1837                                                       blockId,
1838                                                       &data,
1839                                                       &datatype,
1840                                                       frameNumber,
1841                                                       frameTime);
1842     }
1843     else
1844     {
1845         *nAtoms = 1; /* There are not actually any atoms, but it is used for
1846                         allocating memory */
1847         stat    = tng_util_non_particle_data_next_frame_read(input,
1848                                                              blockId,
1849                                                              &data,
1850                                                              &datatype,
1851                                                              frameNumber,
1852                                                              frameTime);
1853     }
1854     if (stat == TNG_CRITICAL)
1855     {
1856         gmx_file("Cannot read next frame of TNG file");
1857     }
1858     if (stat == TNG_FAILURE)
1859     {
1860         *bOK = TRUE;
1861         return FALSE;
1862     }
1863
1864     stat = tng_data_block_num_values_per_frame_get(input, blockId, nValuesPerFrame);
1865     if (stat != TNG_SUCCESS)
1866     {
1867         gmx_file("Cannot read next frame of TNG file");
1868     }
1869     srenew(*values, sizeof(real) * *nValuesPerFrame * *nAtoms);
1870     convert_array_to_real_array(data,
1871                                 *values,
1872                                 getDistanceScaleFactor(gmx_tng_input),
1873                                 *nAtoms,
1874                                 *nValuesPerFrame,
1875                                 datatype);
1876
1877     tng_util_frame_current_compression_get(input, blockId, &codecId, &localPrec);
1878
1879     /* This must be updated if/when more lossy compression methods are added */
1880     if (codecId != TNG_TNG_COMPRESSION)
1881     {
1882         *prec = -1.0;
1883     }
1884     else
1885     {
1886         *prec = localPrec;
1887     }
1888
1889     sfree(data);
1890
1891     *bOK = TRUE;
1892     return TRUE;
1893 #else
1894     GMX_UNUSED_VALUE(gmx_tng_input);
1895     GMX_UNUSED_VALUE(blockId);
1896     GMX_UNUSED_VALUE(values);
1897     GMX_UNUSED_VALUE(frameNumber);
1898     GMX_UNUSED_VALUE(frameTime);
1899     GMX_UNUSED_VALUE(nValuesPerFrame);
1900     GMX_UNUSED_VALUE(nAtoms);
1901     GMX_UNUSED_VALUE(prec);
1902     GMX_UNUSED_VALUE(name);
1903     GMX_UNUSED_VALUE(maxLen);
1904     GMX_UNUSED_VALUE(bOK);
1905     return FALSE;
1906 #endif
1907 }