ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
short_range_cabana.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2025 The ESPResSo project
3 *
4 * This file is part of ESPResSo.
5 *
6 * ESPResSo is free software: you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation, either version 3 of the License, or
9 * (at your option) any later version.
10 *
11 * ESPResSo is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program. If not, see <http://www.gnu.org/licenses/>.
18 */
19
20#pragma once
21
22#include <config/config.hpp>
23
24#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
25
27
28#include "aosoa_pack.hpp"
30#include "forces_cabana.hpp"
31
32#include <Cabana_Core.hpp>
33#include <Cabana_NeighborList.hpp>
34
35#ifdef ESPRESSO_CALIPER
36#include <caliper/cali.h>
37#endif
38
39#include <iterator>
40#include <span>
41#include <utility>
42
43template <class KokkosRangePolicy = Kokkos::RangePolicy<>>
45kokkos_parallel_range_for(auto const &name, auto start, auto end,
46 auto const &kernel) {
47 if (Kokkos::num_threads() > 1) {
48 KokkosRangePolicy policy(start, end);
49 Kokkos::parallel_for(name, policy, kernel);
50 } else {
51 for (auto p_index = start; p_index < end; ++p_index) {
52 kernel(p_index);
53 }
54 }
55}
56
58commit_particle(Particle const &p, auto const index,
59 CellStructure::AoSoA_pack &aosoa, bool const rebuild) {
60 // Always commit: positions, velocities, charges, directors, dipm
61 aosoa.set_vector_at(aosoa.position, index, p.pos());
62#ifdef ESPRESSO_ELECTROSTATICS
63 aosoa.charge(index) = p.q();
64#endif
65#ifdef ESPRESSO_DPD
66 aosoa.set_vector_at(aosoa.velocity, index, p.v());
67#endif
68#if defined(ESPRESSO_GAY_BERNE) or defined(ESPRESSO_DIPOLES)
69 aosoa.set_vector_at(aosoa.director, index,
71#endif
72#ifdef ESPRESSO_DIPOLES
73 aosoa.dipm(index) = p.dipm();
74#endif
75
76 // Only commit on rebuild: id, type
77 if (rebuild) {
78 aosoa.id(index) = p.id();
79 aosoa.type(index) = p.type();
80 }
81
82 // Always update exclusion flags (they can change during simulation)
83#ifdef ESPRESSO_EXCLUSIONS
84 aosoa.set_has_exclusion(index, !p.exclusions().empty());
85#else
86 aosoa.flags(index) = 0;
87#endif
88}
89
91link_cell_kokkos(std::span<Cell *const> cells, BoxGeometry const &box_geo,
92 auto const &verlet_criterion,
93 Kokkos::View<int *> const &id_to_index, int const max_id,
94 auto const &intra_operator, auto const &inter_operator) {
95
96 auto const distance_function = detail::MinimalImageDistance{box_geo};
97
98 // implementation detail: max_id refers to the max local particle id,
99 // but ghost particles from other ranks may have larger particle ids;
100 // -1 is used as a sentinel value for particle ids from other threads
101
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());
109 if (ii >= 0) {
110 // pairs in this cell
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());
115 if (jj >= 0) {
116 intra_operator(ii, jj);
117 }
118 }
119 }
120 }
121 }
122 }
123 }
124 };
125
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());
132 if (ii >= 0) {
133 // pairs with neighboring cells
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());
139 if (jj >= 0) {
140 inter_operator(ii, jj);
141 }
142 }
143 }
144 }
145 }
146 }
147 }
148 }
149 };
150
151 Kokkos::parallel_for("inter", cells.size(), intra_kernel);
152 Kokkos::fence();
153
154 Kokkos::parallel_for("intra", cells.size(), inter_kernel);
155 Kokkos::fence();
156}
157
159update_cabana_state(CellStructure &cell_structure, auto const &verlet_criterion,
160 double const pair_cutoff, auto const integ_switch) {
161#ifdef ESPRESSO_CALIPER
162 CALI_CXX_MARK_FUNCTION;
163#endif
164 using execution_space = Kokkos::DefaultExecutionSpace;
165 using policy_type = Kokkos::RangePolicy<execution_space>;
166 auto const rebuild = cell_structure.prepare_verlet_list_cabana(pair_cutoff);
167 auto const &unique_particles = cell_structure.get_unique_particles();
168 auto const n_part = unique_particles.size();
169 auto const max_id = cell_structure.get_cached_max_local_particle_id();
170 auto &aosoa = cell_structure.get_aosoa();
171
172 if (rebuild) {
173 auto &id_to_index = cell_structure.get_id_to_index();
174
175 // ===================================================
176 // Fill particle storage (full commit)
177 // ===================================================
178#ifdef ESPRESSO_CALIPER
179 CALI_MARK_BEGIN("AoSoA commit full");
180#endif
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);
185 commit_particle(p, index, aosoa, true);
186 id_to_index(p.id()) = index;
187 });
188 Kokkos::fence();
189#ifdef ESPRESSO_CALIPER
190 CALI_MARK_END("AoSoA commit full");
191#endif
192
193 // ===================================================
194 // Get Verlet pairs and fill Verlet list
195 // ===================================================
196 bool rebuild_vl = (integ_switch != INTEG_METHOD_STEEPEST_DESCENT and
197 cell_structure.use_verlet_list);
198#ifdef ESPRESSO_CALIPER
199 CALI_MARK_BEGIN("Verlet list creation");
200#endif
201 cell_structure.rebuild_verlet_list_cabana(
202 [&](std::span<Cell *const> cells, BoxGeometry const &box,
203 CellStructure::ListType &verlet_list) {
205 std::move(cells), box, verlet_criterion, id_to_index, max_id,
206 [&](const int i, const int j) {
207 // intra cell loop
208 verlet_list.addNeighborLB(i, j);
209 },
210 [&](const int i, const int j) {
211 // inter cell loop
212 verlet_list.addNeighbor(i, j);
213 });
214 if (verlet_list.hasOverflow()) {
215 cell_structure.use_verlet_list = false;
217 << "Verlet list overflow detected: neighbor count exceeded "
218 "max_counts. Falling back to the link cell algorithm.";
219 }
220 },
221 rebuild_vl);
222#ifdef ESPRESSO_CALIPER
223 CALI_MARK_END("Verlet list creation");
224#endif
225 } else {
226 // ===================================================
227 // Fill particle storage (partial update)
228 // ===================================================
229#ifdef ESPRESSO_CALIPER
230 CALI_MARK_BEGIN("AoSoA commit partial");
231#endif
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);
236 commit_particle(p, index, aosoa, false);
237 });
238 Kokkos::fence();
239#ifdef ESPRESSO_CALIPER
240 CALI_MARK_END("AoSoA commit partial");
241#endif
242 }
243}
244
245#ifdef ESPRESSO_ELECTROSTATICS
248 using execution_space = Kokkos::DefaultExecutionSpace;
249 using policy_type = Kokkos::RangePolicy<execution_space>;
250 auto const &unique_particles = cell_structure.get_unique_particles();
251 auto const n_part = unique_particles.size();
252 auto &aosoa = cell_structure.get_aosoa();
253
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();
258 });
259}
260#endif
261
262void cabana_short_range(auto const &bond_kernel, auto const &forces_kernel,
263 CellStructure &cell_structure, double pair_cutoff,
264 double bond_cutoff, auto const &verlet_criterion,
265 auto const integ_switch) {
266 using execution_space = Kokkos::DefaultExecutionSpace;
267 assert(cell_structure.get_resort_particles() == Cells::RESORT_NONE);
268
269 if (bond_cutoff >= 0.) {
270#ifdef ESPRESSO_CALIPER
271 CALI_MARK_BEGIN("cabana_bond_loop");
272#endif
273 cell_structure.bond_loop(bond_kernel);
274#ifdef ESPRESSO_CALIPER
275 CALI_MARK_END("cabana_bond_loop");
276#endif
277 }
278
279 // Cabana short range loop
280 if (pair_cutoff > 0.) {
281#ifdef ESPRESSO_CALIPER
282 CALI_MARK_BEGIN("cabana_pair_loop");
283#endif
284 if (integ_switch != INTEG_METHOD_STEEPEST_DESCENT and
285 cell_structure.use_verlet_list) {
286 auto const &verlet_list = cell_structure.get_verlet_list_cabana();
287 Kokkos::RangePolicy<execution_space> policy(
288 std::size_t{0}, cell_structure.get_unique_particles().size());
289 Cabana::neighbor_parallel_for(policy, forces_kernel, verlet_list,
290 Cabana::FirstNeighborsTag(),
291 Cabana::SerialOpTag());
292 } else {
293 cell_structure.cell_list_loop(
294 [&](std::span<Cell *const> cells, BoxGeometry const &box) {
296 std::move(cells), box, verlet_criterion,
297 cell_structure.get_id_to_index(),
298 cell_structure.get_cached_max_local_particle_id(),
299 [&](const int i, const int j) {
300 // intra cell loop
301 forces_kernel(i, j);
302 },
303 [&](const int i, const int j) {
304 // inter cell loop
305 forces_kernel(i, j);
306 });
307 });
308 }
309 Kokkos::fence();
310#ifdef ESPRESSO_CALIPER
311 CALI_MARK_END("cabana_pair_loop");
312#endif
313 }
314}
315
316#endif // ESPRESSO_SHARED_MEMORY_PARALLELISM
#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 & get_id_to_index()
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.
Definition Particle.hpp:450
Utils::compact_vector< int > & exclusions()
Definition Particle.hpp:665
auto const & quat() const
Definition Particle.hpp:532
auto const & q() const
Definition Particle.hpp:593
auto const & v() const
Definition Particle.hpp:488
auto const & type() const
Definition Particle.hpp:473
auto const & dipm() const
Definition Particle.hpp:548
auto const & pos() const
Definition Particle.hpp:486
auto const & id() const
Definition Particle.hpp:469