ohm_electromagnetic_braginskii_m Module

This module solves the electromagnetic and inertial part in Ohm's law

\f[ \beta\frac{\partial}{\partial t}A_\parallel+\frac{\partial}{\partial t}\frac{j_\parallel}{n}=rhs, \f] with Ampere's law \f[ \nabla_\perp^2A_\parallel = -j_\parallel \f]


Uses


Contents


Subroutines

public subroutine ohm_advance(comm_handler, equi, mesh_stag, hsolver_stag, map, penalisation_stag, boundaries, tstep_ohm, nev_stag, nev_adv_stag, tev_adv_stag, aparv, jparv, dohm, apar_adv, jpar_adv, sinfo, res)

Advances Ohm's law in time and solves for parallel elctromagnetic potential and parallel current

Arguments

Type IntentOptional Attributes Name
type(comm_handler_t), intent(in) :: comm_handler

Communicators

class(equilibrium_t), intent(inout) :: equi

Equilibrium

type(mesh_cart_t), intent(in) :: mesh_stag

Mesh (staggered)

class(helmholtz_solver_t), intent(inout) :: hsolver_stag

Elliptic (2D) solver on staggered mesh

type(parallel_map_t), intent(in) :: map

Parallel map

type(penalisation_t), intent(in) :: penalisation_stag

Penalisation (staggered)

type(boundaries_braginskii_t), intent(inout) :: boundaries

Boundary information for the BRAGINSKII model

type(karniadakis_t), intent(inout) :: tstep_ohm

Time-step integrator for Ohm's law

real(kind=GP), intent(in), dimension(mesh_stag%get_n_points()) :: nev_stag

Values of density on staggered grid at time t

real(kind=GP), intent(in), dimension(mesh_stag%get_n_points()) :: nev_adv_stag

Values of density on staggered grid at time t+1

real(kind=GP), intent(in), dimension(mesh_stag%get_n_points()) :: tev_adv_stag

Values of electron temperature on staggered grid at time t+1

real(kind=GP), intent(in), dimension(mesh_stag%get_n_points()) :: aparv

Values of parallel electromagnetic potential at time t

real(kind=GP), intent(in), dimension(mesh_stag%get_n_points()) :: jparv

Values of parallel current at time t

real(kind=GP), intent(in), dimension(mesh_stag%get_n_points()) :: dohm

Right hand side of Ohm's law

real(kind=GP), intent(inout), dimension(mesh_stag%get_n_points()) :: apar_adv

Values of parallel electromagnetic potential at time t+1, on input values for initial guess

real(kind=GP), intent(out), dimension(mesh_stag%get_n_points()) :: jpar_adv

Values of parallel current at time t+1

integer, intent(out) :: sinfo

Info from solver

real(kind=GP), intent(out) :: res

Residual of solution

public subroutine compute_jpar(equi, mesh_stag, boundaries, aparv, jparv)

Computes parallel current from parallel electromagnetic potential \f[ j_\parallel = -\nabla_\perp^2 A_\parallel \f]

Arguments

Type IntentOptional Attributes Name
class(equilibrium_t), intent(inout) :: equi

Equilibrium

type(mesh_cart_t), intent(in) :: mesh_stag

Mesh within poloidal plane

type(boundaries_braginskii_t), intent(inout) :: boundaries
real(kind=GP), intent(in), dimension(mesh_stag%get_n_points()) :: aparv

Values of parallel electromagnetic potential

real(kind=GP), intent(out), dimension(mesh_stag%get_n_points()) :: jparv

Values of parallel current

public subroutine apar_solve(comm_handler, equi, mesh_stag, hsolver_stag, boundaries, psiparv, aparv, sinfo, res, co_inert, co_te, penalisation_stag, tstep_ohm)

Solves for parallel electromagnetic potential if co_inert presents, \f[ \beta A_\parallel -\frac{mu}{n}\nabla_\perp^2A_\parallel= \psi_\parallel \f] otherwise, \f[ -\nabla_\perp^2A_\parallel= \psi_\parallel \f]

Arguments

Type IntentOptional Attributes Name
type(comm_handler_t), intent(in) :: comm_handler

Communicators

class(equilibrium_t), intent(in) :: equi

Equilibrium

type(mesh_cart_t), intent(in) :: mesh_stag

Mesh within poloidal plane

class(helmholtz_solver_t), intent(inout) :: hsolver_stag

Elliptic (2D) solver on staggered mesh

type(boundaries_braginskii_t), intent(inout) :: boundaries

Boundary information for the BRAGINSKII model

real(kind=GP), intent(in), dimension(mesh_stag%get_n_points()) :: psiparv

Values for generalised electromagnetic potential

real(kind=GP), intent(inout), dimension(mesh_stag%get_n_points()) :: aparv

Values for parallel electromagnetic potential, on input initial guess

integer, intent(out) :: sinfo

Info from solver

real(kind=GP), intent(out) :: res

Residual of solution

real(kind=GP), intent(in), optional, dimension(mesh_stag%get_n_points()) :: co_inert

Coefficient in inertial term (= Values of density at time t+1 on staggered grid)

real(kind=GP), intent(in), optional, dimension(mesh_stag%get_n_points()) :: co_te

Coefficient in implicit resistivity (= Values of electron temperature at time t+1 on staggered grid)

type(penalisation_t), intent(in), optional :: penalisation_stag

Penalisation (staggered)

type(karniadakis_t), intent(inout), optional :: tstep_ohm

Time-step integrator for Ohm's law