Merge release-4-6 into master
[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, 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 <gtest/gtest.h>
43
44 #include "gromacs/options/basicoptions.h"
45 #include "gromacs/options/options.h"
46 #include "gromacs/selection/indexutil.h"
47 #include "gromacs/selection/selectioncollection.h"
48 #include "gromacs/selection/selection.h"
49 #include "gromacs/utility/arrayref.h"
50 #include "gromacs/utility/exceptions.h"
51 #include "gromacs/utility/flags.h"
52 #include "gromacs/utility/gmxregex.h"
53 #include "gromacs/utility/stringutil.h"
54
55 #include "testutils/refdata.h"
56 #include "testutils/testasserts.h"
57 #include "testutils/testfilemanager.h"
58 #include "testutils/testoptions.h"
59
60 #include "toputils.h"
61
62 namespace
63 {
64
65 /********************************************************************
66  * Test fixture for selection testing
67  */
68
69 class SelectionCollectionTest : public ::testing::Test
70 {
71     public:
72         static int               s_debugLevel;
73
74         SelectionCollectionTest();
75         ~SelectionCollectionTest();
76
77         void setAtomCount(int natoms)
78         {
79             ASSERT_NO_THROW_GMX(sc_.setTopology(NULL, natoms));
80         }
81         void loadTopology(const char *filename);
82         void setTopology();
83         void loadIndexGroups(const char *filename);
84
85         gmx::test::TopologyManager  topManager_;
86         gmx::SelectionCollection    sc_;
87         gmx::SelectionList          sel_;
88         t_topology                 *top_;
89         t_trxframe                 *frame_;
90         gmx_ana_indexgrps_t        *grps_;
91 };
92
93 int SelectionCollectionTest::s_debugLevel = 0;
94
95 GMX_TEST_OPTIONS(SelectionCollectionTestOptions, options)
96 {
97     options->addOption(gmx::IntegerOption("seldebug")
98                            .store(&SelectionCollectionTest::s_debugLevel)
99                            .description("Set selection debug level"));
100 }
101
102 SelectionCollectionTest::SelectionCollectionTest()
103     : top_(NULL), frame_(NULL), grps_(NULL)
104 {
105     topManager_.requestFrame();
106     sc_.setDebugLevel(s_debugLevel);
107     sc_.setReferencePosType("atom");
108     sc_.setOutputPosType("atom");
109 }
110
111 SelectionCollectionTest::~SelectionCollectionTest()
112 {
113     if (grps_ != NULL)
114     {
115         gmx_ana_indexgrps_free(grps_);
116     }
117 }
118
119 void
120 SelectionCollectionTest::loadTopology(const char *filename)
121 {
122     topManager_.loadTopology(filename);
123     setTopology();
124 }
125
126 void
127 SelectionCollectionTest::setTopology()
128 {
129     top_   = topManager_.topology();
130     frame_ = topManager_.frame();
131
132     ASSERT_NO_THROW_GMX(sc_.setTopology(top_, -1));
133 }
134
135 void
136 SelectionCollectionTest::loadIndexGroups(const char *filename)
137 {
138     GMX_RELEASE_ASSERT(grps_ == NULL,
139                        "External groups can only be loaded once");
140     std::string fullpath =
141         gmx::test::TestFileManager::getInputFilePath(filename);
142     gmx_ana_indexgrps_init(&grps_, NULL, fullpath.c_str());
143     sc_.setIndexGroups(grps_);
144 }
145
146
147 /********************************************************************
148  * Test fixture for selection testing with reference data
149  */
150
151 class SelectionCollectionDataTest : public SelectionCollectionTest
152 {
153     public:
154         enum TestFlag
155         {
156             efTestEvaluation            = 1<<0,
157             efTestPositionAtoms         = 1<<1,
158             efTestPositionCoordinates   = 1<<2,
159             efTestPositionMapping       = 1<<3,
160             efTestPositionMasses        = 1<<4,
161             efTestPositionCharges       = 1<<5,
162             efTestSelectionNames        = 1<<6,
163             efDontTestCompiledAtoms     = 1<<8
164         };
165         typedef gmx::FlagsTemplate<TestFlag> TestFlags;
166
167         SelectionCollectionDataTest()
168             : checker_(data_.rootChecker()), count_(0), framenr_(0)
169         {
170         }
171
172         void setFlags(TestFlags flags) { flags_ = flags; }
173
174         void runParser(const gmx::ConstArrayRef<const char *> &selections);
175         void runCompiler();
176         void runEvaluate();
177         void runEvaluateFinal();
178
179         void runTest(int                                     natoms,
180                      const gmx::ConstArrayRef<const char *> &selections);
181         void runTest(const char                             *filename,
182                      const gmx::ConstArrayRef<const char *> &selections);
183
184     private:
185         static void checkSelection(gmx::test::TestReferenceChecker *checker,
186                                    const gmx::Selection &sel, TestFlags flags);
187
188         void checkCompiled();
189
190         gmx::test::TestReferenceData    data_;
191         gmx::test::TestReferenceChecker checker_;
192         size_t                          count_;
193         int                             framenr_;
194         TestFlags                       flags_;
195 };
196
197
198 void
199 SelectionCollectionDataTest::checkSelection(
200         gmx::test::TestReferenceChecker *checker,
201         const gmx::Selection &sel, TestFlags flags)
202 {
203     using gmx::test::TestReferenceChecker;
204
205     {
206         gmx::ConstArrayRef<int> atoms = sel.atomIndices();
207         checker->checkSequence(atoms.begin(), atoms.end(), "Atoms");
208     }
209     if (flags.test(efTestPositionAtoms)
210         || flags.test(efTestPositionCoordinates)
211         || flags.test(efTestPositionMapping)
212         || flags.test(efTestPositionMasses)
213         || flags.test(efTestPositionCharges))
214     {
215         TestReferenceChecker compound(
216                 checker->checkSequenceCompound("Positions", sel.posCount()));
217         for (int i = 0; i < sel.posCount(); ++i)
218         {
219             TestReferenceChecker          poscompound(compound.checkCompound("Position", NULL));
220             const gmx::SelectionPosition &p = sel.position(i);
221             if (flags.test(efTestPositionAtoms))
222             {
223                 gmx::ConstArrayRef<int> atoms = p.atomIndices();
224                 poscompound.checkSequence(atoms.begin(), atoms.end(), "Atoms");
225             }
226             if (flags.test(efTestPositionCoordinates))
227             {
228                 poscompound.checkVector(p.x(), "Coordinates");
229             }
230             if (flags.test(efTestPositionMapping))
231             {
232                 poscompound.checkInteger(p.refId(), "RefId");
233                 poscompound.checkInteger(p.mappedId(), "MappedId");
234             }
235             if (flags.test(efTestPositionMasses))
236             {
237                 poscompound.checkReal(p.mass(), "Mass");
238             }
239             if (flags.test(efTestPositionCharges))
240             {
241                 poscompound.checkReal(p.charge(), "Charge");
242             }
243         }
244     }
245 }
246
247
248 void
249 SelectionCollectionDataTest::runParser(
250         const gmx::ConstArrayRef<const char *> &selections)
251 {
252     using gmx::test::TestReferenceChecker;
253
254     TestReferenceChecker compound(checker_.checkCompound("ParsedSelections", "Parsed"));
255     size_t               varcount = 0;
256     count_ = 0;
257     for (size_t i = 0; i < selections.size(); ++i)
258     {
259         SCOPED_TRACE(std::string("Parsing selection \"")
260                      + selections[i] + "\"");
261         gmx::SelectionList result;
262         ASSERT_NO_THROW_GMX(result = sc_.parseFromString(selections[i]));
263         sel_.insert(sel_.end(), result.begin(), result.end());
264         if (sel_.size() == count_)
265         {
266             std::string          id = gmx::formatString("Variable%d", static_cast<int>(varcount + 1));
267             TestReferenceChecker varcompound(
268                     compound.checkCompound("ParsedVariable", id.c_str()));
269             varcompound.checkString(selections[i], "Input");
270             ++varcount;
271         }
272         else
273         {
274             std::string          id = gmx::formatString("Selection%d", static_cast<int>(count_ + 1));
275             TestReferenceChecker selcompound(
276                     compound.checkCompound("ParsedSelection", id.c_str()));
277             selcompound.checkString(selections[i], "Input");
278             if (flags_.test(efTestSelectionNames))
279             {
280                 selcompound.checkString(sel_[count_].name(), "Name");
281             }
282             selcompound.checkString(sel_[count_].selectionText(), "Text");
283             selcompound.checkBoolean(sel_[count_].isDynamic(), "Dynamic");
284             ++count_;
285         }
286     }
287 }
288
289
290 void
291 SelectionCollectionDataTest::runCompiler()
292 {
293     ASSERT_NO_THROW_GMX(sc_.compile());
294     ASSERT_EQ(count_, sel_.size());
295     checkCompiled();
296 }
297
298
299 void
300 SelectionCollectionDataTest::checkCompiled()
301 {
302     using gmx::test::TestReferenceChecker;
303     const TestFlags      mask = ~TestFlags(efTestPositionCoordinates);
304
305     TestReferenceChecker compound(checker_.checkCompound("CompiledSelections", "Compiled"));
306     for (size_t i = 0; i < count_; ++i)
307     {
308         SCOPED_TRACE(std::string("Checking selection \"") +
309                      sel_[i].selectionText() + "\"");
310         std::string          id = gmx::formatString("Selection%d", static_cast<int>(i + 1));
311         TestReferenceChecker selcompound(
312                 compound.checkCompound("Selection", id.c_str()));
313         if (flags_.test(efTestSelectionNames))
314         {
315             selcompound.checkString(sel_[i].name(), "Name");
316         }
317         if (!flags_.test(efDontTestCompiledAtoms))
318         {
319             checkSelection(&selcompound, sel_[i], flags_ & mask);
320         }
321     }
322 }
323
324
325 void
326 SelectionCollectionDataTest::runEvaluate()
327 {
328     using gmx::test::TestReferenceChecker;
329
330     ++framenr_;
331     ASSERT_NO_THROW_GMX(sc_.evaluate(frame_, NULL));
332     std::string          frame = gmx::formatString("Frame%d", framenr_);
333     TestReferenceChecker compound(
334             checker_.checkCompound("EvaluatedSelections", frame.c_str()));
335     for (size_t i = 0; i < count_; ++i)
336     {
337         SCOPED_TRACE(std::string("Checking selection \"") +
338                      sel_[i].selectionText() + "\"");
339         std::string          id = gmx::formatString("Selection%d", static_cast<int>(i + 1));
340         TestReferenceChecker selcompound(
341                 compound.checkCompound("Selection", id.c_str()));
342         checkSelection(&selcompound, sel_[i], flags_);
343     }
344 }
345
346
347 void
348 SelectionCollectionDataTest::runEvaluateFinal()
349 {
350     ASSERT_NO_THROW_GMX(sc_.evaluateFinal(framenr_));
351     checkCompiled();
352 }
353
354
355 void
356 SelectionCollectionDataTest::runTest(
357         int natoms, const gmx::ConstArrayRef<const char *> &selections)
358 {
359     ASSERT_NO_FATAL_FAILURE(runParser(selections));
360     ASSERT_NO_FATAL_FAILURE(setAtomCount(natoms));
361     ASSERT_NO_FATAL_FAILURE(runCompiler());
362 }
363
364
365 void
366 SelectionCollectionDataTest::runTest(
367         const char *filename, const gmx::ConstArrayRef<const char *> &selections)
368 {
369     ASSERT_NO_FATAL_FAILURE(runParser(selections));
370     ASSERT_NO_FATAL_FAILURE(loadTopology(filename));
371     ASSERT_NO_FATAL_FAILURE(runCompiler());
372     if (flags_.test(efTestEvaluation))
373     {
374         ASSERT_NO_FATAL_FAILURE(runEvaluate());
375         ASSERT_NO_FATAL_FAILURE(runEvaluateFinal());
376     }
377 }
378
379
380 /********************************************************************
381  * Tests for SelectionCollection functionality without reference data
382  */
383
384 TEST_F(SelectionCollectionTest, HandlesNoSelections)
385 {
386     EXPECT_FALSE(sc_.requiresTopology());
387     EXPECT_NO_THROW_GMX(sc_.compile());
388 }
389
390 TEST_F(SelectionCollectionTest, HandlesVelocityAndForceRequests)
391 {
392     ASSERT_NO_THROW_GMX(sel_ = sc_.parseFromString("atomnr 1 to 10; none"));
393     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
394     ASSERT_EQ(2U, sel_.size());
395     ASSERT_NO_THROW_GMX(sel_[0].setEvaluateVelocities(true));
396     ASSERT_NO_THROW_GMX(sel_[1].setEvaluateVelocities(true));
397     ASSERT_NO_THROW_GMX(sel_[0].setEvaluateForces(true));
398     ASSERT_NO_THROW_GMX(sel_[1].setEvaluateForces(true));
399     ASSERT_NO_THROW_GMX(sc_.compile());
400     EXPECT_TRUE(sel_[0].hasVelocities());
401     EXPECT_TRUE(sel_[1].hasVelocities());
402     EXPECT_TRUE(sel_[0].hasForces());
403     EXPECT_TRUE(sel_[1].hasForces());
404 }
405
406 TEST_F(SelectionCollectionTest, ParsesSelectionsFromFile)
407 {
408     ASSERT_NO_THROW_GMX(sel_ = sc_.parseFromFile(
409                                     gmx::test::TestFileManager::getInputFilePath("selfile.dat")));
410     // These should match the contents of selfile.dat
411     ASSERT_EQ(2U, sel_.size());
412     EXPECT_STREQ("resname RA RB", sel_[0].selectionText());
413     EXPECT_STREQ("resname RB RC", sel_[1].selectionText());
414 }
415
416 TEST_F(SelectionCollectionTest, HandlesInvalidRegularExpressions)
417 {
418     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
419     EXPECT_THROW_GMX({
420                          sc_.parseFromString("resname ~ \"R[A\"");
421                          sc_.compile();
422                      }, gmx::InvalidInputError);
423 }
424
425 TEST_F(SelectionCollectionTest, HandlesUnsupportedRegularExpressions)
426 {
427     if (!gmx::Regex::isSupported())
428     {
429         ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
430         EXPECT_THROW_GMX({
431                              sc_.parseFromString("resname \"R[AD]\"");
432                              sc_.compile();
433                          }, gmx::InvalidInputError);
434     }
435 }
436
437 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue)
438 {
439     EXPECT_THROW_GMX(sc_.parseFromString("mindist from atomnr 1 cutoff"),
440                      gmx::InvalidInputError);
441 }
442
443 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue2)
444 {
445     EXPECT_THROW_GMX(sc_.parseFromString("within 1 of"),
446                      gmx::InvalidInputError);
447 }
448
449 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue3)
450 {
451     EXPECT_THROW_GMX(sc_.parseFromString("within of atomnr 1"),
452                      gmx::InvalidInputError);
453 }
454
455 TEST_F(SelectionCollectionTest, HandlesHelpKeywordInInvalidContext)
456 {
457     EXPECT_THROW_GMX(sc_.parseFromString("resname help"),
458                      gmx::InvalidInputError);
459 }
460
461 // TODO: Tests for more parser errors
462
463 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceParser1)
464 {
465     ASSERT_NO_THROW_GMX(sc_.setIndexGroups(NULL));
466     EXPECT_THROW_GMX(sc_.parseFromString("group \"foo\""), gmx::InconsistentInputError);
467     EXPECT_THROW_GMX(sc_.parseFromString("4"), gmx::InconsistentInputError);
468 }
469
470 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceParser2)
471 {
472     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
473     EXPECT_THROW_GMX(sc_.parseFromString("group \"foo\""), gmx::InconsistentInputError);
474     EXPECT_THROW_GMX(sc_.parseFromString("4"), gmx::InconsistentInputError);
475 }
476
477 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceDelayed1)
478 {
479     ASSERT_NO_THROW_GMX(sc_.parseFromString("group \"foo\""));
480     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
481     EXPECT_THROW_GMX(sc_.setIndexGroups(NULL), gmx::InconsistentInputError);
482     EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
483 }
484
485 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceDelayed2)
486 {
487     ASSERT_NO_THROW_GMX(sc_.parseFromString("group 4; group \"foo\""));
488     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
489     EXPECT_THROW_GMX(loadIndexGroups("simple.ndx"), gmx::InconsistentInputError);
490     EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
491 }
492
493 TEST_F(SelectionCollectionTest, HandlesUnsortedGroupReference)
494 {
495     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
496     EXPECT_THROW_GMX(sc_.parseFromString("group \"GrpUnsorted\""),
497                      gmx::InconsistentInputError);
498     EXPECT_THROW_GMX(sc_.parseFromString("2"), gmx::InconsistentInputError);
499 }
500
501 TEST_F(SelectionCollectionTest, HandlesUnsortedGroupReferenceDelayed)
502 {
503     ASSERT_NO_THROW_GMX(sc_.parseFromString("group 2; group \"GrpUnsorted\""));
504     EXPECT_THROW_GMX(loadIndexGroups("simple.ndx"), gmx::InconsistentInputError);
505     EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
506 }
507
508 TEST_F(SelectionCollectionTest, RecoversFromMissingMoleculeInfo)
509 {
510     ASSERT_NO_THROW_GMX(sc_.parseFromString("molindex 1 to 5"));
511     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
512     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
513 }
514
515 TEST_F(SelectionCollectionTest, RecoversFromMissingAtomTypes)
516 {
517     ASSERT_NO_THROW_GMX(sc_.parseFromString("type CA"));
518     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
519     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
520 }
521
522 TEST_F(SelectionCollectionTest, RecoversFromMissingPDBInfo)
523 {
524     ASSERT_NO_THROW_GMX(sc_.parseFromString("altloc A"));
525     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
526     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
527 }
528
529 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation)
530 {
531     ASSERT_NO_THROW_GMX(sc_.parseFromString("all permute 1 1"));
532     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
533     EXPECT_THROW_GMX(sc_.compile(), gmx::InvalidInputError);
534 }
535
536 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation2)
537 {
538     ASSERT_NO_THROW_GMX(sc_.parseFromString("all permute 3 2 1"));
539     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
540     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
541 }
542
543 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation3)
544 {
545     ASSERT_NO_THROW_GMX(sc_.parseFromString("x < 1.5 permute 3 2 1"));
546     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
547     ASSERT_NO_THROW_GMX(sc_.compile());
548     EXPECT_THROW_GMX(sc_.evaluate(frame_, NULL), gmx::InconsistentInputError);
549 }
550
551 // TODO: Tests for evaluation errors
552
553
554 /********************************************************************
555  * Tests for selection keywords
556  */
557
558 TEST_F(SelectionCollectionDataTest, HandlesAllNone)
559 {
560     static const char * const selections[] = {
561         "all",
562         "none"
563     };
564     runTest(10, selections);
565 }
566
567 TEST_F(SelectionCollectionDataTest, HandlesAtomnr)
568 {
569     static const char * const selections[] = {
570         "atomnr 1 to 3 6 to 8",
571         "atomnr 4 2 5 to 7",
572         "atomnr <= 5"
573     };
574     runTest(10, selections);
575 }
576
577 TEST_F(SelectionCollectionDataTest, HandlesResnr)
578 {
579     static const char * const selections[] = {
580         "resnr 1 2 5",
581         "resid 4 to 3"
582     };
583     runTest("simple.gro", selections);
584 }
585
586 TEST_F(SelectionCollectionDataTest, HandlesResIndex)
587 {
588     static const char * const selections[] = {
589         "resindex 1 4",
590         "residue 1 3"
591     };
592     runTest("simple.pdb", selections);
593 }
594
595 TEST_F(SelectionCollectionDataTest, HandlesMolIndex)
596 {
597     static const char * const selections[] = {
598         "molindex 1 4",
599         "molecule 2 3 5"
600     };
601     ASSERT_NO_FATAL_FAILURE(runParser(selections));
602     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
603     topManager_.initUniformMolecules(3);
604     ASSERT_NO_FATAL_FAILURE(runCompiler());
605 }
606
607 TEST_F(SelectionCollectionDataTest, HandlesAtomname)
608 {
609     static const char * const selections[] = {
610         "name CB",
611         "atomname S1 S2"
612     };
613     runTest("simple.gro", selections);
614 }
615
616 TEST_F(SelectionCollectionDataTest, HandlesPdbAtomname)
617 {
618     static const char * const selections[] = {
619         "name HG21",
620         "name 1HG2",
621         "pdbname HG21 CB",
622         "pdbatomname 1HG2"
623     };
624     runTest("simple.pdb", selections);
625 }
626
627
628 TEST_F(SelectionCollectionDataTest, HandlesAtomtype)
629 {
630     static const char * const selections[] = {
631         "atomtype CA"
632     };
633     ASSERT_NO_FATAL_FAILURE(runParser(selections));
634     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
635     const char *const types[] = { "CA", "SA", "SB" };
636     topManager_.initAtomTypes(types);
637     ASSERT_NO_FATAL_FAILURE(runCompiler());
638 }
639
640 TEST_F(SelectionCollectionDataTest, HandlesChain)
641 {
642     static const char * const selections[] = {
643         "chain A",
644         "chain B"
645     };
646     runTest("simple.pdb", selections);
647 }
648
649 TEST_F(SelectionCollectionDataTest, HandlesMass)
650 {
651     static const char * const selections[] = {
652         "mass > 5"
653     };
654     ASSERT_NO_FATAL_FAILURE(runParser(selections));
655     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
656     for (int i = 0; i < top_->atoms.nr; ++i)
657     {
658         top_->atoms.atom[i].m = 1.0 + i;
659     }
660     ASSERT_NO_FATAL_FAILURE(runCompiler());
661 }
662
663 TEST_F(SelectionCollectionDataTest, HandlesCharge)
664 {
665     static const char * const selections[] = {
666         "charge < 0.5"
667     };
668     ASSERT_NO_FATAL_FAILURE(runParser(selections));
669     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
670     for (int i = 0; i < top_->atoms.nr; ++i)
671     {
672         top_->atoms.atom[i].q = i / 10.0;
673     }
674     ASSERT_NO_FATAL_FAILURE(runCompiler());
675 }
676
677 TEST_F(SelectionCollectionDataTest, HandlesAltLoc)
678 {
679     static const char * const selections[] = {
680         "altloc \" \"",
681         "altloc A"
682     };
683     runTest("simple.pdb", selections);
684 }
685
686 TEST_F(SelectionCollectionDataTest, HandlesInsertCode)
687 {
688     static const char * const selections[] = {
689         "insertcode \" \"",
690         "insertcode A"
691     };
692     runTest("simple.pdb", selections);
693 }
694
695 TEST_F(SelectionCollectionDataTest, HandlesOccupancy)
696 {
697     static const char * const selections[] = {
698         "occupancy 1",
699         "occupancy < .5"
700     };
701     runTest("simple.pdb", selections);
702 }
703
704 TEST_F(SelectionCollectionDataTest, HandlesBeta)
705 {
706     static const char * const selections[] = {
707         "beta 0",
708         "beta >= 0.3"
709     };
710     runTest("simple.pdb", selections);
711 }
712
713 TEST_F(SelectionCollectionDataTest, HandlesResname)
714 {
715     static const char * const selections[] = {
716         "resname RA",
717         "resname RB RC"
718     };
719     runTest("simple.gro", selections);
720 }
721
722 TEST_F(SelectionCollectionDataTest, HandlesCoordinateKeywords)
723 {
724     static const char * const selections[] = {
725         "x < 3",
726         "y >= 3",
727         "x {-1 to 2}"
728     };
729     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
730     runTest("simple.gro", selections);
731 }
732
733
734 TEST_F(SelectionCollectionDataTest, HandlesSameResidue)
735 {
736     static const char * const selections[] = {
737         "same residue as atomnr 1 4 12"
738     };
739     runTest("simple.gro", selections);
740 }
741
742
743 TEST_F(SelectionCollectionDataTest, HandlesSameResidueName)
744 {
745     static const char * const selections[] = {
746         "same resname as atomnr 1 14"
747     };
748     runTest("simple.gro", selections);
749 }
750
751
752 TEST_F(SelectionCollectionDataTest, HandlesPositionKeywords)
753 {
754     static const char * const selections[] = {
755         "cog of resnr 1 3",
756         "res_cog of name CB and resnr 1 3",
757         "whole_res_cog of name CB and resnr 1 3",
758         "part_res_cog of x < 3",
759         "dyn_res_cog of x < 3"
760     };
761     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
762              | efTestPositionAtoms);
763     runTest("simple.gro", selections);
764 }
765
766
767 TEST_F(SelectionCollectionDataTest, HandlesDistanceKeyword)
768 {
769     static const char * const selections[] = {
770         "distance from cog of resnr 1 < 2"
771     };
772     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
773     runTest("simple.gro", selections);
774 }
775
776
777 TEST_F(SelectionCollectionDataTest, HandlesMinDistanceKeyword)
778 {
779     static const char * const selections[] = {
780         "mindistance from resnr 1 < 2"
781     };
782     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
783     runTest("simple.gro", selections);
784 }
785
786
787 TEST_F(SelectionCollectionDataTest, HandlesWithinKeyword)
788 {
789     static const char * const selections[] = {
790         "within 1 of resnr 2"
791     };
792     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
793     runTest("simple.gro", selections);
794 }
795
796
797 TEST_F(SelectionCollectionDataTest, HandlesInSolidAngleKeyword)
798 {
799     // Both of these should evaluate to empty on a correct implementation.
800     static const char * const selections[] = {
801         "resname TP and not insolidangle center cog of resname C span resname R cutoff 20",
802         "resname TN and insolidangle center cog of resname C span resname R cutoff 20"
803     };
804     setFlags(TestFlags() | efDontTestCompiledAtoms | efTestEvaluation);
805     runTest("sphere.gro", selections);
806 }
807
808
809 TEST_F(SelectionCollectionDataTest, HandlesPermuteModifier)
810 {
811     static const char * const selections[] = {
812         "all permute 3 1 2",
813         "res_cog of resnr 1 to 4 permute 2 1",
814         "name CB S1 and res_cog x < 3 permute 2 1"
815     };
816     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
817              | efTestPositionAtoms | efTestPositionMapping);
818     runTest("simple.gro", selections);
819 }
820
821
822 TEST_F(SelectionCollectionDataTest, HandlesPlusModifier)
823 {
824     static const char * const selections[] = {
825         "name S2 plus name S1",
826         "res_cog of resnr 2 plus res_cog of resnr 1 plus res_cog of resnr 3",
827         "name S1 and y < 3 plus res_cog of x < 2.5"
828     };
829     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
830              | efTestPositionAtoms | efTestPositionMapping);
831     runTest("simple.gro", selections);
832 }
833
834
835 TEST_F(SelectionCollectionDataTest, HandlesMergeModifier)
836 {
837     static const char * const selections[] = {
838         "name S2 merge name S1",
839         "resnr 1 2 and name S2 merge resnr 1 2 and name S1 merge res_cog of resnr 1 2",
840         "name S1 and x < 2.5 merge res_cog of x < 2.5"
841     };
842     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
843              | efTestPositionAtoms | efTestPositionMapping);
844     runTest("simple.gro", selections);
845 }
846
847
848 /********************************************************************
849  * Tests for generic selection evaluation
850  */
851
852 TEST_F(SelectionCollectionDataTest, ComputesMassesAndCharges)
853 {
854     static const char * const selections[] = {
855         "name CB",
856         "y > 2",
857         "res_cog of y > 2"
858     };
859     setFlags(TestFlags() | efTestEvaluation | efTestPositionAtoms
860              | efTestPositionMasses | efTestPositionCharges);
861     ASSERT_NO_FATAL_FAILURE(runParser(selections));
862     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
863     for (int i = 0; i < top_->atoms.nr; ++i)
864     {
865         top_->atoms.atom[i].m =   1.0 + i / 100.0;
866         top_->atoms.atom[i].q = -(1.0 + i / 100.0);
867     }
868     ASSERT_NO_FATAL_FAILURE(runCompiler());
869     ASSERT_NO_FATAL_FAILURE(runEvaluate());
870     ASSERT_NO_FATAL_FAILURE(runEvaluateFinal());
871 }
872
873 TEST_F(SelectionCollectionDataTest, ComputesMassesAndChargesWithoutTopology)
874 {
875     static const char * const selections[] = {
876         "atomnr 1 to 3 8 to 9",
877         "y > 2",
878         "cog of (y > 2)"
879     };
880     setFlags(TestFlags() | efTestPositionAtoms
881              | efTestPositionMasses | efTestPositionCharges);
882     runTest(10, selections);
883 }
884
885
886 /********************************************************************
887  * Tests for selection syntactic constructs
888  */
889
890 TEST_F(SelectionCollectionDataTest, HandlesSelectionNames)
891 {
892     static const char * const selections[] = {
893         "\"GroupSelection\" group \"GrpA\"",
894         "\"DynamicSelection\" x < 5",
895         "y < 3"
896     };
897     setFlags(TestFlags() | efTestSelectionNames);
898     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
899     runTest(10, selections);
900 }
901
902 TEST_F(SelectionCollectionDataTest, HandlesIndexGroupsInSelections)
903 {
904     static const char * const selections[] = {
905         "group \"GrpA\"",
906         "GrpB",
907         "1",
908         "group \"GrpB\" and resname RB"
909     };
910     setFlags(TestFlags() | efTestSelectionNames);
911     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
912     runTest("simple.gro", selections);
913 }
914
915 TEST_F(SelectionCollectionDataTest, HandlesIndexGroupsInSelectionsDelayed)
916 {
917     static const char * const selections[] = {
918         "group \"GrpA\"",
919         "GrpB",
920         "1",
921         "group \"GrpB\" and resname RB"
922     };
923     setFlags(TestFlags() | efTestSelectionNames);
924     ASSERT_NO_FATAL_FAILURE(runParser(selections));
925     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
926     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
927     ASSERT_NO_FATAL_FAILURE(runCompiler());
928 }
929
930 TEST_F(SelectionCollectionDataTest, HandlesConstantPositions)
931 {
932     static const char * const selections[] = {
933         "[1, -2, 3.5]"
934     };
935     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
936     runTest("simple.gro", selections);
937 }
938
939
940 TEST_F(SelectionCollectionDataTest, HandlesWithinConstantPositions)
941 {
942     static const char * const selections[] = {
943         "within 1 of [2, 1, 0]"
944     };
945     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
946     runTest("simple.gro", selections);
947 }
948
949
950 TEST_F(SelectionCollectionDataTest, HandlesForcedStringMatchingMode)
951 {
952     static const char * const selections[] = {
953         "name = S1 \"C?\"",
954         "name ? S1 \"C?\""
955     };
956     runTest("simple.gro", selections);
957 }
958
959
960 TEST_F(SelectionCollectionDataTest, HandlesWildcardMatching)
961 {
962     static const char * const selections[] = {
963         "name \"S?\"",
964         "name ? \"S?\""
965     };
966     runTest("simple.gro", selections);
967 }
968
969
970 TEST_F(SelectionCollectionDataTest, HandlesRegexMatching)
971 {
972     static const char * const selections[] = {
973         "resname \"R[BD]\"",
974         "resname ~ \"R[BD]\""
975     };
976     if (gmx::Regex::isSupported())
977     {
978         runTest("simple.gro", selections);
979     }
980 }
981
982
983 TEST_F(SelectionCollectionDataTest, HandlesBasicBoolean)
984 {
985     static const char * const selections[] = {
986         "atomnr 1 to 5 and atomnr 2 to 7",
987         "atomnr 1 to 5 or not atomnr 3 to 8",
988         "not not atomnr 1 to 5 and atomnr 2 to 6 and not not atomnr 3 to 7",
989         "atomnr 1 to 5 and (atomnr 2 to 7 and atomnr 3 to 6)",
990         "x < 5 and atomnr 1 to 5 and y < 3 and atomnr 2 to 4"
991     };
992     runTest(10, selections);
993 }
994
995
996 TEST_F(SelectionCollectionDataTest, HandlesDynamicAtomValuedParameters)
997 {
998     static const char * const selections[] = {
999         "same residue as (atomnr 3 5 13 or y > 5)",
1000         "(resnr 1 3 5 or x > 10) and same residue as (atomnr 3 5 13 or z > 5)"
1001     };
1002     setFlags(TestFlags() | efTestEvaluation);
1003     runTest("simple.gro", selections);
1004 }
1005
1006
1007 TEST_F(SelectionCollectionDataTest, HandlesEmptySelectionWithUnevaluatedExpressions)
1008 {
1009     static const char * const selections[] = {
1010         "none and x > 2",
1011         "none and same resname as resnr 2"
1012     };
1013     runTest("simple.gro", selections);
1014 }
1015
1016
1017 TEST_F(SelectionCollectionDataTest, HandlesNumericComparisons)
1018 {
1019     static const char * const selections[] = {
1020         "x > 2",
1021         "2 < x",
1022         "y > resnr",
1023         "resnr < 2.5",
1024         "2.5 > resnr"
1025     };
1026     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1027     runTest("simple.gro", selections);
1028 }
1029
1030
1031 TEST_F(SelectionCollectionDataTest, HandlesArithmeticExpressions)
1032 {
1033     static const char * const selections[] = {
1034         "x+1 > 3",
1035         "(y-1)^2 <= 1",
1036         "x+--1 > 3",
1037         "-x+-1 < -3"
1038     };
1039     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1040     runTest("simple.gro", selections);
1041 }
1042
1043
1044 TEST_F(SelectionCollectionDataTest, HandlesNumericVariables)
1045 {
1046     static const char * const selections[] = {
1047         "value = x + y",
1048         "value <= 4",
1049         "index = resnr",
1050         "index < 3"
1051     };
1052     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1053     runTest("simple.gro", selections);
1054 }
1055
1056
1057 TEST_F(SelectionCollectionDataTest, HandlesComplexNumericVariables)
1058 {
1059     static const char * const selections[] = {
1060         "value = x + y",
1061         "resname RA and value <= 4",
1062         "resname RA RB and x < 3 and value <= 4",
1063         "index = atomnr",
1064         "resname RA and index < 3",
1065         "resname RB and y < 3 and index < 6"
1066     };
1067     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1068     runTest("simple.gro", selections);
1069 }
1070
1071
1072 TEST_F(SelectionCollectionDataTest, HandlesPositionVariables)
1073 {
1074     static const char * const selections[] = {
1075         "foo = res_cog of resname RA",
1076         "foo",
1077         "within 1 of foo",
1078         "bar = cog of resname RA",
1079         "bar",
1080         "within 1 of bar"
1081     };
1082     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1083     runTest("simple.gro", selections);
1084 }
1085
1086
1087 TEST_F(SelectionCollectionDataTest, HandlesConstantPositionInVariable)
1088 {
1089     static const char * const selections[] = {
1090         "constpos = [1.0, 2.5, 0.5]",
1091         "constpos",
1092         "within 2 of constpos"
1093     };
1094     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
1095              | efTestPositionAtoms);
1096     runTest("simple.gro", selections);
1097 }
1098
1099
1100 TEST_F(SelectionCollectionDataTest, HandlesNumericConstantsInVariables)
1101 {
1102     static const char * const selections[] = {
1103         "constint = 4",
1104         "constreal1 = 0.5",
1105         "constreal2 = 2.7",
1106         "resnr < constint",
1107         "x + constreal1 < constreal2"
1108     };
1109     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1110     runTest("simple.gro", selections);
1111 }
1112
1113
1114 /********************************************************************
1115  * Tests for complex boolean syntax
1116  */
1117
1118 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysis)
1119 {
1120     static const char * const selections[] = {
1121         "atomnr 1 to 5 and atomnr 2 to 7 and x < 2",
1122         "atomnr 1 to 5 and (atomnr 4 to 7 or x < 2)",
1123         "atomnr 1 to 5 and y < 3 and (atomnr 4 to 7 or x < 2)",
1124         "atomnr 1 to 5 and not (atomnr 4 to 7 or x < 2)",
1125         "atomnr 1 to 5 or (atomnr 4 to 6 and (atomnr 5 to 7 or x < 2))"
1126     };
1127     runTest(10, selections);
1128 }
1129
1130
1131 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysisWithVariables)
1132 {
1133     static const char * const selections[] = {
1134         "foo = atomnr 4 to 7 or x < 2",
1135         "atomnr 1 to 4 and foo",
1136         "atomnr 2 to 6 and y < 3 and foo",
1137         "atomnr 6 to 10 and not foo"
1138     };
1139     runTest(10, selections);
1140 }
1141
1142
1143 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysisWithMoreVariables)
1144 {
1145     static const char * const selections[] = {
1146         "foo = atomnr 4 to 7",
1147         "bar = foo and x < 2",
1148         "bar2 = foo and y < 2",
1149         "atomnr 1 to 4 and bar",
1150         "atomnr 2 to 6 and y < 3 and bar2",
1151         "atomnr 6 to 10 and not foo"
1152     };
1153     runTest(10, selections);
1154 }
1155
1156
1157 /********************************************************************
1158  * Tests for complex subexpression cases
1159  *
1160  * These tests use some knowledge of the implementation to trigger different
1161  * paths in the code.
1162  */
1163
1164 TEST_F(SelectionCollectionDataTest, HandlesUnusedVariables)
1165 {
1166     static const char * const selections[] = {
1167         "unused1 = atomnr 1 to 3",
1168         "foo = atomnr 4 to 7",
1169         "atomnr 1 to 6 and foo",
1170         "unused2 = atomnr 3 to 5"
1171     };
1172     runTest(10, selections);
1173 }
1174
1175
1176 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithStaticEvaluationGroups)
1177 {
1178     static const char * const selections[] = {
1179         "foo = atomnr 4 to 7 and x < 2",
1180         "atomnr 1 to 5 and foo",
1181         "atomnr 3 to 7 and foo"
1182     };
1183     runTest(10, selections);
1184 }
1185
1186
1187 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithMixedEvaluationGroups)
1188 {
1189     static const char * const selections[] = {
1190         "foo = atomnr 4 to 7 and x < 2",
1191         "atomnr 1 to 6 and foo",
1192         "within 1 of foo",
1193         "foo"
1194     };
1195     runTest(10, selections);
1196 }
1197
1198
1199 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithMixedEvaluationGroups2)
1200 {
1201     static const char * const selections[] = {
1202         "foo = atomnr 1 to 8 and x < 10",
1203         "atomnr 1 to 5 and y < 10 and foo",
1204         "foo"
1205     };
1206     setFlags(TestFlags() | efTestEvaluation);
1207     runTest("simple.gro", selections);
1208 }
1209
1210
1211 } // namespace