Line data Source code
1 : !-------------------------------------------------------------------------------- 2 : ! Copyright (c) 2016 Peter Grünberg Institut, Forschungszentrum Jülich, Germany 3 : ! This file is part of FLEUR and available as free software under the conditions 4 : ! of the MIT license as expressed in the LICENSE file in more detail. 5 : !-------------------------------------------------------------------------------- 6 : MODULE m_types_hybmpi 7 : #ifdef CPP_MPI 8 : use mpi 9 : #endif 10 : TYPE t_hybmpi 11 : INTEGER :: comm 12 : INTEGER :: rank 13 : INTEGER :: size 14 : contains 15 : procedure :: copy_mpi => t_hybmpi_copy_mpi 16 : procedure :: barrier => t_hybmpi_barrier 17 : procedure :: init => t_hybmpi_init 18 : procedure :: root => t_hybmpi_root 19 : END TYPE t_hybmpi 20 : contains 21 120 : function t_hybmpi_root(mpi_var) result(l_root) 22 : implicit none 23 : class(t_hybmpi), intent(in) :: mpi_var 24 : logical :: l_root 25 : 26 120 : l_root = mpi_var%rank == 0 27 120 : end function t_hybmpi_root 28 : 29 0 : subroutine t_hybmpi_copy_mpi(glob_mpi, mpi_var) 30 : use m_types_mpi 31 : implicit none 32 : class(t_hybmpi), intent(inout) :: glob_mpi 33 : type(t_mpi), intent(in) :: mpi_var 34 : 35 0 : glob_mpi%comm = mpi_var%mpi_comm 36 0 : glob_mpi%size = mpi_var%isize 37 0 : glob_mpi%rank = mpi_var%irank 38 0 : end subroutine 39 : 40 0 : subroutine t_hybmpi_barrier(glob_mpi) 41 : use m_judft 42 : implicit none 43 : class(t_hybmpi), intent(inout) :: glob_mpi 44 : integer :: ierr 45 : #ifdef CPP_MPI 46 0 : call MPI_Barrier(glob_mpi%comm, ierr) 47 0 : if(ierr /= 0) call juDFT_error("barrier failed on process: " // & 48 0 : int2str(glob_mpi%rank)) 49 : #endif 50 0 : end subroutine t_hybmpi_barrier 51 : 52 96 : subroutine t_hybmpi_init(hybmpi, in_comm) 53 : class(t_hybmpi), intent(inout) :: hybmpi 54 : integer, intent(in) :: in_comm 55 : integer :: ierr 56 : 57 48 : hybmpi%comm = in_comm 58 : 59 : #ifdef CPP_MPI 60 48 : call MPI_Comm_size(hybmpi%comm, hybmpi%size, ierr) 61 48 : call MPI_Comm_rank(hybmpi%comm, hybmpi%rank, ierr) 62 : #else 63 : hybmpi%size = 1 64 : hybmpi%rank = 0 65 : #endif 66 48 : end subroutine t_hybmpi_init 67 0 : END MODULE m_types_hybmpi