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