parallel_map_operators_s.f90 Source File


Contents


Source Code

submodule (parallel_map_m) parallel_map_operators_s
    !! Parallel operator routines (mapping, parallel gradient/divergence)
    implicit none

contains

    pure real(GP) module function par_grad_single(self, ucano, ucano_fwd, l_stag)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%grad_stag_cano_bwd%ncol), intent(in) :: ucano
        real(GP), dimension(self%grad_stag_cano_fwd%ncol), intent(in) :: ucano_fwd
        integer, intent(in) :: l_stag
        
        par_grad_single = csr_times_vec_single(self%grad_stag_cano_fwd, ucano_fwd, l_stag) &
                          - csr_times_vec_single(self%grad_stag_cano_bwd, ucano, l_stag)
                              
    end function

    module subroutine par_grad(self, ucano, ucano_fwd, par_grad_u_stag)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%grad_stag_cano_bwd%ncol), intent(in) :: ucano
        real(GP), dimension(self%grad_stag_cano_fwd%ncol), intent(in) :: ucano_fwd
        real(GP), dimension(self%grad_stag_cano_fwd%ndim), intent(out) :: par_grad_u_stag
        
        integer :: l_stag
       
        !$omp do
        do l_stag = 1, self%grad_stag_cano_fwd%ndim
            par_grad_u_stag(l_stag) = self%par_grad_single(ucano, ucano_fwd, l_stag)
        enddo
        !$omp end do
                
    end subroutine
    
    pure real(GP) module function par_divb_single(self, ustag, ustag_bwd, l_cano)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%pdiv_cano_stag_fwd%ncol), intent(in) :: ustag
        real(GP), dimension(self%pdiv_cano_stag_bwd%ncol), intent(in) :: ustag_bwd
        integer, intent(in) :: l_cano
        
        par_divb_single = -csr_times_vec_single(self%pdiv_cano_stag_fwd, ustag, l_cano) &
                          + csr_times_vec_single(self%pdiv_cano_stag_bwd, ustag_bwd, l_cano)
        
    end function
        
    module subroutine par_divb(self, ustag, ustag_bwd, par_pdiv_u_cano)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%pdiv_cano_stag_fwd%ncol), intent(in) :: ustag
        real(GP), dimension(self%pdiv_cano_stag_bwd%ncol), intent(in) :: ustag_bwd
        real(GP), dimension(self%pdiv_cano_stag_fwd%ndim), intent(out) :: par_pdiv_u_cano
         
        integer :: l_cano
     
        !$omp do
        do l_cano = 1, self%pdiv_cano_stag_fwd%ndim
            par_pdiv_u_cano(l_cano) = self%par_divb_single(ustag, ustag_bwd, l_cano)
        enddo
        !$omp end do
        
    end subroutine
    
    pure real(GP) module function lift_cano_to_stag_single(self, ucano, ucano_fwd, l_stag)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_stag_cano_bwd%ncol), intent(in) :: ucano
        real(GP), dimension(self%map_stag_cano_fwd%ncol), intent(in) :: ucano_fwd
        integer, intent(in) :: l_stag
        
        real(GP) :: dfwd, dbwd
        
        dfwd = self%dpar_stag_cano_fwd(l_stag)
        dbwd = self%dpar_stag_cano_bwd(l_stag)
        
        lift_cano_to_stag_single = csr_times_vec_single(self%map_stag_cano_fwd, ucano_fwd, l_stag) * dbwd &
                                   + csr_times_vec_single(self%map_stag_cano_bwd, ucano, l_stag) * dfwd
                                   
        lift_cano_to_stag_single = lift_cano_to_stag_single / (dfwd + dbwd)
        
    end function

    module subroutine lift_cano_to_stag(self, ucano, ucano_fwd, lift_u_stag)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_stag_cano_bwd%ncol), intent(in) :: ucano
        real(GP), dimension(self%map_stag_cano_fwd%ncol), intent(in) :: ucano_fwd
        real(GP), dimension(self%map_stag_cano_fwd%ndim), intent(out) :: lift_u_stag
        
        integer :: l_stag
        
        !$omp do
        do l_stag = 1, self%map_stag_cano_fwd%ndim
            lift_u_stag(l_stag) = self%lift_cano_to_stag_single(ucano, ucano_fwd, l_stag)
        enddo
        !$omp end do
        
    end subroutine
    
    pure real(GP) module function flat_stag_to_cano_single(self, ustag, ustag_bwd, l_cano)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_cano_stag_fwd%ncol), intent(in) :: ustag
        real(GP), dimension(self%map_cano_stag_bwd%ncol), intent(in) :: ustag_bwd
        integer, intent(in) :: l_cano

        real(GP) :: dfwd, dbwd
        
        dfwd = self%dpar_cano_stag_fwd(l_cano)
        dbwd = self%dpar_cano_stag_bwd(l_cano)
        
        
        flat_stag_to_cano_single = csr_times_vec_single(self%map_cano_stag_fwd, ustag, l_cano) * dbwd &
                                   + csr_times_vec_single(self%map_cano_stag_bwd, ustag_bwd, l_cano) * dfwd
                
        flat_stag_to_cano_single = flat_stag_to_cano_single / (dfwd + dbwd)        
        
    end function
        
    module subroutine flat_stag_to_cano(self, ustag, ustag_bwd, flat_u_cano)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_cano_stag_fwd%ncol), intent(in) :: ustag
        real(GP), dimension(self%map_cano_stag_bwd%ncol), intent(in) :: ustag_bwd
        real(GP), dimension(self%map_cano_stag_fwd%ndim), intent(out) :: flat_u_cano
        
        integer :: l_cano
                
        !$omp do
        do l_cano = 1, self%map_cano_stag_fwd%ndim
            flat_u_cano(l_cano) = self%flat_stag_to_cano_single(ustag, ustag_bwd, l_cano)
        enddo
        !$omp end do
        
    end subroutine
    
    pure real(GP) module function upstream_cano_from_cano_fwd_single(self, ucano_fwd, l_cano)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_cano_cano_fwd%ncol), intent(in) :: ucano_fwd
        integer, intent(in) :: l_cano
        
        upstream_cano_from_cano_fwd_single = &
            csr_times_vec_single(self%map_cano_cano_fwd, ucano_fwd, l_cano)
        
    end function
        
    module subroutine upstream_cano_from_cano_fwd(self, ucano_fwd, upstream_cano_fwd)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_cano_cano_fwd%ncol), intent(in) :: ucano_fwd
        real(GP), dimension(self%map_cano_cano_fwd%ndim), intent(out) :: upstream_cano_fwd
        
        integer :: l_cano
                
        !$omp do
        do l_cano = 1, self%map_cano_cano_fwd%ndim
            upstream_cano_fwd(l_cano) = self%upstream_cano_from_cano_fwd_single(ucano_fwd, l_cano)
        enddo
        !$omp end do
        
    end subroutine

    pure real(GP) module function upstream_cano_from_cano_fwd2_single(self, ucano_fwd2, l_cano)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_cano_cano_fwd2%ncol), intent(in) :: ucano_fwd2
        integer, intent(in) :: l_cano

        upstream_cano_from_cano_fwd2_single = &
            csr_times_vec_single(self%map_cano_cano_fwd2, ucano_fwd2, l_cano)

    end function

    module subroutine upstream_cano_from_cano_fwd2(self, ucano_fwd2, upstream_cano_fwd2)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_cano_cano_fwd2%ncol), intent(in) :: ucano_fwd2
        real(GP), dimension(self%map_cano_cano_fwd2%ndim), intent(out) :: upstream_cano_fwd2

        integer :: l_cano

        !$omp do
        do l_cano = 1, self%map_cano_cano_fwd2%ndim
            upstream_cano_fwd2(l_cano) = self%upstream_cano_from_cano_fwd2_single(ucano_fwd2, l_cano)
        enddo
        !$omp end do

    end subroutine

    pure real(GP) module function upstream_cano_from_cano_bwd_single(self, ucano_bwd, l_cano)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_cano_cano_bwd%ncol), intent(in) :: ucano_bwd
        integer, intent(in) :: l_cano
        
        upstream_cano_from_cano_bwd_single = &
            csr_times_vec_single(self%map_cano_cano_bwd, ucano_bwd, l_cano)
        
    end function
        
    module subroutine upstream_cano_from_cano_bwd(self, ucano_bwd, upstream_cano_bwd)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_cano_cano_bwd%ncol), intent(in) :: ucano_bwd
        real(GP), dimension(self%map_cano_cano_bwd%ndim), intent(out) :: upstream_cano_bwd
        
        integer :: l_cano
                
        !$omp do
        do l_cano = 1, self%map_cano_cano_bwd%ndim
            upstream_cano_bwd(l_cano) = self%upstream_cano_from_cano_bwd_single(ucano_bwd, l_cano)
        enddo
        !$omp end do
        
    end subroutine

    pure real(GP) module function upstream_cano_from_cano_bwd2_single(self, ucano_bwd2, l_cano)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_cano_cano_bwd2%ncol), intent(in) :: ucano_bwd2
        integer, intent(in) :: l_cano

        upstream_cano_from_cano_bwd2_single = &
            csr_times_vec_single(self%map_cano_cano_bwd2, ucano_bwd2, l_cano)

    end function

    module subroutine upstream_cano_from_cano_bwd2(self, ucano_bwd2, upstream_cano_bwd2)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_cano_cano_bwd2%ncol), intent(in) :: ucano_bwd2
        real(GP), dimension(self%map_cano_cano_bwd2%ndim), intent(out) :: upstream_cano_bwd2

        integer :: l_cano

        !$omp do
        do l_cano = 1, self%map_cano_cano_bwd2%ndim
            upstream_cano_bwd2(l_cano) = self%upstream_cano_from_cano_bwd2_single(ucano_bwd2, l_cano)
        enddo
        !$omp end do

    end subroutine

    pure real(GP) module function upstream_stag_from_stag_fwd_single(self, ustag_fwd, l_stag)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_stag_stag_fwd%ncol), intent(in) :: ustag_fwd
        integer, intent(in) :: l_stag
        
        upstream_stag_from_stag_fwd_single = &
            csr_times_vec_single(self%map_stag_stag_fwd, ustag_fwd, l_stag)
        
    end function
        
    module subroutine upstream_stag_from_stag_fwd(self, ustag_fwd, upstream_stag_fwd)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_stag_stag_fwd%ncol), intent(in) :: ustag_fwd
        real(GP), dimension(self%map_stag_stag_fwd%ndim), intent(out) :: upstream_stag_fwd
    
        integer :: l_stag
    
        !$omp do
        do l_stag = 1, self%map_stag_stag_fwd%ndim
            upstream_stag_fwd(l_stag) = self%upstream_stag_from_stag_fwd_single(ustag_fwd, l_stag)
        enddo
        !$omp end do
        
    end subroutine
    
    pure real(GP) module function upstream_stag_from_stag_bwd_single(self, ustag_bwd, l_stag)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_stag_stag_bwd%ncol), intent(in) :: ustag_bwd
        integer, intent(in) :: l_stag
        
        upstream_stag_from_stag_bwd_single = &
            csr_times_vec_single(self%map_stag_stag_bwd, ustag_bwd, l_stag)
        
    end function
        
    module subroutine upstream_stag_from_stag_bwd(self, ustag_bwd, upstream_stag_bwd)
        class(parallel_map_t), intent(in) :: self
        real(GP), dimension(self%map_stag_stag_bwd%ncol), intent(in) :: ustag_bwd
        real(GP), dimension(self%map_stag_stag_bwd%ndim), intent(out) :: upstream_stag_bwd
        
        integer :: l_stag
    
        !$omp do
        do l_stag = 1, self%map_stag_stag_bwd%ndim
            upstream_stag_bwd(l_stag) = self%upstream_stag_from_stag_bwd_single(ustag_bwd, l_stag)
        enddo
        !$omp end do
        
    end subroutine
    
end submodule