File: N:\mfix\model\des\mpi_pack_des_mod.f
1
2
3
4
5
6
7
8 MODULE MPI_PACK_DES
9
10 PRIVATE
11 PUBLIC :: DESMPI_PACK_PARCROSS, DESMPI_PACK_GHOSTPAR
12
13 interface pack_dbuf
14 module procedure pack_db0
15 module procedure pack_db1
16 module procedure pack_i0
17 module procedure pack_i1
18 module procedure pack_l0
19 end interface pack_dbuf
20
21 CONTAINS
22
23
24
25
26
27
28
29 SUBROUTINE DESMPI_PACK_GHOSTPAR(PFACE)
30
31
32
33
34 use desmpi, only: iGhostPacketSize
35
36 use discretelement, only: iGHOST_UPDATED
37
38 use desmpi, only: dSENDBUF
39
40 use desmpi, only: iBUFOFFSET
41
42 use run, only: ENERGY_EQ
43
44 use desmpi, only: iNEIGHPROC
45
46 use discretelement, only: DG_PIJKPRV
47
48 use discretelement, only: iGLOBAL_ID
49
50 use discretelement, only: DES_POS_NEW
51
52 use discretelement, only: DES_VEL_NEW
53
54 use discretelement, only: OMEGA_NEW
55
56 use des_thermo, only: DES_T_s
57
58 use discretelement, only: DES_RADIUS
59
60 use discretelement, only: PIJK
61
62 use discretelement, only: MAX_PIP
63
64 use discretelement, only: DES_USR_VAR, DES_USR_VAR_SIZE
65
66 use desgrid, only: dg_ijkconv
67
68 use desmpi, only: isendcnt
69
70 use desmpi, only: dcycl_offset
71
72 use derived_types, only: DG_PIC
73
74 use desmpi, only: iSENDINDICES
75
76
77
78 use constant, only: PI
79
80 use discretelement, only: DIMN
81 use functions, only: is_exiting
82 use functions, only: is_ghost, is_entering_ghost, is_exiting_ghost
83
84 IMPLICIT NONE
85
86
87
88
89 INTEGER, INTENT(IN) :: PFACE
90
91
92
93 integer :: lijk,lindx,ltot_ind,lpicloc,lpar_cnt,lcurpar
94 integer :: lbuf
95
96
97 = 0
98 ltot_ind = isendindices(1,pface)
99 do lindx = 2,ltot_ind+1
100 lijk = isendindices(lindx,pface)
101 do lpicloc =1,dg_pic(lijk)%isize
102 lbuf = lpar_cnt*iGhostPacketSize + ibufoffset
103 lcurpar = dg_pic(lijk)%p(lpicloc)
104
105
106
107 if((is_ghost(lcurpar) .or. &
108 is_entering_ghost(lcurpar) .or. &
109 is_exiting_ghost(lcurpar)) .and. &
110 .not.ighost_updated(lcurpar)) cycle
111
112
113 call pack_dbuf(lbuf,iglobal_id(lcurpar),pface)
114
115 call pack_dbuf(lbuf,dg_ijkconv(lijk,pface, &
116 ineighproc(pface)),pface)
117
118 call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lcurpar),pface, &
119 ineighproc(pface)),pface)
120
121 call pack_dbuf(lbuf,des_radius(lcurpar),pface)
122
123 call pack_dbuf(lbuf,pijk(lcurpar,5),pface)
124
125 call pack_dbuf(lbuf,des_pos_new(lcurpar,:)+ &
126 dcycl_offset(pface,:),pface)
127
128 call pack_dbuf(lbuf,des_vel_new(lcurpar,:),pface)
129
130 call pack_dbuf(lbuf,omega_new(lcurpar,:),pface)
131
132 call pack_dbuf(lbuf,merge(1,0,is_exiting(lcurpar).or.is_exiting_ghost(lcurpar)),pface)
133
134 IF(ENERGY_EQ) &
135 call pack_dbuf(lbuf,des_t_s(lcurpar),pface)
136
137 IF(DES_USR_VAR_SIZE > 0) &
138 call pack_dbuf(lbuf,des_usr_var(:,lcurpar),pface)
139
140 lpar_cnt = lpar_cnt + 1
141 end do
142 end do
143 dsendbuf(1+mod(pface,2))%facebuf(1)=lpar_cnt
144 isendcnt(pface) = lpar_cnt*iGhostPacketSize + ibufoffset
145
146 end subroutine desmpi_pack_ghostpar
147
148
149
150
151
152
153
154 SUBROUTINE DESMPI_PACK_PARCROSS(PFACE)
155
156
157
158
159 use desmpi, only: dSENDBUF
160
161 use desmpi, only: iBUFOFFSET
162
163 use run, only: ENERGY_EQ
164
165 use run, only: ANY_SPECIES_EQ
166
167 use mfix_pic, only: MPPIC
168
169 use discretelement, only: DG_PIJKPRV
170
171 use desmpi, only: iNEIGHPROC
172
173 use mfix_pic, only: DES_STAT_WT
174
175 use discretelement, only: iGLOBAL_ID
176
177 use discretelement, only: DES_POS_NEW, DES_POS_OLD
178
179 use discretelement, only: DES_VEL_NEW, DES_VEL_OLD
180
181 use discretelement, only: OMEGA_NEW, OMEGA_OLD
182
183 use discretelement, only: PARTICLE_ORIENTATION,ORIENTATION
184
185 use discretelement, only: DES_RADIUS, PVOL, RO_SOL, PMASS
186
187 use discretelement, only: DES_ACC_OLD, ROT_ACC_OLD
188
189 use des_rxns, only: DES_X_s
190
191 use des_thermo, only: DES_T_s
192
193 use discretelement, only: FC, TOW
194
195 use discretelement, only: OMOI
196
197 use discretelement, only: PIJK
198
199 use discretelement, only: DO_OLD
200
201 use discretelement, only: PIP, MAX_PIP
202
203 use discretelement, only: iGHOST_CNT
204
205 use discretelement, only: DES_USR_VAR, DES_USR_VAR_SIZE
206
207 use discretelement, only: NEIGHBORS, NEIGHBOR_INDEX, NEIGH_NUM
208
209 use discretelement, only: PFT_NEIGHBOR
210
211 use discretelement, only: DIMN
212
213 use discretelement, only: DES_EXPLICITLY_COUPLED
214
215 use discretelement, only: DRAG_FC
216
217 use des_thermo, only: CONV_Qs, RXNS_Qs
218
219 use desgrid, only: dg_ijkconv, icycoffset
220 use desmpi, only: dcycl_offset, isendcnt
221 use desmpi, only: irecvindices
222
223 use desmpi, only: iParticlePacketSize
224 use desmpi, only: iPairPacketSize
225 use derived_types, only: dg_pic
226
227 use functions
228
229 implicit none
230
231
232
233
234 INTEGER, INTENT(IN) :: PFACE
235
236
237
238 integer :: li, lj, lk
239 integer :: ltot_ind,lindx,cc
240 integer :: lneigh,lijk,&
241 lpicloc,lparcnt,lcurpar
242 integer :: lbuf,num_neighborlists_to_send
243
244 logical, allocatable, dimension(:) :: going_to_send
245
246
247 integer :: num_neighborlists_send_buf_loc
248
249
250
251 = irecvindices(1,pface)
252 lparcnt = 0
253
254 allocate(going_to_send(max_pip))
255 going_to_send(:) = .false.
256
257 do lindx = 2,ltot_ind+1
258 lijk = irecvindices(lindx,pface)
259 do lpicloc = 1,dg_pic(lijk)%isize
260 lcurpar = dg_pic(lijk)%p(lpicloc)
261
262
263 if(is_ghost(lcurpar) .or. &
264 is_entering_ghost(lcurpar) .or. &
265 is_exiting_ghost(lcurpar)) cycle
266
267 going_to_send(lcurpar) = .true.
268 lbuf = lparcnt*iParticlePacketSize + ibufoffset
269
270
271 call pack_dbuf(lbuf,iglobal_id(lcurpar),pface)
272
273 call pack_dbuf(lbuf,dg_ijkconv(lijk,pface, &
274 ineighproc(pface)),pface)
275
276 call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lcurpar),pface, &
277 ineighproc(pface)),pface)
278
279 call pack_dbuf(lbuf,des_radius(lcurpar),pface)
280
281 = pijk(lcurpar,1) + icycoffset(pface,1)
282 call pack_dbuf(lbuf,li,pface)
283
284 = pijk(lcurpar,2) + icycoffset(pface,2)
285 call pack_dbuf(lbuf,lj,pface)
286
287 = pijk(lcurpar,3) + icycoffset(pface,3)
288 call pack_dbuf(lbuf,lk,pface)
289
290 call pack_dbuf(lbuf,funijk_proc(li,lj,lk, &
291 ineighproc(pface)),pface)
292
293 call pack_dbuf(lbuf,pijk(lcurpar,5),pface)
294
295 call pack_dbuf(lbuf, is_entering(lcurpar).or.is_entering_ghost(lcurpar), pface)
296
297 call pack_dbuf(lbuf, is_exiting(lcurpar).or.is_exiting_ghost(lcurpar), pface)
298
299 call pack_dbuf(lbuf,ro_sol(lcurpar),pface)
300
301 call pack_dbuf(lbuf,pvol(lcurpar),pface)
302
303 call pack_dbuf(lbuf,pmass(lcurpar),pface)
304
305 call pack_dbuf(lbuf,omoi(lcurpar),pface)
306
307 call pack_dbuf(lbuf,des_pos_new(lcurpar,:) + &
308 dcycl_offset(pface,:),pface)
309
310 call pack_dbuf(lbuf,des_vel_new(lcurpar,:),pface)
311
312 call pack_dbuf(lbuf,omega_new(lcurpar,:),pface)
313
314 call pack_dbuf(lbuf,fc(lcurpar,:),pface)
315
316 call pack_dbuf(lbuf,tow(lcurpar,:),pface)
317 IF(ENERGY_EQ) THEN
318
319 call pack_dbuf(lbuf,des_t_s(lcurpar),pface)
320
321 call pack_dbuf(lbuf,des_x_s(lcurpar,:),pface)
322 ENDIF
323
324 IF(DES_USR_VAR_SIZE > 0) &
325 call pack_dbuf(lbuf, des_usr_var(:,lcurpar),pface)
326
327 IF(PARTICLE_ORIENTATION) &
328 call pack_dbuf(lbuf,orientation(:,lcurpar),pface)
329
330
331 IF (DO_OLD) THEN
332
333 call pack_dbuf(lbuf,des_pos_old(lcurpar,:) + &
334 dcycl_offset(pface,:),pface)
335
336 call pack_dbuf(lbuf,des_vel_old(lcurpar,:),pface)
337
338 call pack_dbuf(lbuf,omega_old(lcurpar,:),pface)
339
340 call pack_dbuf(lbuf,des_acc_old(lcurpar,:),pface)
341
342 call pack_dbuf(lbuf,rot_acc_old(lcurpar,:),pface)
343 ENDIF
344
345 IF(DES_EXPLICITLY_COUPLED) THEN
346
347 call pack_dbuf(lbuf, drag_fc(lcurpar,:),pface)
348
349 IF(ENERGY_EQ) call pack_dbuf(lbuf,conv_qs(lcurpar),pface)
350
351 IF(ANY_SPECIES_EQ) call pack_dbuf(lbuf, rxns_qs(lcurpar),pface)
352 ENDIF
353
354
355
356 IF (MPPIC) THEN
357
358 call pack_dbuf(lbuf,des_stat_wt(lcurpar),pface)
359 call set_nonexistent(lcurpar)
360 pip = pip - 1
361
362
363
364 ELSE
365 if (is_entering(lcurpar)) then
366 call set_entering_ghost(lcurpar)
367 elseif (is_exiting(lcurpar)) then
368 call set_exiting_ghost(lcurpar)
369 else
370 call set_ghost(lcurpar)
371 endif
372 ighost_cnt = ighost_cnt + 1
373 END IF
374
375
376 (lcurpar,:) = 0.
377 lparcnt = lparcnt + 1
378 end do
379 end do
380
381
382
383
384 = lparcnt*iParticlePacketSize + ibufoffset
385 num_neighborlists_send_buf_loc = lbuf
386 lbuf = lbuf+1
387
388 num_neighborlists_to_send = 0
389 lcurpar = 1
390 do cc = 1, NEIGH_NUM
391 IF (0 .eq. NEIGHBORS(cc)) EXIT
392
393 IF (cc.eq.NEIGHBOR_INDEX(lcurpar)) THEN
394 lcurpar = lcurpar + 1
395 ENDIF
396
397
398 if (.not. going_to_send(lcurpar)) cycle
399
400
401
402 = neighbors(lcurpar)
403 if(is_nonexistent(lneigh)) cycle
404 if(is_exiting(lneigh)) cycle
405
406
407 call pack_dbuf(lbuf,iglobal_id(lcurpar),pface)
408
409 call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lcurpar),pface, &
410 ineighproc(pface)),pface)
411
412 call pack_dbuf(lbuf,iglobal_id(lneigh),pface)
413
414 call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lneigh),pface, &
415 ineighproc(pface)),pface)
416
417 call pack_dbuf(lbuf,PFT_NEIGHBOR(:,CC),pface)
418
419 = num_neighborlists_to_send + 1
420 enddo
421
422
423
424
425 = num_neighborlists_send_buf_loc
426
427 call pack_dbuf(lbuf,num_neighborlists_to_send,pface)
428
429 dsendbuf(1+mod(pface,2))%facebuf(1) = lparcnt
430 isendcnt(pface) = lparcnt*iParticlePacketSize + &
431 num_neighborlists_to_send*iPairPacketSize + ibufoffset + 3
432
433 deallocate(going_to_send)
434
435 RETURN
436 END SUBROUTINE DESMPI_PACK_PARCROSS
437
438
439
440
441 subroutine pack_db0(lbuf,idata,pface)
442 use desmpi, only: dSENDBUF
443 integer, intent(inout) :: lbuf
444 integer, intent(in) :: pface
445 double precision, intent(in) :: idata
446
447 dsendbuf(1+mod(pface,2))%facebuf(lbuf) = idata
448 lbuf = lbuf + 1
449
450 return
451 end subroutine pack_db0
452
453
454
455
456 subroutine pack_db1(lbuf,idata,pface)
457 use desmpi, only: dSENDBUF
458 integer, intent(inout) :: lbuf
459 integer, intent(in) :: pface
460 double precision, intent(in) :: idata(:)
461
462 integer :: lsize
463
464 lsize = size(idata)
465
466 dsendbuf(1+mod(pface,2))%facebuf(lbuf:lbuf+lsize-1) = idata
467 lbuf = lbuf + lsize
468
469 return
470 end subroutine pack_db1
471
472
473
474
475 subroutine pack_i0(lbuf,idata,pface)
476 use desmpi, only: dSENDBUF
477 integer, intent(inout) :: lbuf
478 integer, intent(in) :: pface
479 integer, intent(in) :: idata
480
481 dsendbuf(1+mod(pface,2))%facebuf(lbuf) = idata
482 lbuf = lbuf + 1
483
484 return
485 end subroutine pack_i0
486
487
488
489
490 subroutine pack_i1(lbuf,idata,pface)
491 use desmpi, only: dSENDBUF
492 integer, intent(inout) :: lbuf
493 integer, intent(in) :: pface
494 integer, intent(in) :: idata(:)
495
496 integer :: lsize
497
498 lsize = size(idata)
499
500 dsendbuf(1+mod(pface,2))%facebuf(lbuf:lbuf+lsize-1) = idata
501 lbuf = lbuf + lsize
502
503 return
504 end subroutine pack_i1
505
506
507
508
509 subroutine pack_l0(lbuf,ldata,pface)
510 use desmpi, only: dSENDBUF
511
512 integer, intent(inout) :: lbuf
513 integer, intent(in) :: pface
514 logical, intent(in) :: ldata
515
516 dsendbuf(1+mod(pface,2))%facebuf(lbuf) = merge(1.0, 0.0, ldata)
517 lbuf = lbuf + 1
518
519 return
520 end subroutine pack_l0
521
522 end module mpi_pack_des
523