variable_m.f90 Source File


Contents

Source Code


Source Code

module variable_m
    !! Definition of variable
    use MPI
    use comm_handler_m, only : comm_handler_t
    use netcdf
    use error_handling_grillix_m, only: handle_error_netcdf, &
                                        handle_error, error_info_t
    use status_codes_grillix_m, only : GRILLIX_ERR_ALLOC, GRILLIX_ERR_OTHER
    use precision_grillix_m, only : GP, GP_EPS, GP_NAN 
    use communication_planes_m, only : start_comm_planes, &
                                       finalize_comm_planes
    
    ! from PARALLAX
    use constants_m, only : TWO_PI
    use screen_io_m, only :  get_stdout
    use elementary_functions_m, only : gaussian
    use equilibrium_m, only : equilibrium_t
    use mesh_cart_m, only : mesh_cart_t
    use fieldline_tracer_m, only : trace
    implicit none

    type, public :: variable_t
        !! Datatype for variable_t
        real(GP), allocatable, dimension(:), public :: vals
        !! Values of variable
        real(GP), dimension(:), pointer, public :: hfwd => null()
        !! Pointer to halos_fwd(:,1), 
        !! simplifies call of parallel operators 
        real(GP), dimension(:), pointer, public :: hbwd => null()
        !! Pointer to halos_bwd(:,1),
        !! simplifies call of parallel operators
        integer, private :: ndim
        !! Dimension of data variable
  
        integer, allocatable, dimension(:), private :: ndim_planes
        !! Dimension of data variable for all planes
        integer, private :: nhalo_fwd
        !! Number of halo planes in the forward direction
        integer, private :: nhalo_bwd
        !! Number of halo planes in the backward direction
        real(GP), allocatable, dimension(:,:), private :: halos_fwd
        !! Halo values in forward direction     
        !! e.g. values from plane rank+2 are stored in halo_fwd(:,2)
        !! WILL BE RENAMED to halos_fwd
        real(GP), allocatable, dimension(:,:), private ::  halos_bwd
        !! Halo values in backward direction        
        !! e.g. values from plane rank-2 are stored in halo_bwd(:,2)
        !! WILL BE RENAMED to halos_bwd
        integer, allocatable, dimension(:,:), private :: req_fwd
        !! Communicator request for forward communication
        integer, allocatable, dimension(:,:), private :: req_bwd
        !! Communicator request for backward communication
        
        logical, private :: staggered
        !! Indicates if variable defined on staggered or canonical grid
        !! WILL BE REMOVED
    contains
        procedure, public :: get_ndim   
        procedure, public :: init => init_variable
        procedure, public :: create_halos
        procedure, public :: get_halo
        procedure, public :: start_comm_halos
        procedure, public :: finalize_comm_halos
        procedure, public :: set_blob_toroidal
        procedure, public :: set_blob_aligned        
        procedure, public :: write_netcdf => write_netcdf_variable
        procedure, public :: read_netcdf => read_netcdf_variable
        final :: destructor                        
    end type 
    
    interface

        pure integer module function get_ndim(self)
            !! Returns dimension of variable
            class(variable_t), intent(in) :: self
            !! Instance of the type                   
        end function

        module subroutine init_variable(self, comm_handler, ndim, staggered, &
                                        init_vals)
            !! Initialises variable
            class(variable_t), intent(inout) :: self  
            !! Instance of the type
            type(comm_handler_t), intent(in) :: comm_handler
            !! Communicators
            integer, intent(in) :: ndim
            !! Dimension of variable
            logical, intent(in) :: staggered
            !! Indicates if variable is defined on staggered or canonical grid
            !! WILL BE REMOVED
            real(GP), dimension(ndim), intent(in), optional :: init_vals
            !! Values for initialisation if present, 
            !! otherwise will be initialised with zeros            
        end subroutine
        
        module subroutine create_halos(self, comm_handler, &
                                            nhalo_fwd, nhalo_bwd)
            !! Allocates halo fields
            class(variable_t), intent(inout), target :: self  
            !! Instance of the type
            type(comm_handler_t), intent(in) :: comm_handler
            !! Communicators
            integer, intent(in) :: nhalo_fwd
            !! Number of halo planes in the forward direction
            integer, intent(in) :: nhalo_bwd
            !! Number of halo planes in the backward direction
        end subroutine
        
        module function get_halo(self, comm_handler, ihalo) result(res)
            !! Returns pointer to halo field located at plane=plane+ihalo
            real(GP), dimension(:), pointer :: res
            !! Pointer to halo field
            class(variable_t), intent(inout), target :: self
            !! Instance of the type
            type(comm_handler_t), intent(in) :: comm_handler
            !! Communicators
            integer, intent(in) :: ihalo
            !! Displacement to plane, where information is fetched from
        end function
        
        module subroutine start_comm_halos(self, comm_handler, step)
            !! Starts communication of halo fields
            class(variable_t), intent(inout) :: self  
            !! Instance of the type
            type(comm_handler_t), intent(in) :: comm_handler
            !! Communicators
            integer, intent(in), optional :: step
            !! Only starts communication to get data from plane+step
            !! If not present, communication over all halos is performed
        end subroutine  

        module subroutine finalize_comm_halos(self, comm_handler, step)
            !! Finalizes communication of halos
            class(variable_t), intent(inout) :: self  
            !! Instance of the type
            type(comm_handler_t), intent(in) :: comm_handler
            !! Communicators
            integer, intent(in), optional :: step
            !! Only finalizes communication that gets data from plane+step
            !! If not present, communication over all halos is finalized
        end subroutine
        
        module subroutine set_blob_toroidal(self, mesh, x0, y0, phi0, &
                                            wx, wy, wphi, amp)
            !! Sets Gaussian in (R, phi, Z) for variable data 
            class(variable_t), intent(inout) :: self  
            !! Instance of the type
            type(mesh_cart_t), intent(in) :: mesh
            !! Mesh within poloidal plane (canonical or staggered)
            real(GP), intent(in) :: x0
            !! x-coordinate of center of blob
            real(GP), intent(in) :: y0
            !! y-coordinate of center of blob
            real(GP), intent(in) :: phi0
            !! phi-coordinate of center of blob
            real(GP), intent(in) :: wx
            !! Gaussian width of blob in x-direction 
            real(GP), intent(in) :: wy
            !! Gaussian width of blob in y-direction 
            real(GP), intent(in) :: wphi
            !! Gaussian width of blob in toroidal phi-direction 
            !! A delta peak on first plane is assumed if width is set to zero
            real(GP), intent(in) :: amp
            !! Amplitude of blob  
        end subroutine
        
        module subroutine set_blob_aligned(self, equi, mesh, x0, y0, phi0, &
                                           wx, wy, wpar, amp)
            !! Sets Gaussian in (R, parallel, Z) for variable data,
            !! aligned along magnetic field
            class(variable_t), intent(inout) :: self  
            !! Instance of the type
            class(equilibrium_t), intent(inout) :: equi
            !! Equilibrium
            type(mesh_cart_t), intent(in) :: mesh
            !! Mesh within poloidal plane (canonical or staggered)
            real(GP), intent(in) :: x0
            !! x-coordinate of center of blob
            real(GP), intent(in) :: y0
            !! y-coordinate of center of blob
            real(GP), intent(in) :: phi0
            !! phi-coordinate of center of blob
            real(GP), intent(in) :: wx
            !! Gaussian width of blob in x-direction 
            real(GP), intent(in) :: wy
            !! Gaussian width of blob in y-direction 
            real(GP), intent(in) :: wpar
            !! Gaussian width of blob in toroidal phi-direction 
            !! A delta peak is assumed if width is set to zero
            real(GP), intent(in) :: amp
            !! Amplitude of blob  
        end subroutine
        
        module subroutine write_netcdf_variable(self, comm_handler, fgid, &
                                                varname)
            !! Writes variable to NETCDF file
            class(variable_t), intent(in) :: self
            !! Instance of the type
            type(comm_handler_t), intent(in) :: comm_handler 
            !! MPI communication handler
            integer, intent(in) :: fgid
            !! NETCDF file or group id     
            character(len=*), intent(in) :: varname
            !! Name for variable for identification in NETCDF file            
        end subroutine
        
        module subroutine read_netcdf_variable(self, comm_handler, fgid, &
                                               varname)
            !! Reads variable from NETCDF file
            !! Note: This routine only works yet for axisymmetric conditions,
            !! but not for 3D geometries
            class(variable_t), intent(out) :: self 
            !! Instance of the type
            type(comm_handler_t), intent(in) :: comm_handler 
            !! MPI communication handler
            integer, intent(in) :: fgid
            !! NETCDF file or group id
            character(len=*), intent(in) :: varname
            !! Name for variable for identification in NETCDF file   
        end subroutine

        module subroutine destructor(self)
            !! Frees memory associated with variable
            type(variable_t), intent(inout) :: self
            !! Instance of the type
        end subroutine        

    end interface

end module