ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
forces_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
24#ifdef ESPRESSO_SHARED_MEMORY_PARALLELISM
25
26#include "aosoa_pack.hpp"
27#include "forces_inline.hpp"
28
29#include <utils/Vector.hpp>
30
31#include <Cabana_Core.hpp>
32
33#include <omp.h>
34
35#include <cstddef>
36#include <memory>
37#include <optional>
38#include <variant>
39#include <vector>
40
41#if defined(__GNUG__) or defined(__clang__)
42#define ESPRESSO_ATTR_ALWAYS_INLINE [[gnu::always_inline]]
43#else
44#define ESPRESSO_ATTR_ALWAYS_INLINE
45#endif
46
56 std::vector<Particle *> const &unique_particles;
58#ifdef ESPRESSO_ROTATION
60#endif
61#ifdef ESPRESSO_NPT
64#endif
66#ifdef ESPRESSO_P3M
68#endif
70
72 BondedInteractionsMap const &bonded_ias_,
73 InteractionsNonBonded const &nonbonded_ias_,
77 Coulomb::ShortRangeEnergyKernel::kernel_type const *coulomb_u_kernel_,
78 Coulomb::Solver const &coulomb_,
79 Thermostat::Thermostat const &thermostat_, BoxGeometry const &box_geo_,
80 std::vector<Particle *> const &unique_particles_,
81 CellStructure::ForceType const &local_force_,
82#ifdef ESPRESSO_ROTATION
83 CellStructure::ForceType const &local_torque_,
84#endif
85#ifdef ESPRESSO_NPT
86 Utils::Vector3d *const global_virial_,
87 CellStructure::VirialType const &local_virial_,
88#endif
89 CellStructure::AoSoA_pack const &aosoa_, double system_max_cutoff_)
90 : bonded_ias(bonded_ias_), nonbonded_ias(nonbonded_ias_),
91 coulomb_kernel(coulomb_kernel_), dipoles_kernel(dipoles_kernel_),
92 elc_kernel(elc_kernel_), coulomb_u_kernel(coulomb_u_kernel_),
93 thermostat(thermostat_), box_geo(box_geo_),
94 unique_particles(unique_particles_), local_force(local_force_),
95#ifdef ESPRESSO_ROTATION
96 local_torque(local_torque_),
97#endif
98#ifdef ESPRESSO_NPT
99 global_virial(global_virial_), local_virial(local_virial_),
100#endif
101 aosoa(aosoa_), system_max_cutoff(system_max_cutoff_) {
102#ifdef ESPRESSO_P3M
103 p3m = nullptr;
104 if (auto &solver = coulomb_.impl->solver; solver.has_value()) {
105 if (std::holds_alternative<std::shared_ptr<CoulombP3M>>(*solver)) {
106 p3m = std::get<std::shared_ptr<CoulombP3M>>(*solver).get();
107 }
108 }
109#endif
110 }
111
112 // Helper functions to check if specific algorithms are active
113#ifdef ESPRESSO_GAY_BERNE
114 ESPRESSO_ATTR_ALWAYS_INLINE KOKKOS_INLINE_FUNCTION bool
115 gay_berne_active(double dist, IA_parameters const &ia_params) const {
116 return dist < ia_params.gay_berne.cut;
117 }
118#endif
119
120#ifdef ESPRESSO_NPT
121 ESPRESSO_ATTR_ALWAYS_INLINE KOKKOS_INLINE_FUNCTION bool npt_active() const {
122 return global_virial != nullptr;
123 }
124#endif
125
126#ifdef ESPRESSO_THOLE
127 ESPRESSO_ATTR_ALWAYS_INLINE KOKKOS_INLINE_FUNCTION bool
128 thole_active(IA_parameters const &ia_params) const {
129 return (ia_params.thole.scaling_coeff != 0. and
130 ia_params.thole.q1q2 != 0. and coulomb_kernel != nullptr);
131 }
132#endif
133
134#ifdef ESPRESSO_DIPOLES
135 ESPRESSO_ATTR_ALWAYS_INLINE KOKKOS_INLINE_FUNCTION bool
137 return dipoles_kernel != nullptr;
138 }
139#endif
140
141 ESPRESSO_ATTR_ALWAYS_INLINE KOKKOS_INLINE_FUNCTION void
142 operator()(std::size_t i, std::size_t j) const {
143
144 // calc distance
145 auto const pos1 = aosoa.get_vector_at(aosoa.position, i);
146 auto const pos2 = aosoa.get_vector_at(aosoa.position, j);
147 auto const d = box_geo.get_mi_vector(pos1, pos2);
148 auto const dist = d.norm();
149
150 // Early exit if distance > maximal clobal cutoff
151 if (dist > system_max_cutoff)
152 return;
153
154 auto const &ia_params =
156
157 ParticleForce pf{};
158
159 // Determine which data needs to be loaded based on active algorithms
160#if defined(ESPRESSO_DIPOLES) or defined(ESPRESSO_GAY_BERNE)
161 auto need_directors = false;
162#if defined(ESPRESSO_GAY_BERNE)
163 need_directors |= gay_berne_active(dist, ia_params);
164#endif
165#if defined(ESPRESSO_DIPOLES)
166 need_directors |= dipoles_active();
167#endif
168#endif
169#if defined(ESPRESSO_EXCLUSIONS) or defined(ESPRESSO_THOLE)
170 auto need_particle_pointers = false;
171#if defined(ESPRESSO_EXCLUSIONS)
172 need_particle_pointers |= aosoa.has_exclusion(i) or aosoa.has_exclusion(j);
173#endif
174#if defined(ESPRESSO_THOLE)
175 need_particle_pointers |= thole_active(ia_params);
176#endif
177 Particle const *p1_ptr = nullptr;
178 Particle const *p2_ptr = nullptr;
179 if (need_particle_pointers) {
180 p1_ptr = unique_particles.at(i);
181 p2_ptr = unique_particles.at(j);
182 }
183#endif
184
185 // Load directors only if needed
186#if defined(ESPRESSO_GAY_BERNE) or defined(ESPRESSO_DIPOLES)
187 Utils::Vector3d dir1{}, dir2{};
188 if (need_directors) {
191 }
192#endif
193
194 /***********************************************/
195 /* non-bonded pair potentials */
196 /***********************************************/
197
198 if (dist <= ia_params.max_cut) {
199#ifdef ESPRESSO_EXCLUSIONS
200 bool skip_non_bonded = false;
201 if (aosoa.has_exclusion(i) or aosoa.has_exclusion(j)) {
202 skip_non_bonded = not do_nonbonded(*p1_ptr, *p2_ptr);
203 }
204#else
205 constexpr bool skip_non_bonded = false;
206#endif
207 if (not skip_non_bonded) {
208 pf.f += calc_central_radial_force(ia_params, d, dist);
209
210 // Only call Thole force kernel if active
211#ifdef ESPRESSO_THOLE
212 if (thole_active(ia_params)) {
213 pf.f += thole_pair_force(*p1_ptr, *p2_ptr, ia_params, d, dist,
215 }
216#endif
217 // Only call Gay-Berne force kernel if active
218#ifdef ESPRESSO_GAY_BERNE
219 if (gay_berne_active(dist, ia_params)) {
220 pf += gb_pair_force(dir1, dir2, ia_params, d, dist);
221 }
222#endif
223 } // not skip_non_bonded
224 } // not dist > ia_params.max_cut
225
226 /*********************************************************************/
227 /* everything before this contributes to the virial pressure in NpT, */
228 /* but nothing afterwards, since the contribution to pressure from */
229 /* electrostatic is calculated by energy */
230 /*********************************************************************/
231#ifdef ESPRESSO_NPT
232 Utils::Vector3d virial{};
233 if (npt_active()) {
234 virial = hadamard_product(pf.f, d);
235 }
236#endif // ESPRESSO_NPT
237
238 /***********************************************/
239 /* thermostat */
240 /***********************************************/
241
242 /* The inter dpd force should not be part of the virial */
243#ifdef ESPRESSO_DPD
245 auto const dist2 = dist * dist;
246 auto const vel1 = aosoa.get_vector_at(aosoa.velocity, i);
247 auto const vel2 = aosoa.get_vector_at(aosoa.velocity, j);
248 auto const force =
249 dpd_pair_force(pos1, vel1, aosoa.id(i), pos2, vel2, aosoa.id(j),
250 *thermostat.dpd, box_geo, ia_params, d, dist, dist2);
251 pf += force;
252 }
253#endif // ESPRESSO_DPD
254
255#ifdef ESPRESSO_ELECTROSTATICS
256 Utils::Vector3d f1_asym{};
257 Utils::Vector3d f2_asym{};
258 // real-space electrostatic charge-charge interaction
259 if (coulomb_kernel != nullptr) {
260 if ((aosoa.charge(i) != 0.) and (aosoa.charge(j) != 0.)) {
261 auto const q1q2 = aosoa.charge(i) * aosoa.charge(j);
262#ifdef ESPRESSO_P3M
263 if (p3m) [[likely]] {
264 pf.f += p3m->pair_force(q1q2, d, dist);
265 } else
266#endif
267 {
268 pf.f += (*coulomb_kernel)(q1q2, d, dist);
269 }
270 if (elc_kernel) {
271 (*elc_kernel)(pos1, pos2, f1_asym, f2_asym, q1q2);
272 }
273#ifdef ESPRESSO_NPT
274 if (npt_active()) {
275 virial[0] += (*coulomb_u_kernel)(pos1, pos2, q1q2, d, dist);
276 }
277#endif // ESPRESSO_NPT
278 }
279 }
280#endif // ESPRESSO_ELECTROSTATICS
281
282 // Only call dipole force kernel if active
283#ifdef ESPRESSO_DIPOLES
284 if (dipoles_active()) {
285 auto const d1d2 = aosoa.dipm(i) * aosoa.dipm(j);
286 if (d1d2 != 0.) {
287 pf += (*dipoles_kernel)(d1d2, aosoa.dipm(i) * dir1,
288 aosoa.dipm(j) * dir2, d, dist, dist * dist);
289 }
290 }
291#endif // ESPRESSO_DIPOLES
292
293 auto opf = calc_opposing_force(pf, d);
294#ifdef ESPRESSO_ELECTROSTATICS
295 pf.f += f1_asym;
296 opf.f += f2_asym;
297#endif // ESPRESSO_ELECTROSTATICS
298
299 auto const thread_id = omp_get_thread_num();
300
301 local_force(i, thread_id, 0) += pf.f[0];
302 local_force(i, thread_id, 1) += pf.f[1];
303 local_force(i, thread_id, 2) += pf.f[2];
304#ifdef ESPRESSO_ROTATION
305 local_torque(i, thread_id, 0) += pf.torque[0];
306 local_torque(i, thread_id, 1) += pf.torque[1];
307 local_torque(i, thread_id, 2) += pf.torque[2];
308#endif
309
310 local_force(j, thread_id, 0) += opf.f[0];
311 local_force(j, thread_id, 1) += opf.f[1];
312 local_force(j, thread_id, 2) += opf.f[2];
313#ifdef ESPRESSO_ROTATION
314 local_torque(j, thread_id, 0) += opf.torque[0];
315 local_torque(j, thread_id, 1) += opf.torque[1];
316 local_torque(j, thread_id, 2) += opf.torque[2];
317#endif
318#ifdef ESPRESSO_NPT
319 if (npt_active()) {
320 local_virial(thread_id, 0) += virial[0];
321 local_virial(thread_id, 1) += virial[1];
322 local_virial(thread_id, 2) += virial[2];
323 }
324#endif
325 }
326};
327
328#endif // ESPRESSO_SHARED_MEMORY_PARALLELISM
@ THERMO_DPD
Vector implementation and trait types for boost qvm interoperability.
container for bonded interactions.
ESPRESSO_ATTR_ALWAYS_INLINE Utils::Vector3< T > get_mi_vector(Utils::Vector3< T > const &a, Utils::Vector3< T > const &b) const
Get the minimum-image vector between two coordinates.
auto & get_ia_param(int i, int j)
Get interaction parameters between particle types i and j.
std::shared_ptr< DPDThermostat > dpd
int thermo_switch
Bitmask of currently active thermostats.
Utils::Vector3d dpd_pair_force(Utils::Vector3d const &p1_position, Utils::Vector3d const &p1_velocity, int const &p1_id, Utils::Vector3d const &p2_position, Utils::Vector3d const &p2_velocity, int const &p2_id, DPDThermostat const &dpd, BoxGeometry const &box_geo, IA_parameters const &ia_params, Utils::Vector3d const &d, double dist, double dist2)
Definition dpd.cpp:76
bool do_nonbonded(Particle const &p1, Particle const &p2)
Determine if the non-bonded interactions between p1 and p2 should be calculated.
#define ESPRESSO_ATTR_ALWAYS_INLINE
Force calculation.
ParticleForce calc_opposing_force(ParticleForce const &pf, Utils::Vector3d const &d)
ESPRESSO_ATTR_ALWAYS_INLINE Utils::Vector3d calc_central_radial_force(IA_parameters const &ia_params, Utils::Vector3d const &d, double const dist)
ParticleForce gb_pair_force(Utils::Vector3d const &ui, Utils::Vector3d const &uj, IA_parameters const &ia_params, Utils::Vector3d const &d, double dist)
Calculate Gay-Berne force and torques.
Definition gay_berne.hpp:49
Utils::Vector3d get_vector_at(Kokkos::View< double *[3], array_layout, Kokkos::HostSpace > const &view, std::size_t i) const
bool has_exclusion(std::size_t i) const
PositionViewType position
DirectorViewType director
VelocityViewType velocity
P3M solver.
Definition p3m.hpp:55
Utils::Vector3d pair_force(double q1q2, Utils::Vector3d const &d, double dist) const
Calculate real-space contribution of p3m Coulomb pair forces.
Definition p3m.hpp:145
Solver::ShortRangeEnergyKernel kernel_type
Solver::ShortRangeForceCorrectionsKernel kernel_type
Solver::ShortRangeForceKernel kernel_type
std::unique_ptr< Implementation > impl
Pointer-to-implementation.
Solver::ShortRangeForceKernel kernel_type
ESPRESSO_ATTR_ALWAYS_INLINE KOKKOS_INLINE_FUNCTION bool thole_active(IA_parameters const &ia_params) const
std::vector< Particle * > const & unique_particles
ESPRESSO_ATTR_ALWAYS_INLINE KOKKOS_INLINE_FUNCTION bool gay_berne_active(double dist, IA_parameters const &ia_params) const
Coulomb::ShortRangeEnergyKernel::kernel_type const *const coulomb_u_kernel
Coulomb::ShortRangeForceCorrectionsKernel::kernel_type const * elc_kernel
ESPRESSO_ATTR_ALWAYS_INLINE KOKKOS_INLINE_FUNCTION bool npt_active() const
BoxGeometry const & box_geo
CellStructure::VirialType const & local_virial
ForcesKernel(BondedInteractionsMap const &bonded_ias_, InteractionsNonBonded const &nonbonded_ias_, Coulomb::ShortRangeForceKernel::kernel_type const *coulomb_kernel_, Dipoles::ShortRangeForceKernel::kernel_type const *dipoles_kernel_, Coulomb::ShortRangeForceCorrectionsKernel::kernel_type const *elc_kernel_, Coulomb::ShortRangeEnergyKernel::kernel_type const *coulomb_u_kernel_, Coulomb::Solver const &coulomb_, Thermostat::Thermostat const &thermostat_, BoxGeometry const &box_geo_, std::vector< Particle * > const &unique_particles_, CellStructure::ForceType const &local_force_, CellStructure::ForceType const &local_torque_, Utils::Vector3d *const global_virial_, CellStructure::VirialType const &local_virial_, CellStructure::AoSoA_pack const &aosoa_, double system_max_cutoff_)
ESPRESSO_ATTR_ALWAYS_INLINE KOKKOS_INLINE_FUNCTION bool dipoles_active() const
Dipoles::ShortRangeForceKernel::kernel_type const *const dipoles_kernel
CellStructure::AoSoA_pack const & aosoa
CellStructure::ForceType const & local_torque
double system_max_cutoff
InteractionsNonBonded const & nonbonded_ias
Utils::Vector3d *const global_virial
ESPRESSO_ATTR_ALWAYS_INLINE KOKKOS_INLINE_FUNCTION void operator()(std::size_t i, std::size_t j) const
BondedInteractionsMap const & bonded_ias
CellStructure::ForceType const & local_force
CoulombP3M const * p3m
Coulomb::ShortRangeForceKernel::kernel_type const *const coulomb_kernel
Thermostat::Thermostat const & thermostat
Parameters for non-bonded interactions.
GayBerne_Parameters gay_berne
Force information on a particle.
Definition Particle.hpp:330
Struct holding all information for one particle.
Definition Particle.hpp:435
Utils::Vector3d thole_pair_force(Particle const &p1, Particle const &p2, IA_parameters const &ia_params, Utils::Vector3d const &d, double dist, BondedInteractionsMap const &bonded_ias, Coulomb::ShortRangeForceKernel::kernel_type const *kernel)
Calculate Thole force.
Definition thole.hpp:45