File: RELATIVE:/../../../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     #ifdef MPI
43           integer lerr
44           call mpi_barrier(mpi_comm_world, lerr)
45           call mpi_check( name //':mpi_barrier ', lerr )
46     #endif
47           return
48           end subroutine
49     
50     !------------------------------------------------------------------------
51     ! Subroutine       : des_mpi_irecv_db
52     ! Purpose          : Wrapper class for mpi_irecv
53     !------------------------------------------------------------------------
54           subroutine des_mpi_irecv_db(precvbuf,precvcnt,ptoproc,ptag,preq,perr)
55           implicit none
56     ! dummy variables
57           double precision, dimension(:) :: precvbuf
58           integer :: precvcnt,ptoproc,ptag,preq,perr
59     #ifdef MPI
60           call mpi_irecv(precvbuf,precvcnt,mpi_double_precision,ptoproc,ptag,mpi_comm_world,preq,perr)
61     #endif
62           return
63           end subroutine
64     
65     !------------------------------------------------------------------------
66     ! Subroutine       : des_mpi_isend_db
67     ! Purpose          : Wrapper class for mpi_isend
68     !------------------------------------------------------------------------
69           subroutine des_mpi_isend_db(psendbuf,psendcnt,ptoproc,ptag,preq,perr)
70           implicit none
71     ! dummy variables
72           double precision, dimension(:) :: psendbuf
73           integer :: psendcnt,ptoproc,ptag,preq,perr
74     #ifdef MPI
75           call mpi_isend(psendbuf,psendcnt,mpi_double_precision,ptoproc,ptag,mpi_comm_world,preq,perr)
76     #endif
77           return
78           end subroutine
79     
80     !------------------------------------------------------------------------
81     ! Subroutine       : des_mpi_wait
82     ! Purpose          : Wrapper class for mpi_wait
83     !------------------------------------------------------------------------
84           subroutine des_mpi_wait(preq,perr)
85           implicit none
86     ! dummy variables
87           integer :: preq,perr
88     #ifdef MPI
89     ! local variables
90           integer :: lmpi_status(mpi_status_size)
91           call mpi_wait(preq,lmpi_status,perr)
92     #endif
93           return
94           end subroutine
95     
96     !------------------------------------------------------------------------
97     ! Subroutine       : des_mpi_scatterv_db
98     ! Purpose          : Wrapper class for mpi_scatterv
99     !------------------------------------------------------------------------
100           subroutine des_mpi_scatterv_db(prootbuf,pscattercnts,pdispls, &
101                                pprocbuf,precvcnt,proot,perr )
102           implicit none
103     ! dummy variables
104           double precision, dimension(:):: prootbuf,pprocbuf
105           integer, dimension (:) :: pdispls,pscattercnts
106           integer :: precvcnt,proot,perr
107     
108     #ifdef MPI
109           call mpi_scatterv(prootbuf,pscattercnts,pdispls,mpi_double_precision, &
110                             pprocbuf,precvcnt,mpi_double_precision,proot,mpi_comm_world,perr )
111     #else
112           pprocbuf = prootbuf
113     #endif
114           return
115           end subroutine
116     !------------------------------------------------------------------------
117     ! Subroutine       : des_mpi_scatterv_i
118     ! Purpose          : Wrapper class for mpi_scatterv
119     !------------------------------------------------------------------------
120           subroutine des_mpi_scatterv_i(prootbuf,pscattercnts,pdispls, &
121                                pprocbuf, precvcnt,proot,perr )
122           implicit none
123     ! dummy variables
124           integer, dimension(:) :: prootbuf,pprocbuf
125           integer, dimension (:) :: pdispls,pscattercnts
126           integer :: precvcnt,proot,perr
127     
128     #ifdef MPI
129           call mpi_scatterv(prootbuf,pscattercnts,pdispls,mpi_integer, &
130                             pprocbuf,precvcnt,mpi_integer,proot,mpi_comm_world,perr )
131     #else
132           pprocbuf = prootbuf
133     #endif
134           return
135           end subroutine
136     !------------------------------------------------------------------------
137     ! Subroutine       : des_mpi_gatherv_db
138     ! Purpose          : Wrapper class for mpi_gatherv
139     !------------------------------------------------------------------------
140           subroutine des_mpi_gatherv_db(psendbuf,psendcnt,precvbuf, &
141                                         precvcnts, pdispls,proot,perr )
142           implicit none
143     ! dummy variables
144           double precision, dimension(:) :: psendbuf,precvbuf
145           integer, dimension (:) :: pdispls,precvcnts
146           integer :: psendcnt,proot,perr
147     
148     #ifdef MPI
149           call mpi_gatherv(psendbuf,psendcnt,mpi_double_precision,precvbuf,precvcnts, &
150                            pdispls,mpi_double_precision,proot,mpi_comm_world,perr )
151     #else
152           precvbuf = psendbuf
153     #endif
154           return
155           end subroutine
156     !------------------------------------------------------------------------
157     ! Subroutine       : des_mpi_gatherv_i
158     ! Purpose          : Wrapper class for mpi_gatherv
159     !------------------------------------------------------------------------
160           subroutine des_mpi_gatherv_i(psendbuf,psendcnt,precvbuf, &
161                                         precvcnts, pdispls,proot,perr )
162           implicit none
163     ! dummy variables
164           integer, dimension(:) :: psendbuf,precvbuf
165           integer, dimension (:) :: pdispls,precvcnts
166           integer :: psendcnt,proot,perr
167     
168     #ifdef MPI
169           call mpi_gatherv(psendbuf,psendcnt,mpi_integer,precvbuf,precvcnts, &
170                            pdispls,mpi_integer,proot,mpi_comm_world,perr )
171     #else
172           precvbuf = psendbuf
173     #endif
174           return
175           end subroutine
176     
177     !------------------------------------------------------------------------
178     ! Subroutine       : des_mpi_stop
179     ! Purpose          : Wrapper class for mpi_abort
180     !------------------------------------------------------------------------
181           subroutine des_mpi_stop(myid)
182     
183           USE funits
184           implicit none
185     
186           integer, optional, intent(in) :: myid
187     
188     #ifdef MPI
189     
190           INTEGER :: mylid, ERRORCODE
191     
192           if (.not. present(myid)) then
193              mylid = myPE
194           else
195              mylid = myid
196           endif
197     
198           write(*,100) mylid
199           write(UNIT_LOG,100) mylid
200     
201      100  format(/,'*****************',&
202           '********************************************',/, &
203           '(PE ',I2,') : A fatal error occurred in des routines',/,9X, &
204           '*.LOG file may contain other error messages ',/,'*****************', &
205           '********************************************',/)
206     
207           call MPI_BARRIER(MPI_COMM_WORLD, mpierr)
208           call MPI_ABORT(MPI_COMM_WORLD, ERRORCODE, mpierr)
209           write(*,"('(PE ',I2,') : MPI_ABORT return = ',I2)") &
210           mylid,mpierr
211     
212           call MPI_Finalize(mpierr)
213     
214           STOP 'MPI terminated from des_mpi_stop'
215     #else
216           stop
217     #endif
218           end subroutine
219     
220           end module
221