1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
|
/* Copyright 2019 Luca Fedeli
*
* This file is part of WarpX.
*
* License: BSD-3-Clause-LBNL
*/
#ifndef WARPX_quantum_sync_engine_wrapper_h_
#define WARPX_quantum_sync_engine_wrapper_h_
#include "BreitWheelerEngineWrapper_fwd.H"
#include "QedChiFunctions.H"
#include "QedWrapperCommons.H"
#include "Utils/WarpXConst.H"
#include <AMReX_Extension.H>
#include <AMReX_GpuQualifiers.H>
#include <AMReX_REAL.H>
#include <AMReX_Random.H>
#include <picsar_qed/containers/picsar_array.hpp>
#include <picsar_qed/math/cmath_overloads.hpp>
#include <picsar_qed/math/vec_functions.hpp>
#include <picsar_qed/physics/phys_constants.h>
#include <picsar_qed/physics/quantum_sync/quantum_sync_engine_core.hpp>
#include <picsar_qed/physics/quantum_sync/quantum_sync_engine_tables.hpp>
#include <picsar_qed/physics/unit_conversion.hpp>
#include <cmath>
#include <vector>
namespace amrex { struct RandomEngine; }
// Aliases =============================
using QS_dndt_table_params =
picsar::multi_physics::phys::quantum_sync::
dndt_lookup_table_params<amrex::ParticleReal>;
using QS_dndt_table =
picsar::multi_physics::phys::quantum_sync::
dndt_lookup_table<
amrex::ParticleReal,
PicsarQedVector<amrex::ParticleReal>>;
using QS_dndt_table_view = QS_dndt_table::view_type;
using QS_phot_em_table_params =
picsar::multi_physics::phys::quantum_sync::
photon_emission_lookup_table_params<amrex::ParticleReal>;
using QS_phot_em_table =
picsar::multi_physics::phys::quantum_sync::
photon_emission_lookup_table<
amrex::ParticleReal,
PicsarQedVector<amrex::ParticleReal>>;
using QS_phot_em_table_view = QS_phot_em_table::view_type;
struct PicsarQuantumSyncCtrl
{
QS_dndt_table_params dndt_params;
QS_phot_em_table_params phot_em_params;
};
// Functors ==================================
// These functors allow using the core elementary functions of the library.
// They are generated by a factory class (QuantumSynchrotronEngine, see below).
// They can be included in GPU kernels.
/**
* Functor to initialize the optical depth of leptons for the
* Quantum Synchrotron process
*/
class QuantumSynchrotronGetOpticalDepth
{
public:
/**
* Constructor does nothing because optical depth initialization
* does not require control parameters or lookup tables.
*/
QuantumSynchrotronGetOpticalDepth () = default;
/**
* () operator is just a thin wrapper around a very simple function to
* generate the optical depth. It can be used on GPU.
*/
AMREX_GPU_HOST_DEVICE
AMREX_FORCE_INLINE
amrex::ParticleReal operator() (amrex::RandomEngine const& engine) const noexcept
{
namespace pxr_qs = picsar::multi_physics::phys::quantum_sync;
//A random number in [0,1) should be provided as an argument.
return pxr_qs::get_optical_depth(amrex::Random(engine));
}
};
//____________________________________________
/**
* Functor to evolve the optical depth of leptons due to the
* Quantum Synchrotron process
*/
class QuantumSynchrotronEvolveOpticalDepth
{
public:
/**
* Default constructor: it leaves the functor in a non-initialized state.
*/
QuantumSynchrotronEvolveOpticalDepth () = default;
/**
* Constructor to be used to initialize the functor.
*
* @param[in] table_view a view of a QS_dndt_table lookup table
* @param[in] qs_minimum_chi_part the minimum quantum parameter to evolve the optical depth
*/
QuantumSynchrotronEvolveOpticalDepth (
const QS_dndt_table_view table_view,
const amrex::ParticleReal qs_minimum_chi_part):
m_table_view{table_view}, m_qs_minimum_chi_part{qs_minimum_chi_part}{}
/**
* Evolves the optical depth. It can be used on GPU.
* If the quantum parameter parameter of the particle is
* < qs_minimum_chi_part, the method returns immediately.
*
* @param[in] ux,uy,uz gamma*v components of the lepton.
* @param[in] ex,ey,ez electric field components (SI units)
* @param[in] bx,by,bz magnetic field components (SI units)
* @param[in] dt timestep (SI units)
* @param[in,out] opt_depth optical depth of the particle.
* @return a flag which is 1 if chi_part was out of table.
*/
AMREX_GPU_DEVICE
AMREX_FORCE_INLINE
int operator()(
const amrex::ParticleReal ux, const amrex::ParticleReal uy, const amrex::ParticleReal uz,
const amrex::ParticleReal ex, const amrex::ParticleReal ey, const amrex::ParticleReal ez,
const amrex::ParticleReal bx, const amrex::ParticleReal by, const amrex::ParticleReal bz,
const amrex::Real dt, amrex::ParticleReal& opt_depth) const noexcept
{
using namespace amrex::literals;
namespace pxr_p = picsar::multi_physics::phys;
namespace pxr_qs = picsar::multi_physics::phys::quantum_sync;
constexpr amrex::ParticleReal m_e = PhysConst::m_e;
constexpr amrex::ParticleReal inv_c2 = 1._rt/(PhysConst::c*PhysConst::c);
const amrex::ParticleReal gamma = std::sqrt(1._rt + (ux*ux + uy*uy + uz*uz)*inv_c2);
const auto energy = gamma*m_e*PhysConst::c*PhysConst::c;
const auto chi_part = QedUtils::chi_ele_pos(
m_e*ux, m_e*uy, m_e*uz, ex, ey, ez, bx, by, bz);
if (chi_part < m_qs_minimum_chi_part)
return 0;
const auto is_out = pxr_qs::evolve_optical_depth<
amrex::ParticleReal,
QS_dndt_table_view,
pxr_p::unit_system::SI>(
energy, chi_part, dt, opt_depth, m_table_view);
return is_out;
}
private:
QS_dndt_table_view m_table_view;
amrex::ParticleReal m_qs_minimum_chi_part;
};
/**
* Functor to generate a photon via the Quantum Synchrotron process
* and to update momentum accordingly
*/
class QuantumSynchrotronPhotonEmission
{
public:
/**
* Default constructor: it leaves the functor in a non-initialized state.
*/
QuantumSynchrotronPhotonEmission () = default;
/**
* Constructor acquires pointers to control parameters and
* lookup tables data.
* lookup_table uses non-owning vectors under the hood. So no new data
* allocations should be triggered on GPU
*
* @param[in] table_view a view of a QS_phot_em_table lookup table
*/
QuantumSynchrotronPhotonEmission (
const QS_phot_em_table_view table_view):
m_table_view{table_view}{}
/**
* Generates photons according to Quantum Synchrotron process.
* It can be used on GPU.
*
* @param[in,out] ux,uy,uz gamma*v components of the lepton. They are modified (SI units)
* @param[in] ex,ey,ez electric field components (SI units)
* @param[in] bx,by,bz magnetic field components (SI units)
* @param[out] g_ux,g_uy,g_uz gamma*v components of the generated photon (SI units)
* @param[in] engine random number generator engine
* @return a flag which is 1 if chi_photon was out of table
*/
AMREX_GPU_DEVICE
AMREX_FORCE_INLINE
bool operator()(
amrex::ParticleReal& ux, amrex::ParticleReal& uy,
amrex::ParticleReal& uz,
const amrex::ParticleReal ex, const amrex::ParticleReal ey,
const amrex::ParticleReal ez, const amrex::ParticleReal bx,
const amrex::ParticleReal by, const amrex::ParticleReal bz,
amrex::ParticleReal& g_ux, amrex::ParticleReal& g_uy,
amrex::ParticleReal& g_uz, amrex::RandomEngine const& engine) const noexcept
{
using namespace amrex;
namespace pxr_m = picsar::multi_physics::math;
namespace pxr_p = picsar::multi_physics::phys;
namespace pxr_qs = picsar::multi_physics::phys::quantum_sync;
const auto rand_zero_one_minus_epsi = amrex::Random(engine);
constexpr ParticleReal me = PhysConst::m_e;
constexpr ParticleReal one_over_me = 1._prt/me;
// Particle momentum is stored as gamma * velocity.
// Convert to m * gamma * velocity
auto px = ux*me;
auto py = uy*me;
auto pz = uz*me;
const auto chi_particle = QedUtils::chi_ele_pos(
px, py, pz, ex, ey, ez, bx, by, bz);
auto momentum_particle = pxr_m::vec3<amrex::ParticleReal>{px, py, pz};
auto momentum_photon = pxr_m::vec3<amrex::ParticleReal>();
const auto is_out = pxr_qs::generate_photon_update_momentum<
amrex::ParticleReal,
QS_phot_em_table_view,
pxr_p::unit_system::SI>(
chi_particle, momentum_particle,
rand_zero_one_minus_epsi,
m_table_view,
momentum_photon);
ux = momentum_particle[0]*one_over_me;
uy = momentum_particle[1]*one_over_me;
uz = momentum_particle[2]*one_over_me;
g_ux = momentum_photon[0]*one_over_me;
g_uy = momentum_photon[1]*one_over_me;
g_uz = momentum_photon[2]*one_over_me;
return is_out;
}
private:
QS_phot_em_table_view m_table_view;
};
// Factory class =============================
/**
* Wrapper for the Quantum Synchrotron engine of the PICSAR library
*/
class QuantumSynchrotronEngine
{
public:
/**
* Constructor requires no arguments.
*/
QuantumSynchrotronEngine () = default;
/**
* Builds the functor to initialize the optical depth
*/
QuantumSynchrotronGetOpticalDepth build_optical_depth_functor ();
/**
* Builds the functor to evolve the optical depth
*/
QuantumSynchrotronEvolveOpticalDepth build_evolve_functor ();
/**
* Builds the functor to generate photons
*/
QuantumSynchrotronPhotonEmission build_phot_em_functor ();
/**
* Checks if the optical tables are properly initialized
*/
bool are_lookup_tables_initialized () const;
/**
* Export lookup tables data into a raw binary Vector
*
* @return the data in binary format. The Vector is empty if tables were
* not previously initialized.
*/
std::vector<char> export_lookup_tables_data () const;
/**
* Init lookup tables from raw binary data.
*
* @param[in] raw_data a vector of char
* @param[in] qs_minimum_chi_part minimum chi parameter to evolve the optical depth of a particle.
* @return true if it succeeds, false if it cannot parse raw_data
*/
bool init_lookup_tables_from_raw_data (const std::vector<char>& raw_data,
amrex::ParticleReal qs_minimum_chi_part);
/**
* Init lookup tables using built-in (low resolution) tables
*
* @param[in] qs_minimum_chi_part minimum chi parameter to evolve the optical depth of a particle.
*/
void init_builtin_tables(amrex::ParticleReal qs_minimum_chi_part);
/**
* Computes the lookup tables. It does nothing unless WarpX is compiled with QED_TABLE_GEN=TRUE
*
* @param[in] ctrl control params to generate the tables
* @param[in] qs_minimum_chi_part minimum chi parameter to evolve the optical depth of a particle.
*/
void compute_lookup_tables (PicsarQuantumSyncCtrl ctrl,
amrex::ParticleReal qs_minimum_chi_part);
/**
* gets default values for the control parameters
*
* @return default control params to generate the tables
*/
PicsarQuantumSyncCtrl get_default_ctrl() const;
amrex::ParticleReal get_minimum_chi_part() const;
private:
bool m_lookup_tables_initialized = false;
//Variables to store the minimum chi parameters to enable
//Quantum Synchrotron process
amrex::ParticleReal m_qs_minimum_chi_part;
QS_dndt_table m_dndt_table;
QS_phot_em_table m_phot_em_table;
void init_builtin_dndt_table();
void init_builtin_phot_em_table();
};
//============================================
#endif //WARPX_quantum_sync_engine_wrapper_h_
|