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
80 std::vector<Particle *> const &unique_particles_,
84#endif
88#endif
97#endif
100#endif
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
115#ifdef ESPRESSO_GAY_BERNE
116 return dist < ia_params.gay_berne.cut;
117#else
118 return false;
119#endif
120 }
121
122#ifdef ESPRESSO_NPT
126#endif
127
128#ifdef ESPRESSO_THOLE
131 return (ia_params.thole.scaling_coeff != 0. and
132 ia_params.thole.q1q2 != 0. and coulomb_kernel != nullptr);
133 }
134#endif
135
136#ifdef ESPRESSO_DIPOLES
139 return dipoles_kernel != nullptr;
140 }
141#endif
142
144 operator()(std::size_t i, std::size_t j) const {
145
146 // calc distance
147 auto const pos1 = aosoa.get_vector_at(aosoa.position, i);
148 auto const pos2 = aosoa.get_vector_at(aosoa.position, j);
149 auto const d = box_geo.get_mi_vector(pos1, pos2);
150 auto const dist = d.norm();
151
152 // Early exit if distance > maximal clobal cutoff
154 return;
155
156 auto const &ia_params =
158
160
161 // Determine which data needs to be loaded based on active algorithms
162#if defined(ESPRESSO_DIPOLES) or defined(ESPRESSO_GAY_BERNE)
163 bool const need_directors =
165#endif
166#if defined(ESPRESSO_EXCLUSIONS) or defined(ESPRESSO_THOLE)
170 Particle const *p1_ptr = nullptr;
171 Particle const *p2_ptr = nullptr;
173 p1_ptr = unique_particles.at(i);
175 }
176#endif
177
178 // Load directors only if needed
179#if defined(ESPRESSO_GAY_BERNE) or defined(ESPRESSO_DIPOLES)
181 if (need_directors) {
184 }
185#endif
186
187 /***********************************************/
188 /* non-bonded pair potentials */
189 /***********************************************/
190
191 if (dist <= ia_params.max_cut) {
192#ifdef ESPRESSO_EXCLUSIONS
193 bool skip_non_bonded = false;
196 }
197#else
198 constexpr bool skip_non_bonded = false;
199#endif
200 if (not skip_non_bonded) {
202
203 // Only call Thole force kernel if active
204#ifdef ESPRESSO_THOLE
205 if (thole_active(ia_params)) {
208 }
209#endif
210 // Only call Gay-Berne force kernel if active
211#ifdef ESPRESSO_GAY_BERNE
214 }
215#endif
216 } // not skip_non_bonded
217 } // not dist > ia_params.max_cut
218
219 /*********************************************************************/
220 /* everything before this contributes to the virial pressure in NpT, */
221 /* but nothing afterwards, since the contribution to pressure from */
222 /* electrostatic is calculated by energy */
223 /*********************************************************************/
224#ifdef ESPRESSO_NPT
226 if (npt_active()) {
227 virial = hadamard_product(pf.f, d);
228 }
229#endif // ESPRESSO_NPT
230
231 /***********************************************/
232 /* thermostat */
233 /***********************************************/
234
235 /* The inter dpd force should not be part of the virial */
236#ifdef ESPRESSO_DPD
238 auto const dist2 = dist * dist;
239 auto const vel1 = aosoa.get_vector_at(aosoa.velocity, i);
240 auto const vel2 = aosoa.get_vector_at(aosoa.velocity, j);
241 auto const force =
242 dpd_pair_force(pos1, vel1, aosoa.id(i), pos2, vel2, aosoa.id(j),
243 *thermostat.dpd, box_geo, ia_params, d, dist, dist2);
244 pf += force;
245 }
246#endif // ESPRESSO_DPD
247
248#ifdef ESPRESSO_ELECTROSTATICS
251 // real-space electrostatic charge-charge interaction
252 if (coulomb_kernel != nullptr) {
253 if ((aosoa.charge(i) != 0.) and (aosoa.charge(j) != 0.)) {
254 auto const q1q2 = aosoa.charge(i) * aosoa.charge(j);
255#ifdef ESPRESSO_P3M
256 if (p3m) [[likely]] {
257 pf.f += p3m->pair_force(q1q2, d, dist);
258 } else
259#endif
260 {
261 pf.f += (*coulomb_kernel)(q1q2, d, dist);
262 }
263 if (elc_kernel) {
264 (*elc_kernel)(pos1, pos2, f1_asym, f2_asym, q1q2);
265 }
266#ifdef ESPRESSO_NPT
267 if (npt_active()) {
268 virial[0] += (*coulomb_u_kernel)(pos1, pos2, q1q2, d, dist);
269 }
270#endif // ESPRESSO_NPT
271 }
272 }
273#endif // ESPRESSO_ELECTROSTATICS
274
275 // Only call dipole force kernel if active
276#ifdef ESPRESSO_DIPOLES
277 if (dipoles_active()) {
278 auto const d1d2 = aosoa.dipm(i) * aosoa.dipm(j);
279 if (d1d2 != 0.) {
280 pf += (*dipoles_kernel)(d1d2, aosoa.dipm(i) * dir1,
281 aosoa.dipm(j) * dir2, d, dist, dist * dist);
282 }
283 }
284#endif // ESPRESSO_DIPOLES
285
286 auto opf = calc_opposing_force(pf, d);
287#ifdef ESPRESSO_ELECTROSTATICS
288 pf.f += f1_asym;
289 opf.f += f2_asym;
290#endif // ESPRESSO_ELECTROSTATICS
291
292 auto const thread_id = omp_get_thread_num();
293
294 local_force(i, thread_id, 0) += pf.f[0];
295 local_force(i, thread_id, 1) += pf.f[1];
296 local_force(i, thread_id, 2) += pf.f[2];
297#ifdef ESPRESSO_ROTATION
298 local_torque(i, thread_id, 0) += pf.torque[0];
299 local_torque(i, thread_id, 1) += pf.torque[1];
300 local_torque(i, thread_id, 2) += pf.torque[2];
301#endif
302
303 local_force(j, thread_id, 0) += opf.f[0];
304 local_force(j, thread_id, 1) += opf.f[1];
305 local_force(j, thread_id, 2) += opf.f[2];
306#ifdef ESPRESSO_ROTATION
307 local_torque(j, thread_id, 0) += opf.torque[0];
308 local_torque(j, thread_id, 1) += opf.torque[1];
309 local_torque(j, thread_id, 2) += opf.torque[2];
310#endif
311#ifdef ESPRESSO_NPT
312 if (npt_active()) {
313 local_virial(thread_id, 0) += virial[0];
314 local_virial(thread_id, 1) += virial[1];
315 local_virial(thread_id, 2) += virial[2];
316 }
317#endif
318 }
319};
320
321#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.
cudaStream_t stream[1]
CUDA streams for parallel computing on CPU and GPU.
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
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.
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