File: RELATIVE:/../../../mfix.git/model/des/mpi_funs_des_mod.f
1
2
3
4
5
6
7
8 module mpi_funs_des
9
10 contains
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40 subroutine des_par_exchange()
41
42 use discretelement, only: DO_NSEARCH, DES_PERIODIC_WALLS
43
44 use mfix_pic, only: MPPIC
45 use desmpi, only: iEXCHFLAG
46 use desmpi, only: dSENDBUF, dRECVBUF
47 use discretelement, only: iGHOST_UPDATED
48 use desmpi, only: iSPOT
49
50 use geometry, only: NO_K
51
52
53
54 use mpi_pack_des, only: desmpi_pack_parcross
55 use mpi_unpack_des, only: desmpi_unpack_parcross
56
57 use mpi_pack_des, only: desmpi_pack_ghostpar
58 use mpi_unpack_des, only: desmpi_unpack_ghostpar
59
60 use mpi_comm_des, only: desmpi_sendrecv_init
61 use mpi_comm_des, only: desmpi_sendrecv_wait
62
63 use desgrid, only: desgrid_pic
64 use desmpi_wrapper, only: des_mpi_barrier
65 use compar, only: numpes
66
67 implicit none
68
69
70
71
72 integer :: linter, lface, ii
73
74 integer, save :: lcheckbuf = 0
75
76
77
78 IF (.not.((numPEs>1) .OR. DES_PERIODIC_WALLS)) RETURN
79
80
81
82 if (mod(lcheckbuf,100) .eq. 0) then
83 call desmpi_check_sendrecvbuf
84 lcheckbuf = 0
85 end if
86 lcheckbuf = lcheckbuf + 1
87
88
89 do ii=1, size(dsendbuf)
90 dsendbuf(ii)%facebuf(1) = 0
91 drecvbuf(ii)%facebuf(1) = 0
92 end do
93
94 ispot = 1
95 do linter = merge(2,3,NO_K),1,-1
96 do lface = linter*2-1,linter*2
97 if(.not.iexchflag(lface))cycle
98 call desmpi_pack_parcross(lface)
99 call desmpi_sendrecv_init(lface)
100 end do
101 do lface = linter*2-1,linter*2
102 if(.not.iexchflag(lface)) cycle
103 call desmpi_sendrecv_wait(lface)
104 call desmpi_unpack_parcross(lface)
105 end do
106
107 do lface = linter*2-1,linter*2
108 if(dsendbuf(1+mod(lface,2))%facebuf(1).gt.0 .or. &
109 drecvbuf(1+mod(lface,2))%facebuf(1).gt.0) then
110 call desgrid_pic(plocate=.false.)
111 exit
112 end if
113 end do
114 end do
115 call des_mpi_barrier
116
117
118
119
120 IF(.NOT.MPPIC) THEN
121
122
123 do ii=1, size(dsendbuf)
124 dsendbuf(ii)%facebuf(1) = 0
125 drecvbuf(ii)%facebuf(1) = 0
126 end do
127
128 ighost_updated(:) = .false.
129 ispot = 1
130 do linter = 1,merge(2,3,NO_K)
131 do lface = linter*2-1,linter*2
132 if(.not.iexchflag(lface))cycle
133 call desmpi_pack_ghostpar(lface)
134 call desmpi_sendrecv_init(lface)
135 end do
136 do lface = linter*2-1,linter*2
137 if(.not.iexchflag(lface)) cycle
138 call desmpi_sendrecv_wait(lface)
139 call desmpi_unpack_ghostpar(lface)
140 end do
141
142
143 do lface = linter*2-1,linter*2
144 if(dsendbuf(1+mod(lface,2))%facebuf(1) .gt.0.or.&
145 drecvbuf(1+mod(lface,2))%facebuf(1).gt.0) then
146 call desgrid_pic(plocate=.false.)
147 exit
148 end if
149 end do
150 end do
151 if(do_nsearch) call desmpi_cleanup
152 call des_mpi_barrier
153 ENDIF
154
155
156
157
158
159
160
161 END SUBROUTINE DES_PAR_EXCHANGE
162
163
164
165
166
167
168
169
170
171 SUBROUTINE DESMPI_CHECK_SENDRECVBUF
172
173 use discretelement, only: DIMN, dg_pic
174 use desmpi, only: iMAXBUF
175 use desmpi, only: iBUFOFFSET
176 use desmpi, only: dSENDBUF, dRECVBUF
177 use desmpi, only: iSENDINDICES
178 use desmpi, only: iGhostPacketSize
179
180 use mpi_utility, only: global_all_max
181 use error_manager
182 implicit none
183
184
185
186
187 INTEGER :: lface, lindx, lijk
188
189 INTEGER :: lparcnt
190
191 INTEGER :: lmaxcnt
192
193 INTEGER :: ltot_ind
194
195 INTEGER :: pBUF
196
197 REAL :: lfactor = 1.5
198 DOUBLE PRECISION, PARAMETER :: ONEMBo8 = 131072.0
199
200
201 = 0
202 do lface = 1,6
203 ltot_ind = isendindices(1,lface)
204 lparcnt = 0
205 do lindx = 2,ltot_ind+1
206 lijk = isendindices(lindx,lface)
207 lparcnt = lparcnt + dg_pic(lijk)%isize
208 enddo
209 if(lparcnt.gt.lmaxcnt) lmaxcnt = lparcnt
210 enddo
211
212 call global_all_max(lmaxcnt)
213 if (imaxbuf .lt. lmaxcnt*iGhostPacketSize+ibufoffset) then
214 pbuf = imaxbuf
215 imaxbuf = lmaxcnt*iGhostPacketSize*lfactor
216 do lface = 1,2*dimn
217 if(allocated(dsendbuf(1+mod(lface,2))%facebuf)) then
218 deallocate(dsendbuf(1+mod(lface,2))%facebuf,drecvbuf(1+mod(lface,2))%facebuf)
219 endif
220 allocate(dsendbuf(1+mod(lface,2))%facebuf(imaxbuf),&
221 drecvbuf(1+mod(lface,2))%facebuf(imaxbuf))
222 end do
223
224 WRITE(ERR_MSG, 1000) iMAXBUF/ONEMBo8, &
225 100.0d0+100.0d0*dble(iMAXBUF-pbuf)/dble(pbuf)
226 CALL FLUSH_ERR_MSG(HEADER=.FALSE., FOOTER=.FALSE.)
227
228 1000 FORMAT(/'Resizeing DES MPI buffers: ',F7.1,' MB (+',F5.1'%)')
229
230 endif
231
232 END SUBROUTINE DESMPI_CHECK_SENDRECVBUF
233
234
235
236
237
238
239
240 SUBROUTINE DESMPI_CLEANUP
241
242 use discretelement, only: DIMN
243 use discretelement, only: DES_POS_NEW, DES_POS_OLD
244 use discretelement, only: DES_VEL_NEW, DES_VEL_OLD
245 use discretelement, only: OMEGA_NEW
246 use discretelement, only: PARTICLE_ORIENTATION,ORIENTATION,INIT_ORIENTATION
247 use discretelement, only: FC
248 use discretelement, only: DO_OLD
249 use discretelement, only: PIP
250 use discretelement, only: iGHOST_CNT
251 use discretelement, only: DES_USR_VAR_SIZE, DES_USR_VAR
252 use discretelement, only: dg_pic, pijk
253
254 use run, only: ENERGY_EQ,ANY_SPECIES_EQ
255 use des_thermo, only: DES_T_s_NEW, DES_T_s_OLD
256
257 use des_rxns, only: DES_X_s
258
259 use discretelement, only: iGHOST_UPDATED
260 use functions, only: SET_NONEXISTENT
261 use desmpi, only: iRECVINDICES
262 use desmpi, only: iEXCHFLAG
263
264 use param, only: DIMENSION_N_s
265 use param1, only: ZERO
266
267 implicit none
268
269
270
271 integer ltot_ind,lface,lindx,lijk,lcurpar,lpicloc
272
273 do lface = 1,dimn*2
274 if(.not.iexchflag(lface))cycle
275 ltot_ind = irecvindices(1,lface)
276 do lindx = 2,ltot_ind+1
277 lijk = irecvindices(lindx,lface)
278 do lpicloc =1,dg_pic(lijk)%isize
279 lcurpar = dg_pic(lijk)%p(lpicloc)
280 if(ighost_updated(lcurpar)) cycle
281 pip = pip - 1
282 ighost_cnt = ighost_cnt-1
283 call set_nonexistent(lcurpar)
284 fc(:,lcurpar) = 0.0
285 des_pos_new(:,lcurpar)=0
286 pijk(lcurpar,:) = -10
287 IF (DO_OLD) THEN
288 des_pos_old(:,lcurpar)=0
289 des_vel_old(:,lcurpar)=0
290 ENDIF
291 des_vel_new(:,lcurpar)=0
292 omega_new(:,lcurpar)=0
293
294 IF(PARTICLE_ORIENTATION) ORIENTATION(1:3,lcurpar) = INIT_ORIENTATION
295
296 if(ENERGY_EQ) then
297 des_t_s_new(lcurpar)=0
298 des_t_s_old(lcurpar)=0
299 endif
300
301 if(ANY_SPECIES_EQ)then
302 des_x_s(lcurpar,1:dimension_n_s)= 0
303 endif
304
305 IF(DES_USR_VAR_SIZE > 0)&
306 des_usr_var(:,lcurpar)= 0
307
308 end do
309 end do
310 end do
311 END SUBROUTINE DESMPI_CLEANUP
312
313
314 END MODULE MPI_FUNS_DES
315