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