electron_kinetic_equation

moment_kinetics.electron_kinetic_equation.add_wall_boundary_condition_to_Jacobian!Function
add_wall_boundary_condition_to_Jacobian!(jacobian, phi, pdf, p, 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).

source
moment_kinetics.electron_kinetic_equation.calculate_contribution_from_z_advection!Method

calculates 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

source
moment_kinetics.electron_kinetic_equation.electron_backward_euler!Method
electron_backward_euler!(old_scratch, new_scratch, moments, phi,
    collisions, composition, r, z, vperp, vpa, z_spectral, vperp_spectral,
    vpa_spectral, z_advect, vpa_advect, scratch_dummy, t_params,
    external_source_settings, num_diss_params, nl_solver_params, ir;
    evolve_p=false, ion_dt=nothing) = begin

Take a single backward euler timestep for the electron shape function $g_e$ and parallel pressure $p_{e∥}$.

source
moment_kinetics.electron_kinetic_equation.electron_backward_euler_pseudotimestepping!Method

Update 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_pseudotimestepping!()). t_params.dt[] is adapted according to the iteration counts of the Newton solver.

source
moment_kinetics.electron_kinetic_equation.electron_kinetic_equation_euler_update!Method
electron_kinetic_equation_euler_update!(result_object, f_in, p_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_p=false,
                                        ion_dt=nothing, soft_force_constraints=false,
                                        debug_io=nothing, fields=nothing, r=nothing,
                                        vzeta=nothing, vr=nothing, vz=nothing,
                                        istage=0)

Do a forward-Euler update of the electron kinetic equation.

When evolve_p=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, p_out, f_in, and p_in should have no r-dimension.

result_object should be either a scratch_pdf object (which contains data for all r-indices) or a Tuple containing (pdf_electron, electron_p) fields for r-index ir only. This allows result_object to be (possibly) passed to write_debug_data_to_binary() when result_object is a scratch_pdf.

fields, r, vzeta, vr, vz, and istage are only required when a non-nothing debug_io is passed.

source
moment_kinetics.electron_kinetic_equation.fill_electron_kinetic_equation_Jacobian!Function
fill_electron_kinetic_equation_Jacobian!(jacobian_matrix, f, p, 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_p, include=:all)

Fill a pre-allocated matrix with the Jacobian matrix for electron kinetic equation and (if evolve_p=true) the electron energy equation.

add_identity=false can be passed to skip adding 1's on the diagonal. Doing this would produce the Jacobian for a steady-state solve, rather than a backward-Euler timestep solve. [Note that rows representing boundary points still need the 1 added to their diagonal, as that 1 is used to impose the boundary condition, not to represent the 'new f' in the time derivative term as it is for the non-boundary points.]

source
moment_kinetics.electron_kinetic_equation.fill_electron_kinetic_equation_v_only_Jacobian!Method
fill_electron_kinetic_equation_v_only_Jacobian!()
    jacobian_matrix, f, p, 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_p)

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_p=true) the electron energy equation.

source
moment_kinetics.electron_kinetic_equation.fill_electron_kinetic_equation_z_only_Jacobian_f!Method
fill_electron_kinetic_equation_z_only_Jacobian_f!(
    jacobian_matrix, f, p, 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_p)

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_p=true) the electron energy equation.

source
moment_kinetics.electron_kinetic_equation.fill_electron_kinetic_equation_z_only_Jacobian_p!Method
fill_electron_kinetic_equation_z_only_Jacobian_p!(
    jacobian_matrix, p, 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_p)

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_p=true) the electron energy equation.

source
moment_kinetics.electron_kinetic_equation.get_electron_split_Jacobians!Method
get_electron_split_Jacobians!(ivperp, ivpa, p, 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_p)

Fill a pre-allocated matrix with the Jacobian matrix for electron kinetic equation and (if evolve_p=true) the electron energy equation.

source
moment_kinetics.electron_kinetic_equation.implicit_electron_advance!Method
implicit_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.

source
moment_kinetics.electron_kinetic_equation.update_electron_pdf!Method

updateelectronpdf 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
moments = struct containing electron moments
r = struct containing r-coordinate information
z = struct containing z-coordinate information
vperp = struct containing vperp-coordinate information
vpa = struct containing vpa-coordinate information
z_spectral = struct containing spectral information for the z-coordinate
vperp_spectral = struct containing spectral information for the vperp-coordinate
vpa_spectral = struct containing spectral information for the vpa-coordinate
z_advect = struct containing information for z-advection
vpa_advect = struct containing information for vpa-advection
scratch_dummy = dummy arrays to be used for temporary storage
t_params = parameters, etc. for time-stepping
collisions = parameters for charged species collisions and neutral reactions
composition = parameters describing number and type of species
max_electron_pdf_iterations = maximum number of iterations to use in the solution of the electron kinetic equation
max_electron_sim_time = maximum amount of de-dimensionalised time to use in the solution of the electron kinetic equation
io_electorn = if this is passed, write some output from the electron pseudo-timestepping
initial_time = if this is passed set the initial pseudo-time to this value
residual_tolerance = if this is passed, it overrides the value of `t_params.converged_residual_value`
evolve_p = if set to true, update the electron pressure as part of the solve
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.
solution_method = if this is passed, choose a non-standard method for the solution

OUTPUT: pdf = updated (modified) electron pdf

source
moment_kinetics.electron_kinetic_equation.update_electron_pdf_with_picard_iteration!Method

use 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 p = electron 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

source
moment_kinetics.electron_kinetic_equation.update_electron_pdf_with_shooting_method!Method

updateelectronpdfwithshooting_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
p = electron 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

source
moment_kinetics.electron_kinetic_equation.update_electron_pdf_with_time_advance!Method

updateelectronpdfwithtime_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:
scratch = `scratch_pdf` struct used to store Runge-Kutta stages
pdf = modified electron pdf @ previous time level = (true electron pdf / dens_e) * vth_e
moments = struct containing electron moments
r = struct containing r-coordinate information
z = struct containing z-coordinate information
vperp = struct containing vperp-coordinate information
vpa = struct containing vpa-coordinate information
z_spectral = struct containing spectral information for the z-coordinate
vperp_spectral = struct containing spectral information for the vperp-coordinate
vpa_spectral = struct containing spectral information for the vpa-coordinate
z_advect = struct containing information for z-advection
vpa_advect = struct containing information for vpa-advection
scratch_dummy = dummy arrays to be used for temporary storage
t_params = parameters, etc. for time-stepping
collisions = parameters for charged species collisions and neutral reactions
composition = parameters describing number and type of species
max_electron_pdf_iterations = maximum number of iterations to use in the solution of the electron kinetic equation
max_electron_sim_time = maximum amount of de-dimensionalised time to use in the solution of the electron kinetic equation
io_electorn = if this is passed, write some output from the electron pseudo-timestepping
initial_time = if this is passed set the initial pseudo-time to this value
residual_tolerance = if this is passed, it overrides the value of `t_params.converged_residual_value`
evolve_p = if set to true, update the electron pressure as part of the solve
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

source