1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019,2020, 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 #include "gmxpre.h"
36 
37 #include "gromacs/mdlib/updategroupscog.h"
38 
39 #include <cstdlib>
40 
41 #include <gtest/gtest.h>
42 
43 #include "gromacs/math/functions.h"
44 #include "gromacs/mdlib/updategroups.h"
45 #include "gromacs/topology/idef.h"
46 #include "gromacs/topology/ifunc.h"
47 #include "gromacs/topology/mtop_util.h"
48 #include "gromacs/topology/topology.h"
49 
50 #include "testutils/refdata.h"
51 #include "testutils/testasserts.h"
52 
53 namespace gmx
54 {
55 
56 namespace
57 {
58 
59 //! Database of 51 water atom input positions (taken from spc216.gro) for use as test inputs.
60 gmx::RVec positions[] = { { .130, -.041, -.291 },  { .120, -.056, -.192 },  { .044, -.005, -.327 },
61                           { -.854, -.406, .477 },  { -.900, -.334, .425 },  { -.858, -.386, .575 },
62                           { .351, -.061, .853 },   { .401, -.147, .859 },   { .416, .016, .850 },
63                           { -.067, -.796, .873 },  { -.129, -.811, .797 },  { -.119, -.785, .958 },
64                           { -.635, -.312, -.356 }, { -.629, -.389, -.292 }, { -.687, -.338, -.436 },
65                           { .321, -.919, .242 },   { .403, -.880, .200 },   { .294, -1.001, .193 },
66                           { -.404, .735, .728 },   { -.409, .670, .803 },   { -.324, .794, .741 },
67                           { .461, -.596, -.135 },  { .411, -.595, -.221 },  { .398, -.614, -.059 },
68                           { -.751, -.086, .237 },  { -.811, -.148, .287 },  { -.720, -.130, .152 },
69                           { .202, .285, -.364 },   { .122, .345, -.377 },   { .192, .236, -.278 },
70                           { -.230, -.485, .081 },  { -.262, -.391, .071 },  { -.306, -.548, .069 },
71                           { .464, -.119, .323 },   { .497, -.080, .409 },   { .540, -.126, .258 },
72                           { -.462, .107, .426 },   { -.486, .070, .336 },   { -.363, .123, .430 },
73                           { .249, -.077, -.621 },  { .306, -.142, -.571 },  { .233, -.110, -.714 },
74                           { -.922, -.164, .904 },  { -.842, -.221, .925 },  { -.971, -.204, .827 },
75                           { .382, .700, .480 },    { .427, .610, .477 },    { .288, .689, .513 },
76                           { .781, .264, -.113 },   { .848, .203, -.070 },   { .708, .283, -.048 } };
77 
TEST(UpdateGroupsCog,ComputesCogs)78 TEST(UpdateGroupsCog, ComputesCogs)
79 {
80     const int settleType     = 0;
81     const int atomsPerSettle = NRAL(F_SETTLE);
82     const int numAtoms       = sizeof(positions) / sizeof(positions[0]);
83     const int numMolecules   = gmx::exactDiv(numAtoms, atomsPerSettle);
84 
85     // Set up the topology.
86     gmx_mtop_t mtop;
87 
88     gmx_moltype_t moltype;
89     moltype.atoms.nr         = atomsPerSettle;
90     std::vector<int>& iatoms = moltype.ilist[F_SETTLE].iatoms;
91     iatoms.push_back(settleType);
92     iatoms.push_back(0);
93     iatoms.push_back(1);
94     iatoms.push_back(2);
95     mtop.moltype.push_back(moltype);
96 
97     mtop.molblock.resize(1);
98     mtop.molblock[0].type = 0;
99     mtop.molblock[0].nmol = numMolecules;
100     mtop.natoms           = numAtoms;
101     mtop.finalize();
102 
103     // Set up the SETTLE parameters.
104     const real dOH = 0.1;
105     const real dHH = 0.1633;
106     t_iparams  iparams;
107     iparams.settle.doh = dOH;
108     iparams.settle.dhh = dHH;
109     mtop.ffparams.iparams.push_back(iparams);
110 
111     // Run the test
112     auto updateGroups = makeUpdateGroups(mtop);
113     real temperature  = 300;
114 
115     UpdateGroupsCog updateGroupsCog(mtop, updateGroups, temperature, numAtoms);
116 
117     EXPECT_FLOAT_EQ(updateGroupsCog.maxUpdateGroupRadius(), 0.083887339);
118 
119     std::vector<int> globalAtomIndices(numAtoms);
120     for (int i = 0; i < numAtoms; i++)
121     {
122         globalAtomIndices[i] = i;
123     }
124 
125     // Randomize the atom order
126     for (int i = 0; i < numAtoms; i++)
127     {
128         int a1 = std::rand() % numAtoms;
129         int a2 = std::rand() % numAtoms;
130         if (a1 != a2)
131         {
132             std::swap(globalAtomIndices[a1], globalAtomIndices[a2]);
133             std::swap(positions[a1], positions[a2]);
134         }
135     }
136 
137     updateGroupsCog.addCogs(globalAtomIndices, positions);
138 
139     EXPECT_EQ(updateGroupsCog.numCogs(), numMolecules);
140 
141     std::vector<gmx::RVec> cogPerAtom(numAtoms);
142     for (int i = 0; i < numAtoms; i++)
143     {
144         cogPerAtom[globalAtomIndices[i]] = updateGroupsCog.cogForAtom(i);
145     }
146     gmx::test::TestReferenceData    data;
147     gmx::test::TestReferenceChecker checker(data.rootChecker());
148     checker.checkSequence(cogPerAtom.begin(), cogPerAtom.end(), "cogPerAtom");
149 
150     updateGroupsCog.clear();
151     EXPECT_EQ(updateGroupsCog.numCogs(), 0);
152 }
153 
154 } // namespace
155 } // namespace gmx
156