f57a9de20dc169a23e6d1cb99d5aed74aa61c269
[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, 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(
88             nullptr, &shaked, testData->invmass_.data(), *testData->idef_, testData->ir_, testData->x_,
89             testData->xPrime_, testData->xPrime2_, nullptr, &testData->nrnb_, testData->lambda_,
90             &testData->dHdLambda_, testData->invdt_, testData->v_, testData->computeVirial_,
91             testData->virialScaled_, false, gmx::ConstraintVariable::Positions);
92     EXPECT_TRUE(success) << "Test failed with a false return value in SHAKE.";
93 }
94
95 void LincsConstraintsRunner::applyConstraints(ConstraintsTestData* testData, t_pbc pbc)
96 {
97
98     Lincs* lincsd;
99     int    maxwarn         = 100;
100     int    warncount_lincs = 0;
101     gmx_omp_nthreads_set(emntLINCS, 1);
102
103     // Communication record
104     t_commrec cr;
105     cr.nnodes = 1;
106     cr.dd     = nullptr;
107
108     // Multi-sim record
109     gmx_multisim_t ms{ 1, 0, MPI_COMM_NULL, MPI_COMM_NULL };
110
111     // Make blocka structure for faster LINCS setup
112     std::vector<ListOfLists<int>> at2con_mt;
113     at2con_mt.reserve(testData->mtop_.moltype.size());
114     for (const gmx_moltype_t& moltype : testData->mtop_.moltype)
115     {
116         // This function is in constr.cpp
117         at2con_mt.push_back(make_at2con(moltype, testData->mtop_.ffparams.iparams,
118                                         flexibleConstraintTreatment(EI_DYNAMICS(testData->ir_.eI))));
119     }
120     // Initialize LINCS
121     lincsd = init_lincs(nullptr, testData->mtop_, testData->nflexcon_, at2con_mt, false,
122                         testData->ir_.nLincsIter, testData->ir_.nProjOrder);
123     set_lincs(*testData->idef_, testData->numAtoms_, testData->invmass_.data(), testData->lambda_,
124               EI_DYNAMICS(testData->ir_.eI), &cr, lincsd);
125
126     // Evaluate constraints
127     bool success = constrain_lincs(
128             false, testData->ir_, 0, lincsd, testData->invmass_.data(), &cr, &ms,
129             testData->x_.arrayRefWithPadding(), testData->xPrime_.arrayRefWithPadding(),
130             testData->xPrime2_.arrayRefWithPadding().unpaddedArrayRef(), pbc.box, &pbc,
131             testData->hasMassPerturbed_, testData->lambda_, &testData->dHdLambda_, testData->invdt_,
132             testData->v_.arrayRefWithPadding().unpaddedArrayRef(), testData->computeVirial_,
133             testData->virialScaled_, gmx::ConstraintVariable::Positions, &testData->nrnb_, maxwarn,
134             &warncount_lincs);
135     EXPECT_TRUE(success) << "Test failed with a false return value in LINCS.";
136     EXPECT_EQ(warncount_lincs, 0) << "There were warnings in LINCS.";
137     done_lincs(lincsd);
138 }
139
140 #if !GMX_GPU_CUDA
141 void LincsDeviceConstraintsRunner::applyConstraints(ConstraintsTestData* /* testData */, t_pbc /* pbc */)
142 {
143     GMX_UNUSED_VALUE(testDevice_);
144     FAIL() << "Dummy LINCS CUDA function was called instead of the real one.";
145 }
146 #endif // !GMX_GPU_CUDA
147
148 } // namespace test
149 } // namespace gmx