ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
dipolar_direct_sum_gpu.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2010-2022 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#include "config/config.hpp"
21
22#ifdef ESPRESSO_DIPOLAR_DIRECT_SUM
23
26
27#include "BoxGeometry.hpp"
28#include "communication.hpp"
30#include "system/System.hpp"
31
32static void get_simulation_box(BoxGeometry const &box_geo, float *box,
33 int *per) {
34 for (unsigned int i = 0u; i < 3u; i++) {
35 box[i] = static_cast<float>(box_geo.length()[i]);
36 per[i] = box_geo.periodic(i);
37 }
38}
39
40DipolarDirectSumGpu::DipolarDirectSumGpu(double prefactor, int n_replicas) {
42 this->n_replicas = n_replicas;
43 if (n_replicas < 0) {
44 throw std::domain_error("Parameter 'n_replicas' must be >= 0");
45 }
46}
47
49 auto &gpu_particle_data = get_system().gpu;
50 gpu_particle_data.enable_property(GpuParticleData::prop::force);
51 gpu_particle_data.enable_property(GpuParticleData::prop::torque);
52 gpu_particle_data.enable_property(GpuParticleData::prop::pos);
53 gpu_particle_data.enable_property(GpuParticleData::prop::dip);
54#ifdef ESPRESSO_DIPOLE_FIELD_TRACKING
55 gpu_particle_data.enable_property(GpuParticleData::prop::dip_fld);
56#endif
57}
58
60 auto &system = get_system();
61 auto &gpu = system.gpu;
62 gpu.update();
63 if (this_node != 0) {
64 return;
65 }
66 float box[3];
67 int periodicity[3];
68 get_simulation_box(*system.box_geo, box, periodicity);
69 auto const npart = static_cast<unsigned>(gpu.n_particles());
70 auto const forces_device = gpu.get_particle_forces_device();
71 auto const torques_device = gpu.get_particle_torques_device();
72 auto const positions_device = gpu.get_particle_positions_device();
73 auto const dipoles_device = gpu.get_particle_dipoles_device();
74 float *dipole_fields_device{nullptr};
75#ifdef ESPRESSO_DIPOLE_FIELD_TRACKING
76 dipole_fields_device = gpu.get_particle_dip_fld_device();
77#endif
79 static_cast<float>(prefactor), npart, positions_device, dipoles_device,
80 dipole_fields_device, forces_device, torques_device, box, periodicity,
82}
83
85 auto &system = get_system();
86 auto &gpu = system.gpu;
87 gpu.update();
88 if (this_node != 0) {
89 return;
90 }
91 float box[3];
92 int periodicity[3];
93 get_simulation_box(*system.box_geo, box, periodicity);
94 auto const npart = static_cast<unsigned>(gpu.n_particles());
95 auto const energy_device = &(gpu.get_energy_device()->dipolar);
96 auto const positions_device = gpu.get_particle_positions_device();
97 auto const dipoles_device = gpu.get_particle_dipoles_device();
98 DipolarDirectSum_kernel_wrapper_energy(static_cast<float>(prefactor), npart,
99 positions_device, dipoles_device, box,
100 periodicity, energy_device);
101}
102
103#endif // ESPRESSO_DIPOLAR_DIRECT_SUM
Utils::Vector3d const & length() const
Box length.
constexpr bool periodic(unsigned coord) const
Check periodicity in direction.
double prefactor
Magnetostatics prefactor.
int this_node
The number of this node.
static void get_simulation_box(BoxGeometry const &box_geo, float *box, int *per)
void DipolarDirectSum_kernel_wrapper_energy(float k, unsigned int n, float const *const pos, float const *const dip, float box_l[3], int periodic[3], float *E)
void DipolarDirectSum_kernel_wrapper_force(float k, unsigned int n, float const *const pos, float const *const dip, float *dip_fld, float *f, float *torque, float box_l[3], int periodic[3], int n_replicas)
DipolarDirectSumGpu(double prefactor, int n_replicas)
static constexpr std::size_t force
static constexpr std::size_t dip_fld
static constexpr std::size_t torque
static constexpr std::size_t pos
static constexpr std::size_t dip