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