Merge "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/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             if (flags_.test(efTestSelectionNames))
294             {
295                 selcompound.checkString(sel_[count_].name(), "Name");
296             }
297             selcompound.checkString(sel_[count_].selectionText(), "Text");
298             selcompound.checkBoolean(sel_[count_].isDynamic(), "Dynamic");
299             ++count_;
300         }
301     }
302 }
303
304
305 void
306 SelectionCollectionDataTest::runCompiler()
307 {
308     ASSERT_NO_THROW_GMX(sc_.compile());
309     ASSERT_EQ(count_, sel_.size());
310     checkCompiled();
311 }
312
313
314 void
315 SelectionCollectionDataTest::checkCompiled()
316 {
317     using gmx::test::TestReferenceChecker;
318     const TestFlags      mask = ~TestFlags(efTestPositionCoordinates);
319
320     TestReferenceChecker compound(checker_.checkCompound("CompiledSelections", "Compiled"));
321     for (size_t i = 0; i < count_; ++i)
322     {
323         SCOPED_TRACE(std::string("Checking selection \"") +
324                      sel_[i].selectionText() + "\"");
325         std::string          id = gmx::formatString("Selection%d", static_cast<int>(i + 1));
326         TestReferenceChecker selcompound(
327                 compound.checkCompound("Selection", id.c_str()));
328         if (flags_.test(efTestSelectionNames))
329         {
330             selcompound.checkString(sel_[i].name(), "Name");
331         }
332         if (!flags_.test(efDontTestCompiledAtoms))
333         {
334             checkSelection(&selcompound, sel_[i], flags_ & mask);
335         }
336     }
337 }
338
339
340 void
341 SelectionCollectionDataTest::runEvaluate()
342 {
343     using gmx::test::TestReferenceChecker;
344
345     ++framenr_;
346     ASSERT_NO_THROW_GMX(sc_.evaluate(frame_, NULL));
347     std::string          frame = gmx::formatString("Frame%d", framenr_);
348     TestReferenceChecker compound(
349             checker_.checkCompound("EvaluatedSelections", frame.c_str()));
350     for (size_t i = 0; i < count_; ++i)
351     {
352         SCOPED_TRACE(std::string("Checking selection \"") +
353                      sel_[i].selectionText() + "\"");
354         std::string          id = gmx::formatString("Selection%d", static_cast<int>(i + 1));
355         TestReferenceChecker selcompound(
356                 compound.checkCompound("Selection", id.c_str()));
357         checkSelection(&selcompound, sel_[i], flags_);
358     }
359 }
360
361
362 void
363 SelectionCollectionDataTest::runEvaluateFinal()
364 {
365     ASSERT_NO_THROW_GMX(sc_.evaluateFinal(framenr_));
366     checkCompiled();
367 }
368
369
370 void
371 SelectionCollectionDataTest::runTest(int natoms, const char * const *selections,
372                                      size_t count)
373 {
374     ASSERT_NO_FATAL_FAILURE(runParser(selections, count));
375     ASSERT_NO_FATAL_FAILURE(setAtomCount(natoms));
376     ASSERT_NO_FATAL_FAILURE(runCompiler());
377 }
378
379
380 void
381 SelectionCollectionDataTest::runTest(const char         *filename,
382                                      const char * const *selections,
383                                      size_t              count)
384 {
385     ASSERT_NO_FATAL_FAILURE(runParser(selections, count));
386     ASSERT_NO_FATAL_FAILURE(loadTopology(filename));
387     ASSERT_NO_FATAL_FAILURE(runCompiler());
388     if (flags_.test(efTestEvaluation))
389     {
390         ASSERT_NO_FATAL_FAILURE(runEvaluate());
391         ASSERT_NO_FATAL_FAILURE(runEvaluateFinal());
392     }
393 }
394
395
396 /********************************************************************
397  * Tests for SelectionCollection functionality without reference data
398  */
399
400 TEST_F(SelectionCollectionTest, HandlesNoSelections)
401 {
402     EXPECT_FALSE(sc_.requiresTopology());
403     EXPECT_NO_THROW_GMX(sc_.compile());
404 }
405
406 TEST_F(SelectionCollectionTest, HandlesVelocityAndForceRequests)
407 {
408     ASSERT_NO_THROW_GMX(sel_ = sc_.parseFromString("atomnr 1 to 10; none"));
409     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
410     ASSERT_EQ(2U, sel_.size());
411     ASSERT_NO_THROW_GMX(sel_[0].setEvaluateVelocities(true));
412     ASSERT_NO_THROW_GMX(sel_[1].setEvaluateVelocities(true));
413     ASSERT_NO_THROW_GMX(sel_[0].setEvaluateForces(true));
414     ASSERT_NO_THROW_GMX(sel_[1].setEvaluateForces(true));
415     ASSERT_NO_THROW_GMX(sc_.compile());
416     EXPECT_TRUE(sel_[0].hasVelocities());
417     EXPECT_TRUE(sel_[1].hasVelocities());
418     EXPECT_TRUE(sel_[0].hasForces());
419     EXPECT_TRUE(sel_[1].hasForces());
420 }
421
422 TEST_F(SelectionCollectionTest, ParsesSelectionsFromFile)
423 {
424     ASSERT_NO_THROW_GMX(sel_ = sc_.parseFromFile(
425                                     gmx::test::TestFileManager::getInputFilePath("selfile.dat")));
426     // These should match the contents of selfile.dat
427     ASSERT_EQ(2U, sel_.size());
428     EXPECT_STREQ("resname RA RB", sel_[0].selectionText());
429     EXPECT_STREQ("resname RB RC", sel_[1].selectionText());
430 }
431
432 TEST_F(SelectionCollectionTest, HandlesInvalidRegularExpressions)
433 {
434     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
435     EXPECT_THROW_GMX({
436                          sc_.parseFromString("resname ~ \"R[A\"");
437                          sc_.compile();
438                      }, gmx::InvalidInputError);
439 }
440
441 TEST_F(SelectionCollectionTest, HandlesUnsupportedRegularExpressions)
442 {
443     if (!gmx::Regex::isSupported())
444     {
445         ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
446         EXPECT_THROW_GMX({
447                              sc_.parseFromString("resname \"R[AD]\"");
448                              sc_.compile();
449                          }, gmx::InvalidInputError);
450     }
451 }
452
453 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue)
454 {
455     EXPECT_THROW_GMX(sc_.parseFromString("mindist from atomnr 1 cutoff"),
456                      gmx::InvalidInputError);
457 }
458
459 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue2)
460 {
461     EXPECT_THROW_GMX(sc_.parseFromString("within 1 of"),
462                      gmx::InvalidInputError);
463 }
464
465 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue3)
466 {
467     EXPECT_THROW_GMX(sc_.parseFromString("within of atomnr 1"),
468                      gmx::InvalidInputError);
469 }
470
471 TEST_F(SelectionCollectionTest, HandlesHelpKeywordInInvalidContext)
472 {
473     EXPECT_THROW_GMX(sc_.parseFromString("resname help"),
474                      gmx::InvalidInputError);
475 }
476
477 // TODO: Tests for more parser errors
478
479 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceParser1)
480 {
481     ASSERT_NO_THROW_GMX(sc_.setIndexGroups(NULL));
482     EXPECT_THROW_GMX(sc_.parseFromString("group \"foo\""), gmx::InconsistentInputError);
483     EXPECT_THROW_GMX(sc_.parseFromString("4"), gmx::InconsistentInputError);
484 }
485
486 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceParser2)
487 {
488     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
489     EXPECT_THROW_GMX(sc_.parseFromString("group \"foo\""), gmx::InconsistentInputError);
490     EXPECT_THROW_GMX(sc_.parseFromString("4"), gmx::InconsistentInputError);
491 }
492
493 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceDelayed1)
494 {
495     ASSERT_NO_THROW_GMX(sc_.parseFromString("group \"foo\""));
496     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
497     EXPECT_THROW_GMX(sc_.setIndexGroups(NULL), gmx::InconsistentInputError);
498     EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
499 }
500
501 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceDelayed2)
502 {
503     ASSERT_NO_THROW_GMX(sc_.parseFromString("group 4; group \"foo\""));
504     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
505     EXPECT_THROW_GMX(loadIndexGroups("simple.ndx"), gmx::InconsistentInputError);
506     EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
507 }
508
509 TEST_F(SelectionCollectionTest, HandlesUnsortedGroupReference)
510 {
511     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
512     EXPECT_THROW_GMX(sc_.parseFromString("group \"GrpUnsorted\""),
513                      gmx::InconsistentInputError);
514     EXPECT_THROW_GMX(sc_.parseFromString("2"), gmx::InconsistentInputError);
515 }
516
517 TEST_F(SelectionCollectionTest, HandlesUnsortedGroupReferenceDelayed)
518 {
519     ASSERT_NO_THROW_GMX(sc_.parseFromString("group 2; group \"GrpUnsorted\""));
520     EXPECT_THROW_GMX(loadIndexGroups("simple.ndx"), gmx::InconsistentInputError);
521     EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
522 }
523
524 TEST_F(SelectionCollectionTest, RecoversFromMissingMoleculeInfo)
525 {
526     ASSERT_NO_THROW_GMX(sc_.parseFromString("molindex 1 to 5"));
527     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
528     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
529 }
530
531 TEST_F(SelectionCollectionTest, RecoversFromMissingAtomTypes)
532 {
533     ASSERT_NO_THROW_GMX(sc_.parseFromString("type CA"));
534     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
535     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
536 }
537
538 TEST_F(SelectionCollectionTest, RecoversFromMissingPDBInfo)
539 {
540     ASSERT_NO_THROW_GMX(sc_.parseFromString("altloc A"));
541     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
542     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
543 }
544
545 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation)
546 {
547     ASSERT_NO_THROW_GMX(sc_.parseFromString("all permute 1 1"));
548     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
549     EXPECT_THROW_GMX(sc_.compile(), gmx::InvalidInputError);
550 }
551
552 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation2)
553 {
554     ASSERT_NO_THROW_GMX(sc_.parseFromString("all permute 3 2 1"));
555     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
556     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
557 }
558
559 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation3)
560 {
561     ASSERT_NO_THROW_GMX(sc_.parseFromString("x < 1.5 permute 3 2 1"));
562     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
563     ASSERT_NO_THROW_GMX(sc_.compile());
564     EXPECT_THROW_GMX(sc_.evaluate(frame_, NULL), gmx::InconsistentInputError);
565 }
566
567 // TODO: Tests for evaluation errors
568
569
570 /********************************************************************
571  * Tests for selection keywords
572  */
573
574 TEST_F(SelectionCollectionDataTest, HandlesAllNone)
575 {
576     static const char * const selections[] = {
577         "all",
578         "none"
579     };
580     runTest(10, selections);
581 }
582
583 TEST_F(SelectionCollectionDataTest, HandlesAtomnr)
584 {
585     static const char * const selections[] = {
586         "atomnr 1 to 3 6 to 8",
587         "atomnr 4 2 5 to 7",
588         "atomnr <= 5"
589     };
590     runTest(10, selections);
591 }
592
593 TEST_F(SelectionCollectionDataTest, HandlesResnr)
594 {
595     static const char * const selections[] = {
596         "resnr 1 2 5",
597         "resid 4 to 3"
598     };
599     runTest("simple.gro", selections);
600 }
601
602 TEST_F(SelectionCollectionDataTest, HandlesResIndex)
603 {
604     static const char * const selections[] = {
605         "resindex 1 4",
606         "residue 1 3"
607     };
608     runTest("simple.pdb", selections);
609 }
610
611 TEST_F(SelectionCollectionDataTest, HandlesMolIndex)
612 {
613     static const char * const selections[] = {
614         "molindex 1 4",
615         "molecule 2 3 5"
616     };
617     ASSERT_NO_FATAL_FAILURE(runParser(selections));
618     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
619     topManager_.initUniformMolecules(3);
620     ASSERT_NO_FATAL_FAILURE(runCompiler());
621 }
622
623 TEST_F(SelectionCollectionDataTest, HandlesAtomname)
624 {
625     static const char * const selections[] = {
626         "name CB",
627         "atomname S1 S2"
628     };
629     runTest("simple.gro", selections);
630 }
631
632 TEST_F(SelectionCollectionDataTest, HandlesPdbAtomname)
633 {
634     static const char * const selections[] = {
635         "name HG21",
636         "name 1HG2",
637         "pdbname HG21 CB",
638         "pdbatomname 1HG2"
639     };
640     runTest("simple.pdb", selections);
641 }
642
643
644 TEST_F(SelectionCollectionDataTest, HandlesAtomtype)
645 {
646     static const char * const selections[] = {
647         "atomtype CA"
648     };
649     ASSERT_NO_FATAL_FAILURE(runParser(selections));
650     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
651     const char *const types[] = { "CA", "SA", "SB" };
652     topManager_.initAtomTypes(types);
653     ASSERT_NO_FATAL_FAILURE(runCompiler());
654 }
655
656 TEST_F(SelectionCollectionDataTest, HandlesChain)
657 {
658     static const char * const selections[] = {
659         "chain A",
660         "chain B"
661     };
662     runTest("simple.pdb", selections);
663 }
664
665 TEST_F(SelectionCollectionDataTest, HandlesMass)
666 {
667     static const char * const selections[] = {
668         "mass > 5"
669     };
670     ASSERT_NO_FATAL_FAILURE(runParser(selections));
671     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
672     for (int i = 0; i < top_->atoms.nr; ++i)
673     {
674         top_->atoms.atom[i].m = 1.0 + i;
675     }
676     ASSERT_NO_FATAL_FAILURE(runCompiler());
677 }
678
679 TEST_F(SelectionCollectionDataTest, HandlesCharge)
680 {
681     static const char * const selections[] = {
682         "charge < 0.5"
683     };
684     ASSERT_NO_FATAL_FAILURE(runParser(selections));
685     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
686     for (int i = 0; i < top_->atoms.nr; ++i)
687     {
688         top_->atoms.atom[i].q = i / 10.0;
689     }
690     ASSERT_NO_FATAL_FAILURE(runCompiler());
691 }
692
693 TEST_F(SelectionCollectionDataTest, HandlesAltLoc)
694 {
695     static const char * const selections[] = {
696         "altloc \" \"",
697         "altloc A"
698     };
699     runTest("simple.pdb", selections);
700 }
701
702 TEST_F(SelectionCollectionDataTest, HandlesInsertCode)
703 {
704     static const char * const selections[] = {
705         "insertcode \" \"",
706         "insertcode A"
707     };
708     runTest("simple.pdb", selections);
709 }
710
711 TEST_F(SelectionCollectionDataTest, HandlesOccupancy)
712 {
713     static const char * const selections[] = {
714         "occupancy 1",
715         "occupancy < .5"
716     };
717     runTest("simple.pdb", selections);
718 }
719
720 TEST_F(SelectionCollectionDataTest, HandlesBeta)
721 {
722     static const char * const selections[] = {
723         "beta 0",
724         "beta >= 0.3"
725     };
726     runTest("simple.pdb", selections);
727 }
728
729 TEST_F(SelectionCollectionDataTest, HandlesResname)
730 {
731     static const char * const selections[] = {
732         "resname RA",
733         "resname RB RC"
734     };
735     runTest("simple.gro", selections);
736 }
737
738 TEST_F(SelectionCollectionDataTest, HandlesCoordinateKeywords)
739 {
740     static const char * const selections[] = {
741         "x < 3",
742         "y >= 3",
743         "x {-1 to 2}"
744     };
745     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
746     runTest("simple.gro", selections);
747 }
748
749
750 TEST_F(SelectionCollectionDataTest, HandlesSameResidue)
751 {
752     static const char * const selections[] = {
753         "same residue as atomnr 1 4 12"
754     };
755     runTest("simple.gro", selections);
756 }
757
758
759 TEST_F(SelectionCollectionDataTest, HandlesSameResidueName)
760 {
761     static const char * const selections[] = {
762         "same resname as atomnr 1 14"
763     };
764     runTest("simple.gro", selections);
765 }
766
767
768 TEST_F(SelectionCollectionDataTest, HandlesPositionKeywords)
769 {
770     static const char * const selections[] = {
771         "cog of resnr 1 3",
772         "res_cog of name CB and resnr 1 3",
773         "whole_res_cog of name CB and resnr 1 3",
774         "part_res_cog of x < 3",
775         "dyn_res_cog of x < 3"
776     };
777     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
778              | efTestPositionAtoms);
779     runTest("simple.gro", selections);
780 }
781
782
783 TEST_F(SelectionCollectionDataTest, HandlesDistanceKeyword)
784 {
785     static const char * const selections[] = {
786         "distance from cog of resnr 1 < 2"
787     };
788     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
789     runTest("simple.gro", selections);
790 }
791
792
793 TEST_F(SelectionCollectionDataTest, HandlesMinDistanceKeyword)
794 {
795     static const char * const selections[] = {
796         "mindistance from resnr 1 < 2"
797     };
798     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
799     runTest("simple.gro", selections);
800 }
801
802
803 TEST_F(SelectionCollectionDataTest, HandlesWithinKeyword)
804 {
805     static const char * const selections[] = {
806         "within 1 of resnr 2"
807     };
808     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
809     runTest("simple.gro", selections);
810 }
811
812
813 TEST_F(SelectionCollectionDataTest, HandlesInSolidAngleKeyword)
814 {
815     // Both of these should evaluate to empty on a correct implementation.
816     static const char * const selections[] = {
817         "resname TP and not insolidangle center cog of resname C span resname R cutoff 20",
818         "resname TN and insolidangle center cog of resname C span resname R cutoff 20"
819     };
820     setFlags(TestFlags() | efDontTestCompiledAtoms | efTestEvaluation);
821     runTest("sphere.gro", selections);
822 }
823
824
825 TEST_F(SelectionCollectionDataTest, HandlesPermuteModifier)
826 {
827     static const char * const selections[] = {
828         "all permute 3 1 2",
829         "res_cog of resnr 1 to 4 permute 2 1",
830         "name CB S1 and res_cog x < 3 permute 2 1"
831     };
832     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
833              | efTestPositionAtoms | efTestPositionMapping);
834     runTest("simple.gro", selections);
835 }
836
837
838 TEST_F(SelectionCollectionDataTest, HandlesPlusModifier)
839 {
840     static const char * const selections[] = {
841         "name S2 plus name S1",
842         "res_cog of resnr 2 plus res_cog of resnr 1 plus res_cog of resnr 3",
843         "name S1 and y < 3 plus res_cog of x < 2.5"
844     };
845     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
846              | efTestPositionAtoms | efTestPositionMapping);
847     runTest("simple.gro", selections);
848 }
849
850
851 TEST_F(SelectionCollectionDataTest, HandlesMergeModifier)
852 {
853     static const char * const selections[] = {
854         "name S2 merge name S1",
855         "resnr 1 2 and name S2 merge resnr 1 2 and name S1 merge res_cog of resnr 1 2",
856         "name S1 and x < 2.5 merge res_cog of x < 2.5"
857     };
858     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
859              | efTestPositionAtoms | efTestPositionMapping);
860     runTest("simple.gro", selections);
861 }
862
863
864 /********************************************************************
865  * Tests for generic selection evaluation
866  */
867
868 TEST_F(SelectionCollectionDataTest, ComputesMassesAndCharges)
869 {
870     static const char * const selections[] = {
871         "name CB",
872         "y > 2",
873         "res_cog of y > 2"
874     };
875     setFlags(TestFlags() | efTestEvaluation | efTestPositionAtoms
876              | efTestPositionMasses | efTestPositionCharges);
877     ASSERT_NO_FATAL_FAILURE(runParser(selections));
878     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
879     for (int i = 0; i < top_->atoms.nr; ++i)
880     {
881         top_->atoms.atom[i].m =   1.0 + i / 100.0;
882         top_->atoms.atom[i].q = -(1.0 + i / 100.0);
883     }
884     ASSERT_NO_FATAL_FAILURE(runCompiler());
885     ASSERT_NO_FATAL_FAILURE(runEvaluate());
886     ASSERT_NO_FATAL_FAILURE(runEvaluateFinal());
887 }
888
889 TEST_F(SelectionCollectionDataTest, ComputesMassesAndChargesWithoutTopology)
890 {
891     static const char * const selections[] = {
892         "atomnr 1 to 3 8 to 9",
893         "y > 2",
894         "cog of (y > 2)"
895     };
896     setFlags(TestFlags() | efTestPositionAtoms
897              | efTestPositionMasses | efTestPositionCharges);
898     runTest(10, selections);
899 }
900
901
902 /********************************************************************
903  * Tests for selection syntactic constructs
904  */
905
906 TEST_F(SelectionCollectionDataTest, HandlesSelectionNames)
907 {
908     static const char * const selections[] = {
909         "\"GroupSelection\" group \"GrpA\"",
910         "\"DynamicSelection\" x < 5",
911         "y < 3"
912     };
913     setFlags(TestFlags() | efTestSelectionNames);
914     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
915     runTest(10, selections);
916 }
917
918 TEST_F(SelectionCollectionDataTest, HandlesIndexGroupsInSelections)
919 {
920     static const char * const selections[] = {
921         "group \"GrpA\"",
922         "GrpB",
923         "1",
924         "group \"GrpB\" and resname RB"
925     };
926     setFlags(TestFlags() | efTestSelectionNames);
927     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
928     runTest("simple.gro", selections);
929 }
930
931 TEST_F(SelectionCollectionDataTest, HandlesIndexGroupsInSelectionsDelayed)
932 {
933     static const char * const selections[] = {
934         "group \"GrpA\"",
935         "GrpB",
936         "1",
937         "group \"GrpB\" and resname RB"
938     };
939     setFlags(TestFlags() | efTestSelectionNames);
940     ASSERT_NO_FATAL_FAILURE(runParser(selections));
941     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
942     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
943     ASSERT_NO_FATAL_FAILURE(runCompiler());
944 }
945
946 TEST_F(SelectionCollectionDataTest, HandlesConstantPositions)
947 {
948     static const char * const selections[] = {
949         "[1, -2, 3.5]"
950     };
951     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
952     runTest("simple.gro", selections);
953 }
954
955
956 TEST_F(SelectionCollectionDataTest, HandlesWithinConstantPositions)
957 {
958     static const char * const selections[] = {
959         "within 1 of [2, 1, 0]"
960     };
961     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
962     runTest("simple.gro", selections);
963 }
964
965
966 TEST_F(SelectionCollectionDataTest, HandlesForcedStringMatchingMode)
967 {
968     static const char * const selections[] = {
969         "name = S1 \"C?\"",
970         "name ? S1 \"C?\""
971     };
972     runTest("simple.gro", selections);
973 }
974
975
976 TEST_F(SelectionCollectionDataTest, HandlesWildcardMatching)
977 {
978     static const char * const selections[] = {
979         "name \"S?\"",
980         "name ? \"S?\""
981     };
982     runTest("simple.gro", selections);
983 }
984
985
986 TEST_F(SelectionCollectionDataTest, HandlesRegexMatching)
987 {
988     static const char * const selections[] = {
989         "resname \"R[BD]\"",
990         "resname ~ \"R[BD]\""
991     };
992     if (gmx::Regex::isSupported())
993     {
994         runTest("simple.gro", selections);
995     }
996 }
997
998
999 TEST_F(SelectionCollectionDataTest, HandlesBasicBoolean)
1000 {
1001     static const char * const selections[] = {
1002         "atomnr 1 to 5 and atomnr 2 to 7",
1003         "atomnr 1 to 5 or not atomnr 3 to 8",
1004         "not not atomnr 1 to 5 and atomnr 2 to 6 and not not atomnr 3 to 7",
1005         "atomnr 1 to 5 and (atomnr 2 to 7 and atomnr 3 to 6)",
1006         "x < 5 and atomnr 1 to 5 and y < 3 and atomnr 2 to 4"
1007     };
1008     runTest(10, selections);
1009 }
1010
1011
1012 TEST_F(SelectionCollectionDataTest, HandlesDynamicAtomValuedParameters)
1013 {
1014     static const char * const selections[] = {
1015         "same residue as (atomnr 3 5 13 or y > 5)",
1016         "(resnr 1 3 5 or x > 10) and same residue as (atomnr 3 5 13 or z > 5)"
1017     };
1018     setFlags(TestFlags() | efTestEvaluation);
1019     runTest("simple.gro", selections);
1020 }
1021
1022
1023 TEST_F(SelectionCollectionDataTest, HandlesNumericComparisons)
1024 {
1025     static const char * const selections[] = {
1026         "x > 2",
1027         "2 < x",
1028         "y > resnr",
1029         "resnr < 2.5",
1030         "2.5 > resnr"
1031     };
1032     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1033     runTest("simple.gro", selections);
1034 }
1035
1036
1037 TEST_F(SelectionCollectionDataTest, HandlesArithmeticExpressions)
1038 {
1039     static const char * const selections[] = {
1040         "x+1 > 3",
1041         "(y-1)^2 <= 1",
1042         "x+--1 > 3",
1043         "-x+-1 < -3"
1044     };
1045     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1046     runTest("simple.gro", selections);
1047 }
1048
1049
1050 TEST_F(SelectionCollectionDataTest, HandlesNumericVariables)
1051 {
1052     static const char * const selections[] = {
1053         "value = x + y",
1054         "value <= 4",
1055         "index = resnr",
1056         "index < 3"
1057     };
1058     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1059     runTest("simple.gro", selections);
1060 }
1061
1062
1063 TEST_F(SelectionCollectionDataTest, HandlesComplexNumericVariables)
1064 {
1065     static const char * const selections[] = {
1066         "value = x + y",
1067         "resname RA and value <= 4",
1068         "resname RA RB and x < 3 and value <= 4",
1069         "index = atomnr",
1070         "resname RA and index < 3",
1071         "resname RB and y < 3 and index < 6"
1072     };
1073     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1074     runTest("simple.gro", selections);
1075 }
1076
1077
1078 TEST_F(SelectionCollectionDataTest, HandlesPositionVariables)
1079 {
1080     static const char * const selections[] = {
1081         "foo = res_cog of resname RA",
1082         "foo",
1083         "within 1 of foo",
1084         "bar = cog of resname RA",
1085         "bar",
1086         "within 1 of bar"
1087     };
1088     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1089     runTest("simple.gro", selections);
1090 }
1091
1092
1093 TEST_F(SelectionCollectionDataTest, HandlesConstantPositionInVariable)
1094 {
1095     static const char * const selections[] = {
1096         "constpos = [1.0, 2.5, 0.5]",
1097         "constpos",
1098         "within 2 of constpos"
1099     };
1100     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
1101              | efTestPositionAtoms);
1102     runTest("simple.gro", selections);
1103 }
1104
1105
1106 TEST_F(SelectionCollectionDataTest, HandlesNumericConstantsInVariables)
1107 {
1108     static const char * const selections[] = {
1109         "constint = 4",
1110         "constreal1 = 0.5",
1111         "constreal2 = 2.7",
1112         "resnr < constint",
1113         "x + constreal1 < constreal2"
1114     };
1115     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1116     runTest("simple.gro", selections);
1117 }
1118
1119
1120 /********************************************************************
1121  * Tests for complex boolean syntax
1122  */
1123
1124 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysis)
1125 {
1126     static const char * const selections[] = {
1127         "atomnr 1 to 5 and atomnr 2 to 7 and x < 2",
1128         "atomnr 1 to 5 and (atomnr 4 to 7 or x < 2)",
1129         "atomnr 1 to 5 and y < 3 and (atomnr 4 to 7 or x < 2)",
1130         "atomnr 1 to 5 and not (atomnr 4 to 7 or x < 2)",
1131         "atomnr 1 to 5 or (atomnr 4 to 6 and (atomnr 5 to 7 or x < 2))"
1132     };
1133     runTest(10, selections);
1134 }
1135
1136
1137 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysisWithVariables)
1138 {
1139     static const char * const selections[] = {
1140         "foo = atomnr 4 to 7 or x < 2",
1141         "atomnr 1 to 4 and foo",
1142         "atomnr 2 to 6 and y < 3 and foo",
1143         "atomnr 6 to 10 and not foo"
1144     };
1145     runTest(10, selections);
1146 }
1147
1148
1149 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysisWithMoreVariables)
1150 {
1151     static const char * const selections[] = {
1152         "foo = atomnr 4 to 7",
1153         "bar = foo and x < 2",
1154         "bar2 = foo and y < 2",
1155         "atomnr 1 to 4 and bar",
1156         "atomnr 2 to 6 and y < 3 and bar2",
1157         "atomnr 6 to 10 and not foo"
1158     };
1159     runTest(10, selections);
1160 }
1161
1162
1163 /********************************************************************
1164  * Tests for complex subexpression cases
1165  *
1166  * These tests use some knowledge of the implementation to trigger different
1167  * paths in the code.
1168  */
1169
1170 TEST_F(SelectionCollectionDataTest, HandlesUnusedVariables)
1171 {
1172     static const char * const selections[] = {
1173         "unused1 = atomnr 1 to 3",
1174         "foo = atomnr 4 to 7",
1175         "atomnr 1 to 6 and foo",
1176         "unused2 = atomnr 3 to 5"
1177     };
1178     runTest(10, selections);
1179 }
1180
1181
1182 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithStaticEvaluationGroups)
1183 {
1184     static const char * const selections[] = {
1185         "foo = atomnr 4 to 7 and x < 2",
1186         "atomnr 1 to 5 and foo",
1187         "atomnr 3 to 7 and foo"
1188     };
1189     runTest(10, selections);
1190 }
1191
1192
1193 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithMixedEvaluationGroups)
1194 {
1195     static const char * const selections[] = {
1196         "foo = atomnr 4 to 7 and x < 2",
1197         "atomnr 1 to 6 and foo",
1198         "within 1 of foo",
1199         "foo"
1200     };
1201     runTest(10, selections);
1202 }
1203
1204
1205 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithMixedEvaluationGroups2)
1206 {
1207     static const char * const selections[] = {
1208         "foo = atomnr 1 to 8 and x < 10",
1209         "atomnr 1 to 5 and y < 10 and foo",
1210         "foo"
1211     };
1212     setFlags(TestFlags() | efTestEvaluation);
1213     runTest("simple.gro", selections);
1214 }
1215
1216
1217 } // namespace