File: RELATIVE:/../../../mfix.git/model/des/mpi_unpack_des_mod.f
1
2
3
4
5
6
7
8 MODULE MPI_UNPACK_DES
9
10 PRIVATE
11 PUBLIC :: DESMPI_UNPACK_PARCROSS, DESMPI_UNPACK_GHOSTPAR
12
13 interface unpack_dbuf
14 module procedure unpack_db0
15 module procedure unpack_db1
16 module procedure unpack_i0
17 module procedure unpack_i1
18 module procedure unpack_l0
19 end interface unpack_dbuf
20
21 CONTAINS
22
23
24
25
26
27
28
29 SUBROUTINE DESMPI_UNPACK_GHOSTPAR(pface)
30
31
32
33
34 use desmpi, only: iGhostPacketSize
35
36 use desmpi, only: iSPOT
37
38 use discretelement, only: iGHOST_UPDATED
39
40 use desmpi, only: dRECVBUF
41
42 use desmpi, only: iBUFOFFSET
43
44 use run, only: ENERGY_EQ
45
46 use desgrid, only: DG_IJKSIZE2
47
48 use discretelement, only: DG_PIJK, DG_PIJKPRV
49
50 use discretelement, only: iGLOBAL_ID
51
52 use discretelement, only: DES_POS_NEW, DES_POS_OLD
53
54 use discretelement, only: DES_VEL_NEW, DES_VEL_OLD
55
56 use discretelement, only: OMEGA_NEW, OMEGA_OLD
57
58 use des_thermo, only: DES_T_s_NEW, DES_T_s_OLD
59
60 use discretelement, only: DES_RADIUS, PVOL
61
62 use particle_filter, only: FILTER_SIZE
63
64 use particle_filter, only: FILTER_CELL, FILTER_WEIGHT
65
66 use discretelement, only: PIJK
67
68 use discretelement, only: DO_OLD
69
70 use discretelement, only: PIP
71
72 use discretelement, only: iGHOST_CNT
73
74 use discretelement, only: DES_USR_VAR, DES_USR_VAR_SIZE
75
76 use des_allocate
77
78
79
80 use constant, only: PI
81
82 use discretelement, only: DIMN
83
84 use functions, only: is_nonexistent
85 use functions, only: is_normal, set_normal
86 use functions, only: is_exiting, set_exiting, set_exiting_ghost
87 use functions, only: set_ghost
88
89 IMPLICIT NONE
90
91
92
93
94 INTEGER, INTENT(IN) :: PFACE
95
96
97
98 integer :: lcurpar,lparid,lprvijk,lparijk,lparcnt,ltot_ind
99 integer :: lbuf,llocpar,lnewcnt,lpicloc
100 logical,dimension(:),allocatable :: lfound
101 integer,dimension(:),allocatable :: lnewspot,lnewpic
102 logical :: tmp
103
104
105
106
107
108
109 = drecvbuf(1+mod(pface,2))%facebuf(1)
110 lnewcnt = lparcnt
111 allocate (lfound(lparcnt),lnewspot(lparcnt),lnewpic(dg_ijksize2))
112 lfound(:) = .false.
113 lnewspot(:) =0
114 lnewpic = 0
115
116 do lcurpar = 1,lparcnt
117 lbuf = (lcurpar-1)*iGhostPacketSize+ibufoffset
118
119
120 call unpack_dbuf(lbuf,lparid,pface)
121
122 call unpack_dbuf(lbuf,lparijk,pface)
123
124 call unpack_dbuf(lbuf,lprvijk,pface)
125
126
127
128
129
130 (lcurpar) = locate_par(lparid,lprvijk,llocpar)
131 if (lparijk .ne. lprvijk .and. .not.lfound(lcurpar)) then
132 lfound(lcurpar) = locate_par(lparid,lparijk,llocpar)
133 endif
134
135 if(lfound(lcurpar)) then
136
137 (llocpar) = lparijk
138 dg_pijkprv(llocpar) = lprvijk
139
140
141 call unpack_dbuf(lbuf,des_radius(llocpar),pface)
142
143 call unpack_dbuf(lbuf,pijk(llocpar,5),pface)
144
145 call unpack_dbuf(lbuf,des_pos_new(1:dimn,llocpar),pface)
146
147 call unpack_dbuf(lbuf,des_vel_new(1:dimn,llocpar),pface)
148
149 call unpack_dbuf(lbuf,omega_new(1:3,llocpar),pface)
150
151 call unpack_dbuf(lbuf,tmp,pface)
152 if (tmp) call set_exiting_ghost(llocpar)
153
154 IF(ENERGY_EQ) &
155 call unpack_dbuf(lbuf,des_t_s_new(llocpar),pface)
156
157 IF(DES_USR_VAR_SIZE > 0) &
158 call unpack_dbuf(lbuf,des_usr_var(1:3,llocpar),pface)
159
160 IF(FILTER_SIZE > 0) THEN
161 call unpack_dbuf(lbuf,filter_cell(:,llocpar),pface)
162 call unpack_dbuf(lbuf,filter_weight(:,llocpar),pface)
163 ENDIF
164
165
166
167 (llocpar) = (4.0D0/3.0D0)*PI*DES_RADIUS(llocpar)**3
168
169 (llocpar) = .true.
170 lnewcnt = lnewcnt-1
171
172
173 IF (DO_OLD) THEN
174 des_pos_old(:,llocpar)= des_pos_new(:,llocpar)
175 des_vel_old(:,llocpar)= des_vel_new(:,llocpar)
176 if(ENERGY_EQ)des_t_s_old(llocpar)= des_t_s_new(llocpar)
177 omega_old(:,llocpar)= omega_new(:,llocpar)
178 ENDIF
179
180 else
181 lnewpic(lparijk) = lnewpic(lparijk) + 1
182 endif
183 enddo
184
185
186 if(lnewcnt > 0) then
187 call PARTICLE_GROW(pip+lnewcnt)
188 ighost_cnt = ighost_cnt + lnewcnt
189 pip = pip + lnewcnt
190 do lcurpar = 1,lparcnt
191 if(lfound(lcurpar)) cycle
192 lbuf = (lcurpar-1)*iGhostPacketSize+ibufoffset
193
194
195 call unpack_dbuf(lbuf,lparid,pface)
196
197 call unpack_dbuf(lbuf,lparijk,pface)
198
199 call unpack_dbuf(lbuf,lprvijk,pface)
200
201 do while(.not.is_nonexistent(ispot))
202 ispot = ispot + 1
203 enddo
204
205 call set_ghost(ispot)
206 iglobal_id(ispot) = lparid
207 dg_pijk(ispot) = lparijk
208 dg_pijkprv(ispot) = lprvijk
209
210 call unpack_dbuf(lbuf,des_radius(ispot),pface)
211
212 call unpack_dbuf(lbuf,pijk(ispot,5),pface)
213
214 call unpack_dbuf(lbuf,des_pos_new(1:dimn,ispot),pface)
215
216 call unpack_dbuf(lbuf,des_vel_new(1:dimn,ispot),pface)
217
218 call unpack_dbuf(lbuf,omega_new(1:dimn,ispot),pface)
219
220 call unpack_dbuf(lbuf,tmp,pface)
221 if (tmp) call set_exiting_ghost(ispot)
222
223 IF(ENERGY_EQ) &
224 call unpack_dbuf(lbuf,des_t_s_new(ispot),pface)
225
226 IF(DES_USR_VAR_SIZE > 0)&
227 call unpack_dbuf(lbuf,des_usr_var(:,ispot),pface)
228
229 IF(FILTER_SIZE > 0) THEN
230 call unpack_dbuf(lbuf,filter_cell(:,ispot),pface)
231 call unpack_dbuf(lbuf,filter_weight(:,ispot),pface)
232 ENDIF
233
234 ighost_updated(ispot) = .true.
235 lnewspot(lcurpar) = ispot
236
237 PVOL(ispot) = (4.0D0/3.0D0)*PI*DES_RADIUS(ispot)**3
238
239 IF (DO_OLD) THEN
240 des_pos_old(1:dimn,ispot) = des_pos_new(1:dimn,ispot)
241 des_vel_old(1:dimn,ispot) = des_vel_new(1:dimn,ispot)
242 omega_old(1:3,ispot) = omega_new(1:3,ispot)
243 if(ENERGY_EQ) des_t_s_old(ispot) = des_t_s_new(ispot)
244 ENDIF
245 enddo
246 endif
247
248
249 deallocate (lfound,lnewspot,lnewpic)
250
251 end subroutine desmpi_unpack_ghostpar
252
253
254
255
256
257
258
259 SUBROUTINE DESMPI_UNPACK_PARCROSS(pface)
260
261
262
263
264 use desmpi, only: iParticlePacketSize
265
266 use desmpi, only: iSPOT
267
268 use desmpi, only: dRECVBUF
269
270 use desmpi, only: iBUFOFFSET
271
272 use run, only: ENERGY_EQ
273
274 use run, only: ANY_SPECIES_EQ
275
276 use mfix_pic, only: MPPIC
277
278 use discretelement, only: DG_PIJK, DG_PIJKPRV
279
280 use desmpi, only: iNEIGHPROC
281
282 use mfix_pic, only: DES_STAT_WT
283
284 use discretelement, only: iGLOBAL_ID
285
286 use discretelement, only: DES_POS_NEW, DES_POS_OLD
287
288 use discretelement, only: DES_VEL_NEW, DES_VEL_OLD
289
290 use discretelement, only: OMEGA_NEW, OMEGA_OLD
291
292 use discretelement, only: PARTICLE_ORIENTATION,ORIENTATION
293
294 use discretelement, only: DES_RADIUS, PVOL, RO_SOL, PMASS
295
296 use discretelement, only: DES_ACC_OLD, ROT_ACC_OLD
297
298 use des_rxns, only: DES_X_s
299
300 use des_thermo, only: DES_T_s_NEW, DES_T_s_OLD
301
302 use particle_filter, only: FILTER_WEIGHT
303
304 use discretelement, only: FC, TOW
305
306 use discretelement, only: OMOI
307
308 use discretelement, only: PIJK
309
310 use discretelement, only: DO_OLD
311
312 use discretelement, only: PIP
313
314 use discretelement, only: iGHOST_CNT
315
316 use discretelement, only: DES_EXPLICITLY_COUPLED
317
318 use discretelement, only: DRAG_FC
319
320 use discretelement, only: DES_USR_VAR, DES_USR_VAR_SIZE
321
322 use discretelement, only: PFT_NEIGHBOR
323
324 use discretelement, only: DIMN
325
326 use compar, only: myPE
327
328
329
330 use des_allocate
331 use desmpi_wrapper, only: DES_MPI_STOP
332 use functions, only: IS_NORMAL, IS_NONEXISTENT
333 use functions, only: SET_ENTERING, SET_EXITING, SET_NORMAL
334
335 implicit none
336
337
338
339
340 INTEGER, INTENT(IN) :: PFACE
341
342
343
344 integer :: lcurpar,lparcnt,llocpar,lparid,lparijk,lprvijk
345 integer :: lneigh,lcontactindx,lcontactid,lcontact,&
346 lneighid,lneighijk
347 logical :: lfound
348 integer :: lbuf,lcount
349 logical :: lneighfound
350 integer :: cc,kk,num_neighborlists_sent,nn
351
352 logical :: tmp
353
354
355
356 = drecvbuf(1+mod(pface,2))%facebuf(1)
357
358
359 call PARTICLE_GROW(pip+lparcnt)
360
361 do lcurpar =1,lparcnt
362
363 lfound = .false.
364 lbuf = (lcurpar-1)*iParticlePacketSize + ibufoffset
365
366 call unpack_dbuf(lbuf,lparid,pface)
367
368 call unpack_dbuf(lbuf,lparijk,pface)
369
370 call unpack_dbuf(lbuf,lprvijk,pface)
371
372
373
374
375 IF(MPPIC) THEN
376 DO WHILE(.NOT.IS_NONEXISTENT(ISPOT))
377 ISPOT = ISPOT + 1
378 ENDDO
379 lLOCPAR = iSPOT
380 iGLOBAL_ID(lLOCPAR) = lPARID
381 PIP = PIP + 1
382
383
384
385
386
387 ELSE
388 lFOUND = LOCATE_PAR(lPARID,lPRVIJK,lLOCPAR)
389 IF (.NOT. lFOUND) THEN
390 lFOUND = exten_locate_par(lPARID, lPARIJK, lLOCPAR)
391 IF(.NOT.lFOUND) THEN
392 WRITE(*,1000) iNEIGHPROC(PFACE), MYPE, lPARID
393 CALL DES_MPI_STOP
394 ENDIF
395 ENDIF
396 iGHOST_CNT = iGHOST_CNT - 1
397 ENDIF
398
399 1000 FORMAT(2/1X,72('*'),/1x,'From: DESMPI_UNPACK_PARCROSS: ',/ &
400 ' Error 1000: Unable to match particles crossing processor ', &
401 'boundaries.',/3x,'Source Proc: ',I9,' ---> Destination ', &
402 'Proc: ', I9,/3x,'Global Particle ID: ',I12,/1x,72('*'))
403
404
405 call set_normal(llocpar)
406 dg_pijk(llocpar) = lparijk
407 dg_pijkprv(llocpar) = lprvijk
408
409 call unpack_dbuf(lbuf,des_radius(llocpar),pface)
410
411 call unpack_dbuf(lbuf,pijk(llocpar,:),pface)
412
413 call unpack_dbuf(lbuf,tmp,pface)
414 if (tmp) CALL SET_ENTERING(llocpar)
415
416 call unpack_dbuf(lbuf,tmp,pface)
417 if (tmp) CALL SET_EXITING(llocpar)
418
419 call unpack_dbuf(lbuf,ro_sol(llocpar),pface)
420
421 call unpack_dbuf(lbuf,pvol(llocpar),pface)
422
423 call unpack_dbuf(lbuf,pmass(llocpar),pface)
424
425 call unpack_dbuf(lbuf,omoi(llocpar),pface)
426
427 call unpack_dbuf(lbuf,des_pos_new(:,llocpar),pface)
428
429 call unpack_dbuf(lbuf,des_vel_new(:,llocpar),pface)
430
431 call unpack_dbuf(lbuf,omega_new(:,llocpar),pface)
432
433 call unpack_dbuf(lbuf,fc(:,llocpar),pface)
434
435 call unpack_dbuf(lbuf,tow(:,llocpar),pface)
436
437 IF(ENERGY_EQ) &
438 call unpack_dbuf(lbuf,des_t_s_new(llocpar),pface)
439
440 IF(ANY_SPECIES_EQ) &
441 call unpack_dbuf(lbuf,des_x_s(llocpar,:),pface)
442
443 IF(DES_EXPLICITLY_COUPLED) &
444 call unpack_dbuf(lbuf,drag_fc(:,llocpar),pface)
445
446 IF(DES_USR_VAR_SIZE > 0) &
447 call unpack_dbuf(lbuf,des_usr_var(:,llocpar),pface)
448
449 IF(PARTICLE_ORIENTATION) &
450 call unpack_dbuf(lbuf,orientation(:,llocpar),pface)
451
452
453 IF (DO_OLD) THEN
454
455 call unpack_dbuf(lbuf,des_pos_old(:,llocpar),pface)
456
457 call unpack_dbuf(lbuf,des_vel_old(:,llocpar),pface)
458
459 call unpack_dbuf(lbuf,omega_old(:,llocpar),pface)
460
461 call unpack_dbuf(lbuf,des_acc_old(:,llocpar),pface)
462
463 call unpack_dbuf(lbuf,rot_acc_old(:,llocpar),pface)
464
465 IF(ENERGY_EQ) &
466 call unpack_dbuf(lbuf,des_t_s_old(llocpar),pface)
467 ENDIF
468
469 IF(MPPIC) call unpack_dbuf(lbuf,des_stat_wt(llocpar),pface)
470
471 end do
472
473
474 = lparcnt*iParticlePacketSize + ibufoffset
475 call unpack_dbuf(lbuf,num_neighborlists_sent,pface)
476
477 do nn = 1, num_neighborlists_sent
478
479 call unpack_dbuf(lbuf,lparid,pface)
480
481 call unpack_dbuf(lbuf,lparijk,pface)
482
483
484 if (.not. locate_par(lparid,lparijk,llocpar)) then
485 if (.not. exten_locate_par(lparid,lparijk,llocpar)) then
486 print *,"at buffer location",lbuf," pface = ",pface
487 print *,"COULD NOT FIND PARTICLE ",lparid," IN IJK ",lparijk
488 call des_mpi_stop
489 endif
490 endif
491
492 call unpack_dbuf(lbuf,lneighid,pface)
493
494 call unpack_dbuf(lbuf,lneighijk,pface)
495
496
497 if (.not. locate_par(lneighid,lneighijk,lneigh)) then
498 if (.not. exten_locate_par(lneighid,lparijk,lneigh)) then
499 print *," "
500 print *," "
501 print *," fail on ", myPE
502 print *,"at buffer location",lbuf," pface = ",pface
503 print *,"COULD NOT FIND NEIGHBOR ",lneighid," IN IJK ",lneighijk
504 call des_mpi_stop
505 endif
506 endif
507
508
509
510
511 = add_pair(llocpar,lneigh)
512
513 call unpack_dbuf(lbuf,pft_neighbor(:,cc),pface)
514 enddo
515
516 END SUBROUTINE desmpi_unpack_parcross
517
518
519
520
521
522
523
524
525
526 LOGICAL FUNCTION LOCATE_PAR(pGLOBALID, pIJK, pLOCALNO)
527
528 use discretelement, only: iGLOBAL_ID
529 use desgrid, only: DG_IJKStart2, DG_IJKEnd2
530 use discretelement, only: dg_pic
531
532 implicit none
533
534
535
536
537 INTEGER, INTENT(IN) :: pGlobalID
538
539 INTEGER, INTENT(IN) :: pIJK
540
541 INTEGER, INTENT(OUT) :: pLocalNO
542
543
544
545 INTEGER :: lpicloc, lcurpar
546
547
548 = .false.
549
550
551 if(pIJK < dg_ijkstart2 .or. pIJK > dg_ijkend2) RETURN
552
553
554
555
556 DO lpicloc = 1,dg_pic(pijk)%isize
557 lcurpar = dg_pic(pijk)%p(lpicloc)
558 IF(iGLOBAL_ID(lcurpar) == pGlobalID) THEN
559 plocalno = lcurpar
560 locate_par = .true.
561 RETURN
562 ENDIF
563 ENDDO
564
565 RETURN
566 end function locate_par
567
568
569
570
571
572
573
574
575
576 LOGICAL FUNCTION EXTEN_LOCATE_PAR(pGlobalID, pIJK, pLocalNO)
577
578 use discretelement, only: iGLOBAL_ID, dg_pic
579 use desgrid, only: DG_IJKStart2, DG_IJKEnd2
580 use desgrid, only: dg_Iof_LO, DG_Jof_LO, DG_Kof_LO
581 use geometry, only: NO_K
582
583 use desgrid, only: dg_funijk
584
585 implicit none
586
587
588
589
590 INTEGER, INTENT(IN) :: pGlobalId
591
592 INTEGER, INTENT(IN) :: pIJK
593
594 INTEGER, INTENT(OUT) :: pLocalNo
595
596
597
598
599 INTEGER :: lijk, li, lj, lk, lic, ljc, lkc, lkoffset
600 INTEGER :: lpicloc,lcurpar
601
602 exten_locate_par = .false.
603
604 lic = dg_iof_lo(pijk)
605 ljc = dg_jof_lo(pijk)
606 lkc = dg_kof_lo(pijk)
607 lkoffset = merge(0, 1, NO_K)
608 DO lk = lkc-lkoffset,lkc+lkoffset
609 DO lj = ljc-1,ljc+1
610 DO li = lic-1,lic+1
611 lijk = dg_funijk(li,lj,lk)
612 IF (lijk .lt. dg_ijkstart2 .or. lijk .gt. dg_ijkend2) CYCLE
613 DO lpicloc = 1, dg_pic(lijk)%isize
614 lcurpar = dg_pic(lijk)%p(lpicloc)
615 IF (iglobal_id(lcurpar) .eq. pglobalid) THEN
616 plocalno = lcurpar
617 exten_locate_par = .true.
618 RETURN
619 END IF
620 END DO
621 END DO
622 END DO
623 END DO
624
625 RETURN
626 END FUNCTION EXTEN_LOCATE_PAR
627
628
629
630
631 subroutine unpack_db0(lbuf,idata,pface)
632 use desmpi, only: dRECVBUF
633 integer, intent(inout) :: lbuf
634 integer, intent(in) :: pface
635 double precision, intent(inout) :: idata
636
637 idata = drecvbuf(1+mod(pface,2))%facebuf(lbuf)
638 lbuf = lbuf + 1
639
640 return
641 end subroutine unpack_db0
642
643
644
645
646 subroutine unpack_db1(lbuf,idata,pface)
647 use desmpi, only: dRECVBUF
648 integer, intent(inout) :: lbuf
649 integer, intent(in) :: pface
650 double precision, intent(inout) :: idata(:)
651
652 integer :: lsize
653
654 lsize = size(idata)
655
656 idata = drecvbuf(1+mod(pface,2))%facebuf(lbuf:lbuf+lsize-1)
657 lbuf = lbuf + lsize
658
659 return
660 end subroutine unpack_db1
661
662
663
664
665 subroutine unpack_i0(lbuf,idata,pface)
666 use desmpi, only: dRECVBUF
667 integer, intent(inout) :: lbuf
668 integer, intent(in) :: pface
669 integer, intent(inout) :: idata
670
671 idata = drecvbuf(1+mod(pface,2))%facebuf(lbuf)
672 lbuf = lbuf + 1
673
674 return
675 end subroutine unpack_i0
676
677
678
679
680 subroutine unpack_i1(lbuf,idata,pface)
681 use desmpi, only: dRECVBUF
682 integer, intent(inout) :: lbuf
683 integer, intent(in) :: pface
684 integer, intent(inout) :: idata(:)
685
686 integer :: lsize
687
688 lsize = size(idata)
689
690 idata = drecvbuf(1+mod(pface,2))%facebuf(lbuf:lbuf+lsize-1)
691 lbuf = lbuf + lsize
692
693 return
694 end subroutine unpack_i1
695
696
697
698
699 subroutine unpack_l0(lbuf,idata,pface)
700 use desmpi, only: dRECVBUF
701 integer, intent(inout) :: lbuf
702 integer, intent(in) :: pface
703 logical, intent(inout) :: idata
704
705 idata = merge(.true.,.false.,0.5<drecvbuf(1+mod(pface,2))%facebuf(lbuf))
706 lbuf = lbuf + 1
707
708 return
709 end subroutine unpack_l0
710
711 end module mpi_unpack_des
712