solver_aligned3d_m.f90 Source File


Contents


Source Code

module solver_aligned3d_m
    !! Routines for 3d elliptic solver along magnetic field lines
    use precision_grillix_m, only : GP, GP_EPS
    use error_handling_grillix_m, only : handle_error, error_info_t
    use status_codes_grillix_m, only : GRILLIX_ERR_SOLVER3D
    use comm_handler_m, only : comm_handler_t
    use descriptors_m, only : BND_TYPE_DIRICHLET_ZERO, BND_TYPE_NEUMANN
    use mesh_cart_m, only : mesh_cart_t
    use boundaries_perp_m, only : compute_bndnmn_matrix_row, &
                                  set_perpbnds, extrapolate_ghost_points
    use equilibrium_storage_m, only : equilibrium_storage_t
    use parallel_map_m, only : parallel_map_t
    use inplane_operators_m, only : inplane_operators_t
    use solver3d_m, only : solver3d_t
    use solver3d_factory_m, only : solver3d_factory
    use communication_planes_m, only : start_comm_planes, &
                                       finalize_comm_planes
    use mpi_mapping_auxiliaries_m, only : getdata_fwdbwdplane
    implicit none
    
    type, public :: params_solver_aligned3d_t
        !! Parameters for 3D aligned solver
        character(len=16), public :: krylov_method = 'RGMRES'
        !! Krylov method to be used as iterative solver
        real(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
    end type
    
    public :: solve_aligned3d
    public :: solve_aligned3d_adjoint
    
    interface

        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  
            type(comm_handler_t), intent(in) :: comm_handler
            !! Communicators
            type(mesh_cart_t), intent(in) :: mesh_cano
            !! Canonical mesh
            type(mesh_cart_t), intent(in) :: mesh_stag
            !! Staggered mesh
            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(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(GP), dimension(mesh_cano%get_n_points()), intent(in) :: rhs
            !! Right hand side
            real(GP), dimension(mesh_cano%get_n_points()), intent(inout) :: sol
            !! On input: initial guess, on output :: solution
            integer, intent(out) :: niter
            !! Number of (outer) iterations, negative number if error occured
            real(GP), intent(out) :: res
            !! Pseudo-residual (relative) of solution
            real(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(GP), dimension(mesh_cano%get_n_points()), intent(in), optional :: lambda
            !! Coefficient lambda, default = 1.0
            real(GP), dimension(mesh_cano%get_n_points()), intent(in), optional :: xi
            !! Coefficient xi, default = 1.0
            real(GP), dimension(mesh_stag%get_n_points()), intent(in), optional :: co_stag
            !! Coefficient co, default = 1.0, defined  on staggered mesh
            real(GP), dimension(mesh_cano%get_n_points()), intent(in), optional :: fcx
            !! Coefficient fcx, default = 1.0
            integer, dimension(mesh_cano%get_n_points_boundary()), intent(in), optional :: bnd_descrs
            !! Boundary descriptor for perpendicular boundary points
            real(GP), intent(in), optional :: flutter_fac_cano
            !! Factor (~switch) in front of flutter divergence on canonical mesh
            real(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(GP), dimension(mesh_stag%get_n_points()), intent(in), optional :: apar_stag
            !! Parallel magnetic potential on staggered mesh
            real(GP), dimension(mesh_cano%get_n_points()), intent(in), optional :: apar_cano
            !! Parallel magnetic potential on canonical mesh
            integer, dimension(mesh_stag%get_n_points_boundary()), intent(in), optional :: bnd_descrs_flux
            !! Boundary descriptor for perpendicular boundary points of parallel flux  
        end subroutine
 
        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  
            type(comm_handler_t), intent(in) :: comm_handler
            !! Communicators
            type(mesh_cart_t), intent(in) :: mesh_cano
            !! Canonical mesh
            type(mesh_cart_t), intent(in) :: mesh_stag
            !! Staggered mesh
            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(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(GP), dimension(mesh_stag%get_n_points()), intent(in) :: rhs
            !! Right hand side
            real(GP), dimension(mesh_stag%get_n_points()), intent(inout) :: sol
            !! On input: initial guess, on output :: solution
            integer, intent(out) :: niter
            !! Number of (outer) iterations, negative number if error occured
            real(GP), intent(out) :: res
            !! Pseudo-residual of solution
            real(GP), intent(out) :: true_res
            !! True residual of solution
            integer, intent(out), optional :: cnt_matvec
            !! Counts number that matvec was called (~total iterations)
            real(GP), dimension(mesh_stag%get_n_points()), intent(in), optional :: lambda
            !! Coefficient lambda, default = 1.0
            real(GP), dimension(mesh_stag%get_n_points()), intent(in), optional :: xi
            !! Coefficient xi, default = 1.0
            real(GP), dimension(mesh_cano%get_n_points()), intent(in), optional :: co_cano
            !! Coefficient co, default = 1.0, defined  on staggered mesh
            real(GP), dimension(mesh_stag%get_n_points()), intent(in), optional :: fcx
            !! Coefficient fcx, default = 1.0
            integer, dimension(mesh_stag%get_n_points_boundary()), intent(in), optional :: bnd_descrs
            !! Boundary descriptor for perpendicular boundary points
            real(GP), intent(in), optional :: flutter_fac_cano
            !! Factor (~switch) in front of flutter divergence on canonical mesh
            real(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(GP), dimension(mesh_stag%get_n_points()), intent(in), optional :: apar_stag
            !! Parallel magnetic potential on staggered mesh
            real(GP), dimension(mesh_cano%get_n_points()), intent(in), optional :: apar_cano
            !! Parallel magnetic potential on canonical mesh
            integer, dimension(mesh_cano%get_n_points_boundary()), intent(in), optional :: bnd_descrs_flux
            !! Boundary descriptor for perpendicular boundary points of parallel flux       
        end subroutine
    
    end interface

end module