Clean-up and improvements in C++ tools.
authorTeemu Murtola <teemu.murtola@gmail.com>
Fri, 23 Aug 2013 04:25:39 +0000 (07:25 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Fri, 30 Aug 2013 12:43:03 +0000 (14:43 +0200)
- Make gmx-select accept multiple selections for all output options.
- Make gmx-gangle and gmx-distance accept dynamic selections.
- Remove -dump option from gmx-select and remove commented code related
  to the same in the other tools.

Part of #665.

Change-Id: I5e3c4caf2928a95c2620d969b4f447272ad9b051

15 files changed:
src/gromacs/analysisdata/modules/histogram.cpp
src/gromacs/analysisdata/modules/histogram.h
src/gromacs/trajectoryanalysis/modules/angle.cpp
src/gromacs/trajectoryanalysis/modules/distance.cpp
src/gromacs/trajectoryanalysis/modules/select.cpp
src/gromacs/trajectoryanalysis/modules/select.h
src/gromacs/trajectoryanalysis/tests/angle.cpp
src/gromacs/trajectoryanalysis/tests/distance.cpp
src/gromacs/trajectoryanalysis/tests/refdata/AngleModuleTest_HandlesDynamicSelections.xml [new file with mode: 0644]
src/gromacs/trajectoryanalysis/tests/refdata/DistanceModuleTest_HandlesDynamicSelections.xml [new file with mode: 0644]
src/gromacs/trajectoryanalysis/tests/refdata/SelectModuleTest_BasicTest.xml
src/gromacs/trajectoryanalysis/tests/refdata/SelectModuleTest_HandlesDumpOption.xml [deleted file]
src/gromacs/trajectoryanalysis/tests/refdata/SelectModuleTest_HandlesMaxPDBOutput.xml
src/gromacs/trajectoryanalysis/tests/refdata/SelectModuleTest_HandlesSelectedPDBOutput.xml
src/gromacs/trajectoryanalysis/tests/select.cpp

index 2afc2368b69ed9427b0558c6e6b5c8b4a6410cee..5360430716a577afe4245ad8b814758673fa6523 100644 (file)
@@ -619,7 +619,8 @@ AnalysisDataSimpleHistogramModule::settings() const
 int
 AnalysisDataSimpleHistogramModule::flags() const
 {
-    return efAllowMulticolumn | efAllowMultipoint | efAllowMultipleDataSets;
+    return efAllowMulticolumn | efAllowMultipoint | efAllowMissing
+           | efAllowMultipleDataSets;
 }
 
 
@@ -653,10 +654,13 @@ AnalysisDataSimpleHistogramModule::pointsAdded(const AnalysisDataPointSetRef &po
     frame.selectDataSet(points.dataSetIndex());
     for (int i = 0; i < points.columnCount(); ++i)
     {
-        int bin = settings().findBin(points.y(i));
-        if (bin != -1)
+        if (points.present(i))
         {
-            frame.value(bin) += 1;
+            const int bin = settings().findBin(points.y(i));
+            if (bin != -1)
+            {
+                frame.value(bin) += 1;
+            }
         }
     }
 }
index bbed9f8142dd44c033d9b50ea07ea674cc6f5fd4..58d18bd42abe8ea8b4ed4befaaeef2e14fd27e43 100644 (file)
@@ -338,7 +338,7 @@ class AbstractAverageHistogram : public AbstractAnalysisArrayData
  * Output data contains the same number of frames and data sets as the input
  * data.  Each frame contains the histogram(s) for the points in that frame.
  * Each input data set is processed independently into the corresponding output
- * data set.
+ * data set.  Missing values are ignored.
  * All input columns for a data set are averaged into the same histogram.
  * The number of columns for all data sets equals the number of bins in the
  * histogram.
index edbaec81829b99da8b8725260e078850322fb8c9..b40188071387ef3fc34b47b418e1c3a23187da51 100644 (file)
@@ -135,12 +135,6 @@ Angle::initOptions(Options *options, TrajectoryAnalysisSettings * /*settings*/)
         "set with [TT]-binw[tt].",
         "For [TT]-oav[tt] and [TT]-oh[tt], separate average/histogram is",
         "computed for each selection in [TT]-group1[tt]."
-        /* TODO: Consider if the dump option is necessary and how to best
-         * implement it.
-           "[TT]-od[tt] can be used to dump all the individual angles,",
-           "each on a separate line. This format is better suited for",
-           "further processing, e.g., if angles from multiple runs are needed."
-         */
     };
     static const char *const cGroup1TypeEnum[] =
     { "angle", "dihedral", "vector", "plane" };
@@ -168,15 +162,12 @@ Angle::initOptions(Options *options, TrajectoryAnalysisSettings * /*settings*/)
     options->addOption(DoubleOption("binw").store(&binWidth_)
                            .description("Binwidth for -oh in degrees"));
 
-    // TODO: Consider what is the best way to support dynamic selections.
-    // Again, most of the code already supports it, but it needs to be
-    // considered how should -oall work, and additional checks should be added.
     sel1info_ = options->addOption(SelectionOption("group1")
-                                       .required().onlyStatic().storeVector(&sel1_)
+                                       .required().dynamicMask().storeVector(&sel1_)
                                        .multiValue()
                                        .description("First analysis/vector selection"));
     sel2info_ = options->addOption(SelectionOption("group2")
-                                       .onlyStatic().storeVector(&sel2_)
+                                       .dynamicMask().storeVector(&sel2_)
                                        .multiValue()
                                        .description("Second analysis/vector selection"));
 }
@@ -243,8 +234,8 @@ Angle::checkSelections(const SelectionList &sel1,
 
     for (size_t g = 0; g < sel1.size(); ++g)
     {
-        int na1 = sel1[g].posCount();
-        int na2 = (natoms2_ > 0) ? sel2[g].posCount() : 0;
+        const int na1 = sel1[g].posCount();
+        const int na2 = (natoms2_ > 0) ? sel2[g].posCount() : 0;
         if (natoms1_ > 1 && na1 % natoms1_ != 0)
         {
             GMX_THROW(InconsistentInputError(formatString(
@@ -267,6 +258,30 @@ Angle::checkSelections(const SelectionList &sel1,
             GMX_THROW(InconsistentInputError(
                               "The second group should contain a single position with -g2 sphnorm"));
         }
+        if (sel1[g].isDynamic() || (natoms2_ > 0 && sel2[g].isDynamic()))
+        {
+            for (int i = 0, j = 0; i < na1; i += natoms1_, j += natoms2_)
+            {
+                const bool bSelected = sel1[g].position(i).selected();
+                bool       bOk       = true;
+                for (int k = 1; k < natoms1_ && bOk; ++k)
+                {
+                    bOk = (sel1[g].position(i+k).selected() == bSelected);
+                }
+                for (int k = 1; k < natoms2_ && bOk; ++k)
+                {
+                    bOk = (sel2[g].position(j+k).selected() == bSelected);
+                }
+                if (!bOk)
+                {
+                    std::string message =
+                        formatString("Dynamic selection %d does not select "
+                                     "a consistent set of angles over the frames",
+                                     static_cast<int>(g + 1));
+                    GMX_THROW(InconsistentInputError(message));
+                }
+            }
+        }
     }
 }
 
@@ -422,6 +437,9 @@ Angle::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
         {
             rvec x[4];
             real angle;
+            // checkSelections() ensures that this reflects all the involved
+            // positions.
+            bool bPresent = sel1[g].position(i).selected();
             copy_pos(sel1, natoms1_, g, i, x);
             switch (g1type_[0])
             {
@@ -502,24 +520,7 @@ Angle::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
                 default:
                     GMX_THROW(InternalError("invalid -g1 value"));
             }
-            /* TODO: Should we also calculate distances like g_sgangle?
-             * Could be better to leave that for a separate tool.
-               real dist = 0.0;
-               if (bDumpDist_)
-               {
-                if (pbc)
-                {
-                    rvec dx;
-                    pbc_dx(pbc, c2, c1, dx);
-                    dist = norm(dx);
-                }
-                else
-                {
-                    dist = sqrt(distance2(c1, c2));
-                }
-               }
-             */
-            dh.setPoint(n, angle * RAD2DEG);
+            dh.setPoint(n, angle * RAD2DEG, bPresent);
         }
     }
     dh.finishFrame();
index 9fbadd034b7b7af12455bdf997ab62f75ce95546..de2c2fae1304521f62b3aea089c38691158c327f 100644 (file)
@@ -131,11 +131,8 @@ Distance::initOptions(Options *options, TrajectoryAnalysisSettings * /*settings*
     options->addOption(FileNameOption("oallstat").filetype(eftPlot).outputFile()
                            .store(&fnAllStats_).defaultBasename("diststat")
                            .description("Statistics for individual distances"));
-    // TODO: Consider what is the best way to support dynamic selections.
-    // Again, most of the code already supports it, but it needs to be
-    // considered how should -oall work, and additional checks should be added.
     options->addOption(SelectionOption("select").storeVector(&sel_)
-                           .required().onlyStatic().multiValue()
+                           .required().dynamicMask().multiValue()
                            .description("Position pairs to calculate distances for"));
     // TODO: Extend the histogramming implementation to allow automatic
     // extension of the histograms to cover the data, removing the need for
@@ -167,6 +164,20 @@ void checkSelections(const SelectionList &sel)
                         sel[g].name(), sel[g].posCount());
             GMX_THROW(InconsistentInputError(message));
         }
+        if (sel[g].isDynamic())
+        {
+            for (int i = 0; i < sel[g].posCount(); i += 2)
+            {
+                if (sel[g].position(i).selected() != sel[g].position(i+1).selected())
+                {
+                    std::string message =
+                        formatString("Dynamic selection %d does not select "
+                                     "a consistent set of pairs over the frames",
+                                     static_cast<int>(g + 1));
+                    GMX_THROW(InconsistentInputError(message));
+                }
+            }
+        }
     }
 }
 
@@ -285,8 +296,9 @@ Distance::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
             {
                 rvec_sub(p2.x(), p1.x(), dx);
             }
-            real dist = norm(dx);
-            distHandle.setPoint(n, dist);
+            real dist     = norm(dx);
+            bool bPresent = p1.selected() && p2.selected();
+            distHandle.setPoint(n, dist, bPresent);
             xyzHandle.setPoints(n*3, 3, dx);
         }
     }
index 12f5ee13801d0ef647058505f17973de15d873db..0cee1e54d193c171181ed0869fc2605ad2f62faa 100644 (file)
@@ -45,6 +45,7 @@
 #include <cstring>
 
 #include <algorithm>
+#include <set>
 #include <string>
 #include <vector>
 
@@ -259,11 +260,14 @@ const char Select::shortDescription[] =
 Select::Select()
     : TrajectoryAnalysisModule(name, shortDescription),
       selOpt_(NULL),
-      bDump_(false), bTotNorm_(false), bFracNorm_(false), bResInd_(false),
+      bTotNorm_(false), bFracNorm_(false), bResInd_(false),
       bCumulativeLifetimes_(true), top_(NULL),
       occupancyModule_(new AnalysisDataAverageModule()),
       lifetimeModule_(new AnalysisDataLifetimeModule())
 {
+    mdata_.addModule(occupancyModule_);
+    mdata_.addModule(lifetimeModule_);
+
     registerAnalysisDataset(&sdata_, "size");
     registerAnalysisDataset(&cdata_, "cfrac");
     idata_.setColumnCount(0, 2);
@@ -310,9 +314,7 @@ Select::initOptions(Options *options, TrajectoryAnalysisSettings * /*settings*/)
         "of positions, followed by the atom/residue/molecule numbers.",
         "If more than one selection is specified, the size of the second",
         "group immediately follows the last number of the first group",
-        "and so on. With [TT]-dump[tt], the frame time and the number",
-        "of positions is omitted from the output. In this case, only one",
-        "selection can be given.[PAR]",
+        "and so on.[PAR]",
         "With [TT]-on[tt], the selected atoms are written as a index file",
         "compatible with [TT]make_ndx[tt] and the analyzing tools. Each selection",
         "is written as a selection group and for dynamic selections a",
@@ -328,8 +330,7 @@ Select::initOptions(Options *options, TrajectoryAnalysisSettings * /*settings*/)
         "as a function of time. Each line in the output corresponds to",
         "one frame, and contains either 0/1 for each atom/residue/molecule",
         "possibly selected. 1 stands for the atom/residue/molecule being",
-        "selected for the current frame, 0 for not selected.",
-        "With [TT]-dump[tt], the frame time is omitted from the output.[PAR]",
+        "selected for the current frame, 0 for not selected.[PAR]",
         "With [TT]-of[tt], the occupancy fraction of each position (i.e.,",
         "the fraction of frames where the position is selected) is",
         "printed.[PAR]",
@@ -345,10 +346,8 @@ Select::initOptions(Options *options, TrajectoryAnalysisSettings * /*settings*/)
         "selected positions as a function of the time the position was",
         "continuously selected. [TT]-cumlt[tt] can be used to control whether",
         "subintervals of longer intervals are included in the histogram.[PAR]",
-        "With [TT]-om[tt], [TT]-of[tt], [TT]-olt[tt], and [TT]-ofpdb[tt],",
-        "only one selection can be provided.",
-        "[TT]-om[tt], [TT]-of[tt], and [TT]-olt[tt] only accept dynamic",
-        "selections."
+        "[TT]-om[tt], [TT]-of[tt], and [TT]-olt[tt] only make sense with",
+        "dynamic selections."
     };
 
     options->setDescription(concatenateStrings(desc));
@@ -382,8 +381,6 @@ Select::initOptions(Options *options, TrajectoryAnalysisSettings * /*settings*/)
                                      .required().multiValue()
                                      .description("Selections to analyze"));
 
-    options->addOption(BooleanOption("dump").store(&bDump_)
-                           .description("Do not print the frame time (-om, -oi) or the index size (-oi)"));
     options->addOption(BooleanOption("norm").store(&bTotNorm_)
                            .description("Normalize by total number of positions with -os"));
     options->addOption(BooleanOption("cfnorm").store(&bFracNorm_)
@@ -409,24 +406,12 @@ Select::optionsFinished(Options                     * /*options*/,
         settings->setFlag(TrajectoryAnalysisSettings::efRequireTop);
         settings->setFlag(TrajectoryAnalysisSettings::efUseTopX);
     }
-    if ((!fnIndex_.empty() && bDump_)
-        || !fnMask_.empty() || !fnOccupancy_.empty() || !fnPDB_.empty()
-        || !fnLifetime_.empty())
-    {
-        selOpt_->setValueCount(1);
-    }
 }
 
 void
 Select::initAnalysis(const TrajectoryAnalysisSettings &settings,
                      const TopologyInformation        &top)
 {
-    if (!sel_[0].isDynamic()
-        && (!fnMask_.empty() || !fnOccupancy_.empty() || !fnLifetime_.empty()))
-    {
-        GMX_THROW(InconsistentInputError(
-                          "-om, -of, or -olt are not meaningful with a static selection"));
-    }
     bResInd_ = (resNumberType_ == "index");
 
     for (SelectionList::iterator i = sel_.begin(); i != sel_.end(); ++i)
@@ -473,15 +458,7 @@ Select::initAnalysis(const TrajectoryAnalysisSettings &settings,
         plot->setFileName(fnIndex_);
         plot->setPlainOutput(true);
         plot->setYFormat(4, 0);
-        if (bDump_)
-        {
-            plot->setOmitX(bDump_);
-            idata_.addColumnModule(1, 1, plot);
-        }
-        else
-        {
-            idata_.addModule(plot);
-        }
+        idata_.addModule(plot);
     }
     if (!fnNdx_.empty())
     {
@@ -494,18 +471,17 @@ Select::initAnalysis(const TrajectoryAnalysisSettings &settings,
         idata_.addModule(writer);
     }
 
-    // TODO: Write out the mask and lifetimes for all the selections.
-    mdata_.setColumnCount(0, sel_[0].posCount());
-    mdata_.addModule(occupancyModule_);
-    mdata_.addModule(lifetimeModule_);
+    mdata_.setDataSetCount(sel_.size());
+    for (size_t g = 0; g < sel_.size(); ++g)
+    {
+        mdata_.setColumnCount(g, sel_[g].posCount());
+    }
     lifetimeModule_->setCumulative(bCumulativeLifetimes_);
     if (!fnMask_.empty())
     {
         AnalysisDataPlotModulePointer plot(
                 new AnalysisDataPlotModule(settings.plotSettings()));
         plot->setFileName(fnMask_);
-        plot->setPlainOutput(bDump_);
-        plot->setOmitX(bDump_);
         plot->setTitle("Selection mask");
         plot->setXAxisIsTime();
         plot->setYLabel("Occupancy");
@@ -520,7 +496,7 @@ Select::initAnalysis(const TrajectoryAnalysisSettings &settings,
         plot->setTitle("Fraction of time selection matches");
         plot->setXLabel("Selected position");
         plot->setYLabel("Occupied fraction");
-        occupancyModule_->addColumnModule(0, 1, plot);
+        occupancyModule_->addModule(plot);
     }
     if (!fnLifetime_.empty())
     {
@@ -589,13 +565,17 @@ Select::analyzeFrame(int frnr, const t_trxframe &fr, t_pbc *pbc,
     idh.finishFrame();
 
     mdh.startFrame(frnr, fr.time);
-    for (int i = 0; i < totsize_[0]; ++i)
-    {
-        mdh.setPoint(i, 0);
-    }
-    for (int i = 0; i < sel[0].posCount(); ++i)
+    for (size_t g = 0; g < sel.size(); ++g)
     {
-        mdh.setPoint(sel[0].position(i).refId(), 1);
+        mdh.selectDataSet(g);
+        for (int i = 0; i < totsize_[g]; ++i)
+        {
+            mdh.setPoint(i, 0);
+        }
+        for (int i = 0; i < sel[g].posCount(); ++i)
+        {
+            mdh.setPoint(sel[g].position(i).refId(), 1);
+        }
     }
     mdh.finishFrame();
 }
@@ -628,13 +608,17 @@ Select::writeOutput()
         {
             pdbinfo[i].occup = 0.0;
         }
-        for (int i = 0; i < sel_[0].posCount(); ++i)
+        for (size_t g = 0; g < sel_.size(); ++g)
         {
-            ConstArrayRef<int>                 atomIndices = sel_[0].position(i).atomIndices();
-            ConstArrayRef<int>::const_iterator ai;
-            for (ai = atomIndices.begin(); ai != atomIndices.end(); ++ai)
+            for (int i = 0; i < sel_[g].posCount(); ++i)
             {
-                pdbinfo[*ai].occup = occupancyModule_->average(0, i);
+                ConstArrayRef<int>                 atomIndices
+                    = sel_[g].position(i).atomIndices();
+                ConstArrayRef<int>::const_iterator ai;
+                for (ai = atomIndices.begin(); ai != atomIndices.end(); ++ai)
+                {
+                    pdbinfo[*ai].occup += occupancyModule_->average(g, i);
+                }
             }
         }
 
@@ -654,25 +638,27 @@ Select::writeOutput()
         }
         else if (pdbAtoms_ == "maxsel")
         {
-            ConstArrayRef<int> atomIndices = sel_[0].atomIndices();
-            t_trxstatus       *status      = open_trx(fnPDB_.c_str(), "w");
-            write_trxframe_indexed(status, &fr, atomIndices.size(),
-                                   atomIndices.data(), NULL);
+            std::set<int> atomIndicesSet;
+            for (size_t g = 0; g < sel_.size(); ++g)
+            {
+                ConstArrayRef<int> atomIndices = sel_[g].atomIndices();
+                atomIndicesSet.insert(atomIndices.begin(), atomIndices.end());
+            }
+            std::vector<int>  allAtomIndices(atomIndicesSet.begin(),
+                                             atomIndicesSet.end());
+            t_trxstatus      *status = open_trx(fnPDB_.c_str(), "w");
+            write_trxframe_indexed(status, &fr, allAtomIndices.size(),
+                                   &allAtomIndices[0], NULL);
             close_trx(status);
         }
         else if (pdbAtoms_ == "selected")
         {
             std::vector<int> indices;
-            for (int i = 0; i < sel_[0].posCount(); ++i)
+            for (int i = 0; i < atoms.nr; ++i)
             {
-                if (occupancyModule_->average(0, i) > 0)
+                if (pdbinfo[i].occup > 0.0)
                 {
-                    ConstArrayRef<int>                 atomIndices = sel_[0].position(i).atomIndices();
-                    ConstArrayRef<int>::const_iterator ai;
-                    for (ai = atomIndices.begin(); ai != atomIndices.end(); ++ai)
-                    {
-                        indices.push_back(*ai);
-                    }
+                    indices.push_back(i);
                 }
             }
             t_trxstatus *status = open_trx(fnPDB_.c_str(), "w");
index 22210a33e7488de8d006a4207dc9a6d7503acb76..b25b476c2b1646dc261c81c2936005cc679e64c1 100644 (file)
@@ -93,7 +93,6 @@ class Select : public TrajectoryAnalysisModule
         std::string                         fnOccupancy_;
         std::string                         fnPDB_;
         std::string                         fnLifetime_;
-        bool                                bDump_;
         bool                                bTotNorm_;
         bool                                bFracNorm_;
         bool                                bResInd_;
index b6b7783e4e10a8d3d5968003f5fcb0b3211d3a5e..66f2641566d077cd5149b09f9a01c9f8c1dd6cc5 100644 (file)
@@ -161,4 +161,15 @@ TEST_F(AngleModuleTest, ComputesMultipleAngles)
     runTest(CommandLine::create(cmdline));
 }
 
+TEST_F(AngleModuleTest, HandlesDynamicSelections)
+{
+    const char *const cmdline[] = {
+        "angle",
+        "-g1", "angle", "-group1", "resname RA1 RA2 and name A1 A2 A3 and z < 0.5",
+        "-binw", "60"
+    };
+    setTopology("angle.gro");
+    runTest(CommandLine::create(cmdline));
+}
+
 } // namespace
index 0555114ddff8970a30bc7927ca0dcdd0299ed2db..2ddd62919121b02ed230fb4ac658212d616dec77 100644 (file)
@@ -83,4 +83,15 @@ TEST_F(DistanceModuleTest, ComputesMultipleDistances)
     runTest(CommandLine::create(cmdline));
 }
 
+TEST_F(DistanceModuleTest, HandlesDynamicSelections)
+{
+    const char *const cmdline[] = {
+        "distance",
+        "-select", "atomname S1 S2 and res_cog x < 2.8",
+        "-len", "2", "-binw", "0.5"
+    };
+    setTopology("simple.gro");
+    runTest(CommandLine::create(cmdline));
+}
+
 } // namespace
diff --git a/src/gromacs/trajectoryanalysis/tests/refdata/AngleModuleTest_HandlesDynamicSelections.xml b/src/gromacs/trajectoryanalysis/tests/refdata/AngleModuleTest_HandlesDynamicSelections.xml
new file mode 100644 (file)
index 0000000..bbfadfd
--- /dev/null
@@ -0,0 +1,65 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <String Name="CommandLine">angle -g1 angle -group1 'resname RA1 RA2 and name A1 A2 A3 and z &lt; 0.5' -binw 60</String>
+  <OutputData Name="Data">
+    <AnalysisData Name="angle">
+      <DataFrame Name="Frame0">
+        <Real Name="X">0.000000</Real>
+        <DataValues>
+          <Int Name="Count">2</Int>
+          <DataValue>
+            <Real Name="Value">135.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">45.000000</Real>
+            <Bool Name="Present">false</Bool>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+    </AnalysisData>
+    <AnalysisData Name="average">
+      <DataFrame Name="Frame0">
+        <Real Name="X">0.000000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">135.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+    </AnalysisData>
+    <AnalysisData Name="histogram">
+      <DataFrame Name="Frame0">
+        <Real Name="X">30.000000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame1">
+        <Real Name="X">90.000000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame2">
+        <Real Name="X">150.000000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.016667</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+    </AnalysisData>
+  </OutputData>
+</ReferenceData>
diff --git a/src/gromacs/trajectoryanalysis/tests/refdata/DistanceModuleTest_HandlesDynamicSelections.xml b/src/gromacs/trajectoryanalysis/tests/refdata/DistanceModuleTest_HandlesDynamicSelections.xml
new file mode 100644 (file)
index 0000000..6c879f4
--- /dev/null
@@ -0,0 +1,242 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <String Name="CommandLine">distance -select 'atomname S1 S2 and res_cog x &lt; 2.8' -len 2 -binw 0.5</String>
+  <OutputData Name="Data">
+    <AnalysisData Name="allstats">
+      <DataFrame Name="Frame0">
+        <Real Name="X">0.000000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame1">
+        <Real Name="X">1.000000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame2">
+        <Real Name="X">2.000000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">3.162278</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame3">
+        <Real Name="X">3.000000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame4">
+        <Real Name="X">4.000000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+    </AnalysisData>
+    <AnalysisData Name="average">
+      <DataFrame Name="Frame0">
+        <Real Name="X">0.000000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">1.720759</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+    </AnalysisData>
+    <AnalysisData Name="dist">
+      <DataFrame Name="Frame0">
+        <Real Name="X">0.000000</Real>
+        <DataValues>
+          <Int Name="Count">5</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">3.162278</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Bool Name="Present">false</Bool>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Bool Name="Present">false</Bool>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+    </AnalysisData>
+    <AnalysisData Name="histogram">
+      <DataFrame Name="Frame0">
+        <Real Name="X">0.250000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame1">
+        <Real Name="X">0.750000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame2">
+        <Real Name="X">1.250000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">1.333333</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame3">
+        <Real Name="X">1.750000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame4">
+        <Real Name="X">2.250000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame5">
+        <Real Name="X">2.750000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame6">
+        <Real Name="X">3.250000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.666667</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+      <DataFrame Name="Frame7">
+        <Real Name="X">3.750000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+    </AnalysisData>
+    <AnalysisData Name="stats">
+      <DataFrame Name="Frame0">
+        <Real Name="X">0.000000</Real>
+        <DataValues>
+          <Int Name="Count">1</Int>
+          <DataValue>
+            <Real Name="Value">1.720759</Real>
+            <Real Name="Error">1.248392</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+    </AnalysisData>
+    <AnalysisData Name="xyz">
+      <DataFrame Name="Frame0">
+        <Real Name="X">0.000000</Real>
+        <DataValues>
+          <Int Name="Count">15</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">-3.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+          </DataValue>
+        </DataValues>
+      </DataFrame>
+    </AnalysisData>
+  </OutputData>
+</ReferenceData>
index 617e4789802d820b23a18a54eebd9f13436ba68d..fa89932026cd674cb1c49e1267787cee142af882 100644 (file)
       <DataFrame Name="Frame0">
         <Real Name="X">0.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">8.000000</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">6.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame1">
         <Real Name="X">0.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">3.000000</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">6.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
     </AnalysisData>
         <Real Name="X">0.000000</Real>
         <DataValues>
           <Int Name="Count">15</Int>
+          <Int Name="DataSet">0</Int>
           <DataValue>
             <Real Name="Value">1.000000</Real>
           </DataValue>
             <Real Name="Value">0.000000</Real>
           </DataValue>
         </DataValues>
+        <DataValues>
+          <Int Name="Count">6</Int>
+          <Int Name="DataSet">1</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+        </DataValues>
       </DataFrame>
       <DataFrame Name="Frame1">
         <Real Name="X">0.000000</Real>
         <DataValues>
           <Int Name="Count">15</Int>
+          <Int Name="DataSet">0</Int>
           <DataValue>
             <Real Name="Value">1.000000</Real>
           </DataValue>
             <Real Name="Value">0.000000</Real>
           </DataValue>
         </DataValues>
+        <DataValues>
+          <Int Name="Count">6</Int>
+          <Int Name="DataSet">1</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+          </DataValue>
+        </DataValues>
       </DataFrame>
     </AnalysisData>
     <AnalysisData Name="occupancy">
       <DataFrame Name="Frame0">
         <Real Name="X">1.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
           <DataValue>
             <Real Name="Value">1.000000</Real>
             <Real Name="Error">0.000000</Real>
       <DataFrame Name="Frame1">
         <Real Name="X">2.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame2">
         <Real Name="X">3.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame3">
         <Real Name="X">4.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame4">
         <Real Name="X">5.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame5">
         <Real Name="X">6.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame6">
         <Real Name="X">7.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+            <Bool Name="Present">false</Bool>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame7">
         <Real Name="X">8.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+            <Bool Name="Present">false</Bool>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame8">
         <Real Name="X">9.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">1.000000</Real>
             <Real Name="Error">0.000000</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+            <Bool Name="Present">false</Bool>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame9">
         <Real Name="X">10.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+            <Bool Name="Present">false</Bool>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame10">
         <Real Name="X">11.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
           <DataValue>
             <Real Name="Value">0.000000</Real>
             <Real Name="Error">0.000000</Real>
+            <Bool Name="Present">false</Bool>
           </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame11">
         <Real Name="X">12.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+            <Bool Name="Present">false</Bool>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame12">
         <Real Name="X">13.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">1.000000</Real>
             <Real Name="Error">0.000000</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+            <Bool Name="Present">false</Bool>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame13">
         <Real Name="X">14.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+            <Bool Name="Present">false</Bool>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame14">
         <Real Name="X">15.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
+          <DataValue>
+            <Real Name="Value">0.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
           <DataValue>
             <Real Name="Value">0.000000</Real>
             <Real Name="Error">0.000000</Real>
+            <Bool Name="Present">false</Bool>
           </DataValue>
         </DataValues>
       </DataFrame>
diff --git a/src/gromacs/trajectoryanalysis/tests/refdata/SelectModuleTest_HandlesDumpOption.xml b/src/gromacs/trajectoryanalysis/tests/refdata/SelectModuleTest_HandlesDumpOption.xml
deleted file mode 100644 (file)
index 761e20d..0000000
+++ /dev/null
@@ -1,89 +0,0 @@
-<?xml version="1.0"?>
-<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
-<ReferenceData>
-  <String Name="CommandLine">select -select 'y &lt; 2.5' -dump</String>
-  <OutputData Name="Data">
-    <AnalysisData Name="index">
-      <DataFrame Name="Frame0">
-        <Real Name="X">0.000000</Real>
-        <DataValues>
-          <Int Name="Count">1</Int>
-          <Int Name="FirstColumn">0</Int>
-          <Int Name="LastColumn">0</Int>
-          <DataValue>
-            <Real Name="Value">8.000000</Real>
-          </DataValue>
-        </DataValues>
-        <DataValues>
-          <Int Name="Count">1</Int>
-          <Int Name="FirstColumn">1</Int>
-          <Int Name="LastColumn">1</Int>
-          <DataValue>
-            <Real Name="Value">1.000000</Real>
-          </DataValue>
-        </DataValues>
-        <DataValues>
-          <Int Name="Count">1</Int>
-          <Int Name="FirstColumn">1</Int>
-          <Int Name="LastColumn">1</Int>
-          <DataValue>
-            <Real Name="Value">2.000000</Real>
-          </DataValue>
-        </DataValues>
-        <DataValues>
-          <Int Name="Count">1</Int>
-          <Int Name="FirstColumn">1</Int>
-          <Int Name="LastColumn">1</Int>
-          <DataValue>
-            <Real Name="Value">5.000000</Real>
-          </DataValue>
-        </DataValues>
-        <DataValues>
-          <Int Name="Count">1</Int>
-          <Int Name="FirstColumn">1</Int>
-          <Int Name="LastColumn">1</Int>
-          <DataValue>
-            <Real Name="Value">6.000000</Real>
-          </DataValue>
-        </DataValues>
-        <DataValues>
-          <Int Name="Count">1</Int>
-          <Int Name="FirstColumn">1</Int>
-          <Int Name="LastColumn">1</Int>
-          <DataValue>
-            <Real Name="Value">9.000000</Real>
-          </DataValue>
-        </DataValues>
-        <DataValues>
-          <Int Name="Count">1</Int>
-          <Int Name="FirstColumn">1</Int>
-          <Int Name="LastColumn">1</Int>
-          <DataValue>
-            <Real Name="Value">10.000000</Real>
-          </DataValue>
-        </DataValues>
-        <DataValues>
-          <Int Name="Count">1</Int>
-          <Int Name="FirstColumn">1</Int>
-          <Int Name="LastColumn">1</Int>
-          <DataValue>
-            <Real Name="Value">13.000000</Real>
-          </DataValue>
-        </DataValues>
-        <DataValues>
-          <Int Name="Count">1</Int>
-          <Int Name="FirstColumn">1</Int>
-          <Int Name="LastColumn">1</Int>
-          <DataValue>
-            <Real Name="Value">14.000000</Real>
-          </DataValue>
-        </DataValues>
-      </DataFrame>
-    </AnalysisData>
-  </OutputData>
-  <OutputFiles Name="Files">
-    <String Name="-oi"><![CDATA[
-    1    2    5    6    9   10   13   14
-]]></String>
-  </OutputFiles>
-</ReferenceData>
index f4bf8cf2c7ce4a7f8d256ef25fb15e0ae6590e26..15896efd25d4c5d4216fec3239a2930438e39fc5 100644 (file)
@@ -1,13 +1,17 @@
 <?xml version="1.0"?>
 <?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
 <ReferenceData>
-  <String Name="CommandLine">select -select 'resname RA RD and y &lt; 2.5' -pdbatoms maxsel</String>
+  <String Name="CommandLine">select -select 'resname RA RD and y &lt; 2.5' 'resname RA RB' -pdbatoms maxsel</String>
   <OutputData Name="Data">
     <AnalysisData Name="occupancy">
       <DataFrame Name="Frame0">
         <Real Name="X">1.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
           <DataValue>
             <Real Name="Value">1.000000</Real>
             <Real Name="Error">0.000000</Real>
       <DataFrame Name="Frame1">
         <Real Name="X">2.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame2">
         <Real Name="X">3.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame3">
         <Real Name="X">4.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame4">
         <Real Name="X">5.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame5">
         <Real Name="X">6.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
           <DataValue>
             <Real Name="Value">1.000000</Real>
             <Real Name="Error">0.000000</Real>
       <DataFrame Name="Frame6">
         <Real Name="X">7.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
           <DataValue>
             <Real Name="Value">1.000000</Real>
             <Real Name="Error">0.000000</Real>
       <DataFrame Name="Frame7">
         <Real Name="X">8.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame8">
         <Real Name="X">9.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.000000</Real>
             <Real Name="Error">0.000000</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
     </AnalysisData>
@@ -102,12 +138,15 @@ TITLE     frame t= 0.000
 REMARK    THIS IS A SIMULATION BOX
 CRYST1  100.000  100.000  100.000  90.00  90.00  90.00 P 1           1
 MODEL        1
-ATOM      1  CB   RA A   1      10.000  10.000   0.000  1.00  0.00            
-ATOM      2  S1   RA A   1      10.000  20.000   0.000  0.50  0.00            
-ATOM      3  S2   RA A   1      10.000  30.000   0.000  0.50  0.00            
-ATOM      7  CB   RA B   1      20.000  30.000   0.000  0.50  0.00            
-ATOM      8  S1   RA B   1      20.000  40.000   0.000  0.50  0.30            
-ATOM      9  S2   RA B   1      30.000  10.000   0.000  1.00  0.30            
+ATOM      1  CB   RA A   1      10.000  10.000   0.000  2.00  0.00            
+ATOM      2  S1   RA A   1      10.000  20.000   0.000  1.50  0.00            
+ATOM      3  S2   RA A   1      10.000  30.000   0.000  1.50  0.00            
+ATOM      4  CB   RB A   2      10.000  40.000   0.000  1.00  0.50            
+ATOM      5  S1   RB A   2      20.000  10.000   0.000  1.00  0.50            
+ATOM      6  S2   RB A   2      20.000  20.000   0.000  1.00  0.50            
+ATOM      7  CB   RA B   1      20.000  30.000   0.000  1.50  0.00            
+ATOM      8  S1   RA B   1      20.000  40.000   0.000  1.50  0.30            
+ATOM      9  S2   RA B   1      30.000  10.000   0.000  2.00  0.30            
 ATOM     13  CB   RD B   2      40.000  10.000   0.000  1.00  0.00            
 ATOM     14  S1   RD B   2      40.000  20.000   0.000  0.50  0.00            
 ATOM     15  S2   RD B   2      40.000  30.000   0.000  0.00  0.00            
index 030bd0d310f33a42bdecdc650a0d3914851429aa..c14c5f2dd3c12ab60d41c520e45d3703c4519bdd 100644 (file)
@@ -1,13 +1,17 @@
 <?xml version="1.0"?>
 <?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
 <ReferenceData>
-  <String Name="CommandLine">select -select 'resname RA RD and y &lt; 2.5' -pdbatoms selected</String>
+  <String Name="CommandLine">select -select 'resname RA RD and y &lt; 2.5' 'resname RA RB' -pdbatoms selected</String>
   <OutputData Name="Data">
     <AnalysisData Name="occupancy">
       <DataFrame Name="Frame0">
         <Real Name="X">1.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
           <DataValue>
             <Real Name="Value">1.000000</Real>
             <Real Name="Error">0.000000</Real>
       <DataFrame Name="Frame1">
         <Real Name="X">2.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame2">
         <Real Name="X">3.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame3">
         <Real Name="X">4.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame4">
         <Real Name="X">5.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame5">
         <Real Name="X">6.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
           <DataValue>
             <Real Name="Value">1.000000</Real>
             <Real Name="Error">0.000000</Real>
       <DataFrame Name="Frame6">
         <Real Name="X">7.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
           <DataValue>
             <Real Name="Value">1.000000</Real>
             <Real Name="Error">0.000000</Real>
       <DataFrame Name="Frame7">
         <Real Name="X">8.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.500000</Real>
             <Real Name="Error">0.707107</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
       <DataFrame Name="Frame8">
         <Real Name="X">9.000000</Real>
         <DataValues>
-          <Int Name="Count">1</Int>
+          <Int Name="Count">2</Int>
           <DataValue>
             <Real Name="Value">0.000000</Real>
             <Real Name="Error">0.000000</Real>
           </DataValue>
+          <DataValue>
+            <Real Name="Value">1.000000</Real>
+            <Real Name="Error">0.000000</Real>
+          </DataValue>
         </DataValues>
       </DataFrame>
     </AnalysisData>
@@ -102,12 +138,15 @@ TITLE     frame t= 0.000
 REMARK    THIS IS A SIMULATION BOX
 CRYST1  100.000  100.000  100.000  90.00  90.00  90.00 P 1           1
 MODEL        1
-ATOM      1  CB   RA A   1      10.000  10.000   0.000  1.00  0.00            
-ATOM      2  S1   RA A   1      10.000  20.000   0.000  0.50  0.00            
-ATOM      3  S2   RA A   1      10.000  30.000   0.000  0.50  0.00            
-ATOM      7  CB   RA B   1      20.000  30.000   0.000  0.50  0.00            
-ATOM      8  S1   RA B   1      20.000  40.000   0.000  0.50  0.30            
-ATOM      9  S2   RA B   1      30.000  10.000   0.000  1.00  0.30            
+ATOM      1  CB   RA A   1      10.000  10.000   0.000  2.00  0.00            
+ATOM      2  S1   RA A   1      10.000  20.000   0.000  1.50  0.00            
+ATOM      3  S2   RA A   1      10.000  30.000   0.000  1.50  0.00            
+ATOM      4  CB   RB A   2      10.000  40.000   0.000  1.00  0.50            
+ATOM      5  S1   RB A   2      20.000  10.000   0.000  1.00  0.50            
+ATOM      6  S2   RB A   2      20.000  20.000   0.000  1.00  0.50            
+ATOM      7  CB   RA B   1      20.000  30.000   0.000  1.50  0.00            
+ATOM      8  S1   RA B   1      20.000  40.000   0.000  1.50  0.30            
+ATOM      9  S2   RA B   1      30.000  10.000   0.000  2.00  0.30            
 ATOM     13  CB   RD B   2      40.000  10.000   0.000  1.00  0.00            
 ATOM     14  S1   RD B   2      40.000  20.000   0.000  0.50  0.00            
 TER
index 9c3c3ae551285525b5ed570c6aa282dc2c600a22..8adc76abcedd2d4a7b122f35ddf968223f316275 100644 (file)
@@ -112,7 +112,7 @@ TEST_F(SelectModuleTest, HandlesMaxPDBOutput)
 {
     const char *const cmdline[] = {
         "select",
-        "-select", "resname RA RD and y < 2.5",
+        "-select", "resname RA RD and y < 2.5", "resname RA RB",
         "-pdbatoms", "maxsel"
     };
     setTopology("simple.pdb");
@@ -126,7 +126,7 @@ TEST_F(SelectModuleTest, HandlesSelectedPDBOutput)
 {
     const char *const cmdline[] = {
         "select",
-        "-select", "resname RA RD and y < 2.5",
+        "-select", "resname RA RD and y < 2.5", "resname RA RB",
         "-pdbatoms", "selected"
     };
     setTopology("simple.pdb");
@@ -136,19 +136,6 @@ TEST_F(SelectModuleTest, HandlesSelectedPDBOutput)
     runTest(CommandLine::create(cmdline));
 }
 
-TEST_F(SelectModuleTest, HandlesDumpOption)
-{
-    const char *const cmdline[] = {
-        "select",
-        "-select", "y < 2.5",
-        "-dump"
-    };
-    setTopology("simple.gro");
-    setOutputFile("-oi", "index.dat");
-    includeDataset("index");
-    runTest(CommandLine::create(cmdline));
-}
-
 TEST_F(SelectModuleTest, NormalizesSizes)
 {
     const char *const cmdline[] = {