2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
38 * Implements the Nbnxm class
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_nbnxm
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/nbnxm/atomdata.h"
50 #include "gromacs/timing/wallcycle.h"
52 #include "pairlistsets.h"
53 #include "pairsearch.h"
57 void nbnxn_put_on_grid(nonbonded_verlet_t
*nb_verlet
,
60 const rvec lowerCorner
,
61 const rvec upperCorner
,
62 const gmx::UpdateGroupsCog
*updateGroupsCog
,
67 gmx::ArrayRef
<const gmx::RVec
> x
,
71 nb_verlet
->pairSearch_
->putOnGrid(box
, ddZone
, lowerCorner
, upperCorner
,
72 updateGroupsCog
, atomStart
, atomEnd
, atomDensity
,
73 atinfo
, x
, numAtomsMoved
, move
,
74 nb_verlet
->nbat
.get());
77 /* Calls nbnxn_put_on_grid for all non-local domains */
78 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t
*nbv
,
79 const struct gmx_domdec_zones_t
*zones
,
81 gmx::ArrayRef
<const gmx::RVec
> x
)
83 for (int zone
= 1; zone
< zones
->n
; zone
++)
86 for (int d
= 0; d
< DIM
; d
++)
88 c0
[d
] = zones
->size
[zone
].bb_x0
[d
];
89 c1
[d
] = zones
->size
[zone
].bb_x1
[d
];
92 nbnxn_put_on_grid(nbv
, nullptr,
95 zones
->cg_range
[zone
],
96 zones
->cg_range
[zone
+1],
104 bool nonbonded_verlet_t::isDynamicPruningStepCpu(int64_t step
) const
106 return pairlistSets_
->isDynamicPruningStepCpu(step
);
109 bool nonbonded_verlet_t::isDynamicPruningStepGpu(int64_t step
) const
111 return pairlistSets_
->isDynamicPruningStepGpu(step
);
114 gmx::ArrayRef
<const int> nonbonded_verlet_t::getLocalAtomOrder() const
116 /* Return the atom order for the home cell (index 0) */
117 const Nbnxm::Grid
&grid
= pairSearch_
->gridSet().grids()[0];
119 const int numIndices
= grid
.atomIndexEnd() - grid
.firstAtomInColumn(0);
121 return gmx::constArrayRefFromArray(pairSearch_
->gridSet().atomIndices().data(), numIndices
);
124 void nonbonded_verlet_t::setLocalAtomOrder()
126 pairSearch_
->setLocalAtomOrder();
129 void nonbonded_verlet_t::setAtomProperties(const t_mdatoms
&mdatoms
,
132 nbnxn_atomdata_set(nbat
.get(), pairSearch_
->gridSet(), &mdatoms
, &atinfo
);
135 void nonbonded_verlet_t::setCoordinates(const Nbnxm::AtomLocality locality
,
136 const bool fillLocal
,
137 gmx::ArrayRef
<const gmx::RVec
> x
,
138 gmx_wallcycle
*wcycle
)
140 wallcycle_start(wcycle
, ewcNB_XF_BUF_OPS
);
141 wallcycle_sub_start(wcycle
, ewcsNB_X_BUF_OPS
);
143 nbnxn_atomdata_copy_x_to_nbat_x(pairSearch_
->gridSet(), locality
, fillLocal
,
144 as_rvec_array(x
.data()),
147 wallcycle_sub_stop(wcycle
, ewcsNB_X_BUF_OPS
);
148 wallcycle_stop(wcycle
, ewcNB_XF_BUF_OPS
);
151 gmx::ArrayRef
<const int> nonbonded_verlet_t::getGridIndices() const
153 return pairSearch_
->gridSet().cells();
157 nonbonded_verlet_t::atomdata_add_nbat_f_to_f(const Nbnxm::AtomLocality locality
,
159 gmx_wallcycle
*wcycle
)
161 wallcycle_start(wcycle
, ewcNB_XF_BUF_OPS
);
162 wallcycle_sub_start(wcycle
, ewcsNB_F_BUF_OPS
);
164 reduceForces(nbat
.get(), locality
, pairSearch_
->gridSet(), f
);
166 wallcycle_sub_stop(wcycle
, ewcsNB_F_BUF_OPS
);
167 wallcycle_stop(wcycle
, ewcNB_XF_BUF_OPS
);
170 real
nonbonded_verlet_t::pairlistInnerRadius() const
172 return pairlistSets_
->params().rlistInner
;
175 real
nonbonded_verlet_t::pairlistOuterRadius() const
177 return pairlistSets_
->params().rlistOuter
;
180 void nonbonded_verlet_t::changePairlistRadii(real rlistOuter
,
183 pairlistSets_
->changePairlistRadii(rlistOuter
, rlistInner
);