Write forces and velocities to compressed TNG.
[alexxy/gromacs.git] / src / gromacs / mdlib / mdoutf.cpp
index df1ef1f40989a10366784778814cfa82e52739eb..96b6bd6528f68066777edf7d8e792aa30f6cbbf0 100644 (file)
@@ -112,6 +112,28 @@ gmx_mdoutf_t init_mdoutf(FILE *fplog, int nfile, const t_filenm fnm[],
 
         sprintf(filemode, bAppendFiles ? "a+" : "w+");
 
+        if (EI_DYNAMICS(ir->eI) &&
+            ir->nstxout_compressed > 0)
+        {
+            const char *filename;
+            filename = ftp2fn(efCOMPRESSED, nfile, fnm);
+            switch (fn2ftp(filename))
+            {
+                case efXTC:
+                    of->fp_xtc                  = open_xtc(filename, filemode);
+                    break;
+                case efTNG:
+                    gmx_tng_open(filename, filemode[0], &of->tng_low_prec);
+                    if (filemode[0] == 'w')
+                    {
+                        gmx_tng_prepare_low_prec_writing(of->tng_low_prec, top_global, ir);
+                    }
+                    bCiteTng = TRUE;
+                    break;
+                default:
+                    gmx_incons("Invalid reduced precision file format");
+            }
+        }
         if ((EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
 #ifndef GMX_FAHCORE
             &&
@@ -128,7 +150,14 @@ gmx_mdoutf_t init_mdoutf(FILE *fplog, int nfile, const t_filenm fnm[],
             {
                 case efTRR:
                 case efTRN:
-                    of->fp_trn = gmx_trr_open(filename, filemode);
+                    /* If there is no uncompressed coordinate output and
+                       there is compressed TNG output write forces
+                       and/or velocities to the TNG file instead. */
+                    if (ir->nstxout != 0 || ir->nstxout_compressed == 0 ||
+                        !of->tng_low_prec)
+                    {
+                        of->fp_trn = gmx_trr_open(filename, filemode);
+                    }
                     break;
                 case efTNG:
                     gmx_tng_open(filename, filemode[0], &of->tng);
@@ -142,28 +171,6 @@ gmx_mdoutf_t init_mdoutf(FILE *fplog, int nfile, const t_filenm fnm[],
                     gmx_incons("Invalid full precision file format");
             }
         }
-        if (EI_DYNAMICS(ir->eI) &&
-            ir->nstxout_compressed > 0)
-        {
-            const char *filename;
-            filename = ftp2fn(efCOMPRESSED, nfile, fnm);
-            switch (fn2ftp(filename))
-            {
-                case efXTC:
-                    of->fp_xtc                  = open_xtc(filename, filemode);
-                    break;
-                case efTNG:
-                    gmx_tng_open(filename, filemode[0], &of->tng_low_prec);
-                    if (filemode[0] == 'w')
-                    {
-                        gmx_tng_prepare_low_prec_writing(of->tng_low_prec, top_global, ir);
-                    }
-                    bCiteTng = TRUE;
-                    break;
-                default:
-                    gmx_incons("Invalid reduced precision file format");
-            }
-        }
         if (EI_DYNAMICS(ir->eI) || EI_ENERGY_MINIMIZATION(ir->eI))
         {
             of->fp_ene = open_enx(ftp2fn(efEDR, nfile, fnm), filemode);
@@ -333,12 +340,28 @@ void mdoutf_write_to_trajectory_files(FILE *fplog, t_commrec *cr,
                 }
             }
 
-            gmx_fwrite_tng(of->tng, FALSE, step, t, state_local->lambda[efptFEP],
-                           state_local->box,
-                           top_global->natoms,
-                           (mdof_flags & MDOF_X) ? state_global->x : NULL,
-                           (mdof_flags & MDOF_V) ? global_v : NULL,
-                           (mdof_flags & MDOF_F) ? f_global : NULL);
+            /* If a TNG file is open for uncompressed coordinate output also write
+               velocities and forces to it. */
+            else if (of->tng)
+            {
+                gmx_fwrite_tng(of->tng, FALSE, step, t, state_local->lambda[efptFEP],
+                               state_local->box,
+                               top_global->natoms,
+                               (mdof_flags & MDOF_X) ? state_global->x : NULL,
+                               (mdof_flags & MDOF_V) ? global_v : NULL,
+                               (mdof_flags & MDOF_F) ? f_global : NULL);
+            }
+            /* If only a TNG file is open for compressed coordinate output (no uncompressed
+               coordinate output) also write forces and velocities to it. */
+            else if (of->tng_low_prec)
+            {
+                gmx_fwrite_tng(of->tng_low_prec, FALSE, step, t, state_local->lambda[efptFEP],
+                               state_local->box,
+                               top_global->natoms,
+                               (mdof_flags & MDOF_X) ? state_global->x : NULL,
+                               (mdof_flags & MDOF_V) ? global_v : NULL,
+                               (mdof_flags & MDOF_F) ? f_global : NULL);
+            }
         }
         if (mdof_flags & MDOF_X_COMPRESSED)
         {