Merge branch release-2018
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_trjconv.cpp
index 4556135d889fe30fbd1b1804fb239a443bffab12..37c0eb116217058e3283f36ccb6b7a3fd6aa43d7 100644 (file)
@@ -44,6 +44,7 @@
 
 #include "gromacs/commandline/pargs.h"
 #include "gromacs/commandline/viewit.h"
+#include "gromacs/compat/make_unique.h"
 #include "gromacs/fileio/confio.h"
 #include "gromacs/fileio/g96io.h"
 #include "gromacs/fileio/gmxfio.h"
@@ -77,7 +78,7 @@ enum {
 
 
 static void calc_pbc_cluster(int ecenter, int nrefat, t_topology *top, int ePBC,
-                             rvec x[], int index[], matrix box)
+                             rvec x[], const int index[], matrix box)
 {
     int       m, i, j, j0, j1, jj, ai, aj;
     int       imin, jmin;
@@ -373,7 +374,7 @@ static void put_residue_com_in_box(int unitcell_enum, int ecenter,
                     break;
             }
             rvec_sub(newCom, com, shift);
-            if (norm2(shift))
+            if (norm2(shift) != 0.0f)
             {
                 if (debug)
                 {
@@ -407,7 +408,7 @@ static void put_residue_com_in_box(int unitcell_enum, int ecenter,
     }
 }
 
-static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, int ci[])
+static void center_x(int ecenter, rvec x[], matrix box, int n, int nc, const int ci[])
 {
     int  i, m, ai;
     rvec cmin, cmax, box_center, dx;
@@ -522,7 +523,7 @@ static void do_trunc(const char *fn, real t0)
         {
             fprintf(stderr, "Do you REALLY want to truncate this trajectory (%s) at:\n"
                     "frame %d, time %g, bytes %ld ??? (type YES if so)\n",
-                    fn, j, t, (long int)fpos);
+                    fn, j, t, static_cast<long int>(fpos));
             if (1 != scanf("%s", yesno))
             {
                 gmx_fatal(FARGS, "Error reading user input");
@@ -569,20 +570,21 @@ static void do_trunc(const char *fn, real t0)
  * molecule information will generally be present if the input TNG
  * file was written by a GROMACS tool, this seems like reasonable
  * behaviour. */
-static gmx_mtop_t *read_mtop_for_tng(const char *tps_file,
-                                     const char *input_file,
-                                     const char *output_file)
+static std::unique_ptr<gmx_mtop_t>
+read_mtop_for_tng(const char *tps_file,
+                  const char *input_file,
+                  const char *output_file)
 {
-    gmx_mtop_t *mtop = nullptr;
+    std::unique_ptr<gmx_mtop_t> mtop;
 
     if (fn2bTPX(tps_file) &&
         efTNG != fn2ftp(input_file) &&
         efTNG == fn2ftp(output_file))
     {
         int temp_natoms = -1;
-        snew(mtop, 1);
+        mtop = gmx::compat::make_unique<gmx_mtop_t>();
         read_tpx(tps_file, nullptr, nullptr, &temp_natoms,
-                 nullptr, nullptr, mtop);
+                 nullptr, nullptr, mtop.get());
     }
 
     return mtop;
@@ -871,7 +873,6 @@ int gmx_trjconv(int argc, char *argv[])
     int               m, i, d, frame, outframe, natoms, nout, ncent, newstep = 0, model_nr;
 #define SKIP 10
     t_topology       *top   = nullptr;
-    gmx_mtop_t       *mtop  = nullptr;
     gmx_conect        gc    = nullptr;
     int               ePBC  = -1;
     t_atoms          *atoms = nullptr, useatoms;
@@ -1053,7 +1054,8 @@ int gmx_trjconv(int argc, char *argv[])
             /* Check for number of files disabled, as FOPEN_MAX is not the correct
              * number to check for. In my linux box it is only 16.
              */
-            if (0 && (clust->clust->nr > FOPEN_MAX-4))
+            if (/* DISABLES CODE */ (false))
+            //if (clust->clust->nr > FOPEN_MAX-4)
             {
                 gmx_fatal(FARGS, "Can not open enough (%d) files to write all the"
                           " trajectories.\ntry splitting the index file in %d parts.\n"
@@ -1079,9 +1081,10 @@ int gmx_trjconv(int argc, char *argv[])
         /* skipping */
         if (skip_nr <= 0)
         {
+            gmx_fatal(FARGS, "Argument for -skip (%d) needs to be greater or equal to 1.", skip_nr);
         }
 
-        mtop = read_mtop_for_tng(top_file, in_file, out_file);
+        std::unique_ptr<gmx_mtop_t> mtop = read_mtop_for_tng(top_file, in_file, out_file);
 
         /* Determine whether to read a topology */
         bTPS = (ftp2bSet(efTPS, NFILE, fnm) ||
@@ -1136,7 +1139,7 @@ int gmx_trjconv(int argc, char *argv[])
         if (opt2bSet("-fr", NFILE, fnm))
         {
             printf("Select groups of frame number indices:\n");
-            rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, (int **)&frindex, &frname);
+            rd_index(opt2fn("-fr", NFILE, fnm), 1, &nrfri, &frindex, &frname);
             if (debug)
             {
                 for (i = 0; i < nrfri; i++)
@@ -1361,15 +1364,14 @@ int gmx_trjconv(int argc, char *argv[])
             switch (ftp)
             {
                 case efTNG:
-                    trjtools_gmx_prepare_tng_writing(out_file,
-                                                     filemode[0],
-                                                     trxin,
-                                                     &trxout,
-                                                     nullptr,
-                                                     nout,
-                                                     mtop,
-                                                     index,
-                                                     grpnm);
+                    trxout = trjtools_gmx_prepare_tng_writing(out_file,
+                                                              filemode[0],
+                                                              trxin,
+                                                              nullptr,
+                                                              nout,
+                                                              mtop.get(),
+                                                              index,
+                                                              grpnm);
                     break;
                 case efXTC:
                 case efTRR:
@@ -1441,7 +1443,7 @@ int gmx_trjconv(int argc, char *argv[])
                 if (bSetBox)
                 {
                     /* generate new box */
-                    if (fr.bBox == FALSE)
+                    if (!fr.bBox)
                     {
                         clear_mat(fr.box);
                     }
@@ -1453,7 +1455,7 @@ int gmx_trjconv(int argc, char *argv[])
                         }
                         else
                         {
-                            if (fr.bBox == FALSE)
+                            if (!fr.bBox)
                             {
                                 gmx_fatal(FARGS, "Cannot preserve a box that does not exist.\n");
                             }
@@ -1575,7 +1577,7 @@ int gmx_trjconv(int argc, char *argv[])
                 }
 
                 bWriteFrame =
-                    ( ( !bTDump && !frindex && frame % skip_nr == 0 ) || bDumpFrame );
+                    ( ( !bTDump && (frindex == nullptr) && frame % skip_nr == 0 ) || bDumpFrame );
 
                 if (bWriteFrame && (bDropUnder || bDropOver))
                 {
@@ -1615,8 +1617,7 @@ int gmx_trjconv(int argc, char *argv[])
                     {
                         frout_time = tzero + frame*timestep;
                     }
-                    else
-                    if (bSetTime)
+                    else if (bSetTime)
                     {
                         frout_time += tshift;
                     }
@@ -1851,7 +1852,7 @@ int gmx_trjconv(int argc, char *argv[])
                                 }
                                 if (frout.bStep)
                                 {
-                                    sprintf(stepstr, " step= %" GMX_PRId64, frout.step);
+                                    sprintf(stepstr, " step= %" PRId64, frout.step);
                                 }
                                 else
                                 {
@@ -1978,8 +1979,7 @@ int gmx_trjconv(int argc, char *argv[])
         }
     }
 
-    sfree(mtop);
-    if (top)
+    if (bTPS)
     {
         done_top(top);
         sfree(top);
@@ -1991,7 +1991,6 @@ int gmx_trjconv(int argc, char *argv[])
     sfree(grpnm);
     sfree(index);
     sfree(cindex);
-    done_filenms(NFILE, fnm);
     done_frame(&fr);
 
     do_view(oenv, out_file, nullptr);