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-2026 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
25
26#include "aosoa_pack.hpp"
29#include "forces_cabana.hpp"
30
31#include <Cabana_Core.hpp>
32#include <Cabana_NeighborList.hpp>
33
34#ifdef ESPRESSO_CALIPER
35#include <caliper/cali.h>
36#endif
37
38#include <iterator>
39#include <span>
40#include <utility>
41
42template <class KokkosRangePolicy = Kokkos::RangePolicy<>>
44kokkos_parallel_range_for(auto const &name, auto start, auto end,
45 auto const &kernel) {
46 if (Kokkos::num_threads() > 1) {
47 KokkosRangePolicy policy(start, end);
48 Kokkos::parallel_for(name, policy, kernel);
49 } else {
50 for (auto p_index = start; p_index < end; ++p_index) {
51 kernel(p_index);
52 }
53 }
54}
55
57commit_particle(Particle const &p, auto const index,
58 CellStructure::AoSoA_pack &aosoa, bool const rebuild) {
59 // Always commit: positions, velocities, charges, directors, dipm
60 aosoa.set_vector_at(aosoa.position, index, p.pos());
61#ifdef ESPRESSO_ELECTROSTATICS
62 aosoa.charge(index) = p.q();
63#endif
64 aosoa.set_vector_at(aosoa.velocity, index, p.v());
65#if defined(ESPRESSO_GAY_BERNE) or defined(ESPRESSO_DIPOLES)
66 aosoa.set_vector_at(aosoa.director, index,
68#endif
69#ifdef ESPRESSO_DIPOLES
70 aosoa.dipm(index) = p.dipm();
71#endif
72
73 // Only commit on rebuild: id, type
74 if (rebuild) {
75 aosoa.id(index) = p.id();
76 aosoa.type(index) = p.type();
77 aosoa.set_vector_at(aosoa.image, index, p.image_box());
78#ifdef ESPRESSO_MASS
79 aosoa.mass(index) = p.mass();
80#endif
81 }
82
83 // Always update exclusion flags (they can change during simulation)
84#ifdef ESPRESSO_EXCLUSIONS
85 aosoa.set_has_exclusion(index, !p.exclusions().empty());
86#else
87 aosoa.flags(index) = 0;
88#endif
89}
90
92link_cell_kokkos(std::span<Cell *const> cells, BoxGeometry const &box_geo,
93 auto const &verlet_criterion,
94 Kokkos::View<int *> const &id_to_index, int const max_id,
95 auto const &intra_operator, auto const &inter_operator) {
96
97 auto const distance_function = detail::MinimalImageDistance{box_geo};
98
99 // implementation detail: max_id refers to the max local particle id,
100 // but ghost particles from other ranks may have larger particle ids;
101 // -1 is used as a sentinel value for particle ids from other threads
102
103 auto intra_kernel = [&cells, &distance_function, &verlet_criterion,
104 &id_to_index, &intra_operator, max_id](const int i) {
105 auto &local_particles = cells[i]->particles();
106 for (auto it = local_particles.begin(); it != local_particles.end(); ++it) {
107 auto const &p1 = *it;
108 if (p1.id() <= max_id) {
109 auto const ii = id_to_index(p1.id());
110 if (ii >= 0) {
111 // pairs in this cell
112 for (auto jt = std::next(it); jt != local_particles.end(); ++jt) {
113 if ((*jt).id() <= max_id) {
114 if (verlet_criterion(p1, *jt, distance_function(p1, *jt))) {
115 auto const jj = id_to_index((*jt).id());
116 if (jj >= 0) {
117 intra_operator(ii, jj);
118 }
119 }
120 }
121 }
122 }
123 }
124 }
125 };
126
127 auto inter_kernel = [&cells, &distance_function, &verlet_criterion,
128 &id_to_index, &inter_operator, max_id](const int i) {
129 auto &local_particles = cells[i]->particles();
130 for (auto const &p1 : local_particles) {
131 if (p1.id() <= max_id) {
132 auto const ii = id_to_index(p1.id());
133 if (ii >= 0) {
134 // pairs with neighboring cells
135 for (auto &neighbor : cells[i]->neighbors().red()) {
136 for (auto const &p2 : neighbor->particles()) {
137 if (p2.id() <= max_id) {
138 if (verlet_criterion(p1, p2, distance_function(p1, p2))) {
139 auto const jj = id_to_index(p2.id());
140 if (jj >= 0) {
141 inter_operator(ii, jj);
142 }
143 }
144 }
145 }
146 }
147 }
148 }
149 }
150 };
151
152 Kokkos::parallel_for("intra", cells.size(), intra_kernel);
153 Kokkos::fence();
154
155 Kokkos::parallel_for("inter", cells.size(), inter_kernel);
156 Kokkos::fence();
157}
158
160update_cabana_state(CellStructure &cell_structure, auto const &verlet_criterion,
161 double const pair_cutoff, auto const integ_switch) {
162#ifdef ESPRESSO_CALIPER
163 CALI_CXX_MARK_FUNCTION;
164#endif
165 using execution_space = Kokkos::DefaultExecutionSpace;
166 using policy_type = Kokkos::RangePolicy<execution_space>;
167 auto const rebuild = cell_structure.prepare_verlet_list_cabana(pair_cutoff);
168 auto const &unique_particles = cell_structure.get_unique_particles();
169 auto const n_part = unique_particles.size();
170 auto const max_id = cell_structure.get_cached_max_local_particle_id();
171 auto &aosoa = cell_structure.get_aosoa();
172
173 if (rebuild) {
174 auto &id_to_index = cell_structure.get_id_to_index();
175
176 // ===================================================
177 // Fill particle storage (full commit)
178 // ===================================================
179#ifdef ESPRESSO_CALIPER
180 CALI_MARK_BEGIN("AoSoA commit full");
181#endif
182 int pair_count = 0;
183 int angle_count = 0;
184 int dihedral_count = 0;
185 kokkos_parallel_range_for<policy_type>(
186 "AoSoA write", std::size_t{0}, n_part,
187 [&unique_particles, &aosoa, &id_to_index, &cell_structure, &pair_count,
188 &angle_count, &dihedral_count](int const index) {
189 auto const &p = *unique_particles.at(index);
190 commit_particle(p, index, aosoa, true);
191 id_to_index(p.id()) = index;
192 if (not p.is_ghost()) {
193 cell_structure.update_bond_storage(pair_count, angle_count,
194 dihedral_count, p);
195 }
196 });
197 Kokkos::fence();
198 auto &bs = cell_structure.bond_state();
199 auto &pair_bond_list = bs.pair_list;
200 Kokkos::parallel_for("resolve_pair_bond_indices", pair_count,
201 [&pair_bond_list, &id_to_index](int idx) {
202 for (int col = 0; col < 2; ++col) {
203 pair_bond_list(idx, col) =
204 id_to_index(pair_bond_list(idx, col));
205 }
206 });
207 auto &angle_bond_list = bs.angle_list;
208 Kokkos::parallel_for("resolve_angle_bond_indices", angle_count,
209 [&angle_bond_list, &id_to_index](int idx) {
210 for (int col = 0; col < 3; ++col) {
211 angle_bond_list(idx, col) =
212 id_to_index(angle_bond_list(idx, col));
213 }
214 });
215 auto &dihedral_bond_list = bs.dihedral_list;
216 Kokkos::parallel_for("resolve_dihedral_bond_indices", dihedral_count,
217 [&dihedral_bond_list, &id_to_index](int idx) {
218 for (int col = 0; col < 4; ++col) {
219 dihedral_bond_list(idx, col) =
220 id_to_index(dihedral_bond_list(idx, col));
221 }
222 });
223 Kokkos::fence();
224#ifdef ESPRESSO_CALIPER
225 CALI_MARK_END("AoSoA commit full");
226#endif
227
228 // ===================================================
229 // Get Verlet pairs and fill Verlet list
230 // ===================================================
231 bool rebuild_vl = (integ_switch != INTEG_METHOD_STEEPEST_DESCENT and
232 cell_structure.use_verlet_list);
233#ifdef ESPRESSO_CALIPER
234 CALI_MARK_BEGIN("Verlet list creation");
235#endif
236 cell_structure.rebuild_verlet_list_cabana(
237 [&](std::span<Cell *const> cells, BoxGeometry const &box,
238 CellStructure::ListType &verlet_list) {
240 std::move(cells), box, verlet_criterion, id_to_index, max_id,
241 [&](const int i, const int j) {
242 // intra cell loop
243 verlet_list.addNeighborLB(i, j);
244 },
245 [&](const int i, const int j) {
246 // inter cell loop
247 verlet_list.addNeighbor(i, j);
248 });
249
250 if (verlet_list.hasOverflow()) {
251 cell_structure.use_verlet_list = false;
253 << "Verlet list overflow detected: neighbor count exceeded "
254 "max_counts. Falling back to the link cell algorithm. "
255 "Configured max is "
256 << Cabana::NeighborList<CellStructure::ListType>::maxNeighbor(
257 verlet_list);
258 }
259 },
260 rebuild_vl);
261#ifdef ESPRESSO_CALIPER
262 CALI_MARK_END("Verlet list creation");
263#endif
264 } else {
265 // ===================================================
266 // Fill particle storage (partial update)
267 // ===================================================
268#ifdef ESPRESSO_CALIPER
269 CALI_MARK_BEGIN("AoSoA commit partial");
270#endif
271 kokkos_parallel_range_for<policy_type>(
272 "AoSoA write", std::size_t{0}, n_part,
273 [&unique_particles, &aosoa](int const index) {
274 auto const &p = *unique_particles.at(index);
275 commit_particle(p, index, aosoa, false);
276 });
277 Kokkos::fence();
278#ifdef ESPRESSO_CALIPER
279 CALI_MARK_END("AoSoA commit partial");
280#endif
281 }
282}
283
284#ifdef ESPRESSO_ELECTROSTATICS
287 using execution_space = Kokkos::DefaultExecutionSpace;
288 using policy_type = Kokkos::RangePolicy<execution_space>;
289 auto const &unique_particles = cell_structure.get_unique_particles();
290 auto const n_part = unique_particles.size();
291 auto &aosoa = cell_structure.get_aosoa();
292
293 kokkos_parallel_range_for<policy_type>(
294 "Views update charges", std::size_t{0}, n_part,
295 [&unique_particles, &aosoa](std::size_t const index) {
296 aosoa.charge(index) = unique_particles.at(index)->q();
297 });
298}
299#endif
300
301void cabana_short_range(auto const &pair_bonds_kernel,
302 auto const &angle_bonds_kernel,
303 auto const &dihedral_bonds_kernel,
304 auto const &nonbonded_kernel,
305 CellStructure &cell_structure, double pair_cutoff,
306 double bond_cutoff, auto const &verlet_criterion,
307 auto const integ_switch) {
308 using execution_space = Kokkos::DefaultExecutionSpace;
309 assert(cell_structure.get_resort_particles() == Cells::RESORT_NONE);
310
311 if (bond_cutoff >= 0.) {
312#ifdef ESPRESSO_CALIPER
313 CALI_MARK_BEGIN("cabana_bond_loop");
314#endif
315 auto const n_pair_bonds = cell_structure.get_local_pair_bond_numbers();
316 auto const n_angle_bonds = cell_structure.get_local_angle_bond_numbers();
317 auto const n_dihedral_bonds =
318 cell_structure.get_local_dihedral_bond_numbers();
319 if (n_pair_bonds > 0) {
320 Kokkos::parallel_for( // loop over bonds
321 "for_each_local_pair_bonds", n_pair_bonds, pair_bonds_kernel);
322 Kokkos::fence();
323 }
324 if (n_angle_bonds > 0) {
325 Kokkos::parallel_for( // loop over bonds
326 "for_each_local_angle_bonds", n_angle_bonds, angle_bonds_kernel);
327 Kokkos::fence();
328 }
329 if (n_dihedral_bonds > 0) {
330 Kokkos::parallel_for( // loop over bonds
331 "for_each_local_dihedral_bonds", n_dihedral_bonds,
332 dihedral_bonds_kernel);
333 Kokkos::fence();
334 }
335#ifdef ESPRESSO_CALIPER
336 CALI_MARK_END("cabana_bond_loop");
337#endif
338 }
339
340 // Cabana short range loop
341 if (pair_cutoff > 0.) {
342#ifdef ESPRESSO_CALIPER
343 CALI_MARK_BEGIN("cabana_pair_loop");
344#endif
345 if (integ_switch != INTEG_METHOD_STEEPEST_DESCENT and
346 cell_structure.use_verlet_list) {
347 auto const &verlet_list = cell_structure.get_verlet_list_cabana();
348 Kokkos::RangePolicy<execution_space> policy(
349 std::size_t{0}, cell_structure.get_unique_particles().size());
350 Cabana::neighbor_parallel_for(policy, nonbonded_kernel, verlet_list,
351 Cabana::FirstNeighborsTag(),
352 Cabana::SerialOpTag());
353 } else {
354 cell_structure.cell_list_loop(
355 [&](std::span<Cell *const> cells, BoxGeometry const &box) {
357 std::move(cells), box, verlet_criterion,
358 cell_structure.get_id_to_index(),
359 cell_structure.get_cached_max_local_particle_id(),
360 [&](const int i, const int j) {
361 // intra cell loop
362 nonbonded_kernel(i, j);
363 },
364 [&](const int i, const int j) {
365 // inter cell loop
366 nonbonded_kernel(i, j);
367 });
368 });
369 }
370 Kokkos::fence();
371#ifdef ESPRESSO_CALIPER
372 CALI_MARK_END("cabana_pair_loop");
373#endif
374 }
375}
@ INTEG_METHOD_STEEPEST_DESCENT
#define ESPRESSO_ATTR_ALWAYS_INLINE
Describes a cell structure / cell system.
auto & get_id_to_index()
int get_local_angle_bond_numbers() const
int get_local_pair_bond_numbers() const
auto prepare_verlet_list_cabana(double cutoff)
Reset local properties of the Verlet list.
int get_local_dihedral_bond_numbers() const
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 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)
void cabana_short_range(auto const &pair_bonds_kernel, auto const &angle_bonds_kernel, auto const &dihedral_bonds_kernel, auto const &nonbonded_kernel, CellStructure &cell_structure, double pair_cutoff, double bond_cutoff, auto const &verlet_criterion, auto const integ_switch)
ESPRESSO_ATTR_ALWAYS_INLINE void commit_particle(Particle const &p, auto const index, CellStructure::AoSoA_pack &aosoa, bool const rebuild)
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< T *[N], array_layout, Kokkos::HostSpace > &view, std::size_t i, Utils::Vector< T, N > 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:435
auto const & mass() const
Definition Particle.hpp:492
Utils::compact_vector< int > & exclusions()
Definition Particle.hpp:646
auto const & quat() const
Definition Particle.hpp:517
auto const & q() const
Definition Particle.hpp:578
auto const & v() const
Definition Particle.hpp:473
auto const & image_box() const
Definition Particle.hpp:484
auto const & type() const
Definition Particle.hpp:458
auto const & dipm() const
Definition Particle.hpp:533
auto const & pos() const
Definition Particle.hpp:471
auto const & id() const
Definition Particle.hpp:454