Penning Trap Simulation
Simulate particle behavior inside a Penning Trap
Loading...
Searching...
No Matches
main.cpp
Go to the documentation of this file.
1
13#include <cmath>
14#include <complex>
15#include <fstream>
16#include <omp.h>
17#include <string>
18#include <vector>
19
20#include "PenningTrap.hpp"
21#include "constants.hpp"
22#include "utils.hpp"
23
24#define PARTICLES 100
25#define N 40000
26
27// Particles used for testing
28Particle p1(vec3{20., 0., 20.}, vec3{0., 25., 0.});
29Particle p2(vec3{25., 25., 0.}, vec3{0., 40., 5.});
30
38{
39 double w_0 = T / CA_MASS;
40 double w_z2 = (50. * V / 1000.) / (CA_MASS * 500. * 500.);
41 double w_p = (w_0 + std::sqrt(w_0 * w_0 - 2. * w_z2)) / 2.;
42 double w_n = (w_0 - std::sqrt(w_0 * w_0 - 2. * w_z2)) / 2.;
43 double A_p = (25. + w_n * 20.) / (w_n - w_p);
44 double A_n = -(25. + w_p * 20.) / (w_n - w_p);
45 std::complex<double> f =
46 A_p * std::exp(std::complex<double>(0., -w_p * t))
47 + A_n * std::exp(std::complex<double>(0., -w_n * t));
48 vec3 res{std::real(f), std::imag(f), 20. * std::cos(std::sqrt(w_z2) * t)};
49 return res;
50}
51
55{
56 // Initialize trap with particle 1
57 PenningTrap trap(std::vector<Particle>{p1});
58
59 double time = 50.; // microseconds
60
61 // Simulate and write results to file
62 trap.write_simulation_to_dir("output/simulate_single_particle", time, N,
63 "rk4", false);
64}
65
70{
71 // Initialize traps with particles
72 PenningTrap trap_no_interaction(std::vector<Particle>{p1, p2});
73 PenningTrap trap_with_interaction(std::vector<Particle>{p1, p2});
74
75 double time = 50.; // microseconds
76
77 // Simulate and write results to files
78 trap_no_interaction.write_simulation_to_dir(
79 "output/simulate_2_particles/no_interaction", time, N, "rk4", false);
80 trap_with_interaction.write_simulation_to_dir(
81 "output/simulate_2_particles/with_interaction", time, N);
82}
83
88{
89 double time = 50.; // microseconds
90
91 std::ofstream ofile;
92
93 // Calculate relative error for RK4
94 std::string path = "output/relative_error/RK4/";
95 mkpath(path);
96#pragma omp parallel for
97 for (int i = 0; i < 4; i++) {
98 int steps = 4000 * std::pow(2, i);
99 double dt = time / (double)steps;
100 ofile.open(path + std::to_string(steps) + "_steps.txt");
101 PenningTrap trap(std::vector<Particle>{p1});
102 simulation_t res = trap.simulate(time, steps, "rk4", false);
103 for (int i = 0; i < steps; i++) {
104 ofile << arma::norm(res.r_vecs[0][i]
106 << '\n';
107 }
108 ofile.close();
109 }
110
111 // Calculate relative error for forward Euler
112 path = "output/relative_error/euler/";
113 mkpath(path);
114#pragma omp parallel for
115 for (int i = 0; i < 4; i++) {
116 int steps = 4000 * std::pow(2, i);
117 double dt = time / (double)steps;
118 ofile.open(path + std::to_string(steps) + "_steps.txt");
119 PenningTrap trap(std::vector<Particle>{p1});
120 simulation_t res = trap.simulate(time, steps, "euler", false);
121 for (int i = 0; i < steps; i++) {
122 ofile << arma::norm(res.r_vecs[0][i]
124 << '\n';
125 }
126 ofile.close();
127 }
128}
129
133{
134 PenningTrap trap((unsigned)100);
135
136 double time = 50.; // microseconds
137
138 trap.write_simulation_to_dir("output/simulate_100_particles", time, N,
139 "rk4", false);
140}
141
150{
151 double time = 500.;
152
153 double amplitudes[]{.1, .4, .7};
154
155 double freq_start = .2;
156 double freq_end = 2.5;
157 double freq_increment = .02;
158 size_t freq_iterations =
159 (size_t)((freq_end - freq_start) / freq_increment) + 1;
160
161 double res[4][freq_iterations];
162
163 std::string path = "output/time_dependent_potential/";
164 mkpath(path);
165
166 std::ofstream ofile;
167
168#pragma omp parallel for
169 // Insert frequencies
170 for (size_t i = 0; i < freq_iterations; i++) {
171 res[0][i] = freq_start + freq_increment * i;
172 }
173
174#pragma omp parallel
175 {
176 // Each thread creates a PenningTrap instance and reuses it throughout
177 // the sweep.
178 PenningTrap trap((uint)100);
179#pragma omp for collapse(2)
180 for (size_t i = 0; i < 3; i++) {
181 for (size_t j = 0; j < freq_iterations; j++) {
182 // Reset particles and give new time dependent potential.
183 trap.reinitialize(amplitudes[i], res[0][j]);
184 res[i + 1][j] =
185 trap.fraction_of_particles_left(time, N, "rk4", false);
186 }
187 }
188 }
189
190 // Write results to file
191 ofile.open(path + "wide_sweep.txt");
192 for (size_t i = 0; i < freq_iterations; i++) {
193 ofile << res[0][i] << ',' << res[1][i] << ',' << res[2][i] << ','
194 << res[3][i] << '\n';
195 }
196 ofile.close();
197}
198
207{
208 double time = 500.;
209
210 double amplitudes[]{.1, .4, .7};
211
212 double freq_start = 1.;
213 double freq_end = 1.7;
214 double freq_increment = .002;
215 size_t freq_iterations = (size_t)((freq_end - freq_start) / freq_increment);
216
217 double res[4][freq_iterations];
218
219 std::string path = "output/time_dependent_potential/";
220 mkpath(path);
221
222 std::ofstream ofile;
223
224#pragma omp parallel for
225 // Insert frequencies
226 for (size_t i = 0; i < freq_iterations; i++) {
227 res[0][i] = freq_start + freq_increment * i;
228 }
229
230#pragma omp parallel
231 {
232 // Each thread creates a PenningTrap instance and reuses it throughout
233 // the sweep.
234 PenningTrap trap((uint)100);
235#pragma omp for collapse(2)
236 for (size_t i = 0; i < 3; i++) {
237 for (size_t j = 0; j < freq_iterations; j++) {
238 // Reset particles and give new time dependent potential.
239 trap.reinitialize(amplitudes[i], res[0][j]);
240 res[i + 1][j] =
241 trap.fraction_of_particles_left(time, N, "rk4", false);
242 }
243 }
244 }
245
246 // Write results to file
247 ofile.open(path + "narrow_sweep.txt");
248 for (size_t i = 0; i < freq_iterations; i++) {
249 ofile << res[0][i] << ',' << res[1][i] << ',' << res[2][i] << ','
250 << res[3][i] << '\n';
251 }
252 ofile.close();
253}
254
263{
264 double time = 500.;
265
266 double amplitudes[]{.1, .4, .7};
267
268 double freq_start = 1.;
269 double freq_end = 1.7;
270 double freq_increment = .002;
271 size_t freq_iterations = (size_t)((freq_end - freq_start) / freq_increment);
272
273 double res[4][freq_iterations];
274
275 std::string path = "output/time_dependent_potential/";
276 mkpath(path);
277
278 std::ofstream ofile;
279
280#pragma omp parallel for
281 for (size_t i = 0; i < freq_iterations; i++) {
282 res[0][i] = freq_start + freq_increment * i;
283 }
284
285#pragma omp parallel
286 {
287 // Each thread creates a PenningTrap instance and reuses it throughout
288 // the sweep.
289 PenningTrap trap((uint)100);
290#pragma omp for collapse(2)
291 for (size_t i = 0; i < 3; i++) {
292 for (size_t j = 0; j < freq_iterations; j++) {
293 // Reset particles and give new time dependent potential.
294 trap.reinitialize(amplitudes[i], res[0][j]);
295 res[i + 1][j] = trap.fraction_of_particles_left(time, N);
296 }
297 }
298 }
299
300 // Write results to file
301 ofile.open(path + "narrow_sweep_interactions.txt");
302 for (size_t i = 0; i < freq_iterations; i++) {
303 ofile << res[0][i] << ',' << res[1][i] << ',' << res[2][i] << ','
304 << res[3][i] << '\n';
305 }
306 ofile.close();
307}
308
309int main()
310{
311 double start, end, t1, t2;
312 start = omp_get_wtime();
313
315
317
319
320 t2 = omp_get_wtime();
321
322 std::cout << "Time single and double : " << (t2 - start)
323 << " seconds" << std::endl;
324
325 t1 = omp_get_wtime();
326
328
329 t2 = omp_get_wtime();
330
331 std::cout << "Time 100 particles : " << (t2 - t1)
332 << " seconds" << std::endl;
333
334 t1 = omp_get_wtime();
335
337
338 t2 = omp_get_wtime();
339
340 std::cout << "Time wide sweep : " << (t2 - t1)
341 << " seconds" << std::endl;
342
343 t1 = omp_get_wtime();
344
346
347 t2 = omp_get_wtime();
348
349 std::cout << "Time narrow sweep no interaction : " << (t2 - t1)
350 << " seconds" << std::endl;
351
352 t1 = omp_get_wtime();
353
355
356 t2 = omp_get_wtime();
357
358 std::cout << "Time narrow sweep with interaction: " << (t2 - t1)
359 << " seconds" << std::endl;
360
361 end = omp_get_wtime();
362
363 std::cout << "Time : " << (end - start)
364 << " seconds" << std::endl;
365
366 return 0;
367}
A class for simulating a Penning trap.
A class that holds attributes of a particle.
Definition: Particle.hpp:23
A class that simulates a Penning trap.
Definition: PenningTrap.hpp:32
void reinitialize(double f, double omega_V, double t=0.)
Give all particles new positions and velocities, and change t and V_0.
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.
Library of constants.
#define T
1 Tesla. unit:
Definition: constants.hpp:21
#define CA_MASS
Mass of a single calcium ion. unit: amu.
Definition: constants.hpp:29
#define V
1 Volt. unit:
Definition: constants.hpp:25
void simulate_100_particles()
Simulate 100 particles over 50 .
Definition: main.cpp:132
Particle p1(vec3{20., 0., 20.}, vec3{0., 25., 0.})
Particle 1.
void potential_resonance_narrow_sweep()
Simulate 100 particles over 500 using a time dependent potential.
Definition: main.cpp:206
void potential_resonance_wide_sweep()
Simulate 100 particles over 500 using a time dependent potential.
Definition: main.cpp:149
void simulate_two_particles()
Simulate 2 particles over the period of 50 with and without particle interactions.
Definition: main.cpp:69
vec3 analytical_solution_particle_1(double t)
The analytical solution for particle p1.
Definition: main.cpp:37
void simulate_single_particle()
Simulate a single particle over the period of 50 .
Definition: main.cpp:54
void potential_resonance_narrow_sweep_interaction()
Simulate 100 particles over 500 using a time dependent potential.
Definition: main.cpp:262
void simulate_single_particle_with_different_steps()
Simulate a single particle over 50 using different amount of steps and different methods.
Definition: main.cpp:87
Particle p2(vec3{25., 25., 0.}, vec3{0., 40., 5.})
Particle 2.
Typedef for PenningTrap::simulation return value.
Definition: typedefs.hpp:40
arma::vec::fixed< 3 > vec3
Typedef for a fixed 3d arma vector.
Definition: typedefs.hpp:23
Function prototypes and macros that are useful.
bool mkpath(std::string path, int mode=0777)
Make path given.
Definition: utils.cpp:72