File: N:\mfix\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) == 0) then
83 call desmpi_check_sendrecvbuf(check_global=.true.)
84 lcheckbuf = 0
85 elseif (mod(lcheckbuf,5) == 0) then
86 call desmpi_check_sendrecvbuf(check_global=.false.)
87 end if
88 lcheckbuf = lcheckbuf + 1
89
90
91 do ii=1, size(dsendbuf)
92 dsendbuf(ii)%facebuf(1) = 0
93 drecvbuf(ii)%facebuf(1) = 0
94 end do
95
96 ispot = 1
97 do linter = merge(2,3,NO_K),1,-1
98 do lface = linter*2-1,linter*2
99 if(.not.iexchflag(lface))cycle
100 call desmpi_pack_parcross(lface)
101 call desmpi_sendrecv_init(lface)
102 end do
103 do lface = linter*2-1,linter*2
104 if(.not.iexchflag(lface)) cycle
105 call desmpi_sendrecv_wait(lface)
106 call desmpi_unpack_parcross(lface)
107 end do
108
109 do lface = linter*2-1,linter*2
110 if(dsendbuf(1+mod(lface,2))%facebuf(1).gt.0 .or. &
111 drecvbuf(1+mod(lface,2))%facebuf(1).gt.0) then
112 call desgrid_pic(plocate=.false.)
113 exit
114 end if
115 end do
116 end do
117 call des_mpi_barrier
118
119
120
121
122 IF(.NOT.MPPIC) THEN
123
124
125 do ii=1, size(dsendbuf)
126 dsendbuf(ii)%facebuf(1) = 0
127 drecvbuf(ii)%facebuf(1) = 0
128 end do
129
130 ighost_updated(:) = .false.
131 ispot = 1
132 do linter = 1,merge(2,3,NO_K)
133 do lface = linter*2-1,linter*2
134 if(.not.iexchflag(lface))cycle
135 call desmpi_pack_ghostpar(lface)
136 call desmpi_sendrecv_init(lface)
137 end do
138 do lface = linter*2-1,linter*2
139 if(.not.iexchflag(lface)) cycle
140 call desmpi_sendrecv_wait(lface)
141 call desmpi_unpack_ghostpar(lface)
142 end do
143
144
145 do lface = linter*2-1,linter*2
146 if(dsendbuf(1+mod(lface,2))%facebuf(1) .gt.0.or.&
147 drecvbuf(1+mod(lface,2))%facebuf(1).gt.0) then
148 call desgrid_pic(plocate=.false.)
149 exit
150 end if
151 end do
152 end do
153 call desmpi_cleanup
154 call des_mpi_barrier
155 ENDIF
156
157
158
159
160
161
162
163 END SUBROUTINE DES_PAR_EXCHANGE
164
165
166
167
168
169
170
171
172
173 SUBROUTINE DESMPI_CHECK_SENDRECVBUF(check_global)
174
175 use derived_types, only: dg_pic
176 use discretelement, only: DIMN
177 use desmpi, only: iMAXBUF
178 use desmpi, only: iBUFOFFSET
179 use desmpi, only: dSENDBUF, dRECVBUF
180 use desmpi, only: iSENDINDICES
181 use desmpi, only: iGhostPacketSize
182
183 use mpi_utility, only: global_all_max
184 use error_manager
185 implicit none
186
187 logical, intent(in) :: check_global
188
189
190
191
192 INTEGER :: lface, lindx, lijk
193
194 INTEGER :: lparcnt
195
196 INTEGER :: lmaxcnt
197
198 INTEGER :: ltot_ind
199
200 INTEGER :: pBUF
201
202 REAL :: lfactor = 0.5
203 DOUBLE PRECISION, PARAMETER :: ONEMBo8 = 131072.0
204
205
206 = 0
207 do lface = 1,6
208 ltot_ind = isendindices(1,lface)
209 lparcnt = 0
210 do lindx = 2,ltot_ind+1
211 lijk = isendindices(lindx,lface)
212 lparcnt = lparcnt + dg_pic(lijk)%isize
213 enddo
214 if(lparcnt > lmaxcnt) lmaxcnt = lparcnt
215 enddo
216
217 if(check_global) call global_all_max(lmaxcnt)
218
219 if (imaxbuf < (1.0+0.5*lfactor)*lmaxcnt*iGhostPacketSize) then
220 pbuf = imaxbuf
221 imaxbuf = (1.0+lfactor)*lmaxcnt*iGhostPacketSize
222 do lface = 1,2*dimn
223 if(allocated(dsendbuf(1+mod(lface,2))%facebuf)) then
224 deallocate(dsendbuf(1+mod(lface,2))%&
225 facebuf,drecvbuf(1+mod(lface,2))%facebuf)
226 endif
227 allocate(dsendbuf(1+mod(lface,2))%facebuf(imaxbuf),&
228 drecvbuf(1+mod(lface,2))%facebuf(imaxbuf))
229 end do
230
231 WRITE(ERR_MSG, 1000) iMAXBUF/ONEMBo8, &
232 100.0d0+100.0d0*dble(iMAXBUF-pbuf)/dble(pbuf)
233 CALL FLUSH_ERR_MSG(HEADER=.FALSE., FOOTER=.FALSE.)
234
235 1000 FORMAT(/'Resizeing DES MPI buffers: ',F7.1,' MB (+',F5.1'%)')
236
237 endif
238
239 END SUBROUTINE DESMPI_CHECK_SENDRECVBUF
240
241
242
243
244
245
246
247 SUBROUTINE DESMPI_CLEANUP
248
249 use discretelement, only: DIMN
250 use discretelement, only: DES_POS_NEW, DES_POS_OLD
251 use discretelement, only: DES_VEL_NEW, DES_VEL_OLD
252 use discretelement, only: OMEGA_NEW
253 use discretelement, only: PARTICLE_ORIENTATION
254 use discretelement, only: ORIENTATION,INIT_ORIENTATION
255 use discretelement, only: FC
256 use discretelement, only: DO_OLD
257 use discretelement, only: PIP
258 use discretelement, only: iGHOST_CNT
259 use discretelement, only: DES_USR_VAR_SIZE, DES_USR_VAR
260 use derived_types, only: dg_pic
261 use discretelement, only: pijk
262
263 use run, only: ENERGY_EQ,ANY_SPECIES_EQ
264 use des_thermo, only: DES_T_s
265
266 use des_rxns, only: DES_X_s
267
268 use discretelement, only: iGHOST_UPDATED
269 use functions, only: SET_NONEXISTENT
270 use desmpi, only: iRECVINDICES
271 use desmpi, only: iEXCHFLAG
272
273 use param, only: DIMENSION_N_s
274 use param1, only: ZERO
275
276 implicit none
277
278
279
280 integer ltot_ind,lface,lindx,lijk,lcurpar,lpicloc
281
282 do lface = 1,dimn*2
283 if(.not.iexchflag(lface))cycle
284 ltot_ind = irecvindices(1,lface)
285 do lindx = 2,ltot_ind+1
286 lijk = irecvindices(lindx,lface)
287 do lpicloc =1,dg_pic(lijk)%isize
288 lcurpar = dg_pic(lijk)%p(lpicloc)
289 if(ighost_updated(lcurpar)) cycle
290 pip = pip - 1
291 ighost_cnt = ighost_cnt-1
292 call set_nonexistent(lcurpar)
293 fc(lcurpar,:) = 0.0
294 des_pos_new(lcurpar,:)=0
295 pijk(lcurpar,:) = -10
296 IF (DO_OLD) THEN
297 des_pos_old(lcurpar,:)=0
298 des_vel_old(lcurpar,:)=0
299 ENDIF
300 des_vel_new(lcurpar,:)=0
301 omega_new(lcurpar,:)=0
302
303 IF(PARTICLE_ORIENTATION) &
304 ORIENTATION(1:3,lcurpar) = INIT_ORIENTATION
305
306 if(ENERGY_EQ) des_t_s(lcurpar)=0
307
308 if(ANY_SPECIES_EQ)then
309 des_x_s(lcurpar,1:dimension_n_s)= 0
310 endif
311
312 IF(DES_USR_VAR_SIZE > 0)&
313 des_usr_var(:,lcurpar)= 0
314
315 end do
316 end do
317 end do
318 END SUBROUTINE DESMPI_CLEANUP
319
320
321
322
323
324
325
326
327
328 SUBROUTINE DESMPI_SEND_RECV_FIELD_VARS
329
330
331 use run, only: ENERGY_EQ, ANY_SPECIES_EQ, SOLVE_ROs
332
333 use discretelement, only: DES_EXPLICITLY_COUPLED
334
335 use particle_filter, only: DES_INTERP_ON
336
337 use fldvar, only: EP_g, ROP_g
338
339 use fldvar, only: ROP_s, RO_s
340
341 use des_rxns, only: DES_R_gp, DES_R_gc
342 use des_rxns, only: DES_SUM_R_g, DES_R_PHASE
343
344 use des_thermo, only: CONV_Sc
345 use des_rxns, only: DES_HOR_g
346
347 use sendrecvnode, only: DES_COLLECT_gDATA
348
349 use sendrecv, only: SEND_RECV
350
351 implicit none
352
353
354
355
356 CALL SEND_RECV(EP_G,2)
357 CALL SEND_RECV(ROP_G,2)
358
359 CALL SEND_RECV(ROP_S,2)
360 IF(any(SOLVE_ROs)) CALL SEND_RECV(RO_S,2)
361
362
363 IF(.NOT.DES_EXPLICITLY_COUPLED) THEN
364
365 IF(DES_INTERP_ON) THEN
366 IF(ENERGY_EQ) CALL DES_COLLECT_gDATA(CONV_SC)
367 IF(ANY_SPECIES_EQ) THEN
368 CALL DES_COLLECT_gDATA(DES_R_gp)
369 CALL DES_COLLECT_gDATA(DES_R_gc)
370 CALL DES_COLLECT_gDATA(DES_R_PHASE)
371 CALL DES_COLLECT_gDATA(DES_HOR_g)
372 CALL DES_COLLECT_gDATA(DES_SUM_R_g)
373 ENDIF
374 ENDIF
375
376 IF(ENERGY_EQ) CALL SEND_RECV(CONV_SC, 2)
377 IF(ANY_SPECIES_EQ) THEN
378 CALL SEND_RECV(DES_R_gp, 2)
379 CALL SEND_RECV(DES_R_gc, 2)
380 CALL SEND_RECV(DES_R_PHASE, 2)
381 CALL SEND_RECV(DES_HOR_g, 2)
382 CALL SEND_RECV(DES_SUM_R_g, 2)
383 ENDIF
384 ENDIF
385
386 RETURN
387 END SUBROUTINE DESMPI_SEND_RECV_FIELD_VARS
388
389 END MODULE MPI_FUNS_DES
390