a48c90ff74ed22dca089fde83ab25b8ad31c8021
[alexxy/gromacs.git] / src / gromacs / mdlib / tests / leapfrogtestrunners.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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.
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 Runner for CPU-based implementation of the integrator.
37  *
38  * \author Artem Zhmurov <zhmurov@gmail.com>
39  * \ingroup module_mdlib
40  */
41 #include "gmxpre.h"
42
43 #include "leapfrogtestrunners.h"
44
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/mdlib/gmx_omp_nthreads.h"
47 #include "gromacs/mdlib/update.h"
48 #include "gromacs/utility/arrayref.h"
49
50 #include "testutils/testasserts.h"
51
52 #include "leapfrogtestdata.h"
53
54
55 namespace gmx
56 {
57 namespace test
58 {
59
60 void LeapFrogHostTestRunner::integrate(LeapFrogTestData* testData, int numSteps)
61 {
62     testData->state_.x.resizeWithPadding(testData->numAtoms_);
63     testData->state_.v.resizeWithPadding(testData->numAtoms_);
64     for (int i = 0; i < testData->numAtoms_; i++)
65     {
66         testData->state_.x[i] = testData->x_[i];
67         testData->state_.v[i] = testData->v_[i];
68     }
69
70     gmx_omp_nthreads_set(ModuleMultiThread::Update, 1);
71
72     for (int step = 0; step < numSteps; step++)
73     {
74         testData->update_->update_coords(
75                 testData->inputRecord_,
76                 step,
77                 testData->mdAtoms_.homenr,
78                 testData->mdAtoms_.havePartiallyFrozenAtoms,
79                 testData->mdAtoms_.ptype
80                         ? gmx::arrayRefFromArray(testData->mdAtoms_.ptype, testData->mdAtoms_.nr)
81                         : gmx::ArrayRef<ParticleType>{},
82                 testData->mdAtoms_.invmass
83                         ? gmx::arrayRefFromArray(testData->mdAtoms_.invmass, testData->mdAtoms_.nr)
84                         : gmx::ArrayRef<real>{},
85                 testData->mdAtoms_.invMassPerDim ? gmx::arrayRefFromArray(
86                         testData->mdAtoms_.invMassPerDim, testData->mdAtoms_.nr)
87                                                  : gmx::ArrayRef<rvec>{},
88                 &testData->state_,
89                 testData->f_,
90                 testData->forceCalculationData_,
91                 &testData->kineticEnergyData_,
92                 testData->velocityScalingMatrix_,
93                 etrtNONE,
94                 nullptr,
95                 false);
96         testData->update_->finish_update(testData->inputRecord_,
97                                          testData->mdAtoms_.havePartiallyFrozenAtoms,
98                                          testData->mdAtoms_.homenr,
99                                          &testData->state_,
100                                          nullptr,
101                                          false);
102     }
103     const auto xp = makeArrayRef(*testData->update_->xp()).subArray(0, testData->numAtoms_);
104     for (int i = 0; i < testData->numAtoms_; i++)
105     {
106         for (int d = 0; d < DIM; d++)
107         {
108             testData->x_[i][d]      = testData->state_.x[i][d];
109             testData->v_[i][d]      = testData->state_.v[i][d];
110             testData->xPrime_[i][d] = xp[i][d];
111         }
112     }
113 }
114
115 } // namespace test
116 } // namespace gmx