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