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