Merge release-5-0 into release-5-1
[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 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 #ifdef GMX_USE_TNG
42 #include "tng/tng_io.h"
43 #endif
44
45 #include "gromacs/fileio/gmxfio.h"
46 #include "gromacs/legacyheaders/copyrite.h"
47 #include "gromacs/legacyheaders/types/ifunc.h"
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/basedefinitions.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/programcontext.h"
55 #include "gromacs/utility/sysinfo.h"
56
57 static const char *modeToVerb(char mode)
58 {
59     const char *p;
60     switch (mode)
61     {
62         case 'r':
63             p = "reading";
64             break;
65         case 'w':
66             p = "writing";
67             break;
68         case 'a':
69             p = "appending";
70             break;
71         default:
72             gmx_fatal(FARGS, "Invalid file opening mode %c", mode);
73             p = "";
74             break;
75     }
76     return p;
77 }
78
79 void gmx_tng_open(const char       *filename,
80                   char              mode,
81                   tng_trajectory_t *tng)
82 {
83 #ifdef GMX_USE_TNG
84     /* First check whether we have to make a backup,
85      * only for writing, not for read or append.
86      */
87     if (mode == 'w')
88     {
89         make_backup(filename);
90     }
91
92     /* tng must not be pointing at already allocated memory.
93      * Memory will be allocated by tng_util_trajectory_open() and must
94      * later on be freed by tng_util_trajectory_close(). */
95     if (TNG_SUCCESS != tng_util_trajectory_open(filename, mode, tng))
96     {
97         /* TNG does return more than one degree of error, but there is
98            no use case for GROMACS handling the non-fatal errors
99            gracefully. */
100         gmx_fatal(FARGS,
101                   "%s while opening %s for %s",
102                   gmx_strerror("file"),
103                   filename,
104                   modeToVerb(mode));
105     }
106
107     if (mode == 'w' || mode == 'a')
108     {
109         char hostname[256];
110         gmx_gethostname(hostname, 256);
111         if (mode == 'w')
112         {
113             tng_first_computer_name_set(*tng, hostname);
114         }
115         else
116         {
117             tng_last_computer_name_set(*tng, hostname);
118         }
119
120         char        programInfo[256];
121         const char *precisionString = "";
122 #ifdef GMX_DOUBLE
123         precisionString = " (double precision)";
124 #endif
125         sprintf(programInfo, "%.100s, %.128s%.24s",
126                 gmx::getProgramContext().displayName(),
127                 GromacsVersion(), precisionString);
128         if (mode == 'w')
129         {
130             tng_first_program_name_set(*tng, programInfo);
131         }
132         else
133         {
134             tng_last_program_name_set(*tng, programInfo);
135         }
136
137         char username[256];
138         if (!gmx_getusername(username, 256))
139         {
140             if (mode == 'w')
141             {
142                 tng_first_user_name_set(*tng, username);
143             }
144             else
145             {
146                 tng_last_user_name_set(*tng, username);
147                 tng_file_headers_write(*tng, TNG_USE_HASH);
148             }
149         }
150     }
151 #else
152     gmx_file("GROMACS was compiled without TNG support, cannot handle this file type");
153     GMX_UNUSED_VALUE(filename);
154     GMX_UNUSED_VALUE(mode);
155     GMX_UNUSED_VALUE(tng);
156 #endif
157 }
158
159 void gmx_tng_close(tng_trajectory_t *tng)
160 {
161     /* We have to check that tng is set because
162      * tng_util_trajectory_close wants to return a NULL in it, and
163      * gives a fatal error if it is NULL. */
164 #ifdef GMX_USE_TNG
165     if (tng)
166     {
167         tng_util_trajectory_close(tng);
168     }
169 #else
170     GMX_UNUSED_VALUE(tng);
171 #endif
172 }
173
174 #ifdef GMX_USE_TNG
175 static void addTngMoleculeFromTopology(tng_trajectory_t     tng,
176                                        const char          *moleculeName,
177                                        const t_atoms       *atoms,
178                                        gmx_int64_t          numMolecules,
179                                        tng_molecule_t      *tngMol)
180 {
181     tng_chain_t      tngChain = NULL;
182     tng_residue_t    tngRes   = NULL;
183
184     if (tng_molecule_add(tng, moleculeName, tngMol) != TNG_SUCCESS)
185     {
186         gmx_file("Cannot add molecule to TNG molecular system.");
187     }
188
189     /* FIXME: The TNG atoms should contain mass and atomB info (for free
190      * energy calculations), i.e. in when it's available in TNG (2.0). */
191     for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++)
192     {
193         const t_atom *at = &atoms->atom[atomIndex];
194         /* FIXME: Currently the TNG API can only add atoms belonging to a
195          * residue and chain. Wait for TNG 2.0*/
196         if (atoms->nres > 0)
197         {
198             const t_resinfo *resInfo        = &atoms->resinfo[at->resind];
199             char             chainName[2]   = {resInfo->chainid, 0};
200             tng_atom_t       tngAtom        = NULL;
201             t_atom          *prevAtom;
202
203             if (atomIndex > 0)
204             {
205                 prevAtom = &atoms->atom[atomIndex - 1];
206             }
207             else
208             {
209                 prevAtom = 0;
210             }
211
212             /* If this is the first atom or if the residue changed add the
213              * residue to the TNG molecular system. */
214             if (!prevAtom || resInfo != &atoms->resinfo[prevAtom->resind])
215             {
216                 /* If this is the first atom or if the chain changed add
217                  * the chain to the TNG molecular system. */
218                 if (!prevAtom || resInfo->chainid !=
219                     atoms->resinfo[prevAtom->resind].chainid)
220                 {
221                     tng_molecule_chain_add(tng, *tngMol, chainName,
222                                            &tngChain);
223                 }
224                 /* FIXME: When TNG supports both residue index and residue
225                  * number the latter should be used. Wait for TNG 2.0*/
226                 tng_chain_residue_add(tng, tngChain, *resInfo->name, &tngRes);
227             }
228             tng_residue_atom_add(tng, tngRes, *(atoms->atomname[atomIndex]), *(atoms->atomtype[atomIndex]), &tngAtom);
229         }
230     }
231     tng_molecule_cnt_set(tng, *tngMol, numMolecules);
232 }
233
234 void gmx_tng_add_mtop(tng_trajectory_t  tng,
235                       const gmx_mtop_t *mtop)
236 {
237     int                  i, j;
238     const t_ilist       *ilist;
239     tng_bond_t           tngBond;
240
241     if (!mtop)
242     {
243         /* No topology information available to add. */
244         return;
245     }
246
247     for (int molIndex = 0; molIndex < mtop->nmolblock; molIndex++)
248     {
249         tng_molecule_t       tngMol  = NULL;
250         const gmx_moltype_t *molType =
251             &mtop->moltype[mtop->molblock[molIndex].type];
252
253         /* Add a molecule to the TNG trajectory with the same name as the
254          * current molecule. */
255         addTngMoleculeFromTopology(tng,
256                                    *(molType->name),
257                                    &molType->atoms,
258                                    mtop->molblock[molIndex].nmol,
259                                    &tngMol);
260
261         /* Bonds have to be deduced from interactions (constraints etc). Different
262          * interactions have different sets of parameters. */
263         /* Constraints are specified using two atoms */
264         for (i = 0; i < F_NRE; i++)
265         {
266             if (IS_CHEMBOND(i))
267             {
268                 ilist = &molType->ilist[i];
269                 if (ilist)
270                 {
271                     j = 1;
272                     while (j < ilist->nr)
273                     {
274                         tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
275                         j += 3;
276                     }
277                 }
278             }
279         }
280         /* Settle is described using three atoms */
281         ilist = &molType->ilist[F_SETTLE];
282         if (ilist)
283         {
284             j = 1;
285             while (j < ilist->nr)
286             {
287                 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+1], &tngBond);
288                 tng_molecule_bond_add(tng, tngMol, ilist->iatoms[j], ilist->iatoms[j+2], &tngBond);
289                 j += 4;
290             }
291         }
292     }
293 }
294
295 /*! \libinternal \brief Compute greatest common divisor of n1 and n2
296  * if they are positive.
297  *
298  * If only one of n1 and n2 is positive, then return it.
299  * If neither n1 or n2 is positive, then return -1. */
300 static int
301 greatest_common_divisor_if_positive(int n1, int n2)
302 {
303     if (0 >= n1)
304     {
305         return (0 >= n2) ? -1 : n2;
306     }
307     if (0 >= n2)
308     {
309         return n1;
310     }
311
312     /* We have a non-trivial greatest common divisor to compute. */
313     return gmx_greatest_common_divisor(n1, n2);
314 }
315
316 /* By default try to write 100 frames (of actual output) in each frame set.
317  * This number is the number of outputs of the most frequently written data
318  * type per frame set.
319  * TODO for 5.1: Verify that 100 frames per frame set is efficient for most
320  * setups regarding compression efficiency and compression time. Make this
321  * a hidden command-line option? */
322 const int defaultFramesPerFrameSet = 100;
323
324 /*! \libinternal \brief  Set the number of frames per frame
325  * set according to output intervals.
326  * The default is that 100 frames are written of the data
327  * that is written most often. */
328 static void tng_set_frames_per_frame_set(tng_trajectory_t  tng,
329                                          const gmx_bool    bUseLossyCompression,
330                                          const t_inputrec *ir)
331 {
332     int     gcd = -1;
333
334     /* Set the number of frames per frame set to contain at least
335      * defaultFramesPerFrameSet of the lowest common denominator of
336      * the writing interval of positions and velocities. */
337     /* FIXME after 5.0: consider nstenergy also? */
338     if (bUseLossyCompression)
339     {
340         gcd = ir->nstxout_compressed;
341     }
342     else
343     {
344         gcd = greatest_common_divisor_if_positive(ir->nstxout, ir->nstvout);
345         gcd = greatest_common_divisor_if_positive(gcd, ir->nstfout);
346     }
347     if (0 >= gcd)
348     {
349         return;
350     }
351
352     tng_num_frames_per_frame_set_set(tng, gcd * defaultFramesPerFrameSet);
353 }
354
355 /*! \libinternal \brief Set the data-writing intervals, and number of
356  * frames per frame set */
357 static void set_writing_intervals(tng_trajectory_t  tng,
358                                   const gmx_bool    bUseLossyCompression,
359                                   const t_inputrec *ir)
360 {
361     /* Define pointers to specific writing functions depending on if we
362      * write float or double data */
363     typedef tng_function_status (*set_writing_interval_func_pointer)(tng_trajectory_t,
364                                                                      const gmx_int64_t,
365                                                                      const gmx_int64_t,
366                                                                      const gmx_int64_t,
367                                                                      const char*,
368                                                                      const char,
369                                                                      const char);
370 #ifdef GMX_DOUBLE
371     set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_double_set;
372 #else
373     set_writing_interval_func_pointer set_writing_interval = tng_util_generic_write_interval_set;
374 #endif
375     int  xout, vout, fout;
376 //     int  gcd = -1, lowest = -1;
377     char compression;
378
379     tng_set_frames_per_frame_set(tng, bUseLossyCompression, ir);
380
381     if (bUseLossyCompression)
382     {
383         xout        = ir->nstxout_compressed;
384         vout        = 0;
385         fout        = 0;
386         compression = TNG_TNG_COMPRESSION;
387     }
388     else
389     {
390         xout        = ir->nstxout;
391         vout        = ir->nstvout;
392         fout        = ir->nstfout;
393         compression = TNG_GZIP_COMPRESSION;
394     }
395     if (xout)
396     {
397         set_writing_interval(tng, xout, 3, TNG_TRAJ_POSITIONS,
398                              "POSITIONS", TNG_PARTICLE_BLOCK_DATA,
399                              compression);
400         /* The design of TNG makes it awkward to try to write a box
401          * with multiple periodicities, which might be co-prime. Since
402          * the use cases for the box with a frame consisting only of
403          * velocities seem low, for now we associate box writing with
404          * position writing. */
405         set_writing_interval(tng, xout, 9, TNG_TRAJ_BOX_SHAPE,
406                              "BOX SHAPE", TNG_NON_PARTICLE_BLOCK_DATA,
407                              TNG_GZIP_COMPRESSION);
408         /* TODO: if/when we write energies to TNG also, reconsider how
409          * and when box information is written, because GROMACS
410          * behaviour pre-5.0 was to write the box with every
411          * trajectory frame and every energy frame, and probably
412          * people depend on this. */
413
414         /* TODO: If we need to write lambda values at steps when
415          * positions (or other data) are not also being written, then
416          * code in mdoutf.c will need to match however that is
417          * organized here. */
418         set_writing_interval(tng, xout, 1, TNG_GMX_LAMBDA,
419                              "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
420                              TNG_GZIP_COMPRESSION);
421
422         /* FIXME: gcd and lowest currently not used. */
423 //         gcd = greatest_common_divisor_if_positive(gcd, xout);
424 //         if (lowest < 0 || xout < lowest)
425 //         {
426 //             lowest = xout;
427 //         }
428     }
429     if (vout)
430     {
431         set_writing_interval(tng, ir->nstvout, 3, TNG_TRAJ_VELOCITIES,
432                              "VELOCITIES", TNG_PARTICLE_BLOCK_DATA,
433                              compression);
434
435         /* FIXME: gcd and lowest currently not used. */
436 //         gcd = greatest_common_divisor_if_positive(gcd, vout);
437 //         if (lowest < 0 || vout < lowest)
438 //         {
439 //             lowest = vout;
440 //         }
441     }
442     if (fout)
443     {
444         set_writing_interval(tng, ir->nstfout, 3, TNG_TRAJ_FORCES,
445                              "FORCES", TNG_PARTICLE_BLOCK_DATA,
446                              TNG_GZIP_COMPRESSION);
447
448         /* FIXME: gcd and lowest currently not used. */
449 //         gcd = greatest_common_divisor_if_positive(gcd, fout);
450 //         if (lowest < 0 || fout < lowest)
451 //         {
452 //             lowest = fout;
453 //         }
454     }
455     /* FIXME: See above. gcd interval for lambdas is disabled. */
456 //     if (gcd > 0)
457 //     {
458 //         /* Lambdas written at an interval of the lowest common denominator
459 //          * of other output */
460 //         set_writing_interval(tng, gcd, 1, TNG_GMX_LAMBDA,
461 //                                  "LAMBDAS", TNG_NON_PARTICLE_BLOCK_DATA,
462 //                                  TNG_GZIP_COMPRESSION);
463 //
464 //         if (gcd < lowest / 10)
465 //         {
466 //             gmx_warning("The lowest common denominator of trajectory output is "
467 //                         "every %d step(s), whereas the shortest output interval "
468 //                         "is every %d steps.", gcd, lowest);
469 //         }
470 //     }
471 }
472 #endif
473
474 void gmx_tng_prepare_md_writing(tng_trajectory_t  tng,
475                                 const gmx_mtop_t *mtop,
476                                 const t_inputrec *ir)
477 {
478 #ifdef GMX_USE_TNG
479     gmx_tng_add_mtop(tng, mtop);
480     set_writing_intervals(tng, FALSE, ir);
481     tng_time_per_frame_set(tng, ir->delta_t * PICO);
482 #else
483     GMX_UNUSED_VALUE(tng);
484     GMX_UNUSED_VALUE(mtop);
485     GMX_UNUSED_VALUE(ir);
486 #endif
487 }
488
489 #ifdef GMX_USE_TNG
490 /* Check if all atoms in the molecule system are selected
491  * by a selection group of type specified by gtype. */
492 static gmx_bool all_atoms_selected(const gmx_mtop_t *mtop,
493                                    const int         gtype)
494 {
495     const gmx_moltype_t     *molType;
496     const t_atoms           *atoms;
497
498     /* Iterate through all atoms in the molecule system and
499      * check if they belong to a selection group of the
500      * requested type. */
501     for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
502     {
503         molType = &mtop->moltype[mtop->molblock[molIndex].type];
504
505         atoms = &molType->atoms;
506
507         for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
508         {
509             for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
510             {
511                 if (ggrpnr(&mtop->groups, gtype, i) != 0)
512                 {
513                     return FALSE;
514                 }
515             }
516         }
517     }
518     return TRUE;
519 }
520
521 /* Create TNG molecules which will represent each of the selection
522  * groups for writing. But do that only if there is actually a
523  * specified selection group and it is not the whole system.
524  * TODO: Currently the only selection that is taken into account
525  * is egcCompressedX, but other selections should be added when
526  * e.g. writing energies is implemented.
527  */
528 static void add_selection_groups(tng_trajectory_t  tng,
529                                  const gmx_mtop_t *mtop)
530 {
531     const gmx_moltype_t     *molType;
532     const t_atoms           *atoms;
533     const t_atom            *at;
534     const t_resinfo         *resInfo;
535     const t_ilist           *ilist;
536     int                      nameIndex;
537     int                      atom_offset = 0;
538     tng_molecule_t           mol, iterMol;
539     tng_chain_t              chain;
540     tng_residue_t            res;
541     tng_atom_t               atom;
542     tng_bond_t               tngBond;
543     gmx_int64_t              nMols;
544     char                    *groupName;
545
546     /* TODO: When the TNG molecules block is more flexible TNG selection
547      * groups should not need all atoms specified. It should be possible
548      * just to specify what molecules are selected (and/or which atoms in
549      * the molecule) and how many (if applicable). */
550
551     /* If no atoms are selected we do not need to create a
552      * TNG selection group molecule. */
553     if (mtop->groups.ngrpnr[egcCompressedX] == 0)
554     {
555         return;
556     }
557
558     /* If all atoms are selected we do not have to create a selection
559      * group molecule in the TNG molecule system. */
560     if (all_atoms_selected(mtop, egcCompressedX))
561     {
562         return;
563     }
564
565     /* The name of the TNG molecule containing the selection group is the
566      * same as the name of the selection group. */
567     nameIndex = *mtop->groups.grps[egcCompressedX].nm_ind;
568     groupName = *mtop->groups.grpname[nameIndex];
569
570     tng_molecule_alloc(tng, &mol);
571     tng_molecule_name_set(tng, mol, groupName);
572     tng_molecule_chain_add(tng, mol, "", &chain);
573     for (int molIndex = 0, i = 0; molIndex < mtop->nmoltype; molIndex++)
574     {
575         molType = &mtop->moltype[mtop->molblock[molIndex].type];
576
577         atoms = &molType->atoms;
578
579         for (int j = 0; j < mtop->molblock[molIndex].nmol; j++)
580         {
581             bool bAtomsAdded = FALSE;
582             for (int atomIndex = 0; atomIndex < atoms->nr; atomIndex++, i++)
583             {
584                 char *res_name;
585                 int   res_id;
586
587                 if (ggrpnr(&mtop->groups, egcCompressedX, i) != 0)
588                 {
589                     continue;
590                 }
591                 at = &atoms->atom[atomIndex];
592                 if (atoms->nres > 0)
593                 {
594                     resInfo = &atoms->resinfo[at->resind];
595                     /* FIXME: When TNG supports both residue index and residue
596                      * number the latter should be used. */
597                     res_name = *resInfo->name;
598                     res_id   = at->resind + 1;
599                 }
600                 else
601                 {
602                     res_name = (char *)"";
603                     res_id   = 0;
604                 }
605                 if (tng_chain_residue_find(tng, chain, res_name, res_id, &res)
606                     != TNG_SUCCESS)
607                 {
608                     /* Since there is ONE chain for selection groups do not keep the
609                      * original residue IDs - otherwise there might be conflicts. */
610                     tng_chain_residue_add(tng, chain, res_name, &res);
611                 }
612                 tng_residue_atom_w_id_add(tng, res, *(atoms->atomname[atomIndex]),
613                                           *(atoms->atomtype[atomIndex]),
614                                           atom_offset + atomIndex, &atom);
615                 bAtomsAdded = TRUE;
616             }
617             /* Add bonds. */
618             if (bAtomsAdded)
619             {
620                 for (int k = 0; k < F_NRE; k++)
621                 {
622                     if (IS_CHEMBOND(k))
623                     {
624                         ilist = &molType->ilist[k];
625                         if (ilist)
626                         {
627                             int l = 1;
628                             while (l < ilist->nr)
629                             {
630                                 int atom1, atom2;
631                                 atom1 = ilist->iatoms[l] + atom_offset;
632                                 atom2 = ilist->iatoms[l+1] + atom_offset;
633                                 if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0 &&
634                                     ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
635                                 {
636                                     tng_molecule_bond_add(tng, mol, ilist->iatoms[l],
637                                                           ilist->iatoms[l+1], &tngBond);
638                                 }
639                                 l += 3;
640                             }
641                         }
642                     }
643                 }
644                 /* Settle is described using three atoms */
645                 ilist = &molType->ilist[F_SETTLE];
646                 if (ilist)
647                 {
648                     int l = 1;
649                     while (l < ilist->nr)
650                     {
651                         int atom1, atom2, atom3;
652                         atom1 = ilist->iatoms[l] + atom_offset;
653                         atom2 = ilist->iatoms[l+1] + atom_offset;
654                         atom3 = ilist->iatoms[l+2] + atom_offset;
655                         if (ggrpnr(&mtop->groups, egcCompressedX, atom1) == 0)
656                         {
657                             if (ggrpnr(&mtop->groups, egcCompressedX, atom2) == 0)
658                             {
659                                 tng_molecule_bond_add(tng, mol, atom1,
660                                                       atom2, &tngBond);
661                             }
662                             if (ggrpnr(&mtop->groups, egcCompressedX, atom3) == 0)
663                             {
664                                 tng_molecule_bond_add(tng, mol, atom1,
665                                                       atom3, &tngBond);
666                             }
667                         }
668                         l += 4;
669                     }
670                 }
671             }
672             atom_offset += atoms->nr;
673         }
674     }
675     tng_molecule_existing_add(tng, &mol);
676     tng_molecule_cnt_set(tng, mol, 1);
677     tng_num_molecule_types_get(tng, &nMols);
678     for (gmx_int64_t k = 0; k < nMols; k++)
679     {
680         tng_molecule_of_index_get(tng, k, &iterMol);
681         if (iterMol == mol)
682         {
683             continue;
684         }
685         tng_molecule_cnt_set(tng, iterMol, 0);
686     }
687 }
688 #endif
689
690 void gmx_tng_set_compression_precision(tng_trajectory_t tng,
691                                        real             prec)
692 {
693 #ifdef GMX_USE_TNG
694     tng_compression_precision_set(tng, prec);
695 #else
696     GMX_UNUSED_VALUE(tng);
697     GMX_UNUSED_VALUE(prec);
698 #endif
699 }
700
701 void gmx_tng_prepare_low_prec_writing(tng_trajectory_t  tng,
702                                       const gmx_mtop_t *mtop,
703                                       const t_inputrec *ir)
704 {
705 #ifdef GMX_USE_TNG
706     gmx_tng_add_mtop(tng, mtop);
707     add_selection_groups(tng, mtop);
708     set_writing_intervals(tng, TRUE, ir);
709     tng_time_per_frame_set(tng, ir->delta_t * PICO);
710     gmx_tng_set_compression_precision(tng, ir->x_compression_precision);
711 #else
712     GMX_UNUSED_VALUE(tng);
713     GMX_UNUSED_VALUE(mtop);
714     GMX_UNUSED_VALUE(ir);
715 #endif
716 }
717
718 void gmx_fwrite_tng(tng_trajectory_t tng,
719                     const gmx_bool   bUseLossyCompression,
720                     int              step,
721                     real             elapsedPicoSeconds,
722                     real             lambda,
723                     const rvec      *box,
724                     int              nAtoms,
725                     const rvec      *x,
726                     const rvec      *v,
727                     const rvec      *f)
728 {
729 #ifdef GMX_USE_TNG
730     typedef tng_function_status (*write_data_func_pointer)(tng_trajectory_t,
731                                                            const gmx_int64_t,
732                                                            const double,
733                                                            const real*,
734                                                            const gmx_int64_t,
735                                                            const gmx_int64_t,
736                                                            const char*,
737                                                            const char,
738                                                            const char);
739 #ifdef GMX_DOUBLE
740     static write_data_func_pointer           write_data           = tng_util_generic_with_time_double_write;
741 #else
742     static write_data_func_pointer           write_data           = tng_util_generic_with_time_write;
743 #endif
744     double                                   elapsedSeconds = elapsedPicoSeconds * PICO;
745     gmx_int64_t                              nParticles;
746     char                                     compression;
747
748
749     if (!tng)
750     {
751         /* This function might get called when the type of the
752            compressed trajectory is actually XTC. So we exit and move
753            on. */
754         return;
755     }
756
757     tng_num_particles_get(tng, &nParticles);
758     if (nAtoms != (int)nParticles)
759     {
760         tng_implicit_num_particles_set(tng, nAtoms);
761     }
762
763     if (bUseLossyCompression)
764     {
765         compression = TNG_TNG_COMPRESSION;
766     }
767     else
768     {
769         compression = TNG_GZIP_COMPRESSION;
770     }
771
772     /* The writing is done using write_data, which writes float or double
773      * depending on the GROMACS compilation. */
774     if (x)
775     {
776         GMX_ASSERT(box, "Need a non-NULL box if positions are written");
777
778         if (write_data(tng, step, elapsedSeconds,
779                        reinterpret_cast<const real *>(x),
780                        3, TNG_TRAJ_POSITIONS, "POSITIONS",
781                        TNG_PARTICLE_BLOCK_DATA,
782                        compression) != TNG_SUCCESS)
783         {
784             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
785         }
786         /* TNG-MF1 compression only compresses positions and velocities. Use lossless
787          * compression for box shape regardless of output mode */
788         if (write_data(tng, step, elapsedSeconds,
789                        reinterpret_cast<const real *>(box),
790                        9, TNG_TRAJ_BOX_SHAPE, "BOX SHAPE",
791                        TNG_NON_PARTICLE_BLOCK_DATA,
792                        TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
793         {
794             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
795         }
796     }
797
798     if (v)
799     {
800         if (write_data(tng, step, elapsedSeconds,
801                        reinterpret_cast<const real *>(v),
802                        3, TNG_TRAJ_VELOCITIES, "VELOCITIES",
803                        TNG_PARTICLE_BLOCK_DATA,
804                        compression) != TNG_SUCCESS)
805         {
806             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
807         }
808     }
809
810     if (f)
811     {
812         /* TNG-MF1 compression only compresses positions and velocities. Use lossless
813          * compression for forces regardless of output mode */
814         if (write_data(tng, step, elapsedSeconds,
815                        reinterpret_cast<const real *>(f),
816                        3, TNG_TRAJ_FORCES, "FORCES",
817                        TNG_PARTICLE_BLOCK_DATA,
818                        TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
819         {
820             gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
821         }
822     }
823
824     /* TNG-MF1 compression only compresses positions and velocities. Use lossless
825      * compression for lambdas regardless of output mode */
826     if (write_data(tng, step, elapsedSeconds,
827                    reinterpret_cast<const real *>(&lambda),
828                    1, TNG_GMX_LAMBDA, "LAMBDAS",
829                    TNG_NON_PARTICLE_BLOCK_DATA,
830                    TNG_GZIP_COMPRESSION) != TNG_SUCCESS)
831     {
832         gmx_file("Cannot write TNG trajectory frame; maybe you are out of disk space?");
833     }
834 #else
835     GMX_UNUSED_VALUE(tng);
836     GMX_UNUSED_VALUE(bUseLossyCompression);
837     GMX_UNUSED_VALUE(step);
838     GMX_UNUSED_VALUE(elapsedPicoSeconds);
839     GMX_UNUSED_VALUE(lambda);
840     GMX_UNUSED_VALUE(box);
841     GMX_UNUSED_VALUE(nAtoms);
842     GMX_UNUSED_VALUE(x);
843     GMX_UNUSED_VALUE(v);
844     GMX_UNUSED_VALUE(f);
845 #endif
846 }
847
848 void fflush_tng(tng_trajectory_t tng)
849 {
850 #ifdef GMX_USE_TNG
851     if (!tng)
852     {
853         return;
854     }
855     tng_frame_set_premature_write(tng, TNG_USE_HASH);
856 #else
857     GMX_UNUSED_VALUE(tng);
858 #endif
859 }
860
861 float gmx_tng_get_time_of_final_frame(tng_trajectory_t tng)
862 {
863 #ifdef GMX_USE_TNG
864     gmx_int64_t nFrames;
865     double      time;
866     float       fTime;
867
868     tng_num_frames_get(tng, &nFrames);
869     tng_util_time_of_frame_get(tng, nFrames - 1, &time);
870
871     fTime = time / PICO;
872     return fTime;
873 #else
874     GMX_UNUSED_VALUE(tng);
875     return -1.0;
876 #endif
877 }