Use ObservablesReducer for LINCS RMSD computation
[alexxy/gromacs.git] / src / gromacs / mdlib / tests / constrtestrunners.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019,2020,2021, 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 Function to run SHAKE and LINCS on CPU.
37  *
38  * Functions used in the test to apply constraints on the test data:
39  * CPU-based implementation and a stub for GPU-based implementation.
40  *
41  * \author Artem Zhmurov <zhmurov@gmail.com>
42  * \ingroup module_mdlib
43  */
44
45 #include "gmxpre.h"
46
47 #include "constrtestrunners.h"
48
49 #include "config.h"
50
51 #include <assert.h>
52
53 #include <cmath>
54
55 #include <algorithm>
56 #include <vector>
57
58 #include <gtest/gtest.h>
59
60 #include "gromacs/math/paddedvector.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/math/vectypes.h"
63 #include "gromacs/mdlib/constr.h"
64 #include "gromacs/mdlib/lincs.h"
65 #include "gromacs/mdlib/shake.h"
66 #include "gromacs/mdrunutility/multisim.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/pbcutil/pbc.h"
70 #include "gromacs/topology/idef.h"
71 #include "gromacs/topology/ifunc.h"
72 #include "gromacs/topology/topology.h"
73 #include "gromacs/utility/listoflists.h"
74 #include "gromacs/utility/unique_cptr.h"
75
76 #include "testutils/testasserts.h"
77
78 namespace gmx
79 {
80 namespace test
81 {
82
83 void ShakeConstraintsRunner::applyConstraints(ConstraintsTestData* testData, t_pbc /* pbc */)
84 {
85     shakedata shaked;
86     make_shake_sblock_serial(&shaked, testData->idef_.get(), testData->numAtoms_);
87     bool success = constrain_shake(nullptr,
88                                    &shaked,
89                                    testData->invmass_,
90                                    *testData->idef_,
91                                    testData->ir_,
92                                    testData->x_,
93                                    testData->xPrime_,
94                                    testData->xPrime2_,
95                                    nullptr,
96                                    &testData->nrnb_,
97                                    testData->lambda_,
98                                    &testData->dHdLambda_,
99                                    testData->invdt_,
100                                    testData->v_,
101                                    testData->computeVirial_,
102                                    testData->virialScaled_,
103                                    false,
104                                    gmx::ConstraintVariable::Positions);
105     EXPECT_TRUE(success) << "Test failed with a false return value in SHAKE.";
106 }
107
108 void LincsConstraintsRunner::applyConstraints(ConstraintsTestData* testData, t_pbc pbc)
109 {
110
111     Lincs* lincsd;
112     int    maxwarn         = 100;
113     int    warncount_lincs = 0;
114     gmx_omp_nthreads_set(ModuleMultiThread::Lincs, 1);
115
116     // Communication record
117     t_commrec cr;
118     cr.nnodes = 1;
119     cr.dd     = nullptr;
120
121     // Multi-sim record
122     gmx_multisim_t ms{ 1, 0, MPI_COMM_NULL, MPI_COMM_NULL };
123
124     // Make blocka structure for faster LINCS setup
125     std::vector<ListOfLists<int>> at2con_mt;
126     at2con_mt.reserve(testData->mtop_.moltype.size());
127     for (const gmx_moltype_t& moltype : testData->mtop_.moltype)
128     {
129         // This function is in constr.cpp
130         at2con_mt.push_back(make_at2con(moltype,
131                                         testData->mtop_.ffparams.iparams,
132                                         flexibleConstraintTreatment(EI_DYNAMICS(testData->ir_.eI))));
133     }
134     // Initialize LINCS
135     lincsd = init_lincs(nullptr,
136                         testData->mtop_,
137                         testData->nflexcon_,
138                         at2con_mt,
139                         false,
140                         testData->ir_.nLincsIter,
141                         testData->ir_.nProjOrder,
142                         nullptr);
143     set_lincs(*testData->idef_,
144               testData->numAtoms_,
145               testData->invmass_,
146               testData->lambda_,
147               EI_DYNAMICS(testData->ir_.eI),
148               &cr,
149               lincsd);
150
151     // Evaluate constraints
152     bool success = constrain_lincs(false,
153                                    testData->ir_,
154                                    0,
155                                    lincsd,
156                                    testData->invmass_,
157                                    &cr,
158                                    &ms,
159                                    testData->x_.arrayRefWithPadding(),
160                                    testData->xPrime_.arrayRefWithPadding(),
161                                    testData->xPrime2_.arrayRefWithPadding().unpaddedArrayRef(),
162                                    pbc.box,
163                                    &pbc,
164                                    testData->hasMassPerturbed_,
165                                    testData->lambda_,
166                                    &testData->dHdLambda_,
167                                    testData->invdt_,
168                                    testData->v_.arrayRefWithPadding().unpaddedArrayRef(),
169                                    testData->computeVirial_,
170                                    testData->virialScaled_,
171                                    gmx::ConstraintVariable::Positions,
172                                    &testData->nrnb_,
173                                    maxwarn,
174                                    &warncount_lincs);
175     EXPECT_TRUE(success) << "Test failed with a false return value in LINCS.";
176     EXPECT_EQ(warncount_lincs, 0) << "There were warnings in LINCS.";
177     done_lincs(lincsd);
178 }
179
180 } // namespace test
181 } // namespace gmx