parallel_map_s.f90 Source File


Contents

Source Code


Source Code

submodule (parallel_map_m) parallel_map_s
    implicit none

contains

    module subroutine set_parameters_map(self, comm_handler, filename, &
                                         intorder_in, xorder_in, &
                                         use_fixed_stencil_in, use_gauss_quadrature_in)
        class(parallel_map_t), intent(inout) :: self
        type(comm_handler_t), intent(in) :: comm_handler
        character(*), intent(in), optional :: filename
        integer, intent(in), optional :: intorder_in
        integer, intent(in), optional :: xorder_in
        logical, intent(in), optional :: use_fixed_stencil_in
        logical, intent(in), optional :: use_gauss_quadrature_in
       
        integer :: io_error
        character(len=256) :: io_errmsg 

        integer :: intorder, xorder
        logical :: use_fixed_stencil, use_gauss_quadrature

        namelist / params_map / &
            intorder, xorder, &
            use_fixed_stencil, use_gauss_quadrature

        
        if ( (.not.present(filename)) .and. &
           (present(intorder_in)) .and. & 
           (present(xorder_in)) .and. &
           (present(use_fixed_stencil_in)) .and. &
           (present(use_gauss_quadrature_in)) ) then  

            ! Input via explicit setting                            
            self%intorder = intorder_in
            self%xorder = xorder_in
            self%use_fixed_stencil = use_fixed_stencil_in
            self%use_gauss_quadrature = use_gauss_quadrature_in
            
        elseif ( (present(filename)) .and. &
           (.not.present(intorder_in)) .and. & 
           (.not.present(xorder_in)) .and. &
           (.not.present(use_fixed_stencil_in)) .and. & 
           (.not.present(use_gauss_quadrature_in)) ) then  
            ! Input via file       

            ! Default values
            intorder        = 3
            xorder          = 0
            use_fixed_stencil = .false.
            use_gauss_quadrature = .false.
            
            open(unit = 20, file = filename, status = 'old', action = 'read',&
                iostat = io_error, iomsg = io_errmsg)
            if (io_error /= 0) then
                call handle_error(io_errmsg, GRILLIX_ERR_NAMELIST, __LINE__, __FILE__)
            endif

            read(20, nml = params_map, iostat = io_error, iomsg = io_errmsg)
            if (io_error /= 0) then
                call handle_error(io_errmsg, GRILLIX_ERR_NAMELIST, __LINE__, __FILE__)
            endif          
               
            close(20)

            self%intorder   = intorder
            self%xorder     = xorder
            self%use_fixed_stencil = use_fixed_stencil
            self%use_gauss_quadrature = use_gauss_quadrature
        
        else
            call handle_error('Input to subroutine not consistent', &
                              GRILLIX_ERR_OTHER, __LINE__, __FILE__)
        endif   
         
    end subroutine
    
    module subroutine set_toroidal_location(self, new_iplane, new_phi_cano, new_phi_stag)
        class(parallel_map_t), intent(inout) :: self
        integer, intent(in) :: new_iplane
        real(GP), intent(in) :: new_phi_cano
        real(GP), intent(in) :: new_phi_stag
        
        self%iplane = new_iplane
        self%phi_cano = new_phi_cano
        self%phi_stag = new_phi_stag
        
    end subroutine

    module subroutine display_parallel_map(self)
        class(parallel_map_t), intent(in) :: self
        logical :: axisymmetric
    
        if (allocated(self%dpar_stag_cano_fwd_mem)) then
            axisymmetric = .false.
        else
            axisymmetric = .true.
        endif
    
        if (.not. is_rank_info_writer) then
            return
        endif

        write(get_stdout(),*)'Information on map ------------------------------------'
        write(get_stdout(),101) self%nplanes, self%dphi, self%intorder, self%xorder, &
                                self%use_fixed_stencil, self%use_gauss_quadrature, &
                                axisymmetric
101     FORMAT( 3X,'nplanes              = ',I3 /, &
                3X,'dphi                 = ',ES14.6E3 /, &
                3X,'intorder             = ',I3 /, &
                3X,'xorder               = ',I3 /, &  
                3X,'use_fixed_stencil    = ',L1 /, &  
                3X,'use_gauss_quadrature = ',L1 /, &
                3X,'axisymmetric         = ',L1)
        write(get_stdout(),*)'-------------------------------------------------------'

    end subroutine

end submodule