24#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
32#include <Cabana_Core.hpp>
33#include <Cabana_NeighborList.hpp>
35#ifdef ESPRESSO_CALIPER
36#include <caliper/cali.h>
43template <
class KokkosRangePolicy = Kokkos::RangePolicy<>>
47 if (Kokkos::num_threads() > 1) {
48 KokkosRangePolicy policy(start, end);
49 Kokkos::parallel_for(name, policy, kernel);
51 for (
auto p_index = start; p_index < end; ++p_index) {
62#ifdef ESPRESSO_ELECTROSTATICS
68#if defined(ESPRESSO_GAY_BERNE) or defined(ESPRESSO_DIPOLES)
72#ifdef ESPRESSO_DIPOLES
78 aosoa.
id(index) = p.
id();
83#ifdef ESPRESSO_EXCLUSIONS
86 aosoa.
flags(index) = 0;
92 auto const &verlet_criterion,
94 auto const &intra_operator,
auto const &inter_operator) {
96 auto const distance_function = detail::MinimalImageDistance{box_geo};
102 auto intra_kernel = [&cells, &distance_function, &verlet_criterion,
103 &id_to_index, &intra_operator, max_id](
const int i) {
104 auto &local_particles = cells[i]->particles();
105 for (
auto it = local_particles.begin(); it != local_particles.end(); ++it) {
106 auto const &p1 = *it;
107 if (p1.id() <= max_id) {
108 auto const ii = id_to_index(p1.id());
111 for (
auto jt = std::next(it); jt != local_particles.end(); ++jt) {
112 if ((*jt).id() <= max_id) {
113 if (verlet_criterion(p1, *jt, distance_function(p1, *jt))) {
114 auto const jj = id_to_index((*jt).id());
116 intra_operator(ii, jj);
126 auto inter_kernel = [&cells, &distance_function, &verlet_criterion,
127 &id_to_index, &inter_operator, max_id](
const int i) {
128 auto &local_particles = cells[i]->particles();
129 for (
auto const &p1 : local_particles) {
130 if (p1.id() <= max_id) {
131 auto const ii = id_to_index(p1.id());
134 for (
auto &neighbor : cells[i]->neighbors().red()) {
135 for (
auto const &p2 : neighbor->particles()) {
136 if (p2.id() <= max_id) {
137 if (verlet_criterion(p1, p2, distance_function(p1, p2))) {
138 auto const jj = id_to_index(p2.id());
140 inter_operator(ii, jj);
151 Kokkos::parallel_for(
"inter", cells.size(), intra_kernel);
154 Kokkos::parallel_for(
"intra", cells.size(), inter_kernel);
160 double const pair_cutoff,
auto const integ_switch) {
161#ifdef ESPRESSO_CALIPER
162 CALI_CXX_MARK_FUNCTION;
164 using execution_space = Kokkos::DefaultExecutionSpace;
165 using policy_type = Kokkos::RangePolicy<execution_space>;
168 auto const n_part = unique_particles.size();
170 auto &aosoa = cell_structure.
get_aosoa();
178#ifdef ESPRESSO_CALIPER
179 CALI_MARK_BEGIN(
"AoSoA commit full");
181 kokkos_parallel_range_for<policy_type>(
182 "AoSoA write", std::size_t{0}, n_part,
183 [&unique_particles, &aosoa, &id_to_index](
int const index) {
184 auto const &p = *unique_particles.at(index);
186 id_to_index(p.id()) = index;
189#ifdef ESPRESSO_CALIPER
190 CALI_MARK_END(
"AoSoA commit full");
198#ifdef ESPRESSO_CALIPER
199 CALI_MARK_BEGIN(
"Verlet list creation");
202 [&](std::span<Cell *const> cells,
BoxGeometry const &box,
205 std::move(cells), box, verlet_criterion, id_to_index, max_id,
206 [&](
const int i,
const int j) {
210 [&](
const int i,
const int j) {
217 <<
"Verlet list overflow detected: neighbor count exceeded "
218 "max_counts. Falling back to the link cell algorithm.";
222#ifdef ESPRESSO_CALIPER
223 CALI_MARK_END(
"Verlet list creation");
229#ifdef ESPRESSO_CALIPER
230 CALI_MARK_BEGIN(
"AoSoA commit partial");
232 kokkos_parallel_range_for<policy_type>(
233 "AoSoA write", std::size_t{0}, n_part,
234 [&unique_particles, &aosoa](
int const index) {
235 auto const &p = *unique_particles.at(index);
239#ifdef ESPRESSO_CALIPER
240 CALI_MARK_END(
"AoSoA commit partial");
245#ifdef ESPRESSO_ELECTROSTATICS
248 using execution_space = Kokkos::DefaultExecutionSpace;
249 using policy_type = Kokkos::RangePolicy<execution_space>;
251 auto const n_part = unique_particles.size();
252 auto &aosoa = cell_structure.
get_aosoa();
254 kokkos_parallel_range_for<policy_type>(
255 "Views update charges", std::size_t{0}, n_part,
256 [&unique_particles, &aosoa](std::size_t
const index) {
257 aosoa.charge(index) = unique_particles.at(index)->q();
264 double bond_cutoff,
auto const &verlet_criterion,
265 auto const integ_switch) {
266 using execution_space = Kokkos::DefaultExecutionSpace;
269 if (bond_cutoff >= 0.) {
270#ifdef ESPRESSO_CALIPER
271 CALI_MARK_BEGIN(
"cabana_bond_loop");
274#ifdef ESPRESSO_CALIPER
275 CALI_MARK_END(
"cabana_bond_loop");
280 if (pair_cutoff > 0.) {
281#ifdef ESPRESSO_CALIPER
282 CALI_MARK_BEGIN(
"cabana_pair_loop");
287 Kokkos::RangePolicy<execution_space> policy(
289 Cabana::neighbor_parallel_for(policy, forces_kernel, verlet_list,
290 Cabana::FirstNeighborsTag(),
291 Cabana::SerialOpTag());
294 [&](std::span<Cell *const> cells,
BoxGeometry const &box) {
296 std::move(cells), box, verlet_criterion,
299 [&](
const int i,
const int j) {
303 [&](
const int i,
const int j) {
310#ifdef ESPRESSO_CALIPER
311 CALI_MARK_END(
"cabana_pair_loop");
#define ESPRESSO_ATTR_ALWAYS_INLINE
@ INTEG_METHOD_STEEPEST_DESCENT
Describes a cell structure / cell system.
void rebuild_verlet_list_cabana(auto &&kernel, bool rebuild_verlet_list)
auto prepare_verlet_list_cabana(double cutoff)
Reset local properties of the Verlet list.
int get_cached_max_local_particle_id() const
auto const & get_unique_particles() const
unsigned get_resort_particles() const
Get the currently scheduled resort level.
void bond_loop(BondKernel const &bond_kernel)
Bonded pair loop.
void cell_list_loop(auto &&kernel)
auto const & get_verlet_list_cabana() const
KOKKOS_INLINE_FUNCTION void addNeighborLB(int pid, int nid)
KOKKOS_INLINE_FUNCTION void addNeighbor(int pid, int nid)
KOKKOS_INLINE_FUNCTION bool hasOverflow() const
#define runtimeWarningMsg()
Vector< T, 3 > convert_quaternion_to_director(Quaternion< T > const &quat)
Convert quaternion to director.
ESPRESSO_ATTR_ALWAYS_INLINE void update_cabana_state(CellStructure &cell_structure, auto const &verlet_criterion, double const pair_cutoff, auto const integ_switch)
ESPRESSO_ATTR_ALWAYS_INLINE void link_cell_kokkos(std::span< Cell *const > cells, BoxGeometry const &box_geo, auto const &verlet_criterion, Kokkos::View< int * > const &id_to_index, int const max_id, auto const &intra_operator, auto const &inter_operator)
ESPRESSO_ATTR_ALWAYS_INLINE void commit_particle(Particle const &p, auto const index, CellStructure::AoSoA_pack &aosoa, bool const rebuild)
void cabana_short_range(auto const &bond_kernel, auto const &forces_kernel, CellStructure &cell_structure, double pair_cutoff, double bond_cutoff, auto const &verlet_criterion, auto const integ_switch)
ESPRESSO_ATTR_ALWAYS_INLINE void update_aosoa_charges(CellStructure &cell_structure)
ESPRESSO_ATTR_ALWAYS_INLINE void kokkos_parallel_range_for(auto const &name, auto start, auto end, auto const &kernel)
void set_vector_at(Kokkos::View< double *[3], array_layout, Kokkos::HostSpace > &view, std::size_t i, Utils::Vector3d const &value)
PositionViewType position
void set_has_exclusion(std::size_t i, bool value)
DirectorViewType director
VelocityViewType velocity
Struct holding all information for one particle.
Utils::compact_vector< int > & exclusions()
auto const & quat() const
auto const & type() const
auto const & dipm() const