2 * This file is part of the GROMACS molecular simulation package.
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.
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/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/mdatom.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"
84 * Initialize and apply SHAKE constraints.
86 * \param[in] testData Test data structure.
87 * \param[in] pbc Periodic boundary data.
89 void applyShake(ConstraintsTestData* testData, t_pbc gmx_unused pbc)
91 shakedata* shaked = shake_init();
92 make_shake_sblock_serial(shaked, testData->idef_.get(), testData->md_);
93 bool success = constrain_shake(
94 nullptr, shaked, testData->invmass_.data(), *testData->idef_, testData->ir_,
95 as_rvec_array(testData->x_.data()), as_rvec_array(testData->xPrime_.data()),
96 as_rvec_array(testData->xPrime2_.data()), &testData->nrnb_, testData->md_.lambda,
97 &testData->dHdLambda_, testData->invdt_, as_rvec_array(testData->v_.data()),
98 testData->computeVirial_, testData->virialScaled_, false, gmx::ConstraintVariable::Positions);
99 EXPECT_TRUE(success) << "Test failed with a false return value in SHAKE.";
104 * Initialize and apply LINCS constraints.
106 * \param[in] testData Test data structure.
107 * \param[in] pbc Periodic boundary data.
109 void applyLincs(ConstraintsTestData* testData, t_pbc pbc)
114 int warncount_lincs = 0;
115 gmx_omp_nthreads_set(emntLINCS, 1);
117 // Make blocka structure for faster LINCS setup
118 std::vector<ListOfLists<int>> at2con_mt;
119 at2con_mt.reserve(testData->mtop_.moltype.size());
120 for (const gmx_moltype_t& moltype : testData->mtop_.moltype)
122 // This function is in constr.cpp
123 at2con_mt.push_back(make_at2con(moltype, testData->mtop_.ffparams.iparams,
124 flexibleConstraintTreatment(EI_DYNAMICS(testData->ir_.eI))));
127 lincsd = init_lincs(nullptr, testData->mtop_, testData->nflexcon_, at2con_mt, false,
128 testData->ir_.nLincsIter, testData->ir_.nProjOrder);
129 set_lincs(*testData->idef_, testData->md_, EI_DYNAMICS(testData->ir_.eI), &testData->cr_, lincsd);
131 // Evaluate constraints
132 bool success = constrain_lincs(
133 false, testData->ir_, 0, lincsd, testData->md_, &testData->cr_, &testData->ms_,
134 testData->x_.arrayRefWithPadding(), testData->xPrime_.arrayRefWithPadding(),
135 testData->xPrime2_.arrayRefWithPadding().unpaddedArrayRef(), pbc.box, &pbc,
136 testData->md_.lambda, &testData->dHdLambda_, testData->invdt_,
137 testData->v_.arrayRefWithPadding().unpaddedArrayRef(), testData->computeVirial_,
138 testData->virialScaled_, gmx::ConstraintVariable::Positions, &testData->nrnb_, maxwarn,
140 EXPECT_TRUE(success) << "Test failed with a false return value in LINCS.";
141 EXPECT_EQ(warncount_lincs, 0) << "There were warnings in LINCS.";
145 #if GMX_GPU != GMX_GPU_CUDA
147 * Stub for GPU version of LINCS constraints to satisfy compiler.
149 * \param[in] testData Test data structure.
150 * \param[in] pbc Periodic boundary data.
152 void applyLincsGpu(ConstraintsTestData gmx_unused* testData, t_pbc gmx_unused pbc)
154 FAIL() << "Dummy LINCS CUDA function was called instead of the real one.";