Quasistatic2DWakefield#

class wake_t.physics_models.plasma_wakefields.Quasistatic2DWakefield(density_function, r_max, xi_min, xi_max, n_r, n_xi, ppc=2, dz_fields=None, r_max_plasma=None, parabolic_coefficient=0.0, p_shape='cubic', max_gamma=10, plasma_pusher='rk4', laser=None, laser_evolution=True, laser_envelope_substeps=1, laser_envelope_nxi=None, laser_envelope_nr=None, laser_envelope_use_phase=True)[source]#

This class calculates the plasma wakefields using the gridless quasi-static model in r-z geometry originally developed by P. Baxevanis and G. Stupakov [1].

The model implemented here includes additional features with respect to the original version in [1]. Among them is the support for laser drivers, particle beams (instead of an analytic charge distribution), non-uniform and finite plasma density profiles, as well as an Adams-Bashforth pusher for the plasma particles (in addition the original Runke-Kutta pusher).

As a kinetic quasi-static model in r-z geometry, it computes the plasma response by evolving a 1D radial plasma column from the front to the back of the simulation domain. A special feature of this model is that it does not need a grid in order to compute this evolution, and it allows the fields E and B to be calculated at any radial position in the plasma column.

In the Wake-T implementation, a grid is used only in order to be able to easily interpolate the fields to the beam particles. After evolving the plasma column, E and B are calculated at the locations of the grid nodes. Similarly, the charge density rho and susceptibility chi of the plasma are computed by depositing the charge of the plasma particles on the same grid. This useful for diagnostics and for evolving a laser pulse.

Parameters:
density_functioncallable

Function that returns the density value at the given position z. This parameter is given by the PlasmaStage and does not need to be specified by the user.

r_maxfloat

Maximum radial position up to which plasma wakefield will be calculated.

xi_minfloat

Minimum longitudinal (speed of light frame) position up to which plasma wakefield will be calculated.

xi_maxfloat

Maximum longitudinal (speed of light frame) position up to which plasma wakefield will be calculated.

n_rint

Number of grid elements along r to calculate the wakefields.

n_xiint

Number of grid elements along xi to calculate the wakefields.

ppcint, optional

Number of plasma particles per radial cell. By default ppc=2.

dz_fieldsfloat, optional

Determines how often the plasma wakefields should be updated. For example, if dz_fields=10e-6, the plasma wakefields are only updated every time the simulation window advances by 10 micron. By default dz_fields=xi_max-xi_min, i.e., the length the simulation box.

r_max_plasmafloat, optional

Maximum radial extension of the plasma column. If None, the plasma extends up to the r_max boundary of the simulation box.

parabolic_coefficientfloat or callable, optional

The coefficient for the transverse parabolic density profile. The radial density distribution is calculated as n_r = n_p * (1 + parabolic_coefficient * r**2), where n_p is the local on-axis plasma density. If a float is provided, the same value will be used throwout the stage. Alternatively, a function which returns the value of the coefficient at the given position z (e.g. def func(z)) might also be provided.

p_shapestr, optional

Particle shape to be used for the beam charge deposition. Possible values are 'linear' or 'cubic' (default).

max_gammafloat, optional

Plasma particles whose gamma exceeds max_gamma are considered to violate the quasistatic condition and are put at rest (i.e., gamma=1., pr=pz=0.). By default max_gamma=10.

plasma_pusherstr, optional

The pusher used to evolve the plasma particles. Possible values are 'rk4' (Runge-Kutta 4th order) or 'ab5' (Adams-Bashforth 5th order).

laserLaserPulse, optional

Laser driver of the plasma stage.

laser_evolutionbool, optional

If True (default), the laser pulse is evolved using a laser envelope model. If False, the pulse envelope stays unchanged throughout the computation.

laser_envelope_substepsint, optional

Number of substeps of the laser envelope solver per dz_fields. The time step of the envelope solver is therefore dz_fields / c / laser_envelope_substeps.

laser_envelope_nxi, laser_envelope_nrint, optional

If given, the laser envelope will run in a grid of size (laser_envelope_nxi, laser_envelope_nr) instead of (n_xi, n_r). This allows the laser to run in a finer (or coarser) grid than the plasma wake. It is not necessary to specify both parameters. If one of them is not given, the resolution of the plasma grid with be used for that direction.

laser_envelope_use_phasebool, optional

Determines whether to take into account the terms related to the longitudinal derivative of the complex phase in the envelope solver.

References

[1] (1,2)

P. Baxevanis and G. Stupakov, “Novel fast simulation technique for axisymmetric plasma wakefield acceleration configurations in the blowout regime,” Phys. Rev. Accel. Beams 21, 071301 (2018), https://link.aps.org/doi/10.1103/PhysRevAccelBeams.21.071301

Methods

adjust_dt(t_final)

Autoadjust the time step of the field update.

calculate_field(bunches)

Calculate field using the current properties and given bunches.

evolve_properties(bunches)

Evolve field properties.

gather(x, y, z, t, ex, ey, ez, bx, by, bz, ...)

Gather all field components at the specified locations.

get_openpmd_diagnostics_data(global_time)

Get the data for including the field in the openPMD diagnostics.

initialize_properties(bunches)

Initialize field properties.

update(bunches)

Update field to the next time step.