gyroviscosity_m Module

Computations related with gyroviscosity



Contents


Derived Types

type, public ::  gyroviscosity_t

Holds information and performs computations related with gyroviscosity

Components

Type Visibility Attributes Name Initial
real(kind=GP), public, allocatable, dimension(:) :: eta_flow_v

Values of flow viscosity (on canonical mesh)

Finalizations Procedures

final :: destructor

Type-Bound Procedures

procedure , public :: init => init_gyroviscosity Subroutine
procedure , public :: drive => drive_gyroviscosity Subroutine
procedure , public :: eval_gstress Subroutine
procedure , public :: eval_curvature Subroutine
procedure , public :: eval_pargrad_b32 Subroutine

Functions

public pure function eta_neo_flow(ne_l, ti_l) result(res)

Computes ion flow viscosity with "beyond Braginskii extensions" according to parameters

Arguments

Type IntentOptional Attributes Name
real(kind=GP), intent(in) :: ne_l

Density value

real(kind=GP), intent(in) :: ti_l

Ion temperature value

Return Value real(kind=GP)


Subroutines

public subroutine init_gyroviscosity(self, comm_handler, mesh_cano, mesh_stag)

Initialise gyroviscosity

Arguments

Type IntentOptional Attributes Name
class(gyroviscosity_t), intent(inout) :: self

Instance of type

type(comm_handler_t), intent(in) :: comm_handler

Communicators

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

Mesh (canonical) within poloidal plane

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

Mesh (staggered) within poloidal plane

public subroutine drive_gyroviscosity(self, comm_handler, mesh_cano, mesh_stag, equi_on_stag, ne_v, ti_v, upar_v, qipar_v)

Performs auxiliary precomputations necessary to evaluate guroviscosity - Precomputes viscosities - Communicates uB^3/2 and qiparB^(3/2)

Arguments

Type IntentOptional Attributes Name
class(gyroviscosity_t), intent(inout) :: self

Instance of type

type(comm_handler_t), intent(in) :: comm_handler

Communicators

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

Mesh (canonical) within poloidal plane

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

Mesh (staggered) within poloidal plane

class(equilibrium_storage_t), intent(in) :: equi_on_stag

Equilibrim quantities on staggered mesh

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: ne_v

Values of density

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: ti_v

Values of ion temperature

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

Values of parallel ion velocity

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

Values of parallel ion heat flux, for comoputation of heat-flux related viscosity

public subroutine destructor(self)

Destructor

Arguments

Type IntentOptional Attributes Name
type(gyroviscosity_t), intent(inout) :: self

Instance of type

public subroutine eval_gstress(self, mesh_cano, equi_on_cano, opsinplane_cano, map, ne_v, pot_v, ti_v, apar_fluct_cano_v, gstress, eval_flow_crv, eval_flow_par, eval_heat_crv, eval_heat_par)

Computation of ion stress function G = G_flow + G_heat with "beyond Braginskii extensions" according to parameters

Arguments

Type IntentOptional Attributes Name
class(gyroviscosity_t), intent(in) :: self

Instance of type

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

Mesh (canonical) within poloidal plane

class(equilibrium_storage_t), intent(in) :: equi_on_cano

Equilibrim quantities on canonical mesh

type(inplane_operators_t), intent(in) :: opsinplane_cano

In-plane operators (canonical)

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

Parallel map

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: ne_v

Electron density

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: pot_v

Electrostatic potential

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: ti_v

Ion temperature

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: apar_fluct_cano_v

Fluctuation of apar used for flutter operators It is assumed that communication has already been done to do mapping to canonical mesh

real(kind=GP), intent(out), dimension(mesh_cano%get_n_points()) :: gstress

Values of ion viscous stress function G

logical, intent(in), optional :: eval_flow_crv

If true computes part of G_flow arising from curvature terms (default = T) This feature is useful for technical reasons, when only evaluation of individual parts of G is needed

logical, intent(in), optional :: eval_flow_par

If true computes computes part of G_flow arising from parallel divergence of upar (default = T)

logical, intent(in), optional :: eval_heat_crv

If true computes part of G_heat arising from curvature terms (default = T)

logical, intent(in), optional :: eval_heat_par

If true computes part of G_flow arising from paralle divergence of qipar (default = T)

public subroutine eval_curvature(self, mesh_cano, equi_on_cano, opsinplane_cano, map, ne_v, pot_v, ti_v, apar_fluct_cano_v, crv_gstress)

Evaluates curvature of ion stress function. This routine is present, because there is an issue if one simply evaluated the curvature of gstress obtained via ion_stress_function. Squared curvature operators appear that would imply chequerboarding. Therefore the evaluation of squared curvature operators has to be performed numerically careful and explicit.

Arguments

Type IntentOptional Attributes Name
class(gyroviscosity_t), intent(in) :: self

Instance of type

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

Mesh (canonical) within poloidal plane

class(equilibrium_storage_t), intent(in) :: equi_on_cano

Equilibrim quantities on canonical mesh

type(inplane_operators_t), intent(in) :: opsinplane_cano

In-plane operators (canonical)

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

Parallel map

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: ne_v

Electron density

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: pot_v

Electrostatic potential

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: ti_v

Ion temperature

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: apar_fluct_cano_v

Fluctuation of apar on canonical mesh used for flutter operators

real(kind=GP), intent(out), dimension(mesh_cano%get_n_points_inner()) :: crv_gstress

Values of curvature of ion viscous stress function G

public subroutine eval_pargrad_b32(self, comm_handler, mesh_cano, mesh_stag, equi_on_cano, equi_on_stag, opsinplane_cano, opsinplane_stag, map, ne_v, pot_v, ti_v, apar_fluct_v, apar_fluct_cano_v, pargrad_gstress, eval_flow_crv, eval_flow_par, eval_heat_crv, eval_heat_par)

Evaluates B^(3/2)par_grad(G/B^(3/2)) of ion stress function. Used for the parallel momentum equation

Arguments

Type IntentOptional Attributes Name
class(gyroviscosity_t), intent(in) :: self

Instance of type

type(comm_handler_t), intent(in) :: comm_handler

Communicators

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

Mesh (canonical) within poloidal plane

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

Mesh (staggered) within poloidal plane

class(equilibrium_storage_t), intent(in) :: equi_on_cano

Equilibrim quantities on canonical mesh

class(equilibrium_storage_t), intent(in) :: equi_on_stag

Equilibrim quantities on staggered mesh

type(inplane_operators_t), intent(in) :: opsinplane_cano

In-plane operators (canonical)

type(inplane_operators_t), intent(in) :: opsinplane_stag

In-plane operators (staggered)

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

Parallel map

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: ne_v

Electron density

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: pot_v

Electrostatic potential

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: ti_v

Ion temperature

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

Fluctuation of apar used for flutter operators

real(kind=GP), intent(in), dimension(mesh_cano%get_n_points()) :: apar_fluct_cano_v

Fluctuation of apar on canonical mesh used for flutter operators

real(kind=GP), intent(out), dimension(mesh_stag%get_n_points_inner()) :: pargrad_gstress

Values of ion viscous stress function G

logical, intent(in) :: eval_flow_crv

If true computes part of G_flow arising from curvature terms This feature is useful for technical reasons, when only evaluation of individual parts of G is needed

logical, intent(in) :: eval_flow_par

If true computes computes part of G_flow arising from parallel divergence of upar

logical, intent(in) :: eval_heat_crv

If true computes part of G_heat arising from curvature terms

logical, intent(in) :: eval_heat_par

If true computes computes part of G_flow arising from parallel divergence of qipar