ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
energy_cabana.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 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#include "aosoa_pack.hpp"
25#include "energy_inline.hpp"
27
28#include <utils/Vector.hpp>
29
30#include <Cabana_Core.hpp>
31
32#include <omp.h>
33
34#include <cstddef>
35#include <memory>
36#include <optional>
37#include <variant>
38#include <vector>
39
41 std::size_t n_bonded; // initialized by bonded_ias->get_next_key()
42 std::size_t n_types; // max_seen_particle_type
43 std::size_t off_bonded = 0; // [0, n_bonded)
44 std::size_t off_nb_inter;
45 std::size_t off_nb_intra;
46 std::size_t off_coulomb;
47 std::size_t off_dipolar;
48 std::size_t total;
49
50 EnergyBinLayout(std::size_t n_bonded_, std::size_t n_types_)
51 : n_bonded(n_bonded_), n_types(n_types_) {
52 auto const n_nb = n_types * (n_types + 1) / 2;
57 total = off_dipolar + 1;
58 }
59
60 KOKKOS_INLINE_FUNCTION
61 std::size_t nb_inter_idx(int t1, int t2) const {
62 return off_nb_inter +
63 Utils::lower_triangular(std::max(t1, t2), std::min(t1, t2));
64 }
65
66 KOKKOS_INLINE_FUNCTION
67 std::size_t nb_intra_idx(int t1, int t2) const {
68 return off_nb_intra +
69 Utils::lower_triangular(std::max(t1, t2), std::min(t1, t2));
70 }
71
72 KOKKOS_INLINE_FUNCTION std::size_t dipolar_idx() const { return off_dipolar; }
73 KOKKOS_INLINE_FUNCTION std::size_t coulomb_idx() const { return off_coulomb; }
74 KOKKOS_INLINE_FUNCTION std::size_t bonded_idx(int b) const {
75 return off_bonded + b;
76 }
77};
78
86 std::vector<Particle *> const &unique_particles;
92
94 BondedInteractionsMap const &bonded_ias_,
95 InteractionsNonBonded const &nonbonded_ias_,
96 Coulomb::Solver const &coulomb_,
97 Coulomb::ShortRangeEnergyKernel::kernel_type const *coulomb_u_kernel_,
98 Dipoles::ShortRangeEnergyKernel::kernel_type const *dipoles_u_kernel_,
99 BoxGeometry const &box_geo_,
100 std::vector<Particle *> const &unique_particles_,
102 EnergyBinLayout layout_, CellStructure::AoSoA_pack const &aosoa_,
103 Kokkos::View<int *> mol_id_view_, double system_max_cutoff_)
104 : bonded_ias(bonded_ias_), nonbonded_ias(nonbonded_ias_),
105 coulomb(coulomb_), coulomb_u_kernel(coulomb_u_kernel_),
106 dipoles_u_kernel(dipoles_u_kernel_), box_geo(box_geo_),
107 unique_particles(unique_particles_), local_energy(local_energy_),
108 layout(layout_), aosoa(aosoa_), mol_id_view(std::move(mol_id_view_)),
109 system_max_cutoff(system_max_cutoff_) {}
110
111 KOKKOS_INLINE_FUNCTION
112 void operator()(std::size_t i, std::size_t j) const {
113 auto const pos1 = aosoa.get_vector_at(aosoa.position, i);
114 auto const pos2 = aosoa.get_vector_at(aosoa.position, j);
115 auto const d = box_geo.get_mi_vector(pos1, pos2);
116 auto const dist = d.norm();
117 if (dist > system_max_cutoff)
118 return;
119
120 auto const t1 = aosoa.type(i);
121 auto const t2 = aosoa.type(j);
122 auto const &ia_params = nonbonded_ias.get_ia_param(t1, t2);
123
124 // Determine which data needs to be loaded based on active algorithms
125#if defined(ESPRESSO_GAY_BERNE) or defined(ESPRESSO_DIPOLES) or \
126 defined(ESPRESSO_EXCLUSIONS) or defined(ESPRESSO_THOLE)
127 auto const flag =
128 compute_pair_data_flags(dist, ia_params, coulomb_u_kernel != nullptr,
129 dipoles_u_kernel != nullptr, aosoa, i, j);
130#endif
131
132#if defined(ESPRESSO_EXCLUSIONS) or defined(ESPRESSO_THOLE)
133 Particle const *p1_ptr = nullptr;
134 Particle const *p2_ptr = nullptr;
135 if (flag.need_particle_pointers) {
136 p1_ptr = unique_particles.at(i);
137 p2_ptr = unique_particles.at(j);
138 }
139#endif
140
141 // Load directors only if needed
142#if defined(ESPRESSO_GAY_BERNE) or defined(ESPRESSO_DIPOLES)
143 Utils::Vector3d dir1{}, dir2{};
144 if (flag.need_directors) {
147 }
148#endif
149
150 auto const tid = omp_get_thread_num();
151 double e_nb = 0.0;
152
153 if (dist <= ia_params.max_cut) {
154#ifdef ESPRESSO_EXCLUSIONS
155 bool skip = false;
157 skip = not do_nonbonded(*p1_ptr, *p2_ptr);
158 if (not skip)
159#endif
160 {
161 e_nb += calc_central_radial_energy(ia_params, dist);
162
163 // Only call Thole energy kernel if active
164#ifdef ESPRESSO_THOLE
165 if (thole_active(ia_params, coulomb_u_kernel != nullptr)) {
166 e_nb += thole_pair_energy(*p1_ptr, *p2_ptr, ia_params, d, dist,
168 }
169#endif
170 // Only call Gay-Berne energy kernel if active
171#ifdef ESPRESSO_GAY_BERNE
172 if (gay_berne_active(dist, ia_params)) {
173 e_nb += gb_pair_energy(dir1, dir2, ia_params, d, dist);
174 }
175#endif
176 }
177 }
178 // pick inter vs intra bin like Observable_stat::non_bonded_contribution
179 // does
180 auto const bin = (mol_id_view(i) == mol_id_view(j))
181 ? layout.nb_intra_idx(t1, t2)
182 : layout.nb_inter_idx(t1, t2);
183 local_energy(tid, bin) += e_nb;
184
185#ifdef ESPRESSO_ELECTROSTATICS
186 if (coulomb_u_kernel != nullptr) {
187 auto const q1 = aosoa.charge(i), q2 = aosoa.charge(j);
188 if (q1 != 0. and q2 != 0.) {
189 double const e_c = (*coulomb_u_kernel)(pos1, pos2, q1 * q2, d, dist);
190 local_energy(tid, layout.coulomb_idx()) += e_c;
191 }
192 }
193#endif
194
195#ifdef ESPRESSO_DIPOLES
196 if (dipoles_u_kernel != nullptr) {
197 if (aosoa.dipm(i) != 0. and aosoa.dipm(j) != 0.) {
198 double const e_d = (*dipoles_u_kernel)(
199 aosoa.dipm(i) * dir1, aosoa.dipm(j) * dir2, d, dist, dist * dist);
200 local_energy(tid, layout.dipolar_idx()) += e_d;
201 }
202 }
203#endif
204 }
205};
206
209 EnergyBinLayout const &layout, Observable_stat &obs,
210 BondedInteractionsMap const &bonded_ias, int n_types) {
211 auto const nthreads = local_energy.extent(0);
212 auto host =
213 Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace{}, local_energy);
214
215 auto sum_bin = [&](std::size_t bin) {
216 double s = 0.;
217 for (std::size_t t = 0; t < nthreads; ++t)
218 s += host(t, bin);
219 return s;
220 };
221
222 for (int b = 0; b < int(layout.n_bonded); ++b)
223 obs.bonded_contribution(b)[0] += sum_bin(layout.bonded_idx(b));
224
225 for (int t1 = 0; t1 < n_types; ++t1)
226 for (int t2 = 0; t2 <= t1; ++t2)
227 obs.non_bonded_inter_contribution(t1, t2)[0] +=
228 sum_bin(layout.nb_inter_idx(t1, t2));
229
230 for (int t1 = 0; t1 < n_types; ++t1)
231 for (int t2 = 0; t2 <= t1; ++t2)
232 obs.non_bonded_intra_contribution(t1, t2)[0] +=
233 sum_bin(layout.nb_intra_idx(t1, t2));
234
235 obs.coulomb[0] += sum_bin(layout.coulomb_idx());
236 obs.dipolar[0] += sum_bin(layout.dipolar_idx());
237}
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.
Observable for the pressure and energy.
auto non_bonded_intra_contribution(int type1, int type2) const
Get contribution from a non-bonded intramolecular interaction.
std::span< double > coulomb
Contribution(s) from Coulomb interactions.
auto non_bonded_inter_contribution(int type1, int type2) const
Get contribution from a non-bonded intermolecular interaction.
std::span< double > bonded_contribution(int bond_id) const
Get contribution from a bonded interaction.
std::span< double > dipolar
Contribution(s) from dipolar interactions.
static void reduce_cabana_energy(Kokkos::View< double **, Kokkos::LayoutRight > const &local_energy, EnergyBinLayout const &layout, Observable_stat &obs, BondedInteractionsMap const &bonded_ias, int n_types)
Energy calculation.
double calc_central_radial_energy(IA_parameters const &ia_params, double const dist)
bool do_nonbonded(Particle const &p1, Particle const &p2)
Determine if the non-bonded interactions between p1 and p2 should be calculated.
double gb_pair_energy(Utils::Vector3d const &ui, Utils::Vector3d const &uj, IA_parameters const &ia_params, Utils::Vector3d const &d, double dist)
Calculate Gay-Berne energy.
T lower_triangular(T i, T j)
Linear index into a lower triangular matrix.
Definition index.hpp:83
STL namespace.
KOKKOS_INLINE_FUNCTION PairDataFlags compute_pair_data_flags(double dist, IA_parameters const &ia_params, bool has_coulomb, bool has_dipoles, auto const &aosoa, std::size_t i, std::size_t j)
KOKKOS_INLINE_FUNCTION bool thole_active(IA_parameters const &ia_params, bool has_coulomb_kernel)
KOKKOS_INLINE_FUNCTION bool gay_berne_active(double dist, IA_parameters const &ia_params)
bool has_exclusion(std::size_t i) const
PositionViewType position
Utils::Vector< T, N > get_vector_at(Kokkos::View< T *[N], array_layout, Kokkos::HostSpace > const &view, std::size_t i) const
DirectorViewType director
Solver::ShortRangeEnergyKernel kernel_type
Solver::ShortRangeEnergyKernel kernel_type
std::size_t off_nb_inter
KOKKOS_INLINE_FUNCTION std::size_t nb_intra_idx(int t1, int t2) const
std::size_t n_types
EnergyBinLayout(std::size_t n_bonded_, std::size_t n_types_)
std::size_t off_bonded
std::size_t off_dipolar
KOKKOS_INLINE_FUNCTION std::size_t dipolar_idx() const
std::size_t off_nb_intra
KOKKOS_INLINE_FUNCTION std::size_t bonded_idx(int b) const
KOKKOS_INLINE_FUNCTION std::size_t coulomb_idx() const
std::size_t off_coulomb
std::size_t total
std::size_t n_bonded
KOKKOS_INLINE_FUNCTION std::size_t nb_inter_idx(int t1, int t2) const
CellStructure::AoSoA_pack const & aosoa
std::vector< Particle * > const & unique_particles
BondedInteractionsMap const & bonded_ias
Kokkos::View< double **, Kokkos::LayoutRight > local_energy
EnergyKernel(BondedInteractionsMap const &bonded_ias_, InteractionsNonBonded const &nonbonded_ias_, Coulomb::Solver const &coulomb_, Coulomb::ShortRangeEnergyKernel::kernel_type const *coulomb_u_kernel_, Dipoles::ShortRangeEnergyKernel::kernel_type const *dipoles_u_kernel_, BoxGeometry const &box_geo_, std::vector< Particle * > const &unique_particles_, Kokkos::View< double **, Kokkos::LayoutRight > const &local_energy_, EnergyBinLayout layout_, CellStructure::AoSoA_pack const &aosoa_, Kokkos::View< int * > mol_id_view_, double system_max_cutoff_)
InteractionsNonBonded const & nonbonded_ias
KOKKOS_INLINE_FUNCTION void operator()(std::size_t i, std::size_t j) const
EnergyBinLayout layout
Coulomb::Solver const & coulomb
Dipoles::ShortRangeEnergyKernel::kernel_type const * dipoles_u_kernel
Kokkos::View< int * > mol_id_view
BoxGeometry const & box_geo
Coulomb::ShortRangeEnergyKernel::kernel_type const * coulomb_u_kernel
double system_max_cutoff
Struct holding all information for one particle.
Definition Particle.hpp:435
double thole_pair_energy(Particle const &p1, Particle const &p2, IA_parameters const &ia_params, Utils::Vector3d const &d, double dist, BondedInteractionsMap const &bonded_ias, Coulomb::Solver const &coulomb, Coulomb::ShortRangeEnergyKernel::kernel_type const *kernel)
Calculate Thole energy.
Definition thole.hpp:68