Penning Trap Simulation
Simulate particle behavior inside a Penning Trap
Loading...
Searching...
No Matches
PenningTrap.cpp
Go to the documentation of this file.
1
13#include "PenningTrap.hpp"
14#include "typedefs.hpp"
15#include <algorithm>
16#include <functional>
17#include <sys/types.h>
18#include <vector>
19
20vec3 PenningTrap::v_func(uint i, uint j, double dt)
21{
22 switch (i) {
23 case 0:
24 return .5 * dt * this->k_v[0][j];
25 case 1:
26 return .5 * dt * this->k_v[1][j];
27 case 2:
28 return dt * this->k_v[2][j];
29 case 3:
30 return vec3((dt / 6.)
31 * (this->k_v[0][j] + 2. * this->k_v[1][j]
32 + 2. * this->k_v[2][j] + this->k_v[3][j]));
33 default:
34 std::cout << "Not valid!" << std::endl;
35 abort();
36 }
37}
38
39vec3 PenningTrap::r_func(uint i, uint j, double dt)
40{
41 switch (i) {
42 case 0:
43 return .5 * dt * this->k_r[0][j];
44 case 1:
45 return .5 * dt * this->k_r[1][j];
46 case 2:
47 return dt * this->k_r[2][j];
48 case 3:
49 return vec3((dt / 6.)
50 * (this->k_r[0][j] + 2. * this->k_r[1][j]
51 + 2. * this->k_r[2][j] + this->k_r[3][j]));
52 default:
53 std::cout << "Not valid!" << std::endl;
54 abort();
55 }
56}
57
59{
60 r(2) *= -2.;
61
62 return vec3((this->V_0 * this->perturbation(this->t) / (this->d * this->d))
63 * r);
64}
65
67{
68 return vec3{0., 0., this->B_0};
69}
70
72{
73 // Calculate the difference between the particles' position
74 vec3 res = this->particles[i].r_vec - this->particles[j].r_vec;
75
76 // Get the distance between the particles
77 double norm = arma::norm(res, 2);
78
79 return vec3((this->particles[j].q / (norm * norm * norm)) * res);
80}
81
83{
84 Particle *p = &this->particles[i];
85
86 if (arma::norm(p->r_vec) > this->d) {
87 return vec3{0., 0., 0.};
88 }
89
90 return vec3(p->q
91 * (this->external_E_field(p->r_vec)
92 + arma::cross(p->v_vec, this->external_B_field(p->r_vec))));
93}
94
96{
97 vec3 res;
98
99 for (size_t j = 0; j < this->particles.size(); j++) {
100 if (i != j)
101 res += this->force_on_particle(i, j);
102 }
103
104 return vec3(res * (K_E * this->particles[i].q));
105}
106
108{
109 if (arma::norm(this->particles[i].r_vec) > this->d) {
110 return vec3{0., 0., 0.};
111 }
112 return vec3(this->total_force_external(i) - this->total_force_particles(i));
113}
114
115PenningTrap::PenningTrap(double B_0, double V_0, double d, double t)
116{
117 this->B_0 = B_0;
118 this->V_0 = V_0;
119 this->d = d;
120 this->t = t;
121 this->perturbation = [](double t) { return 1.; };
122}
123
124PenningTrap::PenningTrap(uint i, double B_0, double V_0, double d, double t)
125 : PenningTrap::PenningTrap(B_0, V_0, d)
126{
127 for (size_t j = 0; j < i; j++) {
128 this->particles.push_back(
129 Particle(vec3(vec3().randn() * .1 * this->d),
130 vec3(vec3().randn() * .1 * this->d)));
131 }
132}
133
134PenningTrap::PenningTrap(std::vector<Particle> particles, double B_0,
135 double V_0, double d, double t)
136 : PenningTrap::PenningTrap(B_0, V_0, d)
137{
138 this->particles = particles;
139}
140
141void PenningTrap::set_pertubation(double f, double omega_V)
142{
143 this->perturbation = [f, omega_V](double t) {
144 return 1 + f * std::cos(omega_V * t);
145 };
146}
147
148void PenningTrap::reinitialize(double f, double omega_V, double t)
149{
150 this->t = t;
151 this->set_pertubation(f, omega_V);
152
153 for (size_t i = 0; i < this->particles.size(); i++) {
154 this->particles[i].r_vec = vec3().randn() * .1 * this->d;
155 }
156}
157
159{
160 this->particles.push_back(particle);
161}
162
163void PenningTrap::evolve_RK4(double dt, bool particle_interaction)
164{
165
166 std::vector<Particle> original_particles = this->particles;
167 std::vector<Particle> tmp_particles = this->particles;
168
169 vec3 (PenningTrap::*force)(uint) = particle_interaction
172
173 size_t size = this->particles.size();
174
175 // Allocating takes a long time, so reuse sim_arr if possible
176 if (this->k_v.size() != 4 || this->k_r.size() != 4
177 || this->k_v[0].size() != size || this->k_r[0].size() != size) {
178 this->k_v = sim_arr(4, sim_cols(size));
179 this->k_r = sim_arr(4, sim_cols(size));
180 }
181
182 // Each k_{i+1} is dependent on k_i, so outer loop is not parallelizable
183 for (size_t i = 0; i < 4; i++) {
184 // Inner loop is able to be parallelized
185#pragma omp parallel for
186 for (size_t j = 0; j < size; j++) {
187 this->k_v[i][j] = (this->*force)(j) / this->particles[j].m;
188 this->k_r[i][j] = this->particles[j].v_vec;
189
190 tmp_particles[j].v_vec =
191 original_particles[j].v_vec + this->v_func(i, j, dt);
192 tmp_particles[j].r_vec =
193 original_particles[j].r_vec + this->r_func(i, j, dt);
194 }
195 this->particles = tmp_particles;
196 }
197
198 this->t += dt;
199}
200
201void PenningTrap::evolve_forward_euler(double dt, bool particle_interaction)
202{
203 size_t size = this->particles.size();
204 vec3 force_res[size];
205
206 vec3 (PenningTrap::*force)(uint) = particle_interaction
209
210 // Calculating the force for each particle is independent and therefore
211 // a good candidate for parallel execution
212#pragma omp parallel for
213 for (size_t i = 0; i < size; i++) {
214 force_res[i] = (this->*force)(i);
215 }
216
217 // Updating the particles is also independent, so we can parallelize
218 // this as well
219#pragma omp parallel for
220 for (size_t i = 0; i < size; i++) {
221 this->particles[i].r_vec += dt * this->particles[i].v_vec;
222 this->particles[i].v_vec += dt * force_res[i] / this->particles[i].m;
223 }
224
225 this->t += dt;
226}
227
228simulation_t PenningTrap::simulate(double time, uint steps, std::string method,
229 bool particle_interaction)
230{
231 double dt = time / (double)steps;
232
233 uint size = this->particles.size();
234
235 simulation_t res{sim_arr(size, sim_cols(steps)),
236 sim_arr(size, sim_cols(steps))};
237
238 void (PenningTrap::*func)(double, bool);
239 if (method == "rk4") {
241 } else if (method == "euler") {
243 } else {
244 std::cout << "Not a valid method!" << std::endl;
245 abort();
246 }
247
248 for (size_t j = 0; j < steps; j++) {
249 for (size_t i = 0; i < size; i++) {
250 res.r_vecs[i][j] = this->particles[i].r_vec;
251 res.v_vecs[i][j] = this->particles[i].v_vec;
252 }
253 (this->*func)(dt, particle_interaction);
254 }
255
256 return res;
257}
258
259void PenningTrap::write_simulation_to_dir(std::string path, double time,
260 uint steps, std::string method,
261 bool particle_interaction)
262{
263 if (path.back() != '/') {
264 path += '/';
265 }
266 if (mkpath(path, 0777) != 0) {
267 std::cout << "Failed to make path" << std::endl;
268 abort();
269 }
270
271 simulation_t res =
272 this->simulate(time, steps, method, particle_interaction);
273
274 std::ofstream ofile;
275
276 // Writing each particle to its own file is independent and can be run in
277 // parallel.
278#pragma omp parallel for private(ofile)
279 for (size_t i = 0; i < this->particles.size(); i++) {
280 ofile.open(path + "particle_" + std::to_string(i) + "_r.txt");
281 for (vec3 &vec : res.r_vecs[i]) {
282 ofile << scientific_format(vec(0), 10, 8) << ','
283 << scientific_format(vec(1), 10, 8) << ','
284 << scientific_format(vec(2), 10, 8) << '\n';
285 }
286 ofile.close();
287 ofile.open(path + "particle_" + std::to_string(i) + "_v.txt");
288 for (vec3 &vec : res.v_vecs[i]) {
289 ofile << scientific_format(vec(0), 10, 8) << ','
290 << scientific_format(vec(1), 10, 8) << ','
291 << scientific_format(vec(2), 10, 8) << '\n';
292 }
293 ofile.close();
294 }
295}
296
297double PenningTrap::fraction_of_particles_left(double time, uint steps,
298 std::string method,
299 bool particle_interaction)
300{
301 double dt = time / (double)steps;
302
303 void (PenningTrap::*func)(double, bool);
304 if (method == "rk4") {
306 } else if (method == "euler") {
308 } else {
309 std::cout << "Not a valid method!" << std::endl;
310 abort();
311 }
312
313 for (size_t j = 0; j < steps; j++) {
314 (this->*func)(dt, particle_interaction);
315 }
316
317 int particles_left = 0;
318
319 // A reduction is perfect here
320#pragma omp parallel for reduction(+ : particles_left)
321 for (size_t i = 0; i < this->particles.size(); i++) {
322 if (arma::norm(this->particles[i].r_vec) < this->d) {
323 particles_left++;
324 }
325 }
326
327 return (double)particles_left / (double)this->particles.size();
328}
A class for simulating a Penning trap.
A class that holds attributes of a particle.
Definition: Particle.hpp:23
vec3 r_vec
position
Definition: Particle.hpp:25
vec3 v_vec
velocity
Definition: Particle.hpp:26
double q
Charge.
Definition: Particle.hpp:27
A class that simulates a Penning trap.
Definition: PenningTrap.hpp:32
std::vector< Particle > particles
The particles in the Penning trap.
Definition: PenningTrap.hpp:42
double B_0
Magnetic field strength.
Definition: PenningTrap.hpp:37
vec3 total_force_external(uint i)
Calculate the total external force on a particle.
Definition: PenningTrap.cpp:82
sim_arr k_r
A 2D vector containing all where is the index of a particle.
Definition: PenningTrap.hpp:52
vec3 total_force_particles(uint i)
Calculate the total force on a particle p_i from other particles.
Definition: PenningTrap.cpp:95
vec3 external_B_field(vec3 r)
Calculate B at point r.
Definition: PenningTrap.cpp:66
void evolve_RK4(double dt, bool particle_interaction=true)
Go forward one timestep using the RK4 method.
vec3 v_func(uint i, uint j, double dt)
Helper for evolve_RK4 when calculating values.
Definition: PenningTrap.cpp:20
vec3 external_E_field(vec3 r)
Calculate E at point r.
Definition: PenningTrap.cpp:58
PenningTrap(double B_0=T, double V_0=(25. *V)/1000., double d=500., double t=0.)
Constructor for the PenningTrap class.
double d
Characteristic dimension.
Definition: PenningTrap.hpp:40
void add_particle(Particle particle)
Add a particle to the system.
double V_0
Applied potential.
Definition: PenningTrap.hpp:38
simulation_t simulate(double time, uint steps, std::string method="rk4", bool particle_interaction=true)
Simulate the particle system inside the Penning trap over a certain amount of time.
vec3 force_on_particle(uint i, uint j)
Calculate the force between 2 particles.
Definition: PenningTrap.cpp:71
vec3 r_func(uint i, uint j, double dt)
Helper for evolve_RK4 when calculating values.
Definition: PenningTrap.cpp:39
double t
Current time.
Definition: PenningTrap.hpp:41
vec3 total_force(uint i)
calculate the total force on a particle p_i.
void set_pertubation(double f, double omega_V)
Time dependent perturbation to V_0.
void reinitialize(double f, double omega_V, double t=0.)
Give all particles new positions and velocities, and change t and V_0.
void evolve_forward_euler(double dt, bool particle_interaction=true)
Go forward one timestep using the forward Euler method.
double fraction_of_particles_left(double time, uint steps, std::string method="rk4", bool particle_interaction=true)
Simulate and calculate what fraction of particles are still left inside the Penning trap after the si...
void write_simulation_to_dir(std::string path, double time, uint steps, std::string method="rk4", bool particle_interaction=true)
Simulate and write the displacement of all particles to files.
std::function< double(double)> perturbation
Time-dependent perturbation.
Definition: PenningTrap.hpp:39
sim_arr k_v
A 2D vector containing all where is the index of a particle.
Definition: PenningTrap.hpp:47
#define K_E
Coulomb constant. unit: .
Definition: constants.hpp:17
Typedef for PenningTrap::simulation return value.
Definition: typedefs.hpp:40
Useful typedefs for cleaner code.
arma::vec::fixed< 3 > vec3
Typedef for a fixed 3d arma vector.
Definition: typedefs.hpp:23
std::vector< vec3 > sim_cols
Typedef for the column of the result vector from simulating particles.
Definition: typedefs.hpp:28
std::vector< sim_cols > sim_arr
Typedef for the result of the simulate method.
Definition: typedefs.hpp:36
bool mkpath(std::string path, int mode=0777)
Make path given.
Definition: utils.cpp:72
std::string scientific_format(double d, int width=20, int prec=10)
Turns a double into a string written in scientific format.
Definition: utils.cpp:15