File: /nfs/home/0/users/jenkins/mfix.git/model/des/desmpi_wrapper_mod.f

1     !------------------------------------------------------------------------
2     ! Module           : desmpi_wrapper
3     ! Purpose          : Contains wrapper class for mpi communications- send,recv
4     !                    and scatter, gather
5     ! Author           : Pradeep.G
6     !------------------------------------------------------------------------
7           module desmpi_wrapper
8     
9           use parallel_mpi
10           use mpi_utility
11           use compar
12     
13           interface des_mpi_irecv
14              module procedure des_mpi_irecv_db
15           end interface
16     
17           interface des_mpi_isend
18              module procedure des_mpi_isend_db
19           end interface
20     
21           interface des_mpi_scatterv
22              module procedure des_mpi_scatterv_i
23              module procedure des_mpi_scatterv_db
24           end interface
25     
26           interface des_mpi_gatherv
27              module procedure des_mpi_gatherv_i
28              module procedure des_mpi_gatherv_db
29           end interface
30     
31           contains
32     
33     !------------------------------------------------------------------------
34     ! Subroutine       : des_mpi_barrier
35     ! Purpose          : Wrapper class for barrier
36     !
37     !------------------------------------------------------------------------
38           subroutine des_mpi_barrier
39           implicit none
40     ! local variables
41           character(len=80), parameter :: name = 'des_mpi_barrier'
42           integer lerr
43           call mpi_barrier(mpi_comm_world, lerr)
44           call mpi_check( name //':mpi_barrier ', lerr )
45           return
46           end subroutine
47     
48     !------------------------------------------------------------------------
49     ! Subroutine       : des_mpi_irecv_db
50     ! Purpose          : Wrapper class for mpi_irecv
51     !------------------------------------------------------------------------
52           subroutine des_mpi_irecv_db(precvbuf,precvcnt,ptoproc,ptag,preq,perr)
53           implicit none
54     ! dummy variables
55           double precision, dimension(:) :: precvbuf
56           integer :: precvcnt,ptoproc,ptag,preq,perr
57           call mpi_irecv(precvbuf,precvcnt,mpi_double_precision,ptoproc,ptag,mpi_comm_world,preq,perr)
58           end subroutine
59     
60     !------------------------------------------------------------------------
61     ! Subroutine       : des_mpi_isend_db
62     ! Purpose          : Wrapper class for mpi_isend
63     !------------------------------------------------------------------------
64           subroutine des_mpi_isend_db(psendbuf,psendcnt,ptoproc,ptag,preq,perr)
65           implicit none
66     ! dummy variables
67           double precision, dimension(:) :: psendbuf
68           integer :: psendcnt,ptoproc,ptag,preq,perr
69           call mpi_isend(psendbuf,psendcnt,mpi_double_precision,ptoproc,ptag,mpi_comm_world,preq,perr)
70           end subroutine
71     
72     !------------------------------------------------------------------------
73     ! Subroutine       : des_mpi_wait
74     ! Purpose          : Wrapper class for mpi_wait
75     !------------------------------------------------------------------------
76           subroutine des_mpi_wait(preq,perr)
77           implicit none
78     ! dummy variables
79           integer :: preq,perr
80     ! local variables
81           integer :: lmpi_status(mpi_status_size)
82           call mpi_wait(preq,lmpi_status,perr)
83           end subroutine
84     
85     
86     !------------------------------------------------------------------------
87     ! Subroutine       : des_mpi_scatterv_db
88     ! Purpose          : Wrapper class for mpi_scatterv
89     !------------------------------------------------------------------------
90           subroutine des_mpi_scatterv_db(prootbuf,pscattercnts,pdispls, &
91                                pprocbuf,precvcnt,proot,perr )
92           implicit none
93     ! dummy variables
94           double precision, dimension(:):: prootbuf,pprocbuf
95           integer, dimension (:) :: pdispls,pscattercnts
96           integer :: precvcnt,proot,perr
97     
98           call mpi_scatterv(prootbuf,pscattercnts,pdispls,mpi_double_precision, &
99                             pprocbuf,precvcnt,mpi_double_precision,proot,mpi_comm_world,perr )
100           end subroutine
101     !------------------------------------------------------------------------
102     ! Subroutine       : des_mpi_scatterv_i
103     ! Purpose          : Wrapper class for mpi_scatterv
104     !------------------------------------------------------------------------
105           subroutine des_mpi_scatterv_i(prootbuf,pscattercnts,pdispls, &
106                                pprocbuf, precvcnt,proot,perr )
107           implicit none
108     ! dummy variables
109           integer, dimension(:) :: prootbuf,pprocbuf
110           integer, dimension (:) :: pdispls,pscattercnts
111           integer :: precvcnt,proot,perr
112     
113           call mpi_scatterv(prootbuf,pscattercnts,pdispls,mpi_integer, &
114                             pprocbuf,precvcnt,mpi_integer,proot,mpi_comm_world,perr )
115           end subroutine
116     !------------------------------------------------------------------------
117     ! Subroutine       : des_mpi_gatherv_db
118     ! Purpose          : Wrapper class for mpi_gatherv
119     !------------------------------------------------------------------------
120           subroutine des_mpi_gatherv_db(psendbuf,psendcnt,precvbuf, &
121                                         precvcnts, pdispls,proot,perr )
122           implicit none
123     ! dummy variables
124           double precision, dimension(:) :: psendbuf,precvbuf
125           integer, dimension (:) :: pdispls,precvcnts
126           integer :: psendcnt,proot,perr
127     
128           call mpi_gatherv(psendbuf,psendcnt,mpi_double_precision,precvbuf,precvcnts, &
129                            pdispls,mpi_double_precision,proot,mpi_comm_world,perr )
130           end subroutine
131     !------------------------------------------------------------------------
132     ! Subroutine       : des_mpi_gatherv_i
133     ! Purpose          : Wrapper class for mpi_gatherv
134     !------------------------------------------------------------------------
135           subroutine des_mpi_gatherv_i(psendbuf,psendcnt,precvbuf, &
136                                         precvcnts, pdispls,proot,perr )
137           implicit none
138     ! dummy variables
139           integer, dimension(:) :: psendbuf,precvbuf
140           integer, dimension (:) :: pdispls,precvcnts
141           integer :: psendcnt,proot,perr
142     
143           call mpi_gatherv(psendbuf,psendcnt,mpi_integer,precvbuf,precvcnts, &
144                            pdispls,mpi_integer,proot,mpi_comm_world,perr )
145           end subroutine
146     
147     !------------------------------------------------------------------------
148     ! Subroutine       : des_mpi_stop
149     ! Purpose          : Wrapper class for mpi_abort
150     !------------------------------------------------------------------------
151           subroutine des_mpi_stop(myid)
152     
153           USE funits
154           implicit none
155     
156           integer, optional, intent(in) :: myid
157     
158           INTEGER :: mylid, ERRORCODE
159     
160           if (.not. present(myid)) then
161              mylid = myPE
162           else
163              mylid = myid
164           endif
165     
166           write(*,100) mylid
167           write(UNIT_LOG,100) mylid
168     
169      100  format(/,'*****************',&
170           '********************************************',/, &
171           '(PE ',I2,') : A fatal error occurred in des routines',/,9X, &
172           '*.LOG file may contain other error messages ',/,'*****************', &
173           '********************************************',/)
174     
175           call MPI_BARRIER(MPI_COMM_WORLD, mpierr)
176           call MPI_ABORT(MPI_COMM_WORLD, ERRORCODE, mpierr)
177           write(*,"('(PE ',I2,') : MPI_ABORT return = ',I2)") &
178           mylid,mpierr
179     
180           call MPI_Finalize(mpierr)
181     
182           STOP 'MPI terminated from des_mpi_stop'
183     
184           end subroutine
185     
186           end module
187