710e9ac48c1b83ba40293f68dfb41ca9fba3bf4b
[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, 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/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"
75
76 #include "testutils/testasserts.h"
77
78 namespace gmx
79 {
80 namespace test
81 {
82
83 /*! \brief
84  * Initialize and apply SHAKE constraints.
85  *
86  * \param[in] testData        Test data structure.
87  * \param[in] pbc             Periodic boundary data.
88  */
89 void applyShake(ConstraintsTestData* testData, t_pbc gmx_unused pbc)
90 {
91     shakedata* shaked = shake_init();
92     make_shake_sblock_serial(shaked, &testData->idef_, 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.";
100     done_shake(shaked);
101 }
102
103 /*! \brief
104  * Initialize and apply LINCS constraints.
105  *
106  * \param[in] testData        Test data structure.
107  * \param[in] pbc             Periodic boundary data.
108  */
109 void applyLincs(ConstraintsTestData* testData, t_pbc pbc)
110 {
111
112     Lincs* lincsd;
113     int    maxwarn         = 100;
114     int    warncount_lincs = 0;
115     gmx_omp_nthreads_set(emntLINCS, 1);
116
117     // Make blocka structure for faster LINCS setup
118     std::vector<t_blocka> at2con_mt;
119     at2con_mt.reserve(testData->mtop_.moltype.size());
120     for (const gmx_moltype_t& moltype : testData->mtop_.moltype)
121     {
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))));
125     }
126     // Initialize LINCS
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);
130
131     // Evaluate constraints
132     bool success = constrain_lincs(
133             false, testData->ir_, 0, lincsd, testData->md_, &testData->cr_, &testData->ms_,
134             as_rvec_array(testData->x_.data()), as_rvec_array(testData->xPrime_.data()),
135             as_rvec_array(testData->xPrime2_.data()), pbc.box, &pbc, testData->md_.lambda,
136             &testData->dHdLambda_, testData->invdt_, as_rvec_array(testData->v_.data()),
137             testData->computeVirial_, testData->virialScaled_, gmx::ConstraintVariable::Positions,
138             &testData->nrnb_, maxwarn, &warncount_lincs);
139     EXPECT_TRUE(success) << "Test failed with a false return value in LINCS.";
140     EXPECT_EQ(warncount_lincs, 0) << "There were warnings in LINCS.";
141     for (auto& moltype : at2con_mt)
142     {
143         sfree(moltype.index);
144         sfree(moltype.a);
145     }
146     done_lincs(lincsd);
147 }
148
149 #if GMX_GPU != GMX_GPU_CUDA
150 /*! \brief
151  * Stub for LINCS constraints on CUDA-enabled GPU to satisfy compiler.
152  *
153  * \param[in] testData        Test data structure.
154  * \param[in] pbc             Periodic boundary data.
155  */
156 void applyLincsCuda(ConstraintsTestData gmx_unused* testData, t_pbc gmx_unused pbc)
157 {
158     FAIL() << "Dummy LINCS CUDA function was called instead of the real one.";
159 }
160 #endif
161
162 } // namespace test
163 } // namespace gmx