inplane_operators_m Module

Operators acting within poloidal plane (~perpendicular)


Uses

Used by


Contents


Interfaces

interface

  • public pure module function ddx(self, u, l)

    Evaluates derivative in x direction

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    integer, intent(in) :: l

    Index, for which x-derivative shall be computed

    Return Value real(kind=gp)

interface

  • public pure module function ddy(self, u, l)

    Evaluates derivative in y direction

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    integer, intent(in) :: l

    Index, for which y-derivative shall be computed

    Return Value real(kind=gp)

interface

  • public module function nabla_pol_sqrd(self, u, l, co)

    Evaluates J^(-1) d/dx(J co du/dx u) + d/dy(co du/dy) where J is the Jacobian

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    integer, intent(in) :: l

    Index, for which nabla_pol_sqrd is be computed

    real(kind=GP), intent(in), optional, dimension(self%n_points) :: co

    Coefficient, if not present it is assumed 1

    Return Value real(kind=gp)

interface

  • public pure module function arakawa(self, u, v, l)

    Evaluates Arakawa bracket [u,v] = du/dxdv/dy - du/dydv/dx second order, see [A. Arakawa, J. Comput. Phys. 135, 103 (1997)]

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    real(kind=GP), intent(in), dimension(self%n_points) :: v

    Values of variable v

    integer, intent(in) :: l

    Index, for which Arakawa Bracket is be computed

    Return Value real(kind=gp)

interface

  • public pure module function doubled_arakawa(self, u, v, l)

    Evaluates doubled Arakawa bracket [u,[u,v]], second order

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    real(kind=GP), intent(in), dimension(self%n_points) :: v

    Values of variable v

    integer, intent(in) :: l

    Index, for which double Arakawa bracket is be computed

    Return Value real(kind=gp)

interface

  • public module function curvature_interface(self, u, l)

    Interface for (squared) curvature operator

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    integer, intent(in) :: l

    Index, for which curvature is computed

    Return Value real(kind=gp)

interface

  • public module function curvature_btor_based(self, u, l)

    Curvature operator Computation based on toroidal magnetic field

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    integer, intent(in) :: l

    Index, for which curvature is computed

    Return Value real(kind=gp)

interface

  • public pure module function curvature_ddy_based(self, u, l)

    Curvature operator Simple computation via vertical derivative -2 d/dy Useful for slab/circular geometry

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    integer, intent(in) :: l

    Index, for which curvature is computed

    Return Value real(kind=gp)

interface

  • public module function curvature_sqrd_btor_based(self, u, l)

    Evaluates squared curvature operator Computation based on toroidal magnetic field

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    integer, intent(in) :: l

    Index, for which curvature is computed

    Return Value real(kind=gp)

interface

  • public pure module function curvature_sqrd_ddy_based(self, u, l)

    Evaluates squared curvature operator Simple computation via vertical derivative 4 d^2/dy^2 second order, central finite difference

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    integer, intent(in) :: l

    Index, for which curvature is computed

    Return Value real(kind=gp)

interface

  • public module function polarisation_advection(self, dia_fac, co, adv_pot, potv, tiv, nev, l)

    Evaluates non-linear advection div(co [adv_pot, v]_{x,y}) where v = grad(pot) + diafac*grad(pi) / ne

    Arguments

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

    Instance of the type

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

    Factor in front of diamagnetic term

    real(kind=GP), intent(in), dimension(self%n_points) :: co

    Values of coefficient

    real(kind=GP), intent(in), dimension(self%n_points) :: adv_pot

    Values of the first input in bracket [adv_pot, v]

    real(kind=GP), intent(in), dimension(self%n_points) :: potv

    Values of electrostatic potential

    real(kind=GP), intent(in), dimension(self%n_points) :: tiv

    Values of ion temperature

    real(kind=GP), intent(in), dimension(self%n_points) :: nev

    Values of density

    integer, intent(in) :: l

    Index, for which term is computed

    Return Value real(kind=gp)

interface

  • public pure module function arakawa_loc(self, zeta, psi)

    Evaluates arakawa bracket on local 3x3 stencil without its prefactor 1/(12 hf^2)

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(-1:1,-1:1) :: zeta

    First entry in arakawa bracket given as local 3x3 stencil

    real(kind=GP), intent(in), dimension(-1:1,-1:1) :: psi

    First entry in arakawa bracket given as local 3x3 stencil

    Return Value real(kind=gp)

interface

  • public pure module function doubled_arakawa_loc(self, u_loc, v_loc, jac_loc)

    Evaluates arakawa bracket on local 5x5 stencil If jac_loc not present: [u,[u,v]] If jac_loc present: 1/J[Ju, 1/J[Ju,v]]

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(-2:2,-2:2) :: u_loc

    Values of variable u given as local 5x5 stencil

    real(kind=GP), intent(in), dimension(-2:2,-2:2) :: v_loc

    Values of variable v given as local 5x5 stencil

    real(kind=GP), intent(in), optional, dimension(-2:2,-2:2) :: jac_loc

    Jacobian given as local 5x5 stencil

    Return Value real(kind=gp)

interface

  • public pure module function flutter_grad(self, apar, u, l)

    Evaluates flutter gradient of f using Arakawa bracket co*[apar,u]/btor

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: apar

    Values of variable apar

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    integer, intent(in) :: l

    Index, for which flutter_grad is be computed

    Return Value real(kind=gp)

interface

  • public pure module function flutter_div(self, apar, u, l)

    Evaluates flutter divergence of f using Arakawa bracket co*[apar,u]/btor

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: apar

    Values of variable apar

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    integer, intent(in) :: l

    Index, for which flutter_div is be computed

    Return Value real(kind=gp)

interface

  • public module subroutine init_inplane_operators(self, equi, mesh, equi_on_mesh, norm_scale, norm_flutter)

    Initialises inplane operators

    Arguments

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

    Instance of the type

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

    Equilibrium

    type(mesh_cart_t), intent(inout), target :: mesh

    Mesh

    class(equilibrium_storage_t), intent(inout), target :: equi_on_mesh

    Equilibrim quantities on mesh

    real(kind=GP), intent(in) :: norm_scale
    real(kind=GP), intent(in) :: norm_flutter

interface

  • public module subroutine apply_dissipation_inplane(self, u, ndiss, du, co)

    Applies dissipation of order ndiss: du = div (co grad (nabla_pol_sqrd^{ndiss-1} u) )

    Arguments

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

    Instance of the type

    real(kind=GP), intent(in), dimension(self%n_points) :: u

    Values of variable u

    integer, intent(in) :: ndiss

    Order of dissipation operator

    real(kind=GP), intent(out), dimension(self%n_points_inner) :: du

    Change added with dissipation

    real(kind=GP), intent(in), optional, dimension(self%n_points) :: co

    Coefficient applied in last Laplacian operation

interface

  • public module subroutine destructor(self)

    Frees memory associated with karniadakis

    Arguments

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

    Instance of the type


Derived Types

type, public ::  inplane_operators_t

Datatype for inplane operators

Components

Type Visibility Attributes Name Initial
procedure(curvature_interface), public, pointer :: curvature => null()

Curvature operator

procedure(curvature_interface), public, pointer :: curvature_sqrd => null()

Curvature operator

Finalizations Procedures

final :: destructor

Type-Bound Procedures

procedure , public :: init => init_inplane_operators Interface
procedure , public :: ddx Interface
procedure , public :: ddy Interface
procedure , public :: nabla_pol_sqrd Interface
procedure , public :: arakawa Interface
procedure , public :: doubled_arakawa Interface
procedure , public :: curvature_btor_based Interface
procedure , public :: curvature_ddy_based Interface
procedure , public :: curvature_sqrd_btor_based Interface
procedure , public :: curvature_sqrd_ddy_based Interface
procedure , public :: polarisation_advection Interface
procedure , public :: apply_dissipation_inplane Interface
procedure , public :: flutter_grad Interface
procedure , public :: flutter_div Interface