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