#include <MultiLaser.H>
|
| MultiLaser () |
|
| ~MultiLaser () |
|
void | ReadParameters () |
|
amrex::MultiFab & | getSlices () |
|
const amrex::MultiFab & | getSlices () const |
|
void | InitData (const amrex::BoxArray &slice_ba, const amrex::DistributionMapping &slice_dm, const amrex::Geometry &geom_3D) |
| Allocate laser multifab. More...
|
|
void | InitSliceEnvelope (const int islice, const int comp) |
| Initialize on slice of the 3D laser field. More...
|
|
void | GetEnvelopeFromFileHelper (const amrex::Geometry &gm) |
| Read in a laser from an openPMD file. More...
|
|
template<typename input_type > |
void | GetEnvelopeFromFile (const amrex::Geometry &gm) |
| Read in a laser from an openPMD file. More...
|
|
void | ShiftLaserSlices () |
| Shift 2D slices in zeta. More...
|
|
void | AdvanceSlice (const Fields &fields, amrex::Real dt, int step) |
|
void | AdvanceSliceMG (const Fields &fields, amrex::Real dt, int step) |
|
void | AdvanceSliceFFT (const Fields &fields, amrex::Real dt, int step) |
|
void | InitLaserSlice (const amrex::Geometry &geom, const int islice, const int comp) |
|
void | InSituComputeDiags (int step, amrex::Real time, int islice, const amrex::Geometry &geom3D, int max_step, amrex::Real max_time) |
|
void | InSituWriteToFile (int step, amrex::Real time, const amrex::Geometry &geom3D, int max_step, amrex::Real max_time) |
|
amrex::Real | GetLambda0 () const |
|
◆ SpectralFieldLoc
◆ MultiLaser()
MultiLaser::MultiLaser |
( |
| ) |
|
|
inlineexplicit |
◆ ~MultiLaser()
MultiLaser::~MultiLaser |
( |
| ) |
|
|
inline |
◆ AdvanceSlice()
void MultiLaser::AdvanceSlice |
( |
const Fields & |
fields, |
|
|
amrex::Real |
dt, |
|
|
int |
step |
|
) |
| |
Wrapper function to advance a laser slice by 1 time step.
- Parameters
-
[in] | fields | Field object |
[in] | dt | time step of the simulation |
[in] | step | current iteration. Needed because step 0 needs a specific treatment. |
◆ AdvanceSliceFFT()
void MultiLaser::AdvanceSliceFFT |
( |
const Fields & |
fields, |
|
|
amrex::Real |
dt, |
|
|
int |
step |
|
) |
| |
Advance a laser slice by 1 time step using a FFT solver. The complex phase of the envelope is evaluated on-axis only.
- Parameters
-
[in] | fields | Field object |
[in] | dt | time step of the simulation |
[in] | step | current iteration. Needed because step 0 needs a specific treatment. |
◆ AdvanceSliceMG()
void MultiLaser::AdvanceSliceMG |
( |
const Fields & |
fields, |
|
|
amrex::Real |
dt, |
|
|
int |
step |
|
) |
| |
Advance a laser slice by 1 time step using a multigrid solver. The complex phase of the envelope is evaluated on-axis only, but can be generalized to everywhere.
- Parameters
-
[in] | fields | Field object |
[in] | dt | time step of the simulation |
[in] | step | current iteration. Needed because step 0 needs a specific treatment. |
◆ GetEnvelopeFromFile()
template<typename input_type >
Read in a laser from an openPMD file.
- Parameters
-
[in] | gm | Geometry for level 0 |
◆ GetEnvelopeFromFileHelper()
Read in a laser from an openPMD file.
- Parameters
-
[in] | gm | Geometry for level 0 |
◆ GetLambda0()
amrex::Real MultiLaser::GetLambda0 |
( |
| ) |
const |
|
inline |
Get the central wavelength
◆ getSlices() [1/2]
get function for the 2D slices
◆ getSlices() [2/2]
get function for the 2D slices (const version)
◆ InitData()
Allocate laser multifab.
- Parameters
-
[in] | slice_ba | box array of the slice |
[in] | slice_dm | corresponding distribution mapping |
[in] | geom_3D | 3D Geometry for level 0 |
◆ InitLaserSlice()
Initialize 1 longitudinal slice of the laser, and store it in n00j00 (current time step) and nm1j00 (previous time step).
- Parameters
-
[in] | geom | Geometry object for the slice |
[in] | islice | slice index |
[in] | comp | laser component to initialize |
◆ InitSliceEnvelope()
void MultiLaser::InitSliceEnvelope |
( |
const int |
islice, |
|
|
const int |
comp |
|
) |
| |
Initialize on slice of the 3D laser field.
- Parameters
-
[in] | islice | slice index, referring to the 3D slice |
[in] | comp | laser component to initialize |
◆ InSituComputeDiags()
void MultiLaser::InSituComputeDiags |
( |
int |
step, |
|
|
amrex::Real |
time, |
|
|
int |
islice, |
|
|
const amrex::Geometry & |
geom3D, |
|
|
int |
max_step, |
|
|
amrex::Real |
max_time |
|
) |
| |
Compute in-situ laser diagnostics of current slice, store in member variable
- Parameters
-
[in] | step | current time step |
[in] | time | physical time |
[in] | islice | current slice, on which diags are computed. |
[in] | geom3D | Geometry of the problem |
[in] | max_step | maximum time step of simulation |
[in] | max_time | maximum time of simulation |
◆ InSituWriteToFile()
void MultiLaser::InSituWriteToFile |
( |
int |
step, |
|
|
amrex::Real |
time, |
|
|
const amrex::Geometry & |
geom3D, |
|
|
int |
max_step, |
|
|
amrex::Real |
max_time |
|
) |
| |
Dump in-situ reduced diagnostics to file.
- Parameters
-
[in] | step | current time step |
[in] | time | physical time |
[in] | geom3D | Geometry object for the whole domain |
[in] | max_step | maximum time step of simulation |
[in] | max_time | maximum time of simulation |
◆ ReadParameters()
void MultiLaser::ReadParameters |
( |
| ) |
|
◆ ShiftLaserSlices()
void MultiLaser::ShiftLaserSlices |
( |
| ) |
|
◆ m_all_lasers
◆ m_F_input_file
full 3D laser data stored on the host
◆ m_file_envelope_name
std::string MultiLaser::m_file_envelope_name = "laserEnvelope" |
|
private |
name of the openPMD species in the file
◆ m_file_geometry
std::string MultiLaser::m_file_geometry = "" |
|
private |
Geometry of the laser file, 'rt' or 'xyt'
◆ m_file_num_iteration
int MultiLaser::m_file_num_iteration = 0 |
|
private |
index of the iteration in the openPMD file
◆ m_input_file_path
std::string MultiLaser::m_input_file_path |
|
private |
path to input openPMD file
◆ m_insitu_cdata
All per-slice complex laser properties
◆ m_insitu_file_prefix
std::string MultiLaser::m_insitu_file_prefix = "diags/laser_insitu" |
|
private |
Prefix/path for the output files
◆ m_insitu_ncp
constexpr int MultiLaser::m_insitu_ncp = 1 |
|
staticconstexprprivate |
Number of real complex properties for in-situ per-slice reduced diagnostics.
◆ m_insitu_nrp
constexpr int MultiLaser::m_insitu_nrp = 6 |
|
staticconstexprprivate |
Number of real laser properties for in-situ per-slice reduced diagnostics.
◆ m_insitu_period
int MultiLaser::m_insitu_period {0} |
|
private |
How often the insitu laser diagnostics should be computed and written Default is 0, meaning no output
◆ m_insitu_rdata
All per-slice real laser properties
◆ m_insitu_sum_rdata
Sum of all per-slice real laser properties
◆ m_lambda0
amrex::Real MultiLaser::m_lambda0 {0.} |
|
private |
Laser central wavelength. he central wavelength influences the solver. As long as all the lasers are on the same grid (part of MultiLaser), this must be a property of MultiLaser.
◆ m_laser_from_file
bool MultiLaser::m_laser_from_file = false |
|
private |
if the lasers are initialized from openPMD file
◆ m_laser_geom_3D
◆ m_mg
hpmg solver for the envelope solver
◆ m_MG_average_rhs
bool MultiLaser::m_MG_average_rhs = true |
|
private |
Whether to use time-averaged RHS in envelope solver.
◆ m_MG_tolerance_abs
amrex::Real MultiLaser::m_MG_tolerance_abs = 0. |
|
private |
◆ m_MG_tolerance_rel
amrex::Real MultiLaser::m_MG_tolerance_rel = 1.e-4 |
|
private |
◆ m_MG_verbose
int MultiLaser::m_MG_verbose = 0 |
|
private |
◆ m_names
◆ m_nlasers
int MultiLaser::m_nlasers |
|
private |
◆ m_plan_bkw
FFTW plan for backward C2C transform to solve Complex Poisson equation
◆ m_plan_fwd
FFTW plan for forward C2C transform to solve Complex Poisson equation
◆ m_rhs
Complex FAB to store the RHS in position space
◆ m_rhs_fourier
Complex FAB to store the RHS in Fourier space
◆ m_slice_box
◆ m_slices
Array of N slices required to compute current slice
◆ m_slices_nguards
Number of guard cells for slices MultiFab
◆ m_sol
Complex FAB to store the solution (e.g. laser envelope on current slice)
◆ m_solver_type
std::string MultiLaser::m_solver_type = "multigrid" |
|
private |
◆ m_use_laser
bool MultiLaser::m_use_laser {false} |
whether a laser is used or not
◆ m_use_phase
bool MultiLaser::m_use_phase {true} |
|
private |
The documentation for this class was generated from the following files:
- /home/docs/checkouts/readthedocs.org/user_builds/hipace/checkouts/stable/src/laser/MultiLaser.H
- /home/docs/checkouts/readthedocs.org/user_builds/hipace/checkouts/stable/src/laser/MultiLaser.cpp