Tidy: modernize-use-nullptr
[alexxy/gromacs.git] / src / gromacs / selection / tests / selectioncollection.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  * \brief
37  * Tests selection parsing and compilation.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_selection
41  */
42 #include "gmxpre.h"
43
44 #include "gromacs/selection/selectioncollection.h"
45
46 #include <gtest/gtest.h>
47
48 #include "gromacs/options/basicoptions.h"
49 #include "gromacs/options/ioptionscontainer.h"
50 #include "gromacs/selection/indexutil.h"
51 #include "gromacs/selection/selection.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/trajectory/trajectoryframe.h"
54 #include "gromacs/utility/arrayref.h"
55 #include "gromacs/utility/exceptions.h"
56 #include "gromacs/utility/flags.h"
57 #include "gromacs/utility/gmxregex.h"
58 #include "gromacs/utility/stringutil.h"
59
60 #include "testutils/interactivetest.h"
61 #include "testutils/refdata.h"
62 #include "testutils/testasserts.h"
63 #include "testutils/testfilemanager.h"
64 #include "testutils/testoptions.h"
65
66 #include "toputils.h"
67
68 namespace
69 {
70
71 /********************************************************************
72  * Test fixture for selection testing
73  */
74
75 class SelectionCollectionTest : public ::testing::Test
76 {
77     public:
78         static int               s_debugLevel;
79
80         SelectionCollectionTest();
81         ~SelectionCollectionTest();
82
83         void setAtomCount(int natoms)
84         {
85             ASSERT_NO_THROW_GMX(sc_.setTopology(nullptr, natoms));
86         }
87         void loadTopology(const char *filename);
88         void setTopology();
89         void loadIndexGroups(const char *filename);
90
91         gmx::test::TopologyManager  topManager_;
92         gmx::SelectionCollection    sc_;
93         gmx::SelectionList          sel_;
94         gmx_ana_indexgrps_t        *grps_;
95 };
96
97 int SelectionCollectionTest::s_debugLevel = 0;
98
99 // cond/endcond do not seem to work here with Doxygen 1.8.5 parser.
100 #ifndef DOXYGEN
101 GMX_TEST_OPTIONS(SelectionCollectionTestOptions, options)
102 {
103     options->addOption(gmx::IntegerOption("seldebug")
104                            .store(&SelectionCollectionTest::s_debugLevel)
105                            .description("Set selection debug level"));
106 }
107 #endif
108
109 SelectionCollectionTest::SelectionCollectionTest()
110     : grps_(nullptr)
111 {
112     topManager_.requestFrame();
113     sc_.setDebugLevel(s_debugLevel);
114     sc_.setReferencePosType("atom");
115     sc_.setOutputPosType("atom");
116 }
117
118 SelectionCollectionTest::~SelectionCollectionTest()
119 {
120     if (grps_ != nullptr)
121     {
122         gmx_ana_indexgrps_free(grps_);
123     }
124 }
125
126 void
127 SelectionCollectionTest::loadTopology(const char *filename)
128 {
129     topManager_.loadTopology(filename);
130     setTopology();
131 }
132
133 void
134 SelectionCollectionTest::setTopology()
135 {
136     ASSERT_NO_THROW_GMX(sc_.setTopology(topManager_.topology(), -1));
137 }
138
139 void
140 SelectionCollectionTest::loadIndexGroups(const char *filename)
141 {
142     GMX_RELEASE_ASSERT(grps_ == nullptr,
143                        "External groups can only be loaded once");
144     std::string fullpath =
145         gmx::test::TestFileManager::getInputFilePath(filename);
146     gmx_ana_indexgrps_init(&grps_, nullptr, fullpath.c_str());
147     sc_.setIndexGroups(grps_);
148 }
149
150
151 /********************************************************************
152  * Test fixture for interactive SelectionCollection tests
153  */
154
155 class SelectionCollectionInteractiveTest : public SelectionCollectionTest
156 {
157     public:
158         SelectionCollectionInteractiveTest()
159             : helper_(data_.rootChecker())
160         {
161         }
162
163         void runTest(int count, bool bInteractive,
164                      const gmx::ConstArrayRef<const char *> &input);
165
166         gmx::test::TestReferenceData      data_;
167         gmx::test::InteractiveTestHelper  helper_;
168 };
169
170 void SelectionCollectionInteractiveTest::runTest(
171         int count, bool bInteractive,
172         const gmx::ConstArrayRef<const char *> &inputLines)
173 {
174     helper_.setInputLines(inputLines);
175     // TODO: Check something about the returned selections as well.
176     ASSERT_NO_THROW_GMX(sc_.parseInteractive(
177                                 count, &helper_.inputStream(),
178                                 bInteractive ? &helper_.outputStream() : nullptr,
179                                 "for test context"));
180     helper_.checkSession();
181 }
182
183
184 /********************************************************************
185  * Test fixture for selection testing with reference data
186  */
187
188 class SelectionCollectionDataTest : public SelectionCollectionTest
189 {
190     public:
191         enum TestFlag
192         {
193             efTestEvaluation            = 1<<0,
194             efTestPositionAtoms         = 1<<1,
195             efTestPositionCoordinates   = 1<<2,
196             efTestPositionMapping       = 1<<3,
197             efTestPositionMasses        = 1<<4,
198             efTestPositionCharges       = 1<<5,
199             efTestSelectionNames        = 1<<6,
200             efDontTestCompiledAtoms     = 1<<8
201         };
202         typedef gmx::FlagsTemplate<TestFlag> TestFlags;
203
204         SelectionCollectionDataTest()
205             : checker_(data_.rootChecker()), count_(0), framenr_(0)
206         {
207         }
208
209         void setFlags(TestFlags flags) { flags_ = flags; }
210
211         void runParser(const gmx::ConstArrayRef<const char *> &selections);
212         void runCompiler();
213         void runEvaluate();
214         void runEvaluateFinal();
215
216         void runTest(int                                     natoms,
217                      const gmx::ConstArrayRef<const char *> &selections);
218         void runTest(const char                             *filename,
219                      const gmx::ConstArrayRef<const char *> &selections);
220
221     private:
222         static void checkSelection(gmx::test::TestReferenceChecker *checker,
223                                    const gmx::Selection &sel, TestFlags flags);
224
225         void checkCompiled();
226
227         gmx::test::TestReferenceData    data_;
228         gmx::test::TestReferenceChecker checker_;
229         size_t                          count_;
230         int                             framenr_;
231         TestFlags                       flags_;
232 };
233
234
235 void
236 SelectionCollectionDataTest::checkSelection(
237         gmx::test::TestReferenceChecker *checker,
238         const gmx::Selection &sel, TestFlags flags)
239 {
240     using gmx::test::TestReferenceChecker;
241
242     {
243         gmx::ConstArrayRef<int> atoms = sel.atomIndices();
244         checker->checkSequence(atoms.begin(), atoms.end(), "Atoms");
245     }
246     if (flags.test(efTestPositionAtoms)
247         || flags.test(efTestPositionCoordinates)
248         || flags.test(efTestPositionMapping)
249         || flags.test(efTestPositionMasses)
250         || flags.test(efTestPositionCharges))
251     {
252         TestReferenceChecker compound(
253                 checker->checkSequenceCompound("Positions", sel.posCount()));
254         for (int i = 0; i < sel.posCount(); ++i)
255         {
256             TestReferenceChecker          poscompound(compound.checkCompound("Position", nullptr));
257             const gmx::SelectionPosition &p = sel.position(i);
258             if (flags.test(efTestPositionAtoms))
259             {
260                 gmx::ConstArrayRef<int> atoms = p.atomIndices();
261                 poscompound.checkSequence(atoms.begin(), atoms.end(), "Atoms");
262             }
263             if (flags.test(efTestPositionCoordinates))
264             {
265                 poscompound.checkVector(p.x(), "Coordinates");
266             }
267             if (flags.test(efTestPositionMapping))
268             {
269                 poscompound.checkInteger(p.refId(), "RefId");
270                 poscompound.checkInteger(p.mappedId(), "MappedId");
271             }
272             if (flags.test(efTestPositionMasses))
273             {
274                 poscompound.checkReal(p.mass(), "Mass");
275             }
276             if (flags.test(efTestPositionCharges))
277             {
278                 poscompound.checkReal(p.charge(), "Charge");
279             }
280         }
281     }
282 }
283
284
285 void
286 SelectionCollectionDataTest::runParser(
287         const gmx::ConstArrayRef<const char *> &selections)
288 {
289     using gmx::test::TestReferenceChecker;
290
291     TestReferenceChecker compound(checker_.checkCompound("ParsedSelections", "Parsed"));
292     size_t               varcount = 0;
293     count_ = 0;
294     for (size_t i = 0; i < selections.size(); ++i)
295     {
296         SCOPED_TRACE(std::string("Parsing selection \"")
297                      + selections[i] + "\"");
298         gmx::SelectionList result;
299         ASSERT_NO_THROW_GMX(result = sc_.parseFromString(selections[i]));
300         sel_.insert(sel_.end(), result.begin(), result.end());
301         if (sel_.size() == count_)
302         {
303             std::string          id = gmx::formatString("Variable%d", static_cast<int>(varcount + 1));
304             TestReferenceChecker varcompound(
305                     compound.checkCompound("ParsedVariable", id.c_str()));
306             varcompound.checkString(selections[i], "Input");
307             ++varcount;
308         }
309         else
310         {
311             std::string          id = gmx::formatString("Selection%d", static_cast<int>(count_ + 1));
312             TestReferenceChecker selcompound(
313                     compound.checkCompound("ParsedSelection", id.c_str()));
314             selcompound.checkString(selections[i], "Input");
315             if (flags_.test(efTestSelectionNames))
316             {
317                 selcompound.checkString(sel_[count_].name(), "Name");
318             }
319             selcompound.checkString(sel_[count_].selectionText(), "Text");
320             selcompound.checkBoolean(sel_[count_].isDynamic(), "Dynamic");
321             ++count_;
322         }
323     }
324 }
325
326
327 void
328 SelectionCollectionDataTest::runCompiler()
329 {
330     ASSERT_NO_THROW_GMX(sc_.compile());
331     ASSERT_EQ(count_, sel_.size());
332     checkCompiled();
333 }
334
335
336 void
337 SelectionCollectionDataTest::checkCompiled()
338 {
339     using gmx::test::TestReferenceChecker;
340     const TestFlags      mask = ~TestFlags(efTestPositionCoordinates);
341
342     TestReferenceChecker compound(checker_.checkCompound("CompiledSelections", "Compiled"));
343     for (size_t i = 0; i < count_; ++i)
344     {
345         SCOPED_TRACE(std::string("Checking selection \"") +
346                      sel_[i].selectionText() + "\"");
347         std::string          id = gmx::formatString("Selection%d", static_cast<int>(i + 1));
348         TestReferenceChecker selcompound(
349                 compound.checkCompound("Selection", id.c_str()));
350         if (flags_.test(efTestSelectionNames))
351         {
352             selcompound.checkString(sel_[i].name(), "Name");
353         }
354         if (!flags_.test(efDontTestCompiledAtoms))
355         {
356             checkSelection(&selcompound, sel_[i], flags_ & mask);
357         }
358     }
359 }
360
361
362 void
363 SelectionCollectionDataTest::runEvaluate()
364 {
365     using gmx::test::TestReferenceChecker;
366
367     ++framenr_;
368     ASSERT_NO_THROW_GMX(sc_.evaluate(topManager_.frame(), nullptr));
369     std::string          frame = gmx::formatString("Frame%d", framenr_);
370     TestReferenceChecker compound(
371             checker_.checkCompound("EvaluatedSelections", frame.c_str()));
372     for (size_t i = 0; i < count_; ++i)
373     {
374         SCOPED_TRACE(std::string("Checking selection \"") +
375                      sel_[i].selectionText() + "\"");
376         std::string          id = gmx::formatString("Selection%d", static_cast<int>(i + 1));
377         TestReferenceChecker selcompound(
378                 compound.checkCompound("Selection", id.c_str()));
379         checkSelection(&selcompound, sel_[i], flags_);
380     }
381 }
382
383
384 void
385 SelectionCollectionDataTest::runEvaluateFinal()
386 {
387     ASSERT_NO_THROW_GMX(sc_.evaluateFinal(framenr_));
388     checkCompiled();
389 }
390
391
392 void
393 SelectionCollectionDataTest::runTest(
394         int natoms, const gmx::ConstArrayRef<const char *> &selections)
395 {
396     ASSERT_NO_FATAL_FAILURE(runParser(selections));
397     ASSERT_NO_FATAL_FAILURE(setAtomCount(natoms));
398     ASSERT_NO_FATAL_FAILURE(runCompiler());
399 }
400
401
402 void
403 SelectionCollectionDataTest::runTest(
404         const char *filename, const gmx::ConstArrayRef<const char *> &selections)
405 {
406     ASSERT_NO_FATAL_FAILURE(runParser(selections));
407     ASSERT_NO_FATAL_FAILURE(loadTopology(filename));
408     ASSERT_NO_FATAL_FAILURE(runCompiler());
409     if (flags_.test(efTestEvaluation))
410     {
411         ASSERT_NO_FATAL_FAILURE(runEvaluate());
412         ASSERT_NO_FATAL_FAILURE(runEvaluateFinal());
413     }
414 }
415
416
417 /********************************************************************
418  * Tests for SelectionCollection functionality without reference data
419  */
420
421 TEST_F(SelectionCollectionTest, HandlesNoSelections)
422 {
423     EXPECT_FALSE(sc_.requiredTopologyProperties().hasAny());
424     EXPECT_NO_THROW_GMX(sc_.compile());
425     EXPECT_FALSE(sc_.requiredTopologyProperties().hasAny());
426 }
427
428 TEST_F(SelectionCollectionTest, HandlesNoSelectionsWithDefaultPositionType)
429 {
430     EXPECT_NO_THROW_GMX(sc_.setOutputPosType("res_com"));
431     EXPECT_TRUE(sc_.requiredTopologyProperties().needsTopology);
432     EXPECT_TRUE(sc_.requiredTopologyProperties().needsMasses);
433     EXPECT_NO_THROW_GMX(sc_.setOutputPosType("res_cog"));
434     EXPECT_TRUE(sc_.requiredTopologyProperties().needsTopology);
435     EXPECT_FALSE(sc_.requiredTopologyProperties().needsMasses);
436     ASSERT_NO_THROW_GMX(sc_.parseFromString("atom of atomnr 1 to 10"));
437     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
438     ASSERT_NO_THROW_GMX(sc_.compile());
439     EXPECT_FALSE(sc_.requiredTopologyProperties().hasAny());
440 }
441
442 TEST_F(SelectionCollectionTest, HandlesVelocityAndForceRequests)
443 {
444     ASSERT_NO_THROW_GMX(sel_ = sc_.parseFromString("atomnr 1 to 10; none"));
445     EXPECT_FALSE(sc_.requiredTopologyProperties().hasAny());
446     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
447     ASSERT_EQ(2U, sel_.size());
448     ASSERT_NO_THROW_GMX(sel_[0].setEvaluateVelocities(true));
449     ASSERT_NO_THROW_GMX(sel_[1].setEvaluateVelocities(true));
450     ASSERT_NO_THROW_GMX(sel_[0].setEvaluateForces(true));
451     ASSERT_NO_THROW_GMX(sel_[1].setEvaluateForces(true));
452     EXPECT_FALSE(sc_.requiredTopologyProperties().hasAny());
453     ASSERT_NO_THROW_GMX(sc_.compile());
454     EXPECT_FALSE(sc_.requiredTopologyProperties().hasAny());
455     EXPECT_TRUE(sel_[0].hasVelocities());
456     EXPECT_TRUE(sel_[1].hasVelocities());
457     EXPECT_TRUE(sel_[0].hasForces());
458     EXPECT_TRUE(sel_[1].hasForces());
459 }
460
461 TEST_F(SelectionCollectionTest, HandlesForceRequestForCenterOfGeometry)
462 {
463     ASSERT_NO_THROW_GMX(sel_ = sc_.parseFromString("res_cog of atomnr 1 to 10"));
464     EXPECT_TRUE(sc_.requiredTopologyProperties().needsTopology);
465     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
466     ASSERT_EQ(1U, sel_.size());
467     ASSERT_NO_THROW_GMX(sel_[0].setEvaluateForces(true));
468     // In principle, the code could know here that the masses are required, but
469     // currently it only knows this after compilation.
470     ASSERT_NO_THROW_GMX(sc_.compile());
471     EXPECT_TRUE(sc_.requiredTopologyProperties().needsMasses);
472     EXPECT_TRUE(sel_[0].hasForces());
473 }
474
475 TEST_F(SelectionCollectionTest, ParsesSelectionsFromFile)
476 {
477     ASSERT_NO_THROW_GMX(sel_ = sc_.parseFromFile(
478                                     gmx::test::TestFileManager::getInputFilePath("selfile.dat")));
479     // These should match the contents of selfile.dat
480     ASSERT_EQ(2U, sel_.size());
481     EXPECT_STREQ("resname RA RB", sel_[0].selectionText());
482     EXPECT_STREQ("resname RB RC", sel_[1].selectionText());
483 }
484
485 TEST_F(SelectionCollectionTest, HandlesAtypicalWhitespace)
486 {
487     ASSERT_NO_THROW_GMX(sel_ = sc_.parseFromString("atomnr\n1\r\nto\t10;\vatomnr 3\f to 14\r"));
488     ASSERT_EQ(2U, sel_.size());
489     EXPECT_STREQ("atomnr 1 to 10", sel_[0].selectionText());
490     // TODO: Get rid of the trailing whitespace.
491     EXPECT_STREQ("atomnr 3 to 14 ", sel_[1].selectionText());
492 }
493
494 TEST_F(SelectionCollectionTest, HandlesInvalidRegularExpressions)
495 {
496     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
497     EXPECT_THROW_GMX({
498                          sc_.parseFromString("resname ~ \"R[A\"");
499                          sc_.compile();
500                      }, gmx::InvalidInputError);
501 }
502
503 TEST_F(SelectionCollectionTest, HandlesUnsupportedRegularExpressions)
504 {
505     if (!gmx::Regex::isSupported())
506     {
507         ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
508         EXPECT_THROW_GMX({
509                              sc_.parseFromString("resname \"R[AD]\"");
510                              sc_.compile();
511                          }, gmx::InvalidInputError);
512     }
513 }
514
515 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue)
516 {
517     EXPECT_THROW_GMX(sc_.parseFromString("mindist from atomnr 1 cutoff"),
518                      gmx::InvalidInputError);
519 }
520
521 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue2)
522 {
523     EXPECT_THROW_GMX(sc_.parseFromString("within 1 of"),
524                      gmx::InvalidInputError);
525 }
526
527 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue3)
528 {
529     EXPECT_THROW_GMX(sc_.parseFromString("within of atomnr 1"),
530                      gmx::InvalidInputError);
531 }
532
533 // TODO: Tests for more parser errors
534
535 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceParser1)
536 {
537     ASSERT_NO_THROW_GMX(sc_.setIndexGroups(nullptr));
538     EXPECT_THROW_GMX(sc_.parseFromString("group \"foo\""), gmx::InconsistentInputError);
539     EXPECT_THROW_GMX(sc_.parseFromString("4"), gmx::InconsistentInputError);
540 }
541
542 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceParser2)
543 {
544     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
545     EXPECT_THROW_GMX(sc_.parseFromString("group \"foo\""), gmx::InconsistentInputError);
546     EXPECT_THROW_GMX(sc_.parseFromString("4"), gmx::InconsistentInputError);
547 }
548
549 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceDelayed1)
550 {
551     ASSERT_NO_THROW_GMX(sc_.parseFromString("group \"foo\""));
552     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
553     EXPECT_THROW_GMX(sc_.setIndexGroups(nullptr), gmx::InconsistentInputError);
554     EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
555 }
556
557 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceDelayed2)
558 {
559     ASSERT_NO_THROW_GMX(sc_.parseFromString("group 4; group \"foo\""));
560     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
561     EXPECT_THROW_GMX(loadIndexGroups("simple.ndx"), gmx::InconsistentInputError);
562     EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
563 }
564
565 TEST_F(SelectionCollectionTest, HandlesUnsortedGroupReference)
566 {
567     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
568     EXPECT_THROW_GMX(sc_.parseFromString("atomnr 1 to 3 and group \"GrpUnsorted\""),
569                      gmx::InconsistentInputError);
570     EXPECT_THROW_GMX(sc_.parseFromString("group 2 or atomnr 2 to 5"),
571                      gmx::InconsistentInputError);
572     EXPECT_THROW_GMX(sc_.parseFromString("within 1 of group 2"),
573                      gmx::InconsistentInputError);
574 }
575
576 TEST_F(SelectionCollectionTest, HandlesUnsortedGroupReferenceDelayed)
577 {
578     ASSERT_NO_THROW_GMX(sc_.parseFromString("atomnr 1 to 3 and group \"GrpUnsorted\""));
579     ASSERT_NO_THROW_GMX(sc_.parseFromString("atomnr 1 to 3 and group 2"));
580     EXPECT_THROW_GMX(loadIndexGroups("simple.ndx"), gmx::InconsistentInputError);
581     // TODO: Add a separate check in the selection compiler for a safer API
582     // (makes sense in the future if the compiler needs the information for
583     // other purposes as well).
584     // EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
585 }
586
587 TEST_F(SelectionCollectionTest, HandlesOutOfRangeAtomIndexInGroup)
588 {
589     ASSERT_NO_THROW_GMX(sc_.setTopology(nullptr, 5));
590     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
591     EXPECT_THROW_GMX(sc_.parseFromString("group \"GrpB\""), gmx::InconsistentInputError);
592 }
593
594 TEST_F(SelectionCollectionTest, HandlesOutOfRangeAtomIndexInGroupDelayed)
595 {
596     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
597     ASSERT_NO_THROW_GMX(sc_.parseFromString("group \"GrpB\""));
598     EXPECT_THROW_GMX(sc_.setTopology(nullptr, 5), gmx::InconsistentInputError);
599 }
600
601 TEST_F(SelectionCollectionTest, HandlesOutOfRangeAtomIndexInGroupDelayed2)
602 {
603     ASSERT_NO_THROW_GMX(sc_.setTopology(nullptr, 5));
604     ASSERT_NO_THROW_GMX(sc_.parseFromString("group \"GrpB\""));
605     EXPECT_THROW_GMX(loadIndexGroups("simple.ndx"), gmx::InconsistentInputError);
606 }
607
608 TEST_F(SelectionCollectionTest, RecoversFromMissingMoleculeInfo)
609 {
610     ASSERT_NO_THROW_GMX(sc_.parseFromString("molindex 1 to 5"));
611     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
612     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
613 }
614
615 TEST_F(SelectionCollectionTest, RecoversFromMissingAtomTypes)
616 {
617     ASSERT_NO_THROW_GMX(sc_.parseFromString("type CA"));
618     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
619     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
620 }
621
622 TEST_F(SelectionCollectionTest, RecoversFromMissingPDBInfo)
623 {
624     ASSERT_NO_THROW_GMX(sc_.parseFromString("altloc A"));
625     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
626     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
627 }
628
629 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation)
630 {
631     ASSERT_NO_THROW_GMX(sc_.parseFromString("all permute 1 1"));
632     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
633     EXPECT_THROW_GMX(sc_.compile(), gmx::InvalidInputError);
634 }
635
636 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation2)
637 {
638     ASSERT_NO_THROW_GMX(sc_.parseFromString("all permute 3 2 1"));
639     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
640     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
641 }
642
643 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation3)
644 {
645     ASSERT_NO_THROW_GMX(sc_.parseFromString("x < 1.5 permute 3 2 1"));
646     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
647     ASSERT_NO_THROW_GMX(sc_.compile());
648     EXPECT_THROW_GMX(sc_.evaluate(topManager_.frame(), nullptr), gmx::InconsistentInputError);
649 }
650
651 TEST_F(SelectionCollectionTest, HandlesFramesWithTooSmallAtomSubsets)
652 {
653     ASSERT_NO_THROW_GMX(sc_.parseFromString("atomnr 3 to 10"));
654     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
655     ASSERT_NO_THROW_GMX(sc_.compile());
656     topManager_.frame()->natoms = 8;
657     EXPECT_THROW_GMX(sc_.evaluate(topManager_.frame(), nullptr), gmx::InconsistentInputError);
658 }
659
660 TEST_F(SelectionCollectionTest, HandlesFramesWithTooSmallAtomSubsets2)
661 {
662     const int index[] = { 1, 2, 3, 9 };
663     ASSERT_NO_THROW_GMX(sc_.parseFromString("atomnr 3 4 7 10"));
664     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
665     ASSERT_NO_THROW_GMX(sc_.compile());
666     topManager_.initFrameIndices(index);
667     EXPECT_THROW_GMX(sc_.evaluate(topManager_.frame(), nullptr), gmx::InconsistentInputError);
668 }
669
670 TEST_F(SelectionCollectionTest, HandlesFramesWithTooSmallAtomSubsets3)
671 {
672     const int index[] = { 0, 1, 2, 3, 4, 5, 6, 9, 10, 11 };
673     // Evaluating the positions will require atoms 1-3, 7-12.
674     ASSERT_NO_THROW_GMX(sc_.parseFromString("whole_res_cog of atomnr 2 7 11"));
675     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
676     ASSERT_NO_THROW_GMX(sc_.compile());
677     topManager_.initFrameIndices(index);
678     EXPECT_THROW_GMX(sc_.evaluate(topManager_.frame(), nullptr), gmx::InconsistentInputError);
679 }
680
681 TEST_F(SelectionCollectionTest, HandlesFramesWithTooSmallAtomSubsets4)
682 {
683     ASSERT_NO_THROW_GMX(sc_.parseFromString("mindistance from atomnr 1 to 5 < 2"));
684     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
685     ASSERT_NO_THROW_GMX(sc_.compile());
686     topManager_.frame()->natoms = 10;
687     EXPECT_THROW_GMX(sc_.evaluate(topManager_.frame(), nullptr), gmx::InconsistentInputError);
688 }
689
690 // TODO: Tests for more evaluation errors
691
692 /********************************************************************
693  * Tests for interactive selection input
694  */
695
696 TEST_F(SelectionCollectionInteractiveTest, HandlesBasicInput)
697 {
698     const char *const input[] = {
699         "foo = resname RA",
700         "resname RB",
701         "\"Name\" resname RC"
702     };
703     runTest(-1, true, input);
704 }
705
706 TEST_F(SelectionCollectionInteractiveTest, HandlesContinuation)
707 {
708     const char *const input[] = {
709         "resname RB and \\",
710         "resname RC"
711     };
712     runTest(-1, true, input);
713 }
714
715 TEST_F(SelectionCollectionInteractiveTest, HandlesSingleSelectionInput)
716 {
717     const char *const input[] = {
718         "foo = resname RA",
719         "resname RA"
720     };
721     runTest(1, true, input);
722 }
723
724 TEST_F(SelectionCollectionInteractiveTest, HandlesTwoSelectionInput)
725 {
726     const char *const input[] = {
727         "resname RA",
728         "resname RB"
729     };
730     runTest(2, true, input);
731 }
732
733 TEST_F(SelectionCollectionInteractiveTest, HandlesStatusWithGroups)
734 {
735     const char *const input[] = {
736         "resname RA",
737         ""
738     };
739     loadIndexGroups("simple.ndx");
740     runTest(-1, true, input);
741 }
742
743 TEST_F(SelectionCollectionInteractiveTest, HandlesStatusWithExistingSelections)
744 {
745     const char *const input[] = {
746         "",
747         "bar = resname RC",
748         "resname RA",
749         ""
750     };
751     ASSERT_NO_THROW_GMX(sc_.parseFromString("foo = resname RA"));
752     ASSERT_NO_THROW_GMX(sc_.parseFromString("resname RB"));
753     runTest(-1, true, input);
754 }
755
756 TEST_F(SelectionCollectionInteractiveTest, HandlesSingleSelectionInputStatus)
757 {
758     const char *const input[] = {
759         "foo = resname RA",
760         "",
761         "resname RB"
762     };
763     runTest(1, true, input);
764 }
765
766 TEST_F(SelectionCollectionInteractiveTest, HandlesTwoSelectionInputStatus)
767 {
768     const char *const input[] = {
769         "\"Sel\" resname RA",
770         "",
771         "resname RB"
772     };
773     runTest(2, true, input);
774 }
775
776 TEST_F(SelectionCollectionInteractiveTest, HandlesMultiSelectionInputStatus)
777 {
778     const char *const input[] = {
779         "\"Sel\" resname RA",
780         "\"Sel2\" resname RB",
781         ""
782     };
783     runTest(-1, true, input);
784 }
785
786 TEST_F(SelectionCollectionInteractiveTest, HandlesNoFinalNewline)
787 {
788     // TODO: There is an extra prompt printed after the input is finished; it
789     // would be cleaner not to have it, but it's only a cosmetic issue.
790     const char *const input[] = {
791         "resname RA"
792     };
793     helper_.setLastNewline(false);
794     runTest(-1, true, input);
795 }
796
797 TEST_F(SelectionCollectionInteractiveTest, HandlesEmptySelections)
798 {
799     const char *const input[] = {
800         "resname RA;",
801         "; resname RB;;",
802         " ",
803         ";"
804     };
805     runTest(-1, true, input);
806 }
807
808 TEST_F(SelectionCollectionInteractiveTest, HandlesMultipleSelectionsOnLine)
809 {
810     const char *const input[] = {
811         "resname RA; resname RB and \\",
812         "resname RC"
813     };
814     runTest(2, true, input);
815 }
816
817 TEST_F(SelectionCollectionInteractiveTest, HandlesNoninteractiveInput)
818 {
819     const char *const input[] = {
820         "foo = resname RA",
821         "resname RB",
822         "\"Name\" resname RC"
823     };
824     runTest(-1, false, input);
825 }
826
827 TEST_F(SelectionCollectionInteractiveTest, HandlesSingleSelectionInputNoninteractively)
828 {
829     const char *const input[] = {
830         "foo = resname RA",
831         "resname RA"
832     };
833     runTest(1, false, input);
834 }
835
836
837 /********************************************************************
838  * Tests for selection keywords
839  */
840
841 TEST_F(SelectionCollectionDataTest, HandlesAllNone)
842 {
843     static const char * const selections[] = {
844         "all",
845         "none"
846     };
847     runTest(10, selections);
848 }
849
850 TEST_F(SelectionCollectionDataTest, HandlesAtomnr)
851 {
852     static const char * const selections[] = {
853         "atomnr 1 to 3 6 to 8",
854         "atomnr 4 2 5 to 7",
855         "atomnr <= 5"
856     };
857     runTest(10, selections);
858 }
859
860 TEST_F(SelectionCollectionDataTest, HandlesResnr)
861 {
862     static const char * const selections[] = {
863         "resnr 1 2 5",
864         "resid 4 to 3"
865     };
866     runTest("simple.gro", selections);
867 }
868
869 TEST_F(SelectionCollectionDataTest, HandlesResIndex)
870 {
871     static const char * const selections[] = {
872         "resindex 1 4",
873         "residue 1 3"
874     };
875     runTest("simple.pdb", selections);
876 }
877
878 TEST_F(SelectionCollectionDataTest, HandlesMolIndex)
879 {
880     static const char * const selections[] = {
881         "molindex 1 4",
882         "molecule 2 3 5"
883     };
884     ASSERT_NO_FATAL_FAILURE(runParser(selections));
885     ASSERT_NO_FATAL_FAILURE(topManager_.loadTopology("simple.gro"));
886     topManager_.initUniformMolecules(3);
887     ASSERT_NO_FATAL_FAILURE(setTopology());
888     ASSERT_NO_FATAL_FAILURE(runCompiler());
889 }
890
891 TEST_F(SelectionCollectionDataTest, HandlesAtomname)
892 {
893     static const char * const selections[] = {
894         "name CB",
895         "atomname S1 S2"
896     };
897     runTest("simple.gro", selections);
898 }
899
900 TEST_F(SelectionCollectionDataTest, HandlesPdbAtomname)
901 {
902     static const char * const selections[] = {
903         "name HG21",
904         "name 1HG2",
905         "pdbname HG21 CB",
906         "pdbatomname 1HG2"
907     };
908     runTest("simple.pdb", selections);
909 }
910
911
912 TEST_F(SelectionCollectionDataTest, HandlesAtomtype)
913 {
914     static const char * const selections[] = {
915         "atomtype CA"
916     };
917     ASSERT_NO_FATAL_FAILURE(runParser(selections));
918     ASSERT_NO_FATAL_FAILURE(topManager_.loadTopology("simple.gro"));
919     const char *const types[] = { "CA", "SA", "SB" };
920     topManager_.initAtomTypes(types);
921     ASSERT_NO_FATAL_FAILURE(setTopology());
922     ASSERT_NO_FATAL_FAILURE(runCompiler());
923 }
924
925 TEST_F(SelectionCollectionDataTest, HandlesChain)
926 {
927     static const char * const selections[] = {
928         "chain A",
929         "chain B"
930     };
931     runTest("simple.pdb", selections);
932 }
933
934 TEST_F(SelectionCollectionDataTest, HandlesMass)
935 {
936     static const char * const selections[] = {
937         "mass > 5"
938     };
939     ASSERT_NO_FATAL_FAILURE(runParser(selections));
940     EXPECT_TRUE(sc_.requiredTopologyProperties().needsMasses);
941     ASSERT_NO_FATAL_FAILURE(topManager_.loadTopology("simple.gro"));
942     t_atoms &atoms = topManager_.atoms();
943     for (int i = 0; i < atoms.nr; ++i)
944     {
945         atoms.atom[i].m = 1.0 + i;
946     }
947     atoms.haveMass = TRUE;
948     ASSERT_NO_FATAL_FAILURE(setTopology());
949     ASSERT_NO_FATAL_FAILURE(runCompiler());
950 }
951
952 TEST_F(SelectionCollectionDataTest, HandlesCharge)
953 {
954     static const char * const selections[] = {
955         "charge < 0.5"
956     };
957     ASSERT_NO_FATAL_FAILURE(runParser(selections));
958     ASSERT_NO_FATAL_FAILURE(topManager_.loadTopology("simple.gro"));
959     t_atoms &atoms = topManager_.atoms();
960     for (int i = 0; i < atoms.nr; ++i)
961     {
962         atoms.atom[i].q = i / 10.0;
963     }
964     // Ensure exact representation of 0.5 is used, so that the test is
965     // reproducible.
966     atoms.atom[5].q  = 0.5;
967     atoms.haveCharge = TRUE;
968     ASSERT_NO_FATAL_FAILURE(setTopology());
969     ASSERT_NO_FATAL_FAILURE(runCompiler());
970 }
971
972 TEST_F(SelectionCollectionDataTest, HandlesAltLoc)
973 {
974     static const char * const selections[] = {
975         "altloc \" \"",
976         "altloc A"
977     };
978     runTest("simple.pdb", selections);
979 }
980
981 TEST_F(SelectionCollectionDataTest, HandlesInsertCode)
982 {
983     static const char * const selections[] = {
984         "insertcode \" \"",
985         "insertcode A"
986     };
987     runTest("simple.pdb", selections);
988 }
989
990 TEST_F(SelectionCollectionDataTest, HandlesOccupancy)
991 {
992     static const char * const selections[] = {
993         "occupancy 1",
994         "occupancy < .5"
995     };
996     runTest("simple.pdb", selections);
997 }
998
999 TEST_F(SelectionCollectionDataTest, HandlesBeta)
1000 {
1001     static const char * const selections[] = {
1002         "beta 0",
1003         "beta >= 0.3"
1004     };
1005     runTest("simple.pdb", selections);
1006 }
1007
1008 TEST_F(SelectionCollectionDataTest, HandlesResname)
1009 {
1010     static const char * const selections[] = {
1011         "resname RA",
1012         "resname RB RC"
1013     };
1014     runTest("simple.gro", selections);
1015 }
1016
1017 TEST_F(SelectionCollectionDataTest, HandlesCoordinateKeywords)
1018 {
1019     static const char * const selections[] = {
1020         "x < 3",
1021         "y >= 3",
1022         "x {-1 to 2}"
1023     };
1024     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1025     runTest("simple.gro", selections);
1026 }
1027
1028
1029 TEST_F(SelectionCollectionDataTest, HandlesSameResidue)
1030 {
1031     static const char * const selections[] = {
1032         "same residue as atomnr 1 4 12"
1033     };
1034     runTest("simple.gro", selections);
1035 }
1036
1037
1038 TEST_F(SelectionCollectionDataTest, HandlesSameResidueName)
1039 {
1040     static const char * const selections[] = {
1041         "same resname as atomnr 1 14"
1042     };
1043     runTest("simple.gro", selections);
1044 }
1045
1046
1047 TEST_F(SelectionCollectionDataTest, HandlesPositionKeywords)
1048 {
1049     static const char * const selections[] = {
1050         "cog of resnr 1 3",
1051         "res_cog of name CB and resnr 1 3",
1052         "whole_res_cog of name CB and resnr 1 3",
1053         "part_res_cog of x < 3",
1054         "dyn_res_cog of x < 3"
1055     };
1056     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
1057              | efTestPositionAtoms);
1058     runTest("simple.gro", selections);
1059 }
1060
1061
1062 TEST_F(SelectionCollectionDataTest, HandlesDistanceKeyword)
1063 {
1064     static const char * const selections[] = {
1065         "distance from cog of resnr 1 < 2"
1066     };
1067     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1068     runTest("simple.gro", selections);
1069 }
1070
1071
1072 TEST_F(SelectionCollectionDataTest, HandlesMinDistanceKeyword)
1073 {
1074     static const char * const selections[] = {
1075         "mindistance from resnr 1 < 2"
1076     };
1077     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1078     runTest("simple.gro", selections);
1079 }
1080
1081
1082 TEST_F(SelectionCollectionDataTest, HandlesWithinKeyword)
1083 {
1084     static const char * const selections[] = {
1085         "within 1 of resnr 2"
1086     };
1087     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1088     runTest("simple.gro", selections);
1089 }
1090
1091
1092 TEST_F(SelectionCollectionDataTest, HandlesInSolidAngleKeyword)
1093 {
1094     // Both of these should evaluate to empty on a correct implementation.
1095     static const char * const selections[] = {
1096         "resname TP and not insolidangle center cog of resname C span resname R cutoff 20",
1097         "resname TN and insolidangle center cog of resname C span resname R cutoff 20"
1098     };
1099     setFlags(TestFlags() | efDontTestCompiledAtoms | efTestEvaluation);
1100     runTest("sphere.gro", selections);
1101 }
1102
1103
1104 TEST_F(SelectionCollectionDataTest, HandlesPermuteModifier)
1105 {
1106     static const char * const selections[] = {
1107         "all permute 3 1 2",
1108         "res_cog of resnr 1 to 4 permute 2 1",
1109         "name CB S1 and res_cog x < 3 permute 2 1"
1110     };
1111     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
1112              | efTestPositionAtoms | efTestPositionMapping);
1113     runTest("simple.gro", selections);
1114 }
1115
1116
1117 TEST_F(SelectionCollectionDataTest, HandlesPlusModifier)
1118 {
1119     static const char * const selections[] = {
1120         "name S2 plus name S1",
1121         "res_cog of resnr 2 plus res_cog of resnr 1 plus res_cog of resnr 3",
1122         "name S1 and y < 3 plus res_cog of x < 2.5"
1123     };
1124     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
1125              | efTestPositionAtoms | efTestPositionMapping);
1126     runTest("simple.gro", selections);
1127 }
1128
1129
1130 TEST_F(SelectionCollectionDataTest, HandlesMergeModifier)
1131 {
1132     static const char * const selections[] = {
1133         "name S2 merge name S1",
1134         "resnr 1 2 and name S2 merge resnr 1 2 and name S1 merge res_cog of resnr 1 2",
1135         "name S1 and x < 2.5 merge res_cog of x < 2.5"
1136     };
1137     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
1138              | efTestPositionAtoms | efTestPositionMapping);
1139     runTest("simple.gro", selections);
1140 }
1141
1142
1143 /********************************************************************
1144  * Tests for generic selection evaluation
1145  */
1146
1147 TEST_F(SelectionCollectionDataTest, ComputesMassesAndCharges)
1148 {
1149     static const char * const selections[] = {
1150         "name CB",
1151         "y > 2",
1152         "res_cog of y > 2"
1153     };
1154     setFlags(TestFlags() | efTestEvaluation | efTestPositionAtoms
1155              | efTestPositionMasses | efTestPositionCharges);
1156     ASSERT_NO_FATAL_FAILURE(runParser(selections));
1157     ASSERT_NO_FATAL_FAILURE(topManager_.loadTopology("simple.gro"));
1158     t_atoms &atoms = topManager_.atoms();
1159     for (int i = 0; i < atoms.nr; ++i)
1160     {
1161         atoms.atom[i].m =   1.0 + i / 100.0;
1162         atoms.atom[i].q = -(1.0 + i / 100.0);
1163     }
1164     atoms.haveMass   = TRUE;
1165     atoms.haveCharge = TRUE;
1166     ASSERT_NO_FATAL_FAILURE(setTopology());
1167     ASSERT_NO_FATAL_FAILURE(runCompiler());
1168     ASSERT_NO_FATAL_FAILURE(runEvaluate());
1169     ASSERT_NO_FATAL_FAILURE(runEvaluateFinal());
1170 }
1171
1172 TEST_F(SelectionCollectionDataTest, ComputesMassesAndChargesWithoutTopology)
1173 {
1174     static const char * const selections[] = {
1175         "atomnr 1 to 3 8 to 9",
1176         "y > 2",
1177         "cog of (y > 2)"
1178     };
1179     setFlags(TestFlags() | efTestPositionAtoms
1180              | efTestPositionMasses | efTestPositionCharges);
1181     runTest(10, selections);
1182 }
1183
1184 TEST_F(SelectionCollectionDataTest, HandlesFramesWithAtomSubsets)
1185 {
1186     const int          index[]      = { 0, 1, 2, 3, 4, 5, 9, 10, 11 };
1187     const char * const selections[] = {
1188         "resnr 1 4",
1189         "atomnr 1 2 5 11 and y > 2",
1190         "res_cog of atomnr 2 5 11"
1191     };
1192     setFlags(TestFlags() | efTestEvaluation | efTestPositionAtoms);
1193     ASSERT_NO_FATAL_FAILURE(runParser(selections));
1194     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
1195     ASSERT_NO_FATAL_FAILURE(runCompiler());
1196     topManager_.initFrameIndices(index);
1197     ASSERT_NO_FATAL_FAILURE(runEvaluate());
1198     ASSERT_NO_FATAL_FAILURE(runEvaluateFinal());
1199 }
1200
1201
1202 /********************************************************************
1203  * Tests for selection syntactic constructs
1204  */
1205
1206 TEST_F(SelectionCollectionDataTest, HandlesSelectionNames)
1207 {
1208     static const char * const selections[] = {
1209         "\"GroupSelection\" group \"GrpA\"",
1210         "\"DynamicSelection\" x < 5",
1211         "y < 3"
1212     };
1213     setFlags(TestFlags() | efTestSelectionNames);
1214     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
1215     runTest(10, selections);
1216 }
1217
1218 TEST_F(SelectionCollectionDataTest, HandlesIndexGroupsInSelections)
1219 {
1220     static const char * const selections[] = {
1221         "group \"GrpA\"",
1222         "GrpB",
1223         "1",
1224         // These test that the name of the group is not too eagerly promoted to
1225         // the name of the selection.
1226         "group \"GrpB\" and resname RB",
1227         "group \"GrpA\" permute 5 3 2 1 4",
1228         "group \"GrpA\" plus group \"GrpB\"",
1229         "res_cog of group \"GrpA\""
1230     };
1231     setFlags(TestFlags() | efTestSelectionNames);
1232     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
1233     runTest("simple.gro", selections);
1234 }
1235
1236 TEST_F(SelectionCollectionDataTest, HandlesIndexGroupsInSelectionsDelayed)
1237 {
1238     static const char * const selections[] = {
1239         "group \"GrpA\"",
1240         "GrpB",
1241         "1",
1242         "group \"GrpB\" and resname RB"
1243     };
1244     setFlags(TestFlags() | efTestSelectionNames);
1245     ASSERT_NO_FATAL_FAILURE(runParser(selections));
1246     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
1247     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
1248     ASSERT_NO_FATAL_FAILURE(runCompiler());
1249 }
1250
1251 TEST_F(SelectionCollectionDataTest, HandlesUnsortedIndexGroupsInSelections)
1252 {
1253     static const char * const selections[] = {
1254         "foo = group \"GrpUnsorted\"",
1255         "group \"GrpUnsorted\"",
1256         "GrpUnsorted",
1257         "2",
1258         "res_cog of group \"GrpUnsorted\"",
1259         "group \"GrpUnsorted\" permute 2 1",
1260         "foo"
1261     };
1262     setFlags(TestFlags() | efTestPositionAtoms | efTestPositionMapping
1263              | efTestSelectionNames);
1264     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
1265     runTest("simple.gro", selections);
1266 }
1267
1268 TEST_F(SelectionCollectionDataTest, HandlesUnsortedIndexGroupsInSelectionsDelayed)
1269 {
1270     static const char * const selections[] = {
1271         "foo = group \"GrpUnsorted\"",
1272         "group \"GrpUnsorted\"",
1273         "GrpUnsorted",
1274         "2",
1275         "res_cog of group \"GrpUnsorted\"",
1276         "group \"GrpUnsorted\" permute 2 1",
1277         "foo"
1278     };
1279     ASSERT_NO_FATAL_FAILURE(runParser(selections));
1280     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
1281     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
1282     ASSERT_NO_FATAL_FAILURE(runCompiler());
1283 }
1284
1285
1286 TEST_F(SelectionCollectionDataTest, HandlesConstantPositions)
1287 {
1288     static const char * const selections[] = {
1289         "[1, -2, 3.5]"
1290     };
1291     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
1292              | efTestPositionMapping);
1293     runTest("simple.gro", selections);
1294 }
1295
1296
1297 TEST_F(SelectionCollectionDataTest, HandlesConstantPositionsWithModifiers)
1298 {
1299     static const char * const selections[] = {
1300         "[0, 0, 0] plus [0, 1, 0]"
1301     };
1302     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
1303              | efTestPositionMapping);
1304     runTest("simple.gro", selections);
1305 }
1306
1307
1308 TEST_F(SelectionCollectionDataTest, HandlesWithinConstantPositions)
1309 {
1310     static const char * const selections[] = {
1311         "within 1 of [2, 1, 0]"
1312     };
1313     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1314     runTest("simple.gro", selections);
1315 }
1316
1317
1318 TEST_F(SelectionCollectionDataTest, HandlesOverlappingIntegerRanges)
1319 {
1320     static const char * const selections[] = {
1321         "atomnr 2 to 4 5 to 8",
1322         "atomnr 2 to 5 4 to 7"
1323     };
1324     ASSERT_NO_FATAL_FAILURE(runTest(10, selections));
1325 }
1326
1327
1328 TEST_F(SelectionCollectionDataTest, HandlesOverlappingRealRanges)
1329 {
1330     static const char * const selections[] = {
1331         "charge {-0.35 to -0.05 0.25 to 0.75}",
1332         "charge {0.05 to -0.3 -0.05 to 0.55}"
1333     };
1334     ASSERT_NO_FATAL_FAILURE(runParser(selections));
1335     ASSERT_NO_FATAL_FAILURE(topManager_.loadTopology("simple.gro"));
1336     t_atoms &atoms = topManager_.atoms();
1337     for (int i = 0; i < atoms.nr; ++i)
1338     {
1339         atoms.atom[i].q = i / 10.0 - 0.5;
1340     }
1341     atoms.haveCharge = TRUE;
1342     ASSERT_NO_FATAL_FAILURE(setTopology());
1343     ASSERT_NO_FATAL_FAILURE(runCompiler());
1344 }
1345
1346
1347 TEST_F(SelectionCollectionDataTest, HandlesForcedStringMatchingMode)
1348 {
1349     static const char * const selections[] = {
1350         "name = S1 \"C?\"",
1351         "name ? S1 \"C?\""
1352     };
1353     runTest("simple.gro", selections);
1354 }
1355
1356
1357 TEST_F(SelectionCollectionDataTest, HandlesWildcardMatching)
1358 {
1359     static const char * const selections[] = {
1360         "name \"S?\"",
1361         "name ? \"S?\""
1362     };
1363     runTest("simple.gro", selections);
1364 }
1365
1366
1367 TEST_F(SelectionCollectionDataTest, HandlesRegexMatching)
1368 {
1369     static const char * const selections[] = {
1370         "resname \"R[BD]\"",
1371         "resname ~ \"R[BD]\""
1372     };
1373     if (gmx::Regex::isSupported())
1374     {
1375         runTest("simple.gro", selections);
1376     }
1377 }
1378
1379
1380 TEST_F(SelectionCollectionDataTest, HandlesBasicBoolean)
1381 {
1382     static const char * const selections[] = {
1383         "atomnr 1 to 5 and atomnr 2 to 7",
1384         "atomnr 1 to 5 or not atomnr 3 to 8",
1385         "not not atomnr 1 to 5 and atomnr 2 to 6 and not not atomnr 3 to 7",
1386         "atomnr 1 to 5 and (atomnr 2 to 7 and atomnr 3 to 6)",
1387         "x < 5 and atomnr 1 to 5 and y < 3 and atomnr 2 to 4"
1388     };
1389     runTest(10, selections);
1390 }
1391
1392
1393 TEST_F(SelectionCollectionDataTest, HandlesDynamicAtomValuedParameters)
1394 {
1395     static const char * const selections[] = {
1396         "same residue as (atomnr 3 5 13 or y > 5)",
1397         "(resnr 1 3 5 or x > 10) and same residue as (atomnr 3 5 13 or z > 5)"
1398     };
1399     setFlags(TestFlags() | efTestEvaluation);
1400     runTest("simple.gro", selections);
1401 }
1402
1403
1404 TEST_F(SelectionCollectionDataTest, HandlesEmptySelectionWithUnevaluatedExpressions)
1405 {
1406     static const char * const selections[] = {
1407         "none and x > 2",
1408         "none and same resname as resnr 2"
1409     };
1410     runTest("simple.gro", selections);
1411 }
1412
1413
1414 TEST_F(SelectionCollectionDataTest, HandlesEmptyReferenceForSame)
1415 {
1416     static const char * const selections[] = {
1417         "same residue as none",
1418         "same resname as none"
1419     };
1420     runTest("simple.gro", selections);
1421 }
1422
1423
1424 TEST_F(SelectionCollectionDataTest, HandlesPositionModifiersForKeywords)
1425 {
1426     static const char * const selections[] = {
1427         "res_cog x > 2",
1428         "name CB and res_cog y > 2.5"
1429     };
1430     setFlags(TestFlags() | efTestEvaluation);
1431     runTest("simple.gro", selections);
1432 }
1433
1434
1435 TEST_F(SelectionCollectionDataTest, HandlesPositionModifiersForMethods)
1436 {
1437     static const char * const selections[] = {
1438         "res_cog distance from cog of resnr 1 < 2",
1439         "res_cog within 2 of cog of resnr 1"
1440     };
1441     setFlags(TestFlags() | efTestEvaluation);
1442     runTest("simple.gro", selections);
1443 }
1444
1445
1446 TEST_F(SelectionCollectionDataTest, HandlesKeywordOfPositions)
1447 {
1448     static const char * const selections[] = {
1449         "x < y of cog of resnr 2"
1450     };
1451     setFlags(TestFlags() | efTestEvaluation);
1452     runTest("simple.gro", selections);
1453 }
1454
1455 TEST_F(SelectionCollectionDataTest, HandlesKeywordOfPositionsInArithmetic)
1456 {
1457     static const char * const selections[] = {
1458         "x - y of cog of resnr 2 < 0"
1459     };
1460     setFlags(TestFlags() | efTestEvaluation);
1461     runTest("simple.gro", selections);
1462 }
1463
1464
1465 TEST_F(SelectionCollectionDataTest, HandlesNumericComparisons)
1466 {
1467     static const char * const selections[] = {
1468         "x > 2",
1469         "2 < x",
1470         "y > resnr",
1471         "resnr < 2.5",
1472         "2.5 > resnr"
1473     };
1474     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1475     runTest("simple.gro", selections);
1476 }
1477
1478
1479 TEST_F(SelectionCollectionDataTest, HandlesArithmeticExpressions)
1480 {
1481     static const char * const selections[] = {
1482         "x+1 > 3",
1483         "(y-1)^2 <= 1",
1484         "x+--1 > 3",
1485         "-x+-1 < -3"
1486     };
1487     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1488     runTest("simple.gro", selections);
1489 }
1490
1491
1492 TEST_F(SelectionCollectionDataTest, HandlesNumericVariables)
1493 {
1494     static const char * const selections[] = {
1495         "value = x + y",
1496         "value <= 4",
1497         "index = resnr",
1498         "index < 3"
1499     };
1500     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1501     runTest("simple.gro", selections);
1502 }
1503
1504
1505 TEST_F(SelectionCollectionDataTest, HandlesComplexNumericVariables)
1506 {
1507     static const char * const selections[] = {
1508         "value = x + y",
1509         "resname RA and value <= 4",
1510         "resname RA RB and x < 3 and value <= 4",
1511         "index = atomnr",
1512         "resname RA and index < 3",
1513         "resname RB and y < 3 and index < 6"
1514     };
1515     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1516     runTest("simple.gro", selections);
1517 }
1518
1519
1520 TEST_F(SelectionCollectionDataTest, HandlesPositionVariables)
1521 {
1522     static const char * const selections[] = {
1523         "foo = res_cog of resname RA",
1524         "foo",
1525         "within 1 of foo",
1526         "bar = cog of resname RA",
1527         "bar",
1528         "within 1 of bar"
1529     };
1530     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1531     runTest("simple.gro", selections);
1532 }
1533
1534
1535 TEST_F(SelectionCollectionDataTest, HandlesPositionVariableInModifier)
1536 {
1537     static const char * const selections[] = {
1538         "foo = cog of resnr 1",
1539         "cog of resnr 2 plus foo",
1540         "cog of resnr 3 plus foo"
1541     };
1542     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1543     runTest("simple.gro", selections);
1544 }
1545
1546
1547 TEST_F(SelectionCollectionDataTest, HandlesConstantPositionInVariable)
1548 {
1549     static const char * const selections[] = {
1550         "constpos = [1.0, 2.5, 0.5]",
1551         "constpos",
1552         "within 2 of constpos"
1553     };
1554     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
1555              | efTestPositionAtoms);
1556     runTest("simple.gro", selections);
1557 }
1558
1559
1560 TEST_F(SelectionCollectionDataTest, HandlesNumericConstantsInVariables)
1561 {
1562     static const char * const selections[] = {
1563         "constint = 4",
1564         "constreal1 = 0.5",
1565         "constreal2 = 2.7",
1566         "resnr < constint",
1567         "x + constreal1 < constreal2"
1568     };
1569     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1570     runTest("simple.gro", selections);
1571 }
1572
1573
1574 /********************************************************************
1575  * Tests for complex boolean syntax
1576  */
1577
1578 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysis)
1579 {
1580     static const char * const selections[] = {
1581         "atomnr 1 to 5 and atomnr 2 to 7 and x < 2",
1582         "atomnr 1 to 5 and (atomnr 4 to 7 or x < 2)",
1583         "atomnr 1 to 5 and y < 3 and (atomnr 4 to 7 or x < 2)",
1584         "atomnr 1 to 5 and not (atomnr 4 to 7 or x < 2)",
1585         "atomnr 1 to 5 or (atomnr 4 to 6 and (atomnr 5 to 7 or x < 2))"
1586     };
1587     runTest(10, selections);
1588 }
1589
1590
1591 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysisWithVariables)
1592 {
1593     static const char * const selections[] = {
1594         "foo = atomnr 4 to 7 or x < 2",
1595         "atomnr 1 to 4 and foo",
1596         "atomnr 2 to 6 and y < 3 and foo",
1597         "atomnr 6 to 10 and not foo"
1598     };
1599     runTest(10, selections);
1600 }
1601
1602
1603 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysisWithMoreVariables)
1604 {
1605     static const char * const selections[] = {
1606         "foo = atomnr 4 to 7",
1607         "bar = foo and x < 2",
1608         "bar2 = foo and y < 2",
1609         "atomnr 1 to 4 and bar",
1610         "atomnr 2 to 6 and y < 3 and bar2",
1611         "atomnr 6 to 10 and not foo"
1612     };
1613     runTest(10, selections);
1614 }
1615
1616
1617 /********************************************************************
1618  * Tests for complex subexpression cases
1619  *
1620  * These tests use some knowledge of the implementation to trigger different
1621  * paths in the code.
1622  */
1623
1624 TEST_F(SelectionCollectionDataTest, HandlesUnusedVariables)
1625 {
1626     static const char * const selections[] = {
1627         "unused1 = atomnr 1 to 3",
1628         "foo = atomnr 4 to 7",
1629         "atomnr 1 to 6 and foo",
1630         "unused2 = atomnr 3 to 5"
1631     };
1632     runTest(10, selections);
1633 }
1634
1635
1636 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithStaticEvaluationGroups)
1637 {
1638     static const char * const selections[] = {
1639         "foo = atomnr 4 to 7 and x < 2",
1640         "atomnr 1 to 5 and foo",
1641         "atomnr 3 to 7 and foo"
1642     };
1643     runTest(10, selections);
1644 }
1645
1646
1647 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithMixedEvaluationGroups)
1648 {
1649     static const char * const selections[] = {
1650         "foo = atomnr 4 to 7 and x < 2",
1651         "atomnr 1 to 6 and foo",
1652         "within 1 of foo",
1653         "foo"
1654     };
1655     runTest(10, selections);
1656 }
1657
1658
1659 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithMixedEvaluationGroups2)
1660 {
1661     static const char * const selections[] = {
1662         "foo = atomnr 1 to 8 and x < 10",
1663         "atomnr 1 to 5 and y < 10 and foo",
1664         "foo"
1665     };
1666     setFlags(TestFlags() | efTestEvaluation);
1667     runTest("simple.gro", selections);
1668 }
1669
1670
1671 } // namespace