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 <algorithm>
14#include <functional>
15#include <sys/types.h>
16#include <vector>
17
18#include "PenningTrap.hpp"
19#include "typedefs.hpp"
20
21vec3 PenningTrap::v_func(uint i, uint j, double dt)
22{
23 switch (i) {
24 case 0:
25 return .5 * dt * this->k_v[0][j];
26 case 1:
27 return .5 * dt * this->k_v[1][j];
28 case 2:
29 return dt * this->k_v[2][j];
30 case 3:
31 return vec3((dt / 6.)
32 * (this->k_v[0][j] + 2. * this->k_v[1][j]
33 + 2. * this->k_v[2][j] + this->k_v[3][j]));
34 default:
35 std::cout << "Not valid!" << std::endl;
36 abort();
37 }
38}
39
40vec3 PenningTrap::r_func(uint i, uint j, double dt)
41{
42 switch (i) {
43 case 0:
44 return .5 * dt * this->k_r[0][j];
45 case 1:
46 return .5 * dt * this->k_r[1][j];
47 case 2:
48 return dt * this->k_r[2][j];
49 case 3:
50 return vec3((dt / 6.)
51 * (this->k_r[0][j] + 2. * this->k_r[1][j]
52 + 2. * this->k_r[2][j] + this->k_r[3][j]));
53 default:
54 std::cout << "Not valid!" << std::endl;
55 abort();
56 }
57}
58
60{
61 r(2) *= -2.;
62
63 return vec3(
64 ((this->V_0 * this->perturbation(this->t)) / (this->d * this->d)) * r);
65}
66
68{
69 return vec3{0., 0., this->B_0};
70}
71
73{
74 // Calculate the difference between the particles' position
75 vec3 res = this->particles[i].r_vec - this->particles[j].r_vec;
76
77 // Get the distance between the particles
78 double norm = arma::norm(res, 2);
79
80 return vec3((this->particles[j].q / (norm * norm * norm)) * res);
81}
82
84{
85 Particle *p = &this->particles[i];
86
87 return vec3(p->q
88 * (this->external_E_field(p->r_vec)
89 + arma::cross(p->v_vec, this->external_B_field(p->r_vec))));
90}
91
93{
94 vec3 res;
95
96 for (size_t j = 0; j < this->particles.size(); j++) {
97 if (i != j)
98 res += this->force_on_particle(i, j);
99 }
100
101 return vec3(res * (K_E * this->particles[i].q));
102}
103
105{
106 if (arma::norm(this->particles[i].r_vec) > this->d) {
107 return vec3{0., 0., 0.};
108 }
109 return vec3(this->total_force_external(i) - this->total_force_particles(i));
110}
111
113{
114 if (arma::norm(this->particles[i].r_vec) > this->d) {
115 return vec3{0., 0., 0.};
116 }
117 return vec3(this->total_force_external(i));
118}
119
120PenningTrap::PenningTrap(double B_0, double V_0, double d, double t)
121{
122 this->B_0 = B_0;
123 this->V_0 = V_0;
124 this->d = d;
125 this->t = t;
126 this->perturbation = [](double t) { return 1.; };
127}
128
129PenningTrap::PenningTrap(uint i, double B_0, double V_0, double d, double t)
130 : PenningTrap::PenningTrap(B_0, V_0, d)
131{
132 for (size_t j = 0; j < i; j++) {
133 this->particles.push_back(
134 Particle(vec3(vec3().randn() * .1 * this->d),
135 vec3(vec3().randn() * .1 * this->d)));
136 }
137}
138
139PenningTrap::PenningTrap(std::vector<Particle> particles, double B_0,
140 double V_0, double d, double t)
141 : PenningTrap::PenningTrap(B_0, V_0, d)
142{
143 this->particles = particles;
144}
145
146void PenningTrap::set_pertubation(double f, double omega_V)
147{
148 this->perturbation = [f, omega_V](double t) {
149 return 1 + f * std::cos(omega_V * t);
150 };
151}
152
153void PenningTrap::reinitialize(double f, double omega_V, double t)
154{
155 this->t = t;
156 this->set_pertubation(f, omega_V);
157 Particle *p;
158
159 for (size_t i = 0; i < this->particles.size(); i++) {
160 p = &this->particles[i];
161 p->r_vec = vec3().randn() * .1 * this->d;
162 p->v_vec = vec3().randn() * .1 * this->d;
163 }
164}
165
167{
168 this->particles.push_back(particle);
169}
170
171void PenningTrap::evolve_RK4(double dt, bool particle_interaction)
172{
173 Particle *p;
174
175 // Keep original particles
176 std::vector<Particle> original_particles = this->particles;
177 std::vector<Particle> tmp_particles = this->particles;
178
179 vec3 (PenningTrap::*force)(uint) =
180 particle_interaction ? &PenningTrap::total_force
182
183 size_t size = this->particles.size();
184
185 // Allocating takes a long time, so reuse sim_arr if possible
186 if (this->k_v.size() != 4 || this->k_r.size() != 4
187 || this->k_v[0].size() != size || this->k_r[0].size() != size) {
188 this->k_v = sim_arr(4, sim_cols(size));
189 this->k_r = sim_arr(4, sim_cols(size));
190 }
191
192 // Each k_{i+1} is dependent on k_i, so outer loop is not parallelizable
193 for (size_t i = 0; i < 4; i++) {
194 // Inner loop is able to be parallelized
195#pragma omp parallel for private(p)
196 for (size_t j = 0; j < size; j++) {
197 this->k_v[i][j] = (this->*force)(j) / this->particles[j].m;
198 this->k_r[i][j] = this->particles[j].v_vec;
199
200 p = &tmp_particles[j];
201
202 p->v_vec = original_particles[j].v_vec + this->v_func(i, j, dt);
203 p->r_vec = original_particles[j].r_vec + this->r_func(i, j, dt);
204 }
205 this->particles = tmp_particles;
206 }
207
208 this->t += dt;
209}
210
211void PenningTrap::evolve_forward_euler(double dt, bool particle_interaction)
212{
213 size_t size = this->particles.size();
214 vec3 force_res[size];
215 Particle *p;
216
217 vec3 (PenningTrap::*force)(uint) =
218 particle_interaction ? &PenningTrap::total_force
220
221 // Calculating the force for each particle is independent and therefore
222 // a good candidate for parallel execution
223#pragma omp parallel for
224 for (size_t i = 0; i < size; i++) {
225 force_res[i] = (this->*force)(i);
226 }
227
228 // Updating the particles is also independent, so we can parallelize
229 // this as well
230#pragma omp parallel for
231 for (size_t i = 0; i < size; i++) {
232 p = &this->particles[i];
233 p->r_vec += dt * p->v_vec;
234 p->v_vec += dt * force_res[i] / p->m;
235 }
236
237 this->t += dt;
238}
239
240simulation_t PenningTrap::simulate(double time, uint steps, std::string method,
241 bool particle_interaction)
242{
243 Particle *p;
244 double dt = time / (double)steps;
245
246 uint size = this->particles.size();
247
248 simulation_t res{sim_arr(size, sim_cols(steps)),
249 sim_arr(size, sim_cols(steps))};
250
251 void (PenningTrap::*func)(double, bool);
252 if (method == "rk4") {
254 } else if (method == "euler") {
256 } else {
257 std::cout << "Not a valid method!" << std::endl;
258 abort();
259 }
260
261 for (size_t j = 0; j < steps; j++) {
262 for (size_t i = 0; i < size; i++) {
263 p = &this->particles[i];
264 res.r_vecs[i][j] = p->r_vec;
265 res.v_vecs[i][j] = p->v_vec;
266 }
267 (this->*func)(dt, particle_interaction);
268 }
269
270 return res;
271}
272
273void PenningTrap::write_simulation_to_dir(std::string path, double time,
274 uint steps, std::string method,
275 bool particle_interaction)
276{
277 if (path.back() != '/') {
278 path += '/';
279 }
280 if (mkpath(path, 0777) != 0) {
281 std::cout << "Failed to make path" << std::endl;
282 abort();
283 }
284
285 simulation_t res =
286 this->simulate(time, steps, method, particle_interaction);
287
288 std::ofstream ofile;
289
290 // Writing each particle to its own file is independent and can be run in
291 // parallel.
292#pragma omp parallel for private(ofile)
293 for (size_t i = 0; i < this->particles.size(); i++) {
294 ofile.open(path + "particle_" + std::to_string(i) + "_r.txt");
295 for (vec3 &vec : res.r_vecs[i]) {
296 ofile << scientific_format(vec(0), 10, 8) << ','
297 << scientific_format(vec(1), 10, 8) << ','
298 << scientific_format(vec(2), 10, 8) << '\n';
299 }
300 ofile.close();
301 ofile.open(path + "particle_" + std::to_string(i) + "_v.txt");
302 for (vec3 &vec : res.v_vecs[i]) {
303 ofile << scientific_format(vec(0), 10, 8) << ','
304 << scientific_format(vec(1), 10, 8) << ','
305 << scientific_format(vec(2), 10, 8) << '\n';
306 }
307 ofile.close();
308 }
309}
310
311double PenningTrap::fraction_of_particles_left(double time, uint steps,
312 std::string method,
313 bool particle_interaction)
314{
315 double dt = time / (double)steps;
316
317 void (PenningTrap::*func)(double, bool);
318 if (method == "rk4") {
320 } else if (method == "euler") {
322 } else {
323 std::cout << "Not a valid method!" << std::endl;
324 abort();
325 }
326
327 for (size_t j = 0; j < steps; j++) {
328 (this->*func)(dt, particle_interaction);
329 }
330
331 int particles_left = 0;
332
333 // A reduction is perfect here
334#pragma omp parallel for reduction(+ : particles_left)
335 for (size_t i = 0; i < this->particles.size(); i++) {
336 if (arma::norm(this->particles[i].r_vec) < this->d) {
337 particles_left++;
338 }
339 }
340
341 return (double)particles_left / (double)this->particles.size();
342}
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
double m
Mass.
Definition: Particle.hpp:28
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:83
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:92
vec3 external_B_field(vec3 r)
Calculate B at point r.
Definition: PenningTrap.cpp:67
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:21
vec3 external_E_field(vec3 r)
Calculate E at point r.
Definition: PenningTrap.cpp:59
PenningTrap(double B_0=T, double V_0=(25. *V)/1000., double d=500., double t=0.)
Constructor for the PenningTrap class.
vec3 total_force_no_interaction(uint i)
calculate the total force on a particle p_i without interaction
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:72
vec3 r_func(uint i, uint j, double dt)
Helper for evolve_RK4 when calculating values.
Definition: PenningTrap.cpp:40
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