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