2 * This file is part of the GROMACS molecular simulation package.
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.
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.
36 * \brief Function to run SHAKE and LINCS on CPU.
38 * Functions used in the test to apply constraints on the test data:
39 * CPU-based implementation and a stub for GPU-based implementation.
41 * \author Artem Zhmurov <zhmurov@gmail.com>
42 * \ingroup module_mdlib
47 #include "constrtestrunners.h"
58 #include <gtest/gtest.h>
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"
76 #include "testutils/testasserts.h"
83 void ShakeConstraintsRunner::applyConstraints(ConstraintsTestData* testData, t_pbc /* pbc */)
86 make_shake_sblock_serial(&shaked, testData->idef_.get(), testData->numAtoms_);
87 bool success = constrain_shake(nullptr,
98 &testData->dHdLambda_,
101 testData->computeVirial_,
102 testData->virialScaled_,
104 gmx::ConstraintVariable::Positions);
105 EXPECT_TRUE(success) << "Test failed with a false return value in SHAKE.";
108 void LincsConstraintsRunner::applyConstraints(ConstraintsTestData* testData, t_pbc pbc)
113 int warncount_lincs = 0;
114 gmx_omp_nthreads_set(ModuleMultiThread::Lincs, 1);
116 // Communication record
122 gmx_multisim_t ms{ 1, 0, MPI_COMM_NULL, MPI_COMM_NULL };
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)
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))));
135 lincsd = init_lincs(nullptr,
140 testData->ir_.nLincsIter,
141 testData->ir_.nProjOrder);
142 set_lincs(*testData->idef_,
146 EI_DYNAMICS(testData->ir_.eI),
150 // Evaluate constraints
151 bool success = constrain_lincs(false,
158 testData->x_.arrayRefWithPadding(),
159 testData->xPrime_.arrayRefWithPadding(),
160 testData->xPrime2_.arrayRefWithPadding().unpaddedArrayRef(),
163 testData->hasMassPerturbed_,
165 &testData->dHdLambda_,
167 testData->v_.arrayRefWithPadding().unpaddedArrayRef(),
168 testData->computeVirial_,
169 testData->virialScaled_,
170 gmx::ConstraintVariable::Positions,
174 EXPECT_TRUE(success) << "Test failed with a false return value in LINCS.";
175 EXPECT_EQ(warncount_lincs, 0) << "There were warnings in LINCS.";
180 void LincsDeviceConstraintsRunner::applyConstraints(ConstraintsTestData* /* testData */, t_pbc /* pbc */)
182 GMX_UNUSED_VALUE(testDevice_);
183 FAIL() << "Dummy LINCS CUDA function was called instead of the real one.";
185 #endif // !GMX_GPU_CUDA