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