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