Added tests for gmx clustsize.
authorDavid van der Spoel <spoel@xray.bmc.uu.se>
Tue, 7 Mar 2017 18:39:25 +0000 (19:39 +0100)
committerDavid van der Spoel <davidvanderspoel@gmail.com>
Sun, 12 Mar 2017 09:11:00 +0000 (10:11 +0100)
Some bone headed stuff to allow development of a new
version later. Unfortunately this now touches code
outside the analysis tool due to memory leaks that have to be
fixed.

Part of #2132

Change-Id: I32558e3a7cd7448954cf093208c15ce7014942ff

17 files changed:
src/gromacs/fileio/trxio.cpp
src/gromacs/fileio/trxio.h
src/gromacs/gmxana/gmx_clustsize.cpp
src/gromacs/topology/index.cpp
src/gromacs/trajectory/trajectoryframe.cpp
src/gromacs/trajectory/trajectoryframe.h
src/gromacs/trajectoryanalysis/tests/CMakeLists.txt
src/gromacs/trajectoryanalysis/tests/clustsize.cpp [new file with mode: 0644]
src/gromacs/trajectoryanalysis/tests/clustsize.ndx [new file with mode: 0644]
src/gromacs/trajectoryanalysis/tests/clustsize.pdb [new file with mode: 0644]
src/gromacs/trajectoryanalysis/tests/clustsize.tpr [new file with mode: 0644]
src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_MolDefaultCutoff.xml [new file with mode: 0644]
src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_MolShortCutoff.xml [new file with mode: 0644]
src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_NoMolDefaultCutoff.xml [new file with mode: 0644]
src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_NoMolShortCutoff.xml [new file with mode: 0644]
src/testutils/cmdlinetest.cpp
src/testutils/cmdlinetest.h

index 1815a9b78c8efab4cc3f8e7cdc4e7a56abe90296..16ca193e5c54ccb7f19c27e09d4b6262ea5c040d 100644 (file)
@@ -629,6 +629,10 @@ void close_trx(t_trxstatus *status)
     {
         gmx_fio_close(status->fio);
     }
+    sfree(status->persistent_line);
+#if GMX_USE_PLUGINS
+    sfree(status->vmdplugin);
+#endif
     sfree(status);
 }
 
index 3f9042ff0f626d63ce5e4cf6b17be6b379256fd6..1995e27f73a3cee7bb5effd16657a02684a58e03 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014,2015,2016, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -132,8 +132,9 @@ void write_tng_frame(t_trxstatus        *status,
                      struct t_trxframe  *fr);
 
 void close_trx(t_trxstatus *status);
-/* Close trajectory file as opened with read_first_x, read_frist_frame
+/* Close trajectory file as opened with read_first_x, read_first_frame
  * or open_trx. Identical to close_trj.
+ * Also frees memory in the structure.
  */
 
 t_trxstatus *open_trx(const char *outfile, const char *filemode);
index 8ba0c0ad9488a0c7a4416a157d37e41f12476b94..a88b365a6c9fabed45dc27c41f31f874b0bf452b 100644 (file)
@@ -40,6 +40,7 @@
 
 #include <algorithm>
 
+#include "gromacs/commandline/filenm.h"
 #include "gromacs/commandline/pargs.h"
 #include "gromacs/fileio/matio.h"
 #include "gromacs/fileio/tpxio.h"
@@ -76,7 +77,6 @@ static void clust_size(const char *ndx, const char *trx, const char *xpm,
     t_trxstatus          *status;
     rvec                 *x = nullptr, *v = nullptr, dx;
     t_pbc                 pbc;
-    char                 *gname;
     char                  timebuf[32];
     gmx_bool              bSame, bTPRwarn = TRUE;
     /* Topology stuff */
@@ -148,11 +148,12 @@ static void clust_size(const char *ndx, const char *trx, const char *xpm,
         {
             index[i] = i;
         }
-        gname = gmx_strdup("mols");
     }
     else
     {
+        char *gname;
         rd_index(ndx, 1, &nindex, &index, &gname);
+        sfree(gname);
     }
 
     snew(clust_index, nindex);
@@ -334,6 +335,7 @@ static void clust_size(const char *ndx, const char *trx, const char *xpm,
     }
     while (read_next_frame(oenv, status, &fr));
     close_trx(status);
+    done_frame(&fr);
     xvgrclose(fp);
     xvgrclose(gp);
     xvgrclose(hp);
@@ -426,7 +428,18 @@ static void clust_size(const char *ndx, const char *trx, const char *xpm,
                "Size", n_x, max_size, t_x, t_y, cs_dist, 0, cmid, cmax,
                rlo, rmid, rhi, &nlevels);
     gmx_ffclose(fp);
-
+    if (mtop)
+    {
+        done_mtop(mtop);
+        sfree(mtop);
+    }
+    sfree(t_x);
+    sfree(t_y);
+    for (i = 0; (i < n_x); i++)
+    {
+        sfree(cs_dist[i]);
+    }
+    sfree(cs_dist);
     sfree(clust_index);
     sfree(clust_size);
     sfree(index);
@@ -526,5 +539,8 @@ int gmx_clustsize(int argc, char *argv[])
                bMol, bPBC, fnTPR,
                cutoff, nskip, nlevels, rgblo, rgbhi, ndf, oenv);
 
+    done_filenms(NFILE, fnm);
+    output_env_done(oenv);
+
     return 0;
 }
index 4d148ce2cb13a68e6727cdae1e0474fbf2b4e811..2f5255e34e1372d6afc2662b997219c6c7c11034 100644 (file)
@@ -999,6 +999,14 @@ void rd_index(const char *statfile, int ngrps, int isize[],
     }
     grps = init_index(statfile, &gnames);
     rd_groups(grps, gnames, grpnames, ngrps, isize, index, grpnr);
+    for (int i = 0; i < grps->nr; i++)
+    {
+        sfree(gnames[i]);
+    }
+    sfree(gnames);
+    sfree(grpnr);
+    done_blocka(grps);
+    sfree(grps);
 }
 
 void rd_index_nrs(char *statfile, int ngrps, int isize[],
index c26d39ac6e48e145ae30bbcfed7485d28adc4d61..1abda4ee3f544833c784f97b96b581759054e84f 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2016, by the GROMACS development team, led by
+ * Copyright (c) 2016,2017, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -43,6 +43,7 @@
 #include "gromacs/math/veccompare.h"
 #include "gromacs/topology/atoms.h"
 #include "gromacs/utility/compare.h"
+#include "gromacs/utility/smalloc.h"
 
 void comp_frame(FILE *fp, t_trxframe *fr1, t_trxframe *fr2,
                 gmx_bool bRMSD, real ftol, real abstol)
@@ -92,3 +93,15 @@ void comp_frame(FILE *fp, t_trxframe *fr1, t_trxframe *fr2,
         cmp_rvecs(fp, "box", 3, fr1->box, fr2->box, FALSE, ftol, abstol);
     }
 }
+
+void done_frame(t_trxframe *frame)
+{
+    if (frame->atoms)
+    {
+        done_atom(frame->atoms);
+        sfree(frame->atoms);
+    }
+    sfree(frame->x);
+    sfree(frame->v);
+    sfree(frame->f);
+}
index 9b65402684304e10869a5996a1d4ac28926539de..2c8e34798ae3f92dd3a3183b106e986d35947855 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014,2015,2016, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -86,4 +86,6 @@ typedef struct t_trxframe
 void comp_frame(FILE *fp, t_trxframe *fr1, t_trxframe *fr2,
                 gmx_bool bRMSD, real ftol, real abstol);
 
+void done_frame(t_trxframe *frame);
+
 #endif
index 965330236a953ed53e7f228ee2552bf9693af9e5..db7efdc05c1156fefddee2b2fa35dfbbc6b80fe6 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2010,2012,2013,2014,2015,2016, by the GROMACS development team, led by
+# Copyright (c) 2010,2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
 # Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
 # and including many others, as listed in the AUTHORS file in the
 # top-level source directory and at http://www.gromacs.org.
@@ -34,6 +34,7 @@
 
 gmx_add_unit_test(TrajectoryAnalysisUnitTests trajectoryanalysis-test
                   moduletest.cpp
+                  clustsize.cpp
                   cmdlinerunner.cpp
                   angle.cpp
                   distance.cpp
diff --git a/src/gromacs/trajectoryanalysis/tests/clustsize.cpp b/src/gromacs/trajectoryanalysis/tests/clustsize.cpp
new file mode 100644 (file)
index 0000000..7f853e8
--- /dev/null
@@ -0,0 +1,144 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2017, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \internal \file
+ * \brief
+ * Tests for gmx clustsize
+ *
+ * \todo These will be superseded by tests of the new style analysis
+ * modules.
+ *
+ * \author David van der Spoel <david.vanderspoel@icm.uu.se>
+ * \ingroup module_trajectoryanalysis
+ */
+
+#include "gmxpre.h"
+
+#include <cstring>
+
+#include <string>
+
+#include "gromacs/gmxana/gmx_ana.h"
+#include "gromacs/utility/filestream.h"
+#include "gromacs/utility/path.h"
+
+#include "testutils/cmdlinetest.h"
+#include "testutils/integrationtests.h"
+#include "testutils/refdata.h"
+#include "testutils/testasserts.h"
+#include "testutils/xvgtest.h"
+
+namespace gmx
+{
+
+namespace test
+{
+
+namespace
+{
+
+class ClustsizeTest : public CommandLineTestBase
+{
+    public:
+        ClustsizeTest()
+        {
+            setInputFile("-f", "clustsize.pdb");
+        }
+
+        void runTest(bool        mol,
+                     bool        cutoff)
+        {
+            const char *const command[] = { "clustsize" };
+            CommandLine       args      = CommandLine(command);
+            double            tolerance = 1e-4;
+            test::XvgMatch    xvg;
+
+            setOutputFile("-mc", ".xvg",
+                          xvg.tolerance(gmx::test::relativeToleranceAsFloatingPoint(1, tolerance)));
+            setOutputFile("-nc", ".xvg",
+                          xvg.tolerance(gmx::test::relativeToleranceAsFloatingPoint(1, tolerance)));
+            setOutputFile("-ac", ".xvg",
+                          xvg.tolerance(gmx::test::relativeToleranceAsFloatingPoint(1, tolerance)));
+            setOutputFile("-hc", ".xvg",
+                          xvg.tolerance(gmx::test::relativeToleranceAsFloatingPoint(1, tolerance)));
+
+            if (mol)
+            {
+                setInputFile("-s", "clustsize.tpr");
+                args.addOption("-mol");
+            }
+            else
+            {
+                setInputFile("-n", "clustsize.ndx");
+            }
+            if (cutoff)
+            {
+                args.addOption("-cut", "0.3");
+            }
+            rootChecker().checkString(args.toString(), "CommandLine");
+
+            CommandLine &cmdline = commandLine();
+            cmdline.merge(args);
+
+            ASSERT_EQ(0, gmx_clustsize(cmdline.argc(), cmdline.argv()));
+
+            checkOutputFiles();
+        }
+};
+
+TEST_F(ClustsizeTest, NoMolDefaultCutoff)
+{
+    runTest(false, false);
+}
+
+TEST_F(ClustsizeTest, NoMolShortCutoff)
+{
+    runTest(false, true);
+}
+
+TEST_F(ClustsizeTest, MolDefaultCutoff)
+{
+    runTest(true, false);
+}
+
+TEST_F(ClustsizeTest, MolShortCutoff)
+{
+    runTest(true, true);
+}
+
+} // namespace
+
+} // namespace
+
+} // namespace
diff --git a/src/gromacs/trajectoryanalysis/tests/clustsize.ndx b/src/gromacs/trajectoryanalysis/tests/clustsize.ndx
new file mode 100644 (file)
index 0000000..1fc0f43
--- /dev/null
@@ -0,0 +1,3 @@
+[ SOL ]
+   1    2    3    4    5    6    7    8    9   10   11   12   13   14   15 
+  16   17   18   19   20   21   22   23   24 
diff --git a/src/gromacs/trajectoryanalysis/tests/clustsize.pdb b/src/gromacs/trajectoryanalysis/tests/clustsize.pdb
new file mode 100644 (file)
index 0000000..e8abf0f
--- /dev/null
@@ -0,0 +1,57 @@
+TITLE     216H2O,WATJP01,SPC216,SPC-MODEL,300K,BOX(M)=1.86206NM,WFVG,MAR. 1984
+REMARK    THIS IS A SIMULATION BOX
+CRYST1   18.621   18.621   18.621  90.00  90.00  90.00 P 1           1
+MODEL        1
+ATOM      1  OW  SOL     1       2.300   6.280   1.130  1.00  0.00            
+ATOM      2  HW1 SOL     1       1.370   6.260   1.500  1.00  0.00            
+ATOM      3  HW2 SOL     1       2.310   5.890   0.210  1.00  0.00            
+ATOM     70  OW  SOL    24       2.220   6.430   4.630  1.00  0.00            
+ATOM     71  HW1 SOL    24       0.770   5.550   4.800  1.00  0.00            
+ATOM     72  HW2 SOL    24       1.210   6.970   6.470  1.00  0.00            
+ATOM     73  OW  SOL    25      -0.200  -0.950   3.590  1.00  0.00            
+ATOM     74  HW1 SOL    25       0.340  -1.240   4.390  1.00  0.00            
+ATOM     75  HW2 SOL    25       0.100  -0.050   3.300  1.00  0.00            
+ATOM     76  OW  SOL    26       0.270  -2.660   1.170  1.00  0.00            
+ATOM     77  HW1 SOL    26       0.080  -3.620   1.380  1.00  0.00            
+ATOM     78  HW2 SOL    26      -0.060  -2.080   1.920  1.00  0.00            
+ATOM    637  OW  SOL   213      -5.620   4.530   6.910  1.00  0.00            
+ATOM    638  HW1 SOL   213      -6.210   5.330   6.950  1.00  0.00            
+ATOM    639  HW2 SOL   213      -5.470   4.180   7.840  1.00  0.00            
+ATOM    640  OW  SOL   214      -2.690   2.210   8.820  1.00  0.00            
+ATOM    641  HW1 SOL   214      -3.530   2.200   9.360  1.00  0.00            
+ATOM    642  HW2 SOL   214      -2.670   3.040   8.260  1.00  0.00            
+ATOM    643  OW  SOL   215       0.390  -7.850   3.000  1.00  0.00            
+ATOM    644  HW1 SOL   215       1.380  -7.960   2.910  1.00  0.00            
+ATOM    645  HW2 SOL   215      -0.010  -8.710   3.320  1.00  0.00            
+ATOM    646  OW  SOL   216       8.750  -2.160   3.370  1.00  0.00            
+ATOM    647  HW1 SOL   216       7.980  -2.510   2.830  1.00  0.00            
+ATOM    648  HW2 SOL   216       8.430  -1.450   3.990  1.00  0.00            
+TER
+ENDMDL
+MODEL        2
+ATOM      1  OW  SOL     1       2.300   6.280   1.130  1.00  0.00            
+ATOM      2  HW1 SOL     1       1.370   6.260   1.500  1.00  0.00            
+ATOM      3  HW2 SOL     1       2.310   5.890   0.210  1.00  0.00            
+ATOM     70  OW  SOL    24       1.220   6.430   5.630  1.00  0.00            
+ATOM     71  HW1 SOL    24       0.770   5.550   5.800  1.00  0.00            
+ATOM     72  HW2 SOL    24       1.210   6.970   6.470  1.00  0.00            
+ATOM     73  OW  SOL    25      -0.200  -0.950   3.590  1.00  0.00            
+ATOM     74  HW1 SOL    25       0.340  -1.240   4.390  1.00  0.00            
+ATOM     75  HW2 SOL    25       0.100  -0.050   3.300  1.00  0.00            
+ATOM     76  OW  SOL    26       0.270  -2.660   1.170  1.00  0.00            
+ATOM     77  HW1 SOL    26       0.080  -3.620   1.380  1.00  0.00            
+ATOM     78  HW2 SOL    26      -0.060  -2.080   1.920  1.00  0.00            
+ATOM    637  OW  SOL   213      -5.620   4.530   6.910  1.00  0.00            
+ATOM    638  HW1 SOL   213      -6.210   5.330   6.950  1.00  0.00            
+ATOM    639  HW2 SOL   213      -5.470   4.180   7.840  1.00  0.00            
+ATOM    640  OW  SOL   214      -2.690   2.210   8.820  1.00  0.00            
+ATOM    641  HW1 SOL   214      -3.530   2.200   9.360  1.00  0.00            
+ATOM    642  HW2 SOL   214      -2.670   3.040   8.260  1.00  0.00            
+ATOM    643  OW  SOL   215       0.390  -7.850   3.000  1.00  0.00            
+ATOM    644  HW1 SOL   215       1.380  -7.960   2.910  1.00  0.00            
+ATOM    645  HW2 SOL   215      -0.010  -8.710   3.320  1.00  0.00            
+ATOM    646  OW  SOL   216       8.750  -2.160   3.370  1.00  0.00            
+ATOM    647  HW1 SOL   216       7.980  -2.510   2.830  1.00  0.00            
+ATOM    648  HW2 SOL   216       8.430  -1.450   3.990  1.00  0.00            
+TER
+ENDMDL
diff --git a/src/gromacs/trajectoryanalysis/tests/clustsize.tpr b/src/gromacs/trajectoryanalysis/tests/clustsize.tpr
new file mode 100644 (file)
index 0000000..8eb67f8
Binary files /dev/null and b/src/gromacs/trajectoryanalysis/tests/clustsize.tpr differ
diff --git a/src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_MolDefaultCutoff.xml b/src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_MolDefaultCutoff.xml
new file mode 100644 (file)
index 0000000..8543816
--- /dev/null
@@ -0,0 +1,105 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <String Name="CommandLine">clustsize -mol</String>
+  <OutputFiles Name="Files">
+    <File Name="-mc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Max cluster size"
+xaxis  label "Time (ps)"
+yaxis  label "#molecules"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>2</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>2</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-nc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Number of clusters"
+xaxis  label "Time (ps)"
+yaxis  label "N"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>7</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>7</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-ac">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Average cluster size"
+xaxis  label "Time (ps)"
+yaxis  label "#molecules"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>2.000</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>2.000</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-hc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Cluster size distribution"
+xaxis  label "Cluster size"
+yaxis  label "()"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>0</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>1</Real>
+          <Real>6.000</Real>
+        </Sequence>
+        <Sequence Name="Row2">
+          <Int Name="Length">2</Int>
+          <Real>2</Real>
+          <Real>1.000</Real>
+        </Sequence>
+        <Sequence Name="Row3">
+          <Int Name="Length">2</Int>
+          <Real>3</Real>
+          <Real>0.000</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+  </OutputFiles>
+</ReferenceData>
diff --git a/src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_MolShortCutoff.xml b/src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_MolShortCutoff.xml
new file mode 100644 (file)
index 0000000..018d47d
--- /dev/null
@@ -0,0 +1,105 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <String Name="CommandLine">clustsize -mol -cut 0.3</String>
+  <OutputFiles Name="Files">
+    <File Name="-mc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Max cluster size"
+xaxis  label "Time (ps)"
+yaxis  label "#molecules"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>2</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>2</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-nc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Number of clusters"
+xaxis  label "Time (ps)"
+yaxis  label "N"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>7</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>7</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-ac">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Average cluster size"
+xaxis  label "Time (ps)"
+yaxis  label "#molecules"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>2.000</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>2.000</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-hc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Cluster size distribution"
+xaxis  label "Cluster size"
+yaxis  label "()"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>0</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>1</Real>
+          <Real>6.000</Real>
+        </Sequence>
+        <Sequence Name="Row2">
+          <Int Name="Length">2</Int>
+          <Real>2</Real>
+          <Real>1.000</Real>
+        </Sequence>
+        <Sequence Name="Row3">
+          <Int Name="Length">2</Int>
+          <Real>3</Real>
+          <Real>0.000</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+  </OutputFiles>
+</ReferenceData>
diff --git a/src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_NoMolDefaultCutoff.xml b/src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_NoMolDefaultCutoff.xml
new file mode 100644 (file)
index 0000000..8ae771b
--- /dev/null
@@ -0,0 +1,125 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <String Name="CommandLine">clustsize</String>
+  <OutputFiles Name="Files">
+    <File Name="-mc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Max cluster size"
+xaxis  label "Time (ps)"
+yaxis  label "#molecules"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>6</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>6</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-nc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Number of clusters"
+xaxis  label "Time (ps)"
+yaxis  label "N"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>5</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>6</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-ac">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Average cluster size"
+xaxis  label "Time (ps)"
+yaxis  label "#molecules"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>4.800</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>4.000</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-hc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Cluster size distribution"
+xaxis  label "Cluster size"
+yaxis  label "()"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>0</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>1</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row2">
+          <Int Name="Length">2</Int>
+          <Real>2</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row3">
+          <Int Name="Length">2</Int>
+          <Real>3</Real>
+          <Real>3.000</Real>
+        </Sequence>
+        <Sequence Name="Row4">
+          <Int Name="Length">2</Int>
+          <Real>4</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row5">
+          <Int Name="Length">2</Int>
+          <Real>5</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row6">
+          <Int Name="Length">2</Int>
+          <Real>6</Real>
+          <Real>2.500</Real>
+        </Sequence>
+        <Sequence Name="Row7">
+          <Int Name="Length">2</Int>
+          <Real>7</Real>
+          <Real>0.000</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+  </OutputFiles>
+</ReferenceData>
diff --git a/src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_NoMolShortCutoff.xml b/src/gromacs/trajectoryanalysis/tests/refdata/ClustsizeTest_NoMolShortCutoff.xml
new file mode 100644 (file)
index 0000000..df027f4
--- /dev/null
@@ -0,0 +1,125 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <String Name="CommandLine">clustsize -cut 0.3</String>
+  <OutputFiles Name="Files">
+    <File Name="-mc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Max cluster size"
+xaxis  label "Time (ps)"
+yaxis  label "#molecules"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>6</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>6</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-nc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Number of clusters"
+xaxis  label "Time (ps)"
+yaxis  label "N"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>7</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>7</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-ac">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Average cluster size"
+xaxis  label "Time (ps)"
+yaxis  label "#molecules"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>1.000000e+00</Real>
+          <Real>3.429</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>2.000000e+00</Real>
+          <Real>3.429</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+    <File Name="-hc">
+      <XvgLegend Name="Legend">
+        <String Name="XvgLegend"><![CDATA[
+title "Cluster size distribution"
+xaxis  label "Cluster size"
+yaxis  label "()"
+TYPE xy
+]]></String>
+      </XvgLegend>
+      <XvgData Name="Data">
+        <Sequence Name="Row0">
+          <Int Name="Length">2</Int>
+          <Real>0</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row1">
+          <Int Name="Length">2</Int>
+          <Real>1</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row2">
+          <Int Name="Length">2</Int>
+          <Real>2</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row3">
+          <Int Name="Length">2</Int>
+          <Real>3</Real>
+          <Real>6.000</Real>
+        </Sequence>
+        <Sequence Name="Row4">
+          <Int Name="Length">2</Int>
+          <Real>4</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row5">
+          <Int Name="Length">2</Int>
+          <Real>5</Real>
+          <Real>0.000</Real>
+        </Sequence>
+        <Sequence Name="Row6">
+          <Int Name="Length">2</Int>
+          <Real>6</Real>
+          <Real>1.000</Real>
+        </Sequence>
+        <Sequence Name="Row7">
+          <Int Name="Length">2</Int>
+          <Real>7</Real>
+          <Real>0.000</Real>
+        </Sequence>
+      </XvgData>
+    </File>
+  </OutputFiles>
+</ReferenceData>
index 4cc338fe6c7fdb7eee49486a612018ef099f0df6..1f79ee3486ed62a45596041b045721f91a5ce595 100644 (file)
@@ -173,6 +173,11 @@ std::string value2string(T value)
 
 }       // namespace
 
+void CommandLine::addOption(const char *name)
+{
+    append(name);
+}
+
 void CommandLine::addOption(const char *name, const char *value)
 {
     append(name);
index d4573e27d4e125e384bb663fb2c72eab7342c02c..7bdeed3e136e4b4c61f1b8fef283d7e5d793e855 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012,2013,2014,2015,2016, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -140,6 +140,13 @@ class CommandLine
         void append(const char *arg);
         //! Convenience overload taking a std::string.
         void append(const std::string &arg) { append(arg.c_str()); }
+        /*! \brief
+         * Adds an option to the command line, typically a boolean.
+         *
+         * \param[in] name   Name of the option to append, which
+         *                   should start with "-".
+         */
+        void addOption(const char *name);
         /*! \brief
          * Adds an option-value pair to the command line.
          *