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