solver_aligned3d_m Module

Routines for 3d elliptic solver along magnetic field lines


Uses


Contents


Interfaces

interface

  • public module subroutine solve_aligned3d(comm_handler, equi_on_cano, equi_on_stag, mesh_cano, mesh_stag, map, params_solver, solver_type, rhs, sol, niter, res, true_res, cnt_matvec, lambda, xi, co_stag, fcx, bnd_descrs, flutter_fac_cano, flutter_fac_stag, opsinplane_cano, opsinplane_stag, apar_stag, apar_cano, bnd_descrs_flux)

    Solves an elliptic Helmholtz problem along magnetic field line, optionally along perturbed magnetic field line via magentic flutter \f[ \lambda * u - \xi*\nabla\cdot\left[co\nabla_\parallel u fcx\right] = rhs \f] with coefficients \f$\lambda, \xi, co, fcx $ used e.g. for implicit treatment of parallel heat conduction The problem is stated on the canonical mesh Perpendicular boundary points are set according to bnd_descrs Magnetic flutter is included in parallel gradient and divergence via Possion bracket and arakawa scheme Perpendicular boundary flux_bnd_descrs is specially reqired for dobule arakawa operation

    Arguments

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

    Communicators

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

    Equilibrim quantities on canonical mesh

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

    Equilibrim quantities on staggered mesh

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

    Canonical mesh

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

    Staggered mesh

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

    Parallel map

    type(params_solver_aligned3d_t), intent(in) :: params_solver

    Parameters for solver

    character(len=*), intent(in) :: solver_type

    Type of solver

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

    Right hand side

    real(kind=GP), intent(inout), dimension(mesh_cano%get_n_points()) :: sol

    On input: initial guess, on output :: solution

    integer, intent(out) :: niter

    Number of (outer) iterations, negative number if error occured

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

    Pseudo-residual (relative) of solution

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

    True residual (relative) of solution

    integer, intent(out), optional :: cnt_matvec

    Counts number that matvec was called (~total iterations)

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

    Coefficient lambda, default = 1.0

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

    Coefficient xi, default = 1.0

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

    Coefficient co, default = 1.0, defined on staggered mesh

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

    Coefficient fcx, default = 1.0

    integer, intent(in), optional, dimension(mesh_cano%get_n_points_boundary()) :: bnd_descrs

    Boundary descriptor for perpendicular boundary points

    real(kind=GP), intent(in), optional :: flutter_fac_cano

    Factor (~switch) in front of flutter divergence on canonical mesh

    real(kind=GP), intent(in), optional :: flutter_fac_stag

    Factor (~switch) in front of flutter gradient on staggered mesh

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

    Inplane operators on canonical mesh

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

    Inplane operators on staggered mesh

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

    Parallel magnetic potential on staggered mesh

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

    Parallel magnetic potential on canonical mesh

    integer, intent(in), optional, dimension(mesh_stag%get_n_points_boundary()) :: bnd_descrs_flux

    Boundary descriptor for perpendicular boundary points of parallel flux

interface

  • public module subroutine solve_aligned3d_adjoint(comm_handler, equi_on_cano, equi_on_stag, mesh_cano, mesh_stag, map, params_solver, solver_type, rhs, sol, niter, res, true_res, cnt_matvec, lambda, xi, co_cano, fcx, bnd_descrs, flutter_fac_cano, flutter_fac_stag, opsinplane_cano, opsinplane_stag, apar_stag, apar_cano, bnd_descrs_flux)

    Solves the adjoint elliptic Helmholtz problem along magnetic field line, The problem is stated on the staggered mesh optionally along perturbed magnetic field line via magentic flutter \f[ \lambda * u - \xi*\nabla_parallel\left[co\nabla\cdot\left(\mathbf{b}u fcx\right)\right] = rhs \f] with coefficients \f$\lambda, \xi, co, fcx $ used e.g. for solving the Landau heat flux The problem is stated on the staggered mesh Perpendicular boundary points are set according to bnd_descrs Magnetic flutter is included in parallel gradient and divergence via Possion bracket and arakawa scheme Perpendicular boundary flux_bnd_descrs is specially reqired for dobule arakawa operation

    Arguments

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

    Communicators

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

    Equilibrim quantities on canonical mesh

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

    Equilibrim quantities on staggered mesh

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

    Canonical mesh

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

    Staggered mesh

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

    Parallel map

    type(params_solver_aligned3d_t), intent(in) :: params_solver

    Parameters for solver

    character(len=*), intent(in) :: solver_type

    Type of solver

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

    Right hand side

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

    On input: initial guess, on output :: solution

    integer, intent(out) :: niter

    Number of (outer) iterations, negative number if error occured

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

    Pseudo-residual of solution

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

    True residual of solution

    integer, intent(out), optional :: cnt_matvec

    Counts number that matvec was called (~total iterations)

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

    Coefficient lambda, default = 1.0

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

    Coefficient xi, default = 1.0

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

    Coefficient co, default = 1.0, defined on staggered mesh

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

    Coefficient fcx, default = 1.0

    integer, intent(in), optional, dimension(mesh_stag%get_n_points_boundary()) :: bnd_descrs

    Boundary descriptor for perpendicular boundary points

    real(kind=GP), intent(in), optional :: flutter_fac_cano

    Factor (~switch) in front of flutter divergence on canonical mesh

    real(kind=GP), intent(in), optional :: flutter_fac_stag

    Factor (~switch) in front of flutter gradient on staggered mesh

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

    Inplane operators on canonical mesh

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

    Inplane operators on staggered mesh

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

    Parallel magnetic potential on staggered mesh

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

    Parallel magnetic potential on canonical mesh

    integer, intent(in), optional, dimension(mesh_cano%get_n_points_boundary()) :: bnd_descrs_flux

    Boundary descriptor for perpendicular boundary points of parallel flux


Derived Types

type, public ::  params_solver_aligned3d_t

Parameters for 3D aligned solver

Components

Type Visibility Attributes Name Initial
character(len=16), public :: krylov_method = 'RGMRES'

Krylov method to be used as iterative solver

real(kind=GP), public :: resmax = 1.0E-8_GP

Maximum allowed residual (stopping criterion), default 1E-8

integer, public :: maxiter = 10

Maximum number of (outer) iterations, default = 10

integer, public :: nrestart = 10

Maximum number of iterations before restart (~inner iterations)

character(len=16), public :: precond_type = 'JAC'

Type of preconditioner: 'NONE' (default) 'JAC' (Jacobi)

integer, public :: dbgout = 0

Debug output level