parallel_target_flux_m.f90 Source File


Contents


Source Code

module parallel_target_flux_m
    !! Module related to parallel fluxes onto the target plates
    use MPI
    use netcdf
    use error_handling_grillix_m, only : handle_error_netcdf, &
                                         handle_error, error_info_t
    use status_codes_grillix_m, only : GRILLIX_ERR_OTHER, GRILLIX_ERR_ALLOC
    use precision_grillix_m, only : GP, GP_EPS, MPI_GP
    use comm_handler_m, only : comm_handler_t
    use screen_io_m, only : get_stdout
    use parallelisation_setup_m, only : is_rank_info_writer
    use mesh_cart_m, only : mesh_cart_t
    use parallel_map_m, only : parallel_map_t
    use penalisation_m, only : penalisation_t
    use equilibrium_storage_m, only : equilibrium_storage_t
    use descriptors_braginskii_m, only : BND_TYPE_NEUMANN
    use boundaries_perp_m, only : set_perpbnds, extrapolate_ghost_points
    use variable_m, only : variable_t
    use type_conversion_m, only : logical_to_integer, integer_to_logical
    use mpi_mapping_auxiliaries_m, only : getdata_fwdbwdplane
    implicit none

    type, public :: parallel_target_flux_t
        !! Datatype for computing and storing points which map into
        !! the penalisation volume. 
        !! NOTE: Indices are stored with respect to the INNER mesh
        !! since we require their dirindfun data for computation
        integer, allocatable, dimension(:), public :: inds_fwd
        !! Indices of points on inner mesh which map into pen volume in fwd plane
        real(GP), allocatable, dimension(:), public :: weights_fwd
        !! Penalisation values of points mapped into pen volume in fwd plane
        integer, public :: n_inds_fwd = 0
        !! Number of fwd points to store
        integer, allocatable, dimension(:), public :: inds_bwd
        !! Indices of points on inner mesh which map into pen volume in bwd plane
        real(GP), allocatable, dimension(:), public :: weights_bwd
        !! Penalisation values of points mapped into pen volume in bwd plane
        integer, public :: n_inds_bwd = 0
        !! Number of bwd points to store
        logical, public :: is_staggered
        !! Whether stored points are from the canonical or staggered grid
    contains
        procedure, public :: build => build_target_markers
        procedure, public :: write_netcdf => write_netcdf_target_markers
        procedure, public :: read_netcdf => read_netcdf_target_markers
        procedure, public :: compute => compute_target_flux
        procedure, public :: compute_single => compute_target_flux_single
        final :: destructor
    end type

    interface

        module subroutine build_target_markers(self, comm_handler, mesh, map, penalisation, staggered, &
                                                pen_threshold, dbgout)
            !! Find target points and construct storage
            class(parallel_target_flux_t), intent(inout) :: self
            !! Instance of class
            type(comm_handler_t), intent(in) :: comm_handler
            !! Communicator
            type(mesh_cart_t), intent(in) :: mesh
            !! Mesh
            type(parallel_map_t), intent(in) :: map
            !! Parallel map
            type(penalisation_t), intent(in) :: penalisation
            !! Penalisation   
            logical, intent(in) :: staggered
            !! Whether to construct target points on staggered or canonical grid
            !! Acts on canonical grid if set to .false.
            !! Choice must be consistent with input mesh!
            real(GP), intent(in) :: pen_threshold
            !! Threshold for interpolated penalisation value on mapped point
            !! Only points with values above the threshold are stored as target points
            integer, intent(in), optional :: dbgout
            !! Debug output level 
        end subroutine

        module subroutine write_netcdf_target_markers(self, fgid)
            !! Writes stored indices and weights to netcdf file
            class(parallel_target_flux_t), intent(in) :: self
            !! Instance of class
            integer, intent(in) :: fgid
            !! File or group id
        end subroutine
    
        module subroutine read_netcdf_target_markers(self, fgid)
            !! Reads stored indices and weights from netcdf file
            class(parallel_target_flux_t), intent(inout) :: self
            !! Instance of class
            integer, intent(in) :: fgid
            !! File or group id
        end subroutine

        module subroutine compute_target_flux(self, comm_handler, mesh, penalisation, equi_storage, &
                                              g, flux_fwd, flux_bwd, mode)
            !! Compute parallel fluxes into penalisation region
            !! flux_fwd is defined as [g * dS] integrated over the set of points
            !! which map into the penalisation region in the forward plane.
            !! flux_bwd is defined analogously for penalisation in the backward plane. 
            !! g is a parallel flux, e.g. ne * upar_cano, 
            !! B the magnetic field, dS the boundary surface
            !! Note that input parameters must be consistent with the grid choice
            !! when building the marked points, e.g. only staggered or only canonical
            class(parallel_target_flux_t), intent(in) :: self
            !! Instance of class
            type(comm_handler_t), intent(in) :: comm_handler
            !! Communicator
            type(mesh_cart_t), intent(in) :: mesh
            !! Mesh
            type(penalisation_t), intent(in) :: penalisation
            !! Penalisation
            type(equilibrium_storage_t), intent(in) :: equi_storage
            !! Equilibrium storage
            real(GP), dimension(mesh%get_n_points()), intent(in) :: g
            !! Parallel flux in plane
            real(GP), intent(out) :: flux_fwd
            !! Resulting flux into forward plane penalisation
            real(GP), intent(out) :: flux_bwd
            !! Resulting flux into backward plane penalisation
            character(*), intent(in), optional :: mode
            !! Mode of flux computation, may select between:
            !! 'DEFAULT'  : Fluxes are NOT weighted by local penalisation value
            !! 'WEIGHTED' : Fluxes are weighted by local penalisation value
        end subroutine

        module function compute_target_flux_single(self, mesh, penalisation, equi_storage, &
                                                    i_marker, g, dirind, mode) result(flux_out)
            !! Compute parallel fluxes into penalisation region for a single flux point
            class(parallel_target_flux_t), intent(in) :: self
            !! Instance of class
            type(mesh_cart_t), intent(in) :: mesh
            !! Mesh
            type(penalisation_t), intent(in) :: penalisation
            !! Penalisation
            type(equilibrium_storage_t), intent(in) :: equi_storage
            !! Equilibrium storage
            integer, intent(in) :: i_marker
            !! Index of point in fwd/bwd marker array
            real(GP), dimension(mesh%get_n_points()), intent(in) :: g
            !! Parallel flux in plane
            character(*), intent(in) :: dirind
            !! Field direction indication, either:
            !! 'fwd' : Compute forward flux into penalisation
            !! 'bwd' : Compute backward flux into penalisation
            character(*), intent(in), optional :: mode
            !! Mode of flux computation, may select between:
            !! 'DEFAULT'  : Fluxes are NOT weighted by local penalisation value
            !! 'WEIGHTED' : Fluxes are weighted by local penalisation value

            real(GP) :: flux_out
            !! Resulting flux into into penalisation
        end function

        module subroutine destructor(self)
            !! Destructor for target point storage
            type(parallel_target_flux_t), intent(inout) :: self
            !! Instance of type
        end subroutine

    end interface

end module