2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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.
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.
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.
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.
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.
37 * Tests the position mapping engine.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 #include <gtest/gtest.h>
46 #include "gromacs/legacyheaders/smalloc.h"
47 #include "gromacs/legacyheaders/typedefs.h"
48 #include "gromacs/legacyheaders/vec.h"
50 #include "gromacs/selection/indexutil.h"
51 #include "gromacs/selection/poscalc.h"
52 #include "gromacs/selection/position.h"
54 #include "testutils/refdata.h"
61 /********************************************************************
62 * PositionCalculationTest
65 class PositionCalculationTest : public ::testing::Test
68 PositionCalculationTest();
69 ~PositionCalculationTest();
71 void generateCoordinates();
73 gmx_ana_poscalc_t *createCalculation(e_poscalc_t type, int flags);
74 void setMaximumGroup(gmx_ana_poscalc_t *pc,
75 int count, const int atoms[]);
76 gmx_ana_pos_t *initPositions(gmx_ana_poscalc_t *pc, const char *name);
78 void checkInitialized();
79 void updateAndCheck(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p,
80 int count, const int atoms[],
81 gmx::test::TestReferenceChecker *checker,
84 void testSingleStatic(e_poscalc_t type, int flags, bool bExpectTop,
85 int atomCount, const int atoms[]);
86 void testSingleDynamic(e_poscalc_t type, int flags, bool bExpectTop,
87 int initCount, const int initAtoms[],
88 int evalCount, const int evalAtoms[]);
91 void setMaximumGroup(gmx_ana_poscalc_t *pc, const int (&atoms)[count])
93 setMaximumGroup(pc, count, atoms);
96 void updateAndCheck(gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p,
97 const int (&atoms)[count],
98 gmx::test::TestReferenceChecker *checker,
101 updateAndCheck(pc, p, count, atoms, checker, name);
103 template <int atomCount>
104 void testSingleStatic(e_poscalc_t type, int flags, bool bExpectTop,
105 const int (&atoms)[atomCount])
107 testSingleStatic(type, flags, bExpectTop, atomCount, atoms);
109 template <int initCount, int evalCount>
110 void testSingleDynamic(e_poscalc_t type, int flags, bool bExpectTop,
111 const int (&initAtoms)[initCount],
112 const int (&evalAtoms)[evalCount])
114 testSingleDynamic(type, flags, bExpectTop,
115 initCount, initAtoms, evalCount, evalAtoms);
118 gmx::test::TestReferenceData data_;
119 gmx::test::TestReferenceChecker checker_;
120 gmx::test::TopologyManager topManager_;
121 gmx::PositionCalculationCollection pcc_;
126 PositionTest() : pos(NULL), pc(NULL), name(NULL) {}
127 PositionTest(gmx_ana_pos_t *pos, gmx_ana_poscalc_t *pc,
129 : pos(pos), pc(pc), name(name)
134 gmx_ana_poscalc_t *pc;
138 typedef std::vector<PositionTest> PositionTestList;
140 void setTopologyIfRequired();
141 void checkPositions(gmx::test::TestReferenceChecker *checker,
142 const char *name, gmx_ana_pos_t *p,
145 std::vector<gmx_ana_poscalc_t *> pcList_;
146 PositionTestList posList_;
150 PositionCalculationTest::PositionCalculationTest()
151 : checker_(data_.rootChecker()), bTopSet_(false)
153 topManager_.requestFrame();
156 PositionCalculationTest::~PositionCalculationTest()
158 std::vector<gmx_ana_poscalc_t *>::reverse_iterator pci;
159 for (pci = pcList_.rbegin(); pci != pcList_.rend(); ++pci)
161 gmx_ana_poscalc_free(*pci);
164 PositionTestList::iterator pi;
165 for (pi = posList_.begin(); pi != posList_.end(); ++pi)
167 gmx_ana_pos_free(pi->pos);
171 void PositionCalculationTest::generateCoordinates()
173 t_topology *top = topManager_.topology();
174 t_trxframe *frame = topManager_.frame();
175 for (int i = 0; i < top->atoms.nr; ++i)
178 frame->x[i][YY] = top->atoms.atom[i].resind;
179 frame->x[i][ZZ] = 0.0;
182 copy_rvec(frame->x[i], frame->v[i]);
183 frame->v[i][ZZ] = 1.0;
187 copy_rvec(frame->x[i], frame->f[i]);
188 frame->f[i][ZZ] = -1.0;
194 PositionCalculationTest::createCalculation(e_poscalc_t type, int flags)
196 pcList_.reserve(pcList_.size() + 1);
197 pcList_.push_back(pcc_.createCalculation(type, flags));
198 return pcList_.back();
201 void PositionCalculationTest::setMaximumGroup(gmx_ana_poscalc_t *pc,
202 int count, const int atoms[])
204 setTopologyIfRequired();
207 g.index = const_cast<int *>(atoms);
208 gmx_ana_poscalc_set_maxindex(pc, &g);
212 PositionCalculationTest::initPositions(gmx_ana_poscalc_t *pc, const char *name)
214 posList_.reserve(posList_.size() + 1);
217 posList_.push_back(PositionTest(p, pc, name));
218 gmx_ana_poscalc_init_pos(pc, p);
222 void PositionCalculationTest::checkInitialized()
224 gmx::test::TestReferenceChecker compound(
225 checker_.checkCompound("InitializedPositions", NULL));
226 PositionTestList::const_iterator pi;
227 for (pi = posList_.begin(); pi != posList_.end(); ++pi)
229 checkPositions(&compound, pi->name, pi->pos, false);
233 void PositionCalculationTest::updateAndCheck(
234 gmx_ana_poscalc_t *pc, gmx_ana_pos_t *p, int count, const int atoms[],
235 gmx::test::TestReferenceChecker *checker, const char *name)
239 g.index = const_cast<int *>(atoms);
240 gmx_ana_poscalc_update(pc, p, &g, topManager_.frame(), NULL);
241 checkPositions(checker, name, p, true);
244 void PositionCalculationTest::testSingleStatic(
245 e_poscalc_t type, int flags, bool bExpectTop,
246 int atomCount, const int atoms[])
248 t_trxframe *frame = topManager_.frame();
251 flags |= POS_VELOCITIES;
257 gmx_ana_poscalc_t *pc = createCalculation(type, flags);
258 EXPECT_EQ(bExpectTop, gmx_ana_poscalc_requires_top(pc));
259 setMaximumGroup(pc, atomCount, atoms);
260 gmx_ana_pos_t *p = initPositions(pc, NULL);
263 pcc_.initEvaluation();
265 generateCoordinates();
266 gmx::test::TestReferenceChecker frameCompound(
267 checker_.checkCompound("EvaluatedPositions", "Frame0"));
268 updateAndCheck(pc, p, atomCount, atoms, &frameCompound, NULL);
272 void PositionCalculationTest::testSingleDynamic(
273 e_poscalc_t type, int flags, bool bExpectTop,
274 int initCount, const int initAtoms[],
275 int evalCount, const int evalAtoms[])
277 gmx_ana_poscalc_t *pc = createCalculation(type, flags | POS_DYNAMIC);
278 EXPECT_EQ(bExpectTop, gmx_ana_poscalc_requires_top(pc));
279 setMaximumGroup(pc, initCount, initAtoms);
280 gmx_ana_pos_t *p = initPositions(pc, NULL);
283 pcc_.initEvaluation();
285 generateCoordinates();
286 gmx::test::TestReferenceChecker frameCompound(
287 checker_.checkCompound("EvaluatedPositions", "Frame0"));
288 updateAndCheck(pc, p, evalCount, evalAtoms, &frameCompound, NULL);
292 void PositionCalculationTest::setTopologyIfRequired()
298 std::vector<gmx_ana_poscalc_t *>::const_iterator pci;
299 for (pci = pcList_.begin(); pci != pcList_.end(); ++pci)
301 if (gmx_ana_poscalc_requires_top(*pci))
304 pcc_.setTopology(topManager_.topology());
310 void PositionCalculationTest::checkPositions(
311 gmx::test::TestReferenceChecker *checker,
312 const char *name, gmx_ana_pos_t *p, bool bCoordinates)
314 EXPECT_EQ(p->nr, p->m.nr);
315 EXPECT_EQ(p->nr, p->m.mapb.nr);
316 gmx::test::TestReferenceChecker compound(
317 checker->checkCompound("Positions", name));
318 compound.checkInteger(p->nr, "Count");
319 const char *type = "???";
322 case INDEX_UNKNOWN: type = "unknown"; break;
323 case INDEX_ATOM: type = "atoms"; break;
324 case INDEX_RES: type = "residues"; break;
325 case INDEX_MOL: type = "molecules"; break;
326 case INDEX_ALL: type = "single"; break;
328 compound.checkString(type, "Type");
329 compound.checkSequenceArray(p->nr + 1, p->m.mapb.index, "Block");
330 for (int i = 0; i < p->nr; ++i)
332 gmx::test::TestReferenceChecker posCompound(
333 compound.checkCompound("Position", NULL));
334 posCompound.checkSequence(&p->m.mapb.a[p->m.mapb.index[i]],
335 &p->m.mapb.a[p->m.mapb.index[i+1]],
337 posCompound.checkInteger(p->m.refid[i], "RefId");
340 posCompound.checkVector(p->x[i], "Coordinates");
342 if (bCoordinates && p->v != NULL)
344 posCompound.checkVector(p->v[i], "Velocity");
346 if (bCoordinates && p->f != NULL)
348 posCompound.checkVector(p->f[i], "Force");
350 int originalIdIndex = (p->m.refid[i] != -1 ? p->m.refid[i] : i);
351 EXPECT_EQ(p->m.orgid[originalIdIndex], p->m.mapid[i]);
355 /********************************************************************
359 TEST_F(PositionCalculationTest, ComputesAtomPositions)
361 const int group[] = { 0, 1, 2, 3 };
362 topManager_.requestVelocities();
363 topManager_.requestForces();
364 topManager_.initAtoms(4);
365 testSingleStatic(POS_ATOM, 0, false, group);
368 TEST_F(PositionCalculationTest, ComputesResidueCOGPositions)
370 const int group[] = { 0, 1, 2, 3, 4, 8 };
371 topManager_.requestVelocities();
372 topManager_.requestForces();
373 topManager_.initAtoms(9);
374 topManager_.initUniformResidues(3);
375 testSingleStatic(POS_RES, 0, true, group);
378 TEST_F(PositionCalculationTest, ComputesResidueCOMPositions)
380 const int group[] = { 0, 1, 2, 3, 4, 8 };
381 topManager_.requestVelocities();
382 topManager_.requestForces();
383 topManager_.initAtoms(9);
384 topManager_.initUniformResidues(3);
385 testSingleStatic(POS_RES, POS_MASS, true, group);
388 TEST_F(PositionCalculationTest, ComputesGroupCOGPositions)
390 const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8 };
391 topManager_.requestVelocities();
392 topManager_.requestForces();
393 topManager_.initAtoms(9);
394 // Topology (masses) is requires for computing the force
395 testSingleStatic(POS_ALL, 0, true, group);
398 TEST_F(PositionCalculationTest, ComputesGroupCOMPositions)
400 const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8 };
401 topManager_.requestVelocities();
402 topManager_.requestForces();
403 topManager_.initAtoms(9);
404 testSingleStatic(POS_ALL, POS_MASS, true, group);
407 TEST_F(PositionCalculationTest, ComputesPositionsWithCompleteWhole)
409 const int group[] = { 0, 1, 2, 3, 4, 8 };
410 topManager_.initAtoms(9);
411 topManager_.initUniformResidues(3);
412 testSingleStatic(POS_RES, POS_COMPLWHOLE, true, group);
415 TEST_F(PositionCalculationTest, ComputesPositionsWithCompleteMax)
417 const int maxGroup[] = { 0, 1, 4, 5, 6, 8 };
418 const int evalGroup[] = { 0, 1, 5, 6 };
419 topManager_.initAtoms(9);
420 topManager_.initUniformResidues(3);
421 testSingleDynamic(POS_RES, POS_COMPLMAX, true, maxGroup, evalGroup);
424 TEST_F(PositionCalculationTest, ComputesPositionMask)
426 const int maxGroup[] = { 0, 1, 2, 3, 4, 5 };
427 const int evalGroup[] = { 1, 2, 4 };
428 topManager_.initAtoms(6);
429 testSingleDynamic(POS_ATOM, POS_MASKONLY, false, maxGroup, evalGroup);
432 // TODO: Check for POS_ALL_PBC
434 TEST_F(PositionCalculationTest, HandlesIdenticalStaticCalculations)
436 const int group[] = { 0, 1, 4, 5, 6, 7 };
437 topManager_.initAtoms(9);
438 topManager_.initUniformResidues(3);
440 gmx_ana_poscalc_t *pc1 = createCalculation(POS_RES, 0);
441 gmx_ana_poscalc_t *pc2 = createCalculation(POS_RES, 0);
442 gmx_ana_poscalc_t *pc3 = createCalculation(POS_RES, 0);
443 setMaximumGroup(pc1, group);
444 setMaximumGroup(pc2, group);
445 setMaximumGroup(pc3, group);
446 gmx_ana_pos_t *p1 = initPositions(pc1, "Positions");
447 gmx_ana_pos_t *p2 = initPositions(pc2, "Positions");
448 gmx_ana_pos_t *p3 = initPositions(pc3, "Positions");
451 pcc_.initEvaluation();
453 generateCoordinates();
454 gmx::test::TestReferenceChecker frameCompound(
455 checker_.checkCompound("EvaluatedPositions", "Frame0"));
456 updateAndCheck(pc1, p1, group, &frameCompound, "Positions");
457 updateAndCheck(pc2, p2, group, &frameCompound, "Positions");
458 updateAndCheck(pc3, p3, group, &frameCompound, "Positions");
462 TEST_F(PositionCalculationTest, HandlesOverlappingStaticCalculations)
464 const int group1[] = { 0, 1, 4, 5 };
465 const int group2[] = { 4, 5, 7, 8 };
466 topManager_.initAtoms(9);
467 topManager_.initUniformResidues(3);
469 gmx_ana_poscalc_t *pc1 = createCalculation(POS_RES, 0);
470 gmx_ana_poscalc_t *pc2 = createCalculation(POS_RES, 0);
471 setMaximumGroup(pc1, group1);
472 setMaximumGroup(pc2, group2);
473 gmx_ana_pos_t *p1 = initPositions(pc1, "P1");
474 gmx_ana_pos_t *p2 = initPositions(pc2, "P2");
477 pcc_.initEvaluation();
479 generateCoordinates();
480 gmx::test::TestReferenceChecker frameCompound(
481 checker_.checkCompound("EvaluatedPositions", "Frame0"));
482 updateAndCheck(pc1, p1, group1, &frameCompound, "P1");
483 updateAndCheck(pc2, p2, group2, &frameCompound, "P2");
487 // TODO: Check for handling of more multiple calculation cases