2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019, 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/block.h"
71 #include "gromacs/topology/idef.h"
72 #include "gromacs/topology/ifunc.h"
73 #include "gromacs/topology/topology.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_
, testData
->md_
);
93 bool success
= constrain_shake(
96 testData
->invmass_
.data(),
99 as_rvec_array(testData
->x_
.data()),
100 as_rvec_array(testData
->xPrime_
.data()),
101 as_rvec_array(testData
->xPrime2_
.data()),
103 testData
->md_
.lambda
,
104 &testData
->dHdLambda_
,
106 as_rvec_array(testData
->v_
.data()),
107 testData
->computeVirial_
,
108 testData
->virialScaled_
,
110 gmx::ConstraintVariable::Positions
);
111 EXPECT_TRUE(success
) << "Test failed with a false return value in SHAKE.";
116 * Initialize and apply LINCS constraints.
118 * \param[in] testData Test data structure.
119 * \param[in] pbc Periodic boundary data.
121 void applyLincs(ConstraintsTestData
*testData
, t_pbc pbc
)
126 int warncount_lincs
= 0;
127 gmx_omp_nthreads_set(emntLINCS
, 1);
129 // Make blocka structure for faster LINCS setup
130 std::vector
<t_blocka
> at2con_mt
;
131 at2con_mt
.reserve(testData
->mtop_
.moltype
.size());
132 for (const gmx_moltype_t
&moltype
: testData
->mtop_
.moltype
)
134 // This function is in constr.cpp
135 at2con_mt
.push_back(make_at2con(moltype
,
136 testData
->mtop_
.ffparams
.iparams
,
137 flexibleConstraintTreatment(EI_DYNAMICS(testData
->ir_
.eI
))));
140 lincsd
= init_lincs(nullptr, testData
->mtop_
,
141 testData
->nflexcon_
, at2con_mt
,
143 testData
->ir_
.nLincsIter
, testData
->ir_
.nProjOrder
);
144 set_lincs(testData
->idef_
, testData
->md_
, EI_DYNAMICS(testData
->ir_
.eI
), &testData
->cr_
, lincsd
);
146 // Evaluate constraints
147 bool success
= constrain_lincs(false, testData
->ir_
, 0, lincsd
, testData
->md_
,
150 as_rvec_array(testData
->x_
.data()),
151 as_rvec_array(testData
->xPrime_
.data()),
152 as_rvec_array(testData
->xPrime2_
.data()),
154 testData
->md_
.lambda
, &testData
->dHdLambda_
,
156 as_rvec_array(testData
->v_
.data()),
157 testData
->computeVirial_
, testData
->virialScaled_
,
158 gmx::ConstraintVariable::Positions
,
160 maxwarn
, &warncount_lincs
);
161 EXPECT_TRUE(success
) << "Test failed with a false return value in LINCS.";
162 EXPECT_EQ(warncount_lincs
, 0) << "There were warnings in LINCS.";
163 for (unsigned int i
= 0; i
< testData
->mtop_
.moltype
.size(); i
++)
165 sfree(at2con_mt
.at(i
).index
);
166 sfree(at2con_mt
.at(i
).a
);
171 #if GMX_GPU != GMX_GPU_CUDA
173 * Stub for LINCS constraints on CUDA-enabled GPU to satisfy compiler.
175 * \param[in] testData Test data structure.
176 * \param[in] pbc Periodic boundary data.
178 void applyLincsCuda(ConstraintsTestData gmx_unused
*testData
, t_pbc gmx_unused pbc
)
180 FAIL() << "Dummy LINCS CUDA function was called instead of the real one.";