e4c10b132a402120f45c34dcf20e87cba5608007
[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,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  * \brief
37  * Tests selection parsing and compilation.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_selection
41  */
42 #include <gtest/gtest.h>
43
44 #include "gromacs/options/basicoptions.h"
45 #include "gromacs/options/options.h"
46 #include "gromacs/selection/indexutil.h"
47 #include "gromacs/selection/selectioncollection.h"
48 #include "gromacs/selection/selection.h"
49 #include "gromacs/utility/arrayref.h"
50 #include "gromacs/utility/exceptions.h"
51 #include "gromacs/utility/flags.h"
52 #include "gromacs/utility/gmxregex.h"
53 #include "gromacs/utility/stringutil.h"
54
55 #include "testutils/refdata.h"
56 #include "testutils/testasserts.h"
57 #include "testutils/testfilemanager.h"
58 #include "testutils/testoptions.h"
59
60 #include "toputils.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 int               s_debugLevel;
73
74         SelectionCollectionTest();
75         ~SelectionCollectionTest();
76
77         void setAtomCount(int natoms)
78         {
79             ASSERT_NO_THROW_GMX(sc_.setTopology(NULL, natoms));
80         }
81         void loadTopology(const char *filename);
82         void setTopology();
83         void loadIndexGroups(const char *filename);
84
85         gmx::test::TopologyManager  topManager_;
86         gmx::SelectionCollection    sc_;
87         gmx::SelectionList          sel_;
88         t_topology                 *top_;
89         t_trxframe                 *frame_;
90         gmx_ana_indexgrps_t        *grps_;
91 };
92
93 int SelectionCollectionTest::s_debugLevel = 0;
94
95 // \cond/\endcond do not seem to work here with Doxygen 1.8.5 parser.
96 #ifndef DOXYGEN
97 GMX_TEST_OPTIONS(SelectionCollectionTestOptions, options)
98 {
99     options->addOption(gmx::IntegerOption("seldebug")
100                            .store(&SelectionCollectionTest::s_debugLevel)
101                            .description("Set selection debug level"));
102 }
103 #endif
104
105 SelectionCollectionTest::SelectionCollectionTest()
106     : top_(NULL), frame_(NULL), grps_(NULL)
107 {
108     topManager_.requestFrame();
109     sc_.setDebugLevel(s_debugLevel);
110     sc_.setReferencePosType("atom");
111     sc_.setOutputPosType("atom");
112 }
113
114 SelectionCollectionTest::~SelectionCollectionTest()
115 {
116     if (grps_ != NULL)
117     {
118         gmx_ana_indexgrps_free(grps_);
119     }
120 }
121
122 void
123 SelectionCollectionTest::loadTopology(const char *filename)
124 {
125     topManager_.loadTopology(filename);
126     setTopology();
127 }
128
129 void
130 SelectionCollectionTest::setTopology()
131 {
132     top_   = topManager_.topology();
133     frame_ = topManager_.frame();
134
135     ASSERT_NO_THROW_GMX(sc_.setTopology(top_, -1));
136 }
137
138 void
139 SelectionCollectionTest::loadIndexGroups(const char *filename)
140 {
141     GMX_RELEASE_ASSERT(grps_ == NULL,
142                        "External groups can only be loaded once");
143     std::string fullpath =
144         gmx::test::TestFileManager::getInputFilePath(filename);
145     gmx_ana_indexgrps_init(&grps_, NULL, fullpath.c_str());
146     sc_.setIndexGroups(grps_);
147 }
148
149
150 /********************************************************************
151  * Test fixture for selection testing with reference data
152  */
153
154 class SelectionCollectionDataTest : public SelectionCollectionTest
155 {
156     public:
157         enum TestFlag
158         {
159             efTestEvaluation            = 1<<0,
160             efTestPositionAtoms         = 1<<1,
161             efTestPositionCoordinates   = 1<<2,
162             efTestPositionMapping       = 1<<3,
163             efTestPositionMasses        = 1<<4,
164             efTestPositionCharges       = 1<<5,
165             efTestSelectionNames        = 1<<6,
166             efDontTestCompiledAtoms     = 1<<8
167         };
168         typedef gmx::FlagsTemplate<TestFlag> TestFlags;
169
170         SelectionCollectionDataTest()
171             : checker_(data_.rootChecker()), count_(0), framenr_(0)
172         {
173         }
174
175         void setFlags(TestFlags flags) { flags_ = flags; }
176
177         void runParser(const gmx::ConstArrayRef<const char *> &selections);
178         void runCompiler();
179         void runEvaluate();
180         void runEvaluateFinal();
181
182         void runTest(int                                     natoms,
183                      const gmx::ConstArrayRef<const char *> &selections);
184         void runTest(const char                             *filename,
185                      const gmx::ConstArrayRef<const char *> &selections);
186
187     private:
188         static void checkSelection(gmx::test::TestReferenceChecker *checker,
189                                    const gmx::Selection &sel, TestFlags flags);
190
191         void checkCompiled();
192
193         gmx::test::TestReferenceData    data_;
194         gmx::test::TestReferenceChecker checker_;
195         size_t                          count_;
196         int                             framenr_;
197         TestFlags                       flags_;
198 };
199
200
201 void
202 SelectionCollectionDataTest::checkSelection(
203         gmx::test::TestReferenceChecker *checker,
204         const gmx::Selection &sel, TestFlags flags)
205 {
206     using gmx::test::TestReferenceChecker;
207
208     {
209         gmx::ConstArrayRef<int> atoms = sel.atomIndices();
210         checker->checkSequence(atoms.begin(), atoms.end(), "Atoms");
211     }
212     if (flags.test(efTestPositionAtoms)
213         || flags.test(efTestPositionCoordinates)
214         || flags.test(efTestPositionMapping)
215         || flags.test(efTestPositionMasses)
216         || flags.test(efTestPositionCharges))
217     {
218         TestReferenceChecker compound(
219                 checker->checkSequenceCompound("Positions", sel.posCount()));
220         for (int i = 0; i < sel.posCount(); ++i)
221         {
222             TestReferenceChecker          poscompound(compound.checkCompound("Position", NULL));
223             const gmx::SelectionPosition &p = sel.position(i);
224             if (flags.test(efTestPositionAtoms))
225             {
226                 gmx::ConstArrayRef<int> atoms = p.atomIndices();
227                 poscompound.checkSequence(atoms.begin(), atoms.end(), "Atoms");
228             }
229             if (flags.test(efTestPositionCoordinates))
230             {
231                 poscompound.checkVector(p.x(), "Coordinates");
232             }
233             if (flags.test(efTestPositionMapping))
234             {
235                 poscompound.checkInteger(p.refId(), "RefId");
236                 poscompound.checkInteger(p.mappedId(), "MappedId");
237             }
238             if (flags.test(efTestPositionMasses))
239             {
240                 poscompound.checkReal(p.mass(), "Mass");
241             }
242             if (flags.test(efTestPositionCharges))
243             {
244                 poscompound.checkReal(p.charge(), "Charge");
245             }
246         }
247     }
248 }
249
250
251 void
252 SelectionCollectionDataTest::runParser(
253         const gmx::ConstArrayRef<const char *> &selections)
254 {
255     using gmx::test::TestReferenceChecker;
256
257     TestReferenceChecker compound(checker_.checkCompound("ParsedSelections", "Parsed"));
258     size_t               varcount = 0;
259     count_ = 0;
260     for (size_t i = 0; i < selections.size(); ++i)
261     {
262         SCOPED_TRACE(std::string("Parsing selection \"")
263                      + selections[i] + "\"");
264         gmx::SelectionList result;
265         ASSERT_NO_THROW_GMX(result = sc_.parseFromString(selections[i]));
266         sel_.insert(sel_.end(), result.begin(), result.end());
267         if (sel_.size() == count_)
268         {
269             std::string          id = gmx::formatString("Variable%d", static_cast<int>(varcount + 1));
270             TestReferenceChecker varcompound(
271                     compound.checkCompound("ParsedVariable", id.c_str()));
272             varcompound.checkString(selections[i], "Input");
273             ++varcount;
274         }
275         else
276         {
277             std::string          id = gmx::formatString("Selection%d", static_cast<int>(count_ + 1));
278             TestReferenceChecker selcompound(
279                     compound.checkCompound("ParsedSelection", id.c_str()));
280             selcompound.checkString(selections[i], "Input");
281             if (flags_.test(efTestSelectionNames))
282             {
283                 selcompound.checkString(sel_[count_].name(), "Name");
284             }
285             selcompound.checkString(sel_[count_].selectionText(), "Text");
286             selcompound.checkBoolean(sel_[count_].isDynamic(), "Dynamic");
287             ++count_;
288         }
289     }
290 }
291
292
293 void
294 SelectionCollectionDataTest::runCompiler()
295 {
296     ASSERT_NO_THROW_GMX(sc_.compile());
297     ASSERT_EQ(count_, sel_.size());
298     checkCompiled();
299 }
300
301
302 void
303 SelectionCollectionDataTest::checkCompiled()
304 {
305     using gmx::test::TestReferenceChecker;
306     const TestFlags      mask = ~TestFlags(efTestPositionCoordinates);
307
308     TestReferenceChecker compound(checker_.checkCompound("CompiledSelections", "Compiled"));
309     for (size_t i = 0; i < count_; ++i)
310     {
311         SCOPED_TRACE(std::string("Checking selection \"") +
312                      sel_[i].selectionText() + "\"");
313         std::string          id = gmx::formatString("Selection%d", static_cast<int>(i + 1));
314         TestReferenceChecker selcompound(
315                 compound.checkCompound("Selection", id.c_str()));
316         if (flags_.test(efTestSelectionNames))
317         {
318             selcompound.checkString(sel_[i].name(), "Name");
319         }
320         if (!flags_.test(efDontTestCompiledAtoms))
321         {
322             checkSelection(&selcompound, sel_[i], flags_ & mask);
323         }
324     }
325 }
326
327
328 void
329 SelectionCollectionDataTest::runEvaluate()
330 {
331     using gmx::test::TestReferenceChecker;
332
333     ++framenr_;
334     ASSERT_NO_THROW_GMX(sc_.evaluate(frame_, NULL));
335     std::string          frame = gmx::formatString("Frame%d", framenr_);
336     TestReferenceChecker compound(
337             checker_.checkCompound("EvaluatedSelections", frame.c_str()));
338     for (size_t i = 0; i < count_; ++i)
339     {
340         SCOPED_TRACE(std::string("Checking selection \"") +
341                      sel_[i].selectionText() + "\"");
342         std::string          id = gmx::formatString("Selection%d", static_cast<int>(i + 1));
343         TestReferenceChecker selcompound(
344                 compound.checkCompound("Selection", id.c_str()));
345         checkSelection(&selcompound, sel_[i], flags_);
346     }
347 }
348
349
350 void
351 SelectionCollectionDataTest::runEvaluateFinal()
352 {
353     ASSERT_NO_THROW_GMX(sc_.evaluateFinal(framenr_));
354     checkCompiled();
355 }
356
357
358 void
359 SelectionCollectionDataTest::runTest(
360         int natoms, const gmx::ConstArrayRef<const char *> &selections)
361 {
362     ASSERT_NO_FATAL_FAILURE(runParser(selections));
363     ASSERT_NO_FATAL_FAILURE(setAtomCount(natoms));
364     ASSERT_NO_FATAL_FAILURE(runCompiler());
365 }
366
367
368 void
369 SelectionCollectionDataTest::runTest(
370         const char *filename, const gmx::ConstArrayRef<const char *> &selections)
371 {
372     ASSERT_NO_FATAL_FAILURE(runParser(selections));
373     ASSERT_NO_FATAL_FAILURE(loadTopology(filename));
374     ASSERT_NO_FATAL_FAILURE(runCompiler());
375     if (flags_.test(efTestEvaluation))
376     {
377         ASSERT_NO_FATAL_FAILURE(runEvaluate());
378         ASSERT_NO_FATAL_FAILURE(runEvaluateFinal());
379     }
380 }
381
382
383 /********************************************************************
384  * Tests for SelectionCollection functionality without reference data
385  */
386
387 TEST_F(SelectionCollectionTest, HandlesNoSelections)
388 {
389     EXPECT_FALSE(sc_.requiresTopology());
390     EXPECT_NO_THROW_GMX(sc_.compile());
391 }
392
393 TEST_F(SelectionCollectionTest, HandlesVelocityAndForceRequests)
394 {
395     ASSERT_NO_THROW_GMX(sel_ = sc_.parseFromString("atomnr 1 to 10; none"));
396     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
397     ASSERT_EQ(2U, sel_.size());
398     ASSERT_NO_THROW_GMX(sel_[0].setEvaluateVelocities(true));
399     ASSERT_NO_THROW_GMX(sel_[1].setEvaluateVelocities(true));
400     ASSERT_NO_THROW_GMX(sel_[0].setEvaluateForces(true));
401     ASSERT_NO_THROW_GMX(sel_[1].setEvaluateForces(true));
402     ASSERT_NO_THROW_GMX(sc_.compile());
403     EXPECT_TRUE(sel_[0].hasVelocities());
404     EXPECT_TRUE(sel_[1].hasVelocities());
405     EXPECT_TRUE(sel_[0].hasForces());
406     EXPECT_TRUE(sel_[1].hasForces());
407 }
408
409 TEST_F(SelectionCollectionTest, ParsesSelectionsFromFile)
410 {
411     ASSERT_NO_THROW_GMX(sel_ = sc_.parseFromFile(
412                                     gmx::test::TestFileManager::getInputFilePath("selfile.dat")));
413     // These should match the contents of selfile.dat
414     ASSERT_EQ(2U, sel_.size());
415     EXPECT_STREQ("resname RA RB", sel_[0].selectionText());
416     EXPECT_STREQ("resname RB RC", sel_[1].selectionText());
417 }
418
419 TEST_F(SelectionCollectionTest, HandlesInvalidRegularExpressions)
420 {
421     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
422     EXPECT_THROW_GMX({
423                          sc_.parseFromString("resname ~ \"R[A\"");
424                          sc_.compile();
425                      }, gmx::InvalidInputError);
426 }
427
428 TEST_F(SelectionCollectionTest, HandlesUnsupportedRegularExpressions)
429 {
430     if (!gmx::Regex::isSupported())
431     {
432         ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
433         EXPECT_THROW_GMX({
434                              sc_.parseFromString("resname \"R[AD]\"");
435                              sc_.compile();
436                          }, gmx::InvalidInputError);
437     }
438 }
439
440 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue)
441 {
442     EXPECT_THROW_GMX(sc_.parseFromString("mindist from atomnr 1 cutoff"),
443                      gmx::InvalidInputError);
444 }
445
446 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue2)
447 {
448     EXPECT_THROW_GMX(sc_.parseFromString("within 1 of"),
449                      gmx::InvalidInputError);
450 }
451
452 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue3)
453 {
454     EXPECT_THROW_GMX(sc_.parseFromString("within of atomnr 1"),
455                      gmx::InvalidInputError);
456 }
457
458 // TODO: Tests for more parser errors
459
460 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceParser1)
461 {
462     ASSERT_NO_THROW_GMX(sc_.setIndexGroups(NULL));
463     EXPECT_THROW_GMX(sc_.parseFromString("group \"foo\""), gmx::InconsistentInputError);
464     EXPECT_THROW_GMX(sc_.parseFromString("4"), gmx::InconsistentInputError);
465 }
466
467 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceParser2)
468 {
469     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
470     EXPECT_THROW_GMX(sc_.parseFromString("group \"foo\""), gmx::InconsistentInputError);
471     EXPECT_THROW_GMX(sc_.parseFromString("4"), gmx::InconsistentInputError);
472 }
473
474 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceDelayed1)
475 {
476     ASSERT_NO_THROW_GMX(sc_.parseFromString("group \"foo\""));
477     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
478     EXPECT_THROW_GMX(sc_.setIndexGroups(NULL), gmx::InconsistentInputError);
479     EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
480 }
481
482 TEST_F(SelectionCollectionTest, HandlesUnknownGroupReferenceDelayed2)
483 {
484     ASSERT_NO_THROW_GMX(sc_.parseFromString("group 4; group \"foo\""));
485     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
486     EXPECT_THROW_GMX(loadIndexGroups("simple.ndx"), gmx::InconsistentInputError);
487     EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
488 }
489
490 TEST_F(SelectionCollectionTest, HandlesUnsortedGroupReference)
491 {
492     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
493     EXPECT_THROW_GMX(sc_.parseFromString("atomnr 1 to 3 and group \"GrpUnsorted\""),
494                      gmx::InconsistentInputError);
495     EXPECT_THROW_GMX(sc_.parseFromString("group 2 or atomnr 2 to 5"),
496                      gmx::InconsistentInputError);
497     EXPECT_THROW_GMX(sc_.parseFromString("within 1 of group 2"),
498                      gmx::InconsistentInputError);
499 }
500
501 TEST_F(SelectionCollectionTest, HandlesUnsortedGroupReferenceDelayed)
502 {
503     ASSERT_NO_THROW_GMX(sc_.parseFromString("atomnr 1 to 3 and group \"GrpUnsorted\""));
504     ASSERT_NO_THROW_GMX(sc_.parseFromString("atomnr 1 to 3 and group 2"));
505     EXPECT_THROW_GMX(loadIndexGroups("simple.ndx"), gmx::InconsistentInputError);
506     // TODO: Add a separate check in the selection compiler for a safer API
507     // (makes sense in the future if the compiler needs the information for
508     // other purposes as well).
509     // EXPECT_THROW_GMX(sc_.compile(), gmx::APIError);
510 }
511
512 TEST_F(SelectionCollectionTest, RecoversFromMissingMoleculeInfo)
513 {
514     ASSERT_NO_THROW_GMX(sc_.parseFromString("molindex 1 to 5"));
515     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
516     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
517 }
518
519 TEST_F(SelectionCollectionTest, RecoversFromMissingAtomTypes)
520 {
521     ASSERT_NO_THROW_GMX(sc_.parseFromString("type CA"));
522     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
523     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
524 }
525
526 TEST_F(SelectionCollectionTest, RecoversFromMissingPDBInfo)
527 {
528     ASSERT_NO_THROW_GMX(sc_.parseFromString("altloc A"));
529     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
530     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
531 }
532
533 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation)
534 {
535     ASSERT_NO_THROW_GMX(sc_.parseFromString("all permute 1 1"));
536     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
537     EXPECT_THROW_GMX(sc_.compile(), gmx::InvalidInputError);
538 }
539
540 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation2)
541 {
542     ASSERT_NO_THROW_GMX(sc_.parseFromString("all permute 3 2 1"));
543     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
544     EXPECT_THROW_GMX(sc_.compile(), gmx::InconsistentInputError);
545 }
546
547 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation3)
548 {
549     ASSERT_NO_THROW_GMX(sc_.parseFromString("x < 1.5 permute 3 2 1"));
550     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
551     ASSERT_NO_THROW_GMX(sc_.compile());
552     EXPECT_THROW_GMX(sc_.evaluate(frame_, NULL), gmx::InconsistentInputError);
553 }
554
555 // TODO: Tests for evaluation errors
556
557
558 /********************************************************************
559  * Tests for selection keywords
560  */
561
562 TEST_F(SelectionCollectionDataTest, HandlesAllNone)
563 {
564     static const char * const selections[] = {
565         "all",
566         "none"
567     };
568     runTest(10, selections);
569 }
570
571 TEST_F(SelectionCollectionDataTest, HandlesAtomnr)
572 {
573     static const char * const selections[] = {
574         "atomnr 1 to 3 6 to 8",
575         "atomnr 4 2 5 to 7",
576         "atomnr <= 5"
577     };
578     runTest(10, selections);
579 }
580
581 TEST_F(SelectionCollectionDataTest, HandlesResnr)
582 {
583     static const char * const selections[] = {
584         "resnr 1 2 5",
585         "resid 4 to 3"
586     };
587     runTest("simple.gro", selections);
588 }
589
590 TEST_F(SelectionCollectionDataTest, HandlesResIndex)
591 {
592     static const char * const selections[] = {
593         "resindex 1 4",
594         "residue 1 3"
595     };
596     runTest("simple.pdb", selections);
597 }
598
599 TEST_F(SelectionCollectionDataTest, HandlesMolIndex)
600 {
601     static const char * const selections[] = {
602         "molindex 1 4",
603         "molecule 2 3 5"
604     };
605     ASSERT_NO_FATAL_FAILURE(runParser(selections));
606     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
607     topManager_.initUniformMolecules(3);
608     ASSERT_NO_FATAL_FAILURE(runCompiler());
609 }
610
611 TEST_F(SelectionCollectionDataTest, HandlesAtomname)
612 {
613     static const char * const selections[] = {
614         "name CB",
615         "atomname S1 S2"
616     };
617     runTest("simple.gro", selections);
618 }
619
620 TEST_F(SelectionCollectionDataTest, HandlesPdbAtomname)
621 {
622     static const char * const selections[] = {
623         "name HG21",
624         "name 1HG2",
625         "pdbname HG21 CB",
626         "pdbatomname 1HG2"
627     };
628     runTest("simple.pdb", selections);
629 }
630
631
632 TEST_F(SelectionCollectionDataTest, HandlesAtomtype)
633 {
634     static const char * const selections[] = {
635         "atomtype CA"
636     };
637     ASSERT_NO_FATAL_FAILURE(runParser(selections));
638     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
639     const char *const types[] = { "CA", "SA", "SB" };
640     topManager_.initAtomTypes(types);
641     ASSERT_NO_FATAL_FAILURE(runCompiler());
642 }
643
644 TEST_F(SelectionCollectionDataTest, HandlesChain)
645 {
646     static const char * const selections[] = {
647         "chain A",
648         "chain B"
649     };
650     runTest("simple.pdb", selections);
651 }
652
653 TEST_F(SelectionCollectionDataTest, HandlesMass)
654 {
655     static const char * const selections[] = {
656         "mass > 5"
657     };
658     ASSERT_NO_FATAL_FAILURE(runParser(selections));
659     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
660     for (int i = 0; i < top_->atoms.nr; ++i)
661     {
662         top_->atoms.atom[i].m = 1.0 + i;
663     }
664     ASSERT_NO_FATAL_FAILURE(runCompiler());
665 }
666
667 TEST_F(SelectionCollectionDataTest, HandlesCharge)
668 {
669     static const char * const selections[] = {
670         "charge < 0.5"
671     };
672     ASSERT_NO_FATAL_FAILURE(runParser(selections));
673     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
674     for (int i = 0; i < top_->atoms.nr; ++i)
675     {
676         top_->atoms.atom[i].q = i / 10.0;
677     }
678     ASSERT_NO_FATAL_FAILURE(runCompiler());
679 }
680
681 TEST_F(SelectionCollectionDataTest, HandlesAltLoc)
682 {
683     static const char * const selections[] = {
684         "altloc \" \"",
685         "altloc A"
686     };
687     runTest("simple.pdb", selections);
688 }
689
690 TEST_F(SelectionCollectionDataTest, HandlesInsertCode)
691 {
692     static const char * const selections[] = {
693         "insertcode \" \"",
694         "insertcode A"
695     };
696     runTest("simple.pdb", selections);
697 }
698
699 TEST_F(SelectionCollectionDataTest, HandlesOccupancy)
700 {
701     static const char * const selections[] = {
702         "occupancy 1",
703         "occupancy < .5"
704     };
705     runTest("simple.pdb", selections);
706 }
707
708 TEST_F(SelectionCollectionDataTest, HandlesBeta)
709 {
710     static const char * const selections[] = {
711         "beta 0",
712         "beta >= 0.3"
713     };
714     runTest("simple.pdb", selections);
715 }
716
717 TEST_F(SelectionCollectionDataTest, HandlesResname)
718 {
719     static const char * const selections[] = {
720         "resname RA",
721         "resname RB RC"
722     };
723     runTest("simple.gro", selections);
724 }
725
726 TEST_F(SelectionCollectionDataTest, HandlesCoordinateKeywords)
727 {
728     static const char * const selections[] = {
729         "x < 3",
730         "y >= 3",
731         "x {-1 to 2}"
732     };
733     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
734     runTest("simple.gro", selections);
735 }
736
737
738 TEST_F(SelectionCollectionDataTest, HandlesSameResidue)
739 {
740     static const char * const selections[] = {
741         "same residue as atomnr 1 4 12"
742     };
743     runTest("simple.gro", selections);
744 }
745
746
747 TEST_F(SelectionCollectionDataTest, HandlesSameResidueName)
748 {
749     static const char * const selections[] = {
750         "same resname as atomnr 1 14"
751     };
752     runTest("simple.gro", selections);
753 }
754
755
756 TEST_F(SelectionCollectionDataTest, HandlesPositionKeywords)
757 {
758     static const char * const selections[] = {
759         "cog of resnr 1 3",
760         "res_cog of name CB and resnr 1 3",
761         "whole_res_cog of name CB and resnr 1 3",
762         "part_res_cog of x < 3",
763         "dyn_res_cog of x < 3"
764     };
765     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
766              | efTestPositionAtoms);
767     runTest("simple.gro", selections);
768 }
769
770
771 TEST_F(SelectionCollectionDataTest, HandlesDistanceKeyword)
772 {
773     static const char * const selections[] = {
774         "distance from cog of resnr 1 < 2"
775     };
776     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
777     runTest("simple.gro", selections);
778 }
779
780
781 TEST_F(SelectionCollectionDataTest, HandlesMinDistanceKeyword)
782 {
783     static const char * const selections[] = {
784         "mindistance from resnr 1 < 2"
785     };
786     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
787     runTest("simple.gro", selections);
788 }
789
790
791 TEST_F(SelectionCollectionDataTest, HandlesWithinKeyword)
792 {
793     static const char * const selections[] = {
794         "within 1 of resnr 2"
795     };
796     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
797     runTest("simple.gro", selections);
798 }
799
800
801 TEST_F(SelectionCollectionDataTest, HandlesInSolidAngleKeyword)
802 {
803     // Both of these should evaluate to empty on a correct implementation.
804     static const char * const selections[] = {
805         "resname TP and not insolidangle center cog of resname C span resname R cutoff 20",
806         "resname TN and insolidangle center cog of resname C span resname R cutoff 20"
807     };
808     setFlags(TestFlags() | efDontTestCompiledAtoms | efTestEvaluation);
809     runTest("sphere.gro", selections);
810 }
811
812
813 TEST_F(SelectionCollectionDataTest, HandlesPermuteModifier)
814 {
815     static const char * const selections[] = {
816         "all permute 3 1 2",
817         "res_cog of resnr 1 to 4 permute 2 1",
818         "name CB S1 and res_cog x < 3 permute 2 1"
819     };
820     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
821              | efTestPositionAtoms | efTestPositionMapping);
822     runTest("simple.gro", selections);
823 }
824
825
826 TEST_F(SelectionCollectionDataTest, HandlesPlusModifier)
827 {
828     static const char * const selections[] = {
829         "name S2 plus name S1",
830         "res_cog of resnr 2 plus res_cog of resnr 1 plus res_cog of resnr 3",
831         "name S1 and y < 3 plus res_cog of x < 2.5"
832     };
833     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
834              | efTestPositionAtoms | efTestPositionMapping);
835     runTest("simple.gro", selections);
836 }
837
838
839 TEST_F(SelectionCollectionDataTest, HandlesMergeModifier)
840 {
841     static const char * const selections[] = {
842         "name S2 merge name S1",
843         "resnr 1 2 and name S2 merge resnr 1 2 and name S1 merge res_cog of resnr 1 2",
844         "name S1 and x < 2.5 merge res_cog of x < 2.5"
845     };
846     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
847              | efTestPositionAtoms | efTestPositionMapping);
848     runTest("simple.gro", selections);
849 }
850
851
852 /********************************************************************
853  * Tests for generic selection evaluation
854  */
855
856 TEST_F(SelectionCollectionDataTest, ComputesMassesAndCharges)
857 {
858     static const char * const selections[] = {
859         "name CB",
860         "y > 2",
861         "res_cog of y > 2"
862     };
863     setFlags(TestFlags() | efTestEvaluation | efTestPositionAtoms
864              | efTestPositionMasses | efTestPositionCharges);
865     ASSERT_NO_FATAL_FAILURE(runParser(selections));
866     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
867     for (int i = 0; i < top_->atoms.nr; ++i)
868     {
869         top_->atoms.atom[i].m =   1.0 + i / 100.0;
870         top_->atoms.atom[i].q = -(1.0 + i / 100.0);
871     }
872     ASSERT_NO_FATAL_FAILURE(runCompiler());
873     ASSERT_NO_FATAL_FAILURE(runEvaluate());
874     ASSERT_NO_FATAL_FAILURE(runEvaluateFinal());
875 }
876
877 TEST_F(SelectionCollectionDataTest, ComputesMassesAndChargesWithoutTopology)
878 {
879     static const char * const selections[] = {
880         "atomnr 1 to 3 8 to 9",
881         "y > 2",
882         "cog of (y > 2)"
883     };
884     setFlags(TestFlags() | efTestPositionAtoms
885              | efTestPositionMasses | efTestPositionCharges);
886     runTest(10, selections);
887 }
888
889
890 /********************************************************************
891  * Tests for selection syntactic constructs
892  */
893
894 TEST_F(SelectionCollectionDataTest, HandlesSelectionNames)
895 {
896     static const char * const selections[] = {
897         "\"GroupSelection\" group \"GrpA\"",
898         "\"DynamicSelection\" x < 5",
899         "y < 3"
900     };
901     setFlags(TestFlags() | efTestSelectionNames);
902     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
903     runTest(10, selections);
904 }
905
906 TEST_F(SelectionCollectionDataTest, HandlesIndexGroupsInSelections)
907 {
908     static const char * const selections[] = {
909         "group \"GrpA\"",
910         "GrpB",
911         "1",
912         "group \"GrpB\" and resname RB"
913     };
914     setFlags(TestFlags() | efTestSelectionNames);
915     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
916     runTest("simple.gro", selections);
917 }
918
919 TEST_F(SelectionCollectionDataTest, HandlesIndexGroupsInSelectionsDelayed)
920 {
921     static const char * const selections[] = {
922         "group \"GrpA\"",
923         "GrpB",
924         "1",
925         "group \"GrpB\" and resname RB"
926     };
927     setFlags(TestFlags() | efTestSelectionNames);
928     ASSERT_NO_FATAL_FAILURE(runParser(selections));
929     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
930     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
931     ASSERT_NO_FATAL_FAILURE(runCompiler());
932 }
933
934 TEST_F(SelectionCollectionDataTest, HandlesUnsortedIndexGroupsInSelections)
935 {
936     static const char * const selections[] = {
937         "foo = group \"GrpUnsorted\"",
938         "group \"GrpUnsorted\"",
939         "GrpUnsorted",
940         "2",
941         "res_cog of group \"GrpUnsorted\"",
942         "group \"GrpUnsorted\" permute 2 1",
943         "foo"
944     };
945     setFlags(TestFlags() | efTestPositionAtoms | efTestPositionMapping
946              | efTestSelectionNames);
947     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
948     runTest("simple.gro", selections);
949 }
950
951 TEST_F(SelectionCollectionDataTest, HandlesUnsortedIndexGroupsInSelectionsDelayed)
952 {
953     static const char * const selections[] = {
954         "foo = group \"GrpUnsorted\"",
955         "group \"GrpUnsorted\"",
956         "GrpUnsorted",
957         "2",
958         "res_cog of group \"GrpUnsorted\"",
959         "group \"GrpUnsorted\" permute 2 1",
960         "foo"
961     };
962     ASSERT_NO_FATAL_FAILURE(runParser(selections));
963     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
964     ASSERT_NO_THROW_GMX(loadIndexGroups("simple.ndx"));
965     ASSERT_NO_FATAL_FAILURE(runCompiler());
966 }
967
968 TEST_F(SelectionCollectionDataTest, HandlesConstantPositions)
969 {
970     static const char * const selections[] = {
971         "[1, -2, 3.5]"
972     };
973     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
974     runTest("simple.gro", selections);
975 }
976
977
978 TEST_F(SelectionCollectionDataTest, HandlesWithinConstantPositions)
979 {
980     static const char * const selections[] = {
981         "within 1 of [2, 1, 0]"
982     };
983     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
984     runTest("simple.gro", selections);
985 }
986
987
988 TEST_F(SelectionCollectionDataTest, HandlesForcedStringMatchingMode)
989 {
990     static const char * const selections[] = {
991         "name = S1 \"C?\"",
992         "name ? S1 \"C?\""
993     };
994     runTest("simple.gro", selections);
995 }
996
997
998 TEST_F(SelectionCollectionDataTest, HandlesWildcardMatching)
999 {
1000     static const char * const selections[] = {
1001         "name \"S?\"",
1002         "name ? \"S?\""
1003     };
1004     runTest("simple.gro", selections);
1005 }
1006
1007
1008 TEST_F(SelectionCollectionDataTest, HandlesRegexMatching)
1009 {
1010     static const char * const selections[] = {
1011         "resname \"R[BD]\"",
1012         "resname ~ \"R[BD]\""
1013     };
1014     if (gmx::Regex::isSupported())
1015     {
1016         runTest("simple.gro", selections);
1017     }
1018 }
1019
1020
1021 TEST_F(SelectionCollectionDataTest, HandlesBasicBoolean)
1022 {
1023     static const char * const selections[] = {
1024         "atomnr 1 to 5 and atomnr 2 to 7",
1025         "atomnr 1 to 5 or not atomnr 3 to 8",
1026         "not not atomnr 1 to 5 and atomnr 2 to 6 and not not atomnr 3 to 7",
1027         "atomnr 1 to 5 and (atomnr 2 to 7 and atomnr 3 to 6)",
1028         "x < 5 and atomnr 1 to 5 and y < 3 and atomnr 2 to 4"
1029     };
1030     runTest(10, selections);
1031 }
1032
1033
1034 TEST_F(SelectionCollectionDataTest, HandlesDynamicAtomValuedParameters)
1035 {
1036     static const char * const selections[] = {
1037         "same residue as (atomnr 3 5 13 or y > 5)",
1038         "(resnr 1 3 5 or x > 10) and same residue as (atomnr 3 5 13 or z > 5)"
1039     };
1040     setFlags(TestFlags() | efTestEvaluation);
1041     runTest("simple.gro", selections);
1042 }
1043
1044
1045 TEST_F(SelectionCollectionDataTest, HandlesEmptySelectionWithUnevaluatedExpressions)
1046 {
1047     static const char * const selections[] = {
1048         "none and x > 2",
1049         "none and same resname as resnr 2"
1050     };
1051     runTest("simple.gro", selections);
1052 }
1053
1054
1055 TEST_F(SelectionCollectionDataTest, HandlesNumericComparisons)
1056 {
1057     static const char * const selections[] = {
1058         "x > 2",
1059         "2 < x",
1060         "y > resnr",
1061         "resnr < 2.5",
1062         "2.5 > resnr"
1063     };
1064     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1065     runTest("simple.gro", selections);
1066 }
1067
1068
1069 TEST_F(SelectionCollectionDataTest, HandlesArithmeticExpressions)
1070 {
1071     static const char * const selections[] = {
1072         "x+1 > 3",
1073         "(y-1)^2 <= 1",
1074         "x+--1 > 3",
1075         "-x+-1 < -3"
1076     };
1077     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1078     runTest("simple.gro", selections);
1079 }
1080
1081
1082 TEST_F(SelectionCollectionDataTest, HandlesNumericVariables)
1083 {
1084     static const char * const selections[] = {
1085         "value = x + y",
1086         "value <= 4",
1087         "index = resnr",
1088         "index < 3"
1089     };
1090     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1091     runTest("simple.gro", selections);
1092 }
1093
1094
1095 TEST_F(SelectionCollectionDataTest, HandlesComplexNumericVariables)
1096 {
1097     static const char * const selections[] = {
1098         "value = x + y",
1099         "resname RA and value <= 4",
1100         "resname RA RB and x < 3 and value <= 4",
1101         "index = atomnr",
1102         "resname RA and index < 3",
1103         "resname RB and y < 3 and index < 6"
1104     };
1105     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1106     runTest("simple.gro", selections);
1107 }
1108
1109
1110 TEST_F(SelectionCollectionDataTest, HandlesPositionVariables)
1111 {
1112     static const char * const selections[] = {
1113         "foo = res_cog of resname RA",
1114         "foo",
1115         "within 1 of foo",
1116         "bar = cog of resname RA",
1117         "bar",
1118         "within 1 of bar"
1119     };
1120     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1121     runTest("simple.gro", selections);
1122 }
1123
1124
1125 TEST_F(SelectionCollectionDataTest, HandlesConstantPositionInVariable)
1126 {
1127     static const char * const selections[] = {
1128         "constpos = [1.0, 2.5, 0.5]",
1129         "constpos",
1130         "within 2 of constpos"
1131     };
1132     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates
1133              | efTestPositionAtoms);
1134     runTest("simple.gro", selections);
1135 }
1136
1137
1138 TEST_F(SelectionCollectionDataTest, HandlesNumericConstantsInVariables)
1139 {
1140     static const char * const selections[] = {
1141         "constint = 4",
1142         "constreal1 = 0.5",
1143         "constreal2 = 2.7",
1144         "resnr < constint",
1145         "x + constreal1 < constreal2"
1146     };
1147     setFlags(TestFlags() | efTestEvaluation | efTestPositionCoordinates);
1148     runTest("simple.gro", selections);
1149 }
1150
1151
1152 /********************************************************************
1153  * Tests for complex boolean syntax
1154  */
1155
1156 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysis)
1157 {
1158     static const char * const selections[] = {
1159         "atomnr 1 to 5 and atomnr 2 to 7 and x < 2",
1160         "atomnr 1 to 5 and (atomnr 4 to 7 or x < 2)",
1161         "atomnr 1 to 5 and y < 3 and (atomnr 4 to 7 or x < 2)",
1162         "atomnr 1 to 5 and not (atomnr 4 to 7 or x < 2)",
1163         "atomnr 1 to 5 or (atomnr 4 to 6 and (atomnr 5 to 7 or x < 2))"
1164     };
1165     runTest(10, selections);
1166 }
1167
1168
1169 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysisWithVariables)
1170 {
1171     static const char * const selections[] = {
1172         "foo = atomnr 4 to 7 or x < 2",
1173         "atomnr 1 to 4 and foo",
1174         "atomnr 2 to 6 and y < 3 and foo",
1175         "atomnr 6 to 10 and not foo"
1176     };
1177     runTest(10, selections);
1178 }
1179
1180
1181 TEST_F(SelectionCollectionDataTest, HandlesBooleanStaticAnalysisWithMoreVariables)
1182 {
1183     static const char * const selections[] = {
1184         "foo = atomnr 4 to 7",
1185         "bar = foo and x < 2",
1186         "bar2 = foo and y < 2",
1187         "atomnr 1 to 4 and bar",
1188         "atomnr 2 to 6 and y < 3 and bar2",
1189         "atomnr 6 to 10 and not foo"
1190     };
1191     runTest(10, selections);
1192 }
1193
1194
1195 /********************************************************************
1196  * Tests for complex subexpression cases
1197  *
1198  * These tests use some knowledge of the implementation to trigger different
1199  * paths in the code.
1200  */
1201
1202 TEST_F(SelectionCollectionDataTest, HandlesUnusedVariables)
1203 {
1204     static const char * const selections[] = {
1205         "unused1 = atomnr 1 to 3",
1206         "foo = atomnr 4 to 7",
1207         "atomnr 1 to 6 and foo",
1208         "unused2 = atomnr 3 to 5"
1209     };
1210     runTest(10, selections);
1211 }
1212
1213
1214 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithStaticEvaluationGroups)
1215 {
1216     static const char * const selections[] = {
1217         "foo = atomnr 4 to 7 and x < 2",
1218         "atomnr 1 to 5 and foo",
1219         "atomnr 3 to 7 and foo"
1220     };
1221     runTest(10, selections);
1222 }
1223
1224
1225 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithMixedEvaluationGroups)
1226 {
1227     static const char * const selections[] = {
1228         "foo = atomnr 4 to 7 and x < 2",
1229         "atomnr 1 to 6 and foo",
1230         "within 1 of foo",
1231         "foo"
1232     };
1233     runTest(10, selections);
1234 }
1235
1236
1237 TEST_F(SelectionCollectionDataTest, HandlesVariablesWithMixedEvaluationGroups2)
1238 {
1239     static const char * const selections[] = {
1240         "foo = atomnr 1 to 8 and x < 10",
1241         "atomnr 1 to 5 and y < 10 and foo",
1242         "foo"
1243     };
1244     setFlags(TestFlags() | efTestEvaluation);
1245     runTest("simple.gro", selections);
1246 }
1247
1248
1249 } // namespace