electron_kinetic_equation
moment_kinetics.electron_kinetic_equation.add_wall_boundary_condition_to_Jacobian!
— Functionadd_wall_boundary_condition_to_Jacobian!(jacobian, phi, pdf, ppar, vthe, upar, z,
vperp, vpa, vperp_spectral, vpa_spectral,
vpa_adv, moments, vpa_diffusion, me_over_mi,
ir)
All the contributions that we add in this function have to be added with a -'ve sign so that they combine with the 1 on the diagonal of the preconditioner matrix to make rows corresponding to the boundary points which define constraint equations imposing the boundary condition on those entries of δg (when the right-hand-side is set to zero).
moment_kinetics.electron_kinetic_equation.calculate_contribution_from_z_advection!
— Methodcalculates the contribution to the residual of the electron kinetic equation from the z advection term: residual = zdot * d(pdf)/dz + wpadot * d(pdf)/dwpa - pdf * prefactor INPUTS: zadvectionterm = dummy array to be filled with the contribution to the residual from the z advection term pdf = modified electron pdf used in the kinetic equation = (true electron pdf / dense) * vthe vthe = electron thermal speed z = z grid vpa = vparallel grid zspectral = spectral representation of the z grid scratchdummy = dummy array to be used for temporary storage OUTPUT: zadvectionterm = updated contribution to the residual from the z advection term
moment_kinetics.electron_kinetic_equation.electron_adaptive_timestep_update!
— Methodelectron_adaptive_timestep_update!(scratch, t, t_params, moments, phi, z_advect,
vpa_advect, composition, r, z, vperp, vpa,
vperp_spectral, vpa_spectral,
external_source_settings, num_diss_params;
evolve_ppar=false)
Check the error estimate for the embedded RK method and adjust the timestep if appropriate.
moment_kinetics.electron_kinetic_equation.electron_backward_euler!
— MethodUpdate the electron distribution function using backward-Euler for an artifical time advance of the electron kinetic equation until a steady-state solution is reached.
Note that this function does not use the runge_kutta
timestep functionality. t_params.previous_dt[]
is used to store the (adaptively updated) initial timestep of the pseudotimestepping loop (initial value of t_params.dt[]
within electron_backward_euler!()
). t_params.dt[]
is adapted according to the iteration counts of the Newton solver.
moment_kinetics.electron_kinetic_equation.electron_kinetic_equation_euler_update!
— Methodelectron_kinetic_equation_euler_update!(f_out, ppar_out, f_in, ppar_in, moments,
z, vperp, vpa, z_spectral, vpa_spectral,
z_advect, vpa_advect, scratch_dummy,
collisions, composition,
external_source_settings,
num_diss_params, t_params, ir;
evolve_ppar=false, ion_dt=nothing)
Do a forward-Euler update of the electron kinetic equation.
When evolve_ppar=true
is passed, also updates the electron parallel pressure.
Note that this function operates on a single point in r
, given by ir
, and f_out
, ppar_out
, f_in
, and ppar_in
should have no r-dimension.
moment_kinetics.electron_kinetic_equation.fill_electron_kinetic_equation_Jacobian!
— Functionfill_electron_kinetic_equation_Jacobian!(jacobian_matrix, f, ppar, moments,
collisions, composition, z, vperp, vpa,
z_spectral, vperp_specral,
vpa_spectral, z_advect, vpa_advect,
scratch_dummy, external_source_settings,
num_diss_params, t_params, ion_dt,
ir, evolve_ppar, include=:all)
Fill a pre-allocated matrix with the Jacobian matrix for electron kinetic equation and (if evolve_ppar=true
) the electron energy equation.
moment_kinetics.electron_kinetic_equation.fill_electron_kinetic_equation_v_only_Jacobian!
— Methodfill_electron_kinetic_equation_v_only_Jacobian!()
jacobian_matrix, f, ppar, dpdf_dz, dpdf_dvpa, z_speed, moments, zeroth_moment,
first_moment, second_moment, third_moment, dthird_moment_dz, phi, collisions,
composition, z, vperp, vpa, z_spectral, vperp_spectral, vpa_spectral, z_advect,
vpa_advect, scratch_dummy, external_source_settings, num_diss_params, t_params,
ion_dt, ir, iz, evolve_ppar)
Fill a pre-allocated matrix with the Jacobian matrix for a velocity-space solve part of the ADI method for electron kinetic equation and (if evolve_ppar=true
) the electron energy equation.
moment_kinetics.electron_kinetic_equation.fill_electron_kinetic_equation_z_only_Jacobian_f!
— Methodfill_electron_kinetic_equation_z_only_Jacobian_f!(
jacobian_matrix, f, ppar, dpdf_dz, dpdf_dvpa, z_speed, moments, zeroth_moment,
first_moment, second_moment, third_moment, dthird_moment_dz, collisions,
composition, z, vperp, vpa, z_spectral, vperp_spectral, vpa_spectral, z_advect,
vpa_advect, scratch_dummy, external_source_settings, num_diss_params, t_params,
ion_dt, ir, ivperp, ivpa, evolve_ppar)
Fill a pre-allocated matrix with the Jacobian matrix for a z-direction solve part of the ADI method for electron kinetic equation and (if evolve_ppar=true
) the electron energy equation.
moment_kinetics.electron_kinetic_equation.fill_electron_kinetic_equation_z_only_Jacobian_ppar!
— Methodfill_electron_kinetic_equation_z_only_Jacobian_ppar!(
jacobian_matrix, ppar, moments, zeroth_moment, first_moment, second_moment,
third_moment, dthird_moment_dz, collisions, composition, z, vperp, vpa,
z_spectral, vperp_spectral, vpa_spectral, z_advect, vpa_advect, scratch_dummy,
external_source_settings, num_diss_params, t_params, ion_dt, ir, evolve_ppar)
Fill a pre-allocated matrix with the Jacobian matrix for a z-direction solve part of the ADI method for electron kinetic equation and (if evolve_ppar=true
) the electron energy equation.
moment_kinetics.electron_kinetic_equation.get_electron_split_Jacobians!
— Methodget_electron_split_Jacobians!(ivperp, ivpa, ppar, moments, collisions, composition,
z, vperp, vpa, z_spectral, vperp_spectral, vpa_spectral,
z_advect, vpa_advect, scratch_dummy,
external_source_settings, num_diss_params, t_params,
ion_dt, ir, evolve_ppar
Fill a pre-allocated matrix with the Jacobian matrix for electron kinetic equation and (if evolve_ppar=true
) the electron energy equation.
moment_kinetics.electron_kinetic_equation.implicit_electron_advance!
— Methodimplicit_electron_advance!()
Do an implicit solve which finds: the steady-state electron shape function $g_e$; the backward-Euler advanced electron pressure which is updated using $g_e$ at the new time-level.
The r-dimension is not parallelised. For 1D runs this makes no difference. In 2D it might or might not be necessary. If r-dimension parallelisation is needed, it would need some work. The simplest option would be a non-parallelised outer loop over r, with each nonlinear solve being parallelised over {z,vperp,vpa}. More efficient might be to add an equivalent to the 'anyv' parallelisation used for the collision operator (e.g. 'anyzv'?) to allow the outer r-loop to be parallelised.
moment_kinetics.electron_kinetic_equation.update_electron_pdf!
— Methodupdateelectronpdf is a function that uses the electron kinetic equation to solve for the updated electron pdf
The electron kinetic equation is: zdot * d(pdf)/dz + wpadot * d(pdf)/dwpa = pdf * pre_factor
INPUTS:
scratch = `scratch_pdf` struct used to store Runge-Kutta stages
pdf = modified electron pdf @ previous time level = (true electron pdf / dens_e) * vth_e
dens = electron density
vthe = electron thermal speed
ppar = electron parallel pressure
ddens_dz = z-derivative of the electron density
dppar_dz = z-derivative of the electron parallel pressure
dqpar_dz = z-derivative of the electron parallel heat flux
dvth_dz = z-derivative of the electron thermal speed
z = struct containing z-coordinate information
vpa = struct containing vpa-coordinate information
z_spectral = struct containing spectral information for the z-coordinate
vpa_spectral = struct containing spectral information for the vpa-coordinate
scratch_dummy = dummy arrays to be used for temporary storage
dt = time step size
max_electron_pdf_iterations = maximum number of iterations to use in the solution of the electron kinetic equation
ion_dt = if this is passed, the electron pressure is evolved in a form that results in
a backward-Euler update on the ion timestep (ion_dt) once the electron
pseudo-timestepping reaches steady state.
OUTPUT: pdf = updated (modified) electron pdf
moment_kinetics.electron_kinetic_equation.update_electron_pdf_with_picard_iteration!
— Methoduse Picard iteration to solve the electron kinetic equation
The electron kinetic equation is: zdot * d(pdf)/dz + wpadot * d(pdf)/dwpa = pdf * pre_factor Picard iteration uses the previous iteration of the electron pdf to calculate the next iteration: zdot * d(pdf^{i+1})/dz + wpadot^{i} * d(pdf^{i})/dwpa = pdf^{i} * prefactor^{i}
INPUTS: pdf = modified electron pdf @ previous time level = (true electron pdf / dense) * vthe dens = electron density vthe = electron thermal speed ppar = electron parallel pressure ddensdz = z-derivative of the electron density dppardz = z-derivative of the electron parallel pressure dqpardz = z-derivative of the electron parallel heat flux dvthdz = z-derivative of the electron thermal speed z = struct containing z-coordinate information vpa = struct containing vpa-coordinate information zspectral = struct containing spectral information for the z-coordinate vpaspectral = struct containing spectral information for the vpa-coordinate scratchdummy = dummy arrays to be used for temporary storage maxelectronpdfiterations = maximum number of iterations to use in the solution of the electron kinetic equation OUTPUT: pdf = updated (modified) electron pdf
moment_kinetics.electron_kinetic_equation.update_electron_pdf_with_shooting_method!
— Methodupdateelectronpdfwithshooting_method is a function that using a shooting method to solve for the electron pdf
The electron kinetic equation is: zdot * d(pdf)/dz + wpadot * d(pdf)/dwpa = pdf * prefactor The shooting method is 'explicit' in z, solving zdoti * (pdf{i+1} - pdf{i})/dz{i} + wpadot{i} * d(pdf{i})/dwpa = pdf{i} * prefactor{i}
INPUTS:
pdf = modified electron pdf @ previous time level = (true electron pdf / dens_e) * vth_e
dens = electron density
vthe = electron thermal speed
ppar = electron parallel pressure
ddens_dz = z-derivative of the electron density
dppar_dz = z-derivative of the electron parallel pressure
dqpar_dz = z-derivative of the electron parallel heat flux
dvth_dz = z-derivative of the electron thermal speed
z = struct containing z-coordinate information
vpa = struct containing vpa-coordinate information
vpa_spectral = struct containing spectral information for the vpa-coordinate
scratch_dummy = dummy arrays to be used for temporary storage
OUTPUT: pdf = updated (modified) electron pdf
moment_kinetics.electron_kinetic_equation.update_electron_pdf_with_time_advance!
— Methodupdateelectronpdfwithtime_advance is a function that introduces an artifical time derivative to advance the electron kinetic equation until a steady-state solution is reached.
The electron kinetic equation is: zdot * d(pdf)/dz + wpadot * d(pdf)/dwpa = pdf * pre_factor
INPUTS:
pdf = modified electron pdf @ previous time level = (true electron pdf / dens_e) * vth_e
dens = electron density
vthe = electron thermal speed
ppar = electron parallel pressure
ddens_dz = z-derivative of the electron density
dppar_dz = z-derivative of the electron parallel pressure
dqpar_dz = z-derivative of the electron parallel heat flux
dvth_dz = z-derivative of the electron thermal speed
z = struct containing z-coordinate information
vpa = struct containing vpa-coordinate information
z_spectral = struct containing spectral information for the z-coordinate
vpa_spectral = struct containing spectral information for the vpa-coordinate
scratch_dummy = dummy arrays to be used for temporary storage
max_electron_pdf_iterations = maximum number of iterations to use in the solution of the electron kinetic equation
io_electron = info struct for binary file I/O
initial_time = initial value for the (pseudo-)time
ion_dt = if this is passed, the electron pressure is evolved in a form that results in
a backward-Euler update on the ion timestep (ion_dt) once the electron
pseudo-timestepping reaches steady state.
OUTPUT: pdf = updated (modified) electron pdf