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