MFIX  2016-1
mpi_unpack_des_mod.f
Go to the documentation of this file.
1 !----------------------------------------------------------------------!
2 ! Module: MPI_UNPACK_DES !
3 ! Author: Pradeep Gopalakrishnan !
4 ! !
5 ! Purpose: Contains routines for unpacking real and ghost particles !
6 ! from the MPI recv buffers. !
7 !----------------------------------------------------------------------!
9 
10  PRIVATE
12 
13  interface unpack_dbuf
14  module procedure unpack_db0 ! real scalars
15  module procedure unpack_db1 ! real arrays
16  module procedure unpack_i0 ! integer scalars
17  module procedure unpack_i1 ! integer arrays
18  module procedure unpack_l0 ! logical scalars
19  end interface unpack_dbuf
20 
21  CONTAINS
22 
23 !----------------------------------------------------------------------!
24 ! Subroutine: DESMPI_UNPACK_GHOSTPAR !
25 ! Author: Pradeep Gopalakrishnan !
26 ! !
27 ! Purpose: Unpacks ghost particle from the recv buffer. !
28 !----------------------------------------------------------------------!
29  SUBROUTINE desmpi_unpack_ghostpar(pface)
30 
31 ! Global Variables:
32 !---------------------------------------------------------------------//
33 ! Size of Particle data packet
34  use desmpi, only: ighostpacketsize
35 ! Index of last particle added to this process.
36  use desmpi, only: ispot
37 ! Flag indicating that the ghost particle was updated
38  use discretelement, only: ighost_updated
39 ! The MPI receive buffer
40  use desmpi, only: drecvbuf
41 ! Buffer offset
42  use desmpi, only: ibufoffset
43 ! Runtime flag for solving the energy equations
44  use run, only: energy_eq
45 ! Dimensions of DES grid
46  use desgrid, only: dg_ijksize2
47 ! DES grid cell containing each particle: current/previous
48  use discretelement, only: dg_pijk, dg_pijkprv
49 ! The global ID for each particle
50  use discretelement, only: iglobal_id
51 ! Particle positions: current/previous
52  use discretelement, only: des_pos_new, des_pos_old
53 ! Particle tangential velocities: current/previous
54  use discretelement, only: des_vel_new, des_vel_old
55 ! Particle rotational velocities: current/previous
56  use discretelement, only: omega_new, omega_old
57 ! Particle tempertures
58  use des_thermo, only: des_t_s
59 ! Particle radius, volume
60  use discretelement, only: des_radius, pvol
61 ! Map to fluid grid cells and solids phase (I,J,K,IJK,M)
62  use discretelement, only: pijk
63 ! Flag to send/recv old (previous) values
64  use discretelement, only: do_old
65 ! Number of particles on the process (max particle array size)
66  use discretelement, only: pip
67 ! Number of ghost particles on the current process
68  use discretelement, only: ighost_cnt
69 ! User-defined variables for each particle.
70  use discretelement, only: des_usr_var, des_usr_var_size
71 
72  use des_allocate
73 
74 ! Global Constants:
75 !---------------------------------------------------------------------//
76  use constant, only: pi
77 ! Dimension of particle spatial arrays.
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 ! Dummy arguments:
89 !---------------------------------------------------------------------//
90 ! Processor boundary being packed (Top/Bottom/North/South/East/West)
91  INTEGER, INTENT(IN) :: PFACE
92 
93 ! Local variables
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 ! unpack the particles:
103 ! if it already exists update the position
104 ! if not and do_nsearch is true then add to the particle array
105 
106  lparcnt = 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 ! 1) Global ID
117  call unpack_dbuf(lbuf,lparid,pface)
118 ! 2) DES Grid IJK
119  call unpack_dbuf(lbuf,lparijk,pface)
120 ! 3) DES Grid IJK - Previous
121  call unpack_dbuf(lbuf,lprvijk,pface)
122 
123 ! Determine if this particle already exists on this process as a
124 ! ghost particle. If so, (lfound), then the current infomration is
125 ! updated on the current process. Otherwise (.NOT.lfound) a new
126 ! ghost particle is created on this process.
127  lfound(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 ! Store the local variables
134  dg_pijk(llocpar) = lparijk
135  dg_pijkprv(llocpar) = lprvijk
136 
137 ! 4) Radious
138  call unpack_dbuf(lbuf,des_radius(llocpar),pface)
139 ! 5) Phase index
140  call unpack_dbuf(lbuf,pijk(llocpar,5),pface)
141 ! 6) Position
142  call unpack_dbuf(lbuf,des_pos_new(llocpar,1:dimn),pface)
143 ! 7) Translational Velocity
144  call unpack_dbuf(lbuf,des_vel_new(llocpar,1:dimn),pface)
145 ! 8) Rotational Velocity
146  call unpack_dbuf(lbuf,omega_new(llocpar,1:3),pface)
147 ! 9) Exiting particle flag
148  call unpack_dbuf(lbuf,tmp,pface)
149  if (tmp) call set_exiting_ghost(llocpar)
150 ! 10) Temperature
151  IF(energy_eq) &
152  call unpack_dbuf(lbuf,des_t_s(llocpar),pface)
153 ! 11) User Variables
154  IF(des_usr_var_size > 0) &
155  call unpack_dbuf(lbuf,des_usr_var(:,llocpar),pface)
156 
157 ! Calculate the volume of the ghost particle.
158  pvol(llocpar) = (4.0d0/3.0d0)*pi*des_radius(llocpar)**3
159 ! Flag that the ghost particle was updated.
160  ighost_updated(llocpar) = .true.
161  lnewcnt = lnewcnt-1
162 
163 ! Copy the current value to the previous value if needed.
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 ! iAdd new ghost particles
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 ! 1) Global particle ID
186  call unpack_dbuf(lbuf,lparid,pface)
187 ! 2) DES grid IJK
188  call unpack_dbuf(lbuf,lparijk,pface)
189 ! 3) DES grid IJK - Previous
190  call unpack_dbuf(lbuf,lprvijk,pface)
191 ! Locate the first open space in the particle array.
192  do while(.not.is_nonexistent(ispot))
193  ispot = ispot + 1
194  enddo
195 ! Set the flags for the ghost particle and store the local variables.
196  call set_ghost(ispot)
197  iglobal_id(ispot) = lparid
198  dg_pijk(ispot) = lparijk
199  dg_pijkprv(ispot) = lprvijk
200 ! 4) Radius
201  call unpack_dbuf(lbuf,des_radius(ispot),pface)
202 ! 5) Phase index
203  call unpack_dbuf(lbuf,pijk(ispot,5),pface)
204 ! 6) Position
205  call unpack_dbuf(lbuf,des_pos_new(ispot,1:dimn),pface)
206 ! 7) Translational velocity
207  call unpack_dbuf(lbuf,des_vel_new(ispot,1:dimn),pface)
208 ! 8) Rotational velocity
209  call unpack_dbuf(lbuf,omega_new(ispot,1:dimn),pface)
210 ! 9) Exiting particle flag
211  call unpack_dbuf(lbuf,tmp,pface)
212  if (tmp) call set_exiting_ghost(ispot)
213 ! 10) Temperature.
214  IF(energy_eq) &
215  call unpack_dbuf(lbuf,des_t_s(ispot),pface)
216 ! 11) User varaible
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 !deallocate temporary variablies
234  deallocate (lfound,lnewspot,lnewpic)
235 
236  end subroutine desmpi_unpack_ghostpar
237 
238 !----------------------------------------------------------------------!
239 ! Subroutine: DESMPI_UNPACK_PARCROSS !
240 ! Author: Pradeep Gopalakrishnan !
241 ! !
242 ! Purpose: Unpacks real particle from the recv buffer. !
243 !----------------------------------------------------------------------!
244  SUBROUTINE desmpi_unpack_parcross(pface)
246 ! Global Variables:
247 !---------------------------------------------------------------------//
248 ! Size of ghost particle data packet
249  use desmpi, only: iparticlepacketsize
250 ! Index of last particle added to this process.
251  use desmpi, only: ispot
252 ! The MPI receive buffer
253  use desmpi, only: drecvbuf
254 ! Buffer offset
255  use desmpi, only: ibufoffset
256 ! Runtime flag for solving the energy equations
257  use run, only: energy_eq
258 ! Runtime flag for solving species equations
259  use run, only: any_species_eq
260 ! Runtime flag for MPPIC solids
261  use mfix_pic, only: mppic
262 ! DES grid cell containing each particle: current/previous
263  use discretelement, only: dg_pijk, dg_pijkprv
264 ! The neighbor processor's rank
265  use desmpi, only: ineighproc
266 ! The statistical weight of each particle.
267  use mfix_pic, only: des_stat_wt
268 ! The global ID for each particle
269  use discretelement, only: iglobal_id
270 ! Particle positions: current/previous
271  use discretelement, only: des_pos_new, des_pos_old
272 ! Particle tangential velocities: current/previous
273  use discretelement, only: des_vel_new, des_vel_old
274 ! Particle rotational velocities: current/previous
275  use discretelement, only: omega_new, omega_old
276 !Particle orientation
277  use discretelement, only: particle_orientation,orientation
278 ! Particle radius, volume, density, mass
279  use discretelement, only: des_radius, pvol, ro_sol, pmass
280 ! Previous value for particle acceleration (tangential/rotational)
281  use discretelement, only: des_acc_old, rot_acc_old
282 ! Particle species composition
283  use des_rxns, only: des_x_s
284 ! Particle tempertures.
285  use des_thermo, only: des_t_s
286 ! Force arrays acting on the particle
287  use discretelement, only: fc, tow
288 ! One of the moment of inertia
289  use discretelement, only: omoi
290 ! Map to fluid grid cells and solids phase (I,J,K,IJK,M)
291  use discretelement, only: pijk
292 ! Flag to send/recv old (previous) values
293  use discretelement, only: do_old
294 ! Number of particles on the process (max particle array size)
295  use discretelement, only: pip
296 ! Number of ghost particles on the current process
297  use discretelement, only: ighost_cnt
298 ! Flag indicating the the fluid-particle drag is explicitly coupled.
299  use discretelement, only: des_explicitly_coupled
300 ! Explicit fluid-particle drag force
301  use discretelement, only: drag_fc
302 ! Explict convection and HOR
303  use des_thermo, only: conv_qs, rxns_qs
304 ! User-defined variables for each particle.
305  use discretelement, only: des_usr_var, des_usr_var_size
306 ! Neighbor collision history information
307  use discretelement, only: pft_neighbor
308 ! Dimension of particle spatial arrays.
309  use discretelement, only: dimn
310 ! The ID of the current process
311  use compar, only: mype
312 
313 ! Module Procedures:
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 ! Dummy arguments:
324 !---------------------------------------------------------------------//
325 ! Processor boundary being packed (Top/Bottom/North/South/East/West)
326  INTEGER, INTENT(IN) :: PFACE
327 
328 ! Local variables
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 ! loop through particles and locate them and make changes
342  lparcnt = drecvbuf(1+mod(pface,2))%facebuf(1)
343 
344 ! if mppic make sure enough space available
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 ! 1) Global ID
353  call unpack_dbuf(lbuf,lparid,pface)
354 ! 2) DES Grid IJK
355  call unpack_dbuf(lbuf,lparijk,pface)
356 ! 3) DES grid IJK - previous
357  call unpack_dbuf(lbuf,lprvijk,pface)
358 
359 ! PIC particles are always 'new' to the receiving process. Find the
360 ! first available array position and store the global ID. Increment
361 ! the PIP counter to include the new particle.
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 ! A DEM particle should already exist on the current processor as a
371 ! ghost particle. Match the sent particle to the local ghost particle
372 ! by matching the global IDs. Decrement the iGHOST_CNT counter to
373 ! account for the switch from ghost to real particle.
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 ! convert the local particle from ghost to existing and update its position
392  call set_normal(llocpar)
393  dg_pijk(llocpar) = lparijk
394  dg_pijkprv(llocpar) = lprvijk
395 ! 4) Radius
396  call unpack_dbuf(lbuf,des_radius(llocpar),pface)
397 ! 5-9) Fluid cell I,J,K,IJK, and solids phase index
398  call unpack_dbuf(lbuf,pijk(llocpar,:),pface)
399 ! 10) Entering particle flag.
400  call unpack_dbuf(lbuf,tmp,pface)
401  if (tmp) CALL set_entering(llocpar)
402 ! 11) Exiting particle flag.
403  call unpack_dbuf(lbuf,tmp,pface)
404  if (tmp) CALL set_exiting(llocpar)
405 ! 12) Density
406  call unpack_dbuf(lbuf,ro_sol(llocpar),pface)
407 ! 13) Volume
408  call unpack_dbuf(lbuf,pvol(llocpar),pface)
409 ! 14) Mass
410  call unpack_dbuf(lbuf,pmass(llocpar),pface)
411 ! 15) 1/Moment of Inertia
412  call unpack_dbuf(lbuf,omoi(llocpar),pface)
413 ! 16) Position with cyclic shift
414  call unpack_dbuf(lbuf,des_pos_new(llocpar,:),pface)
415 ! 17) Translational velocity
416  call unpack_dbuf(lbuf,des_vel_new(llocpar,:),pface)
417 ! 18) Rotational velocity
418  call unpack_dbuf(lbuf,omega_new(llocpar,:),pface)
419 ! 19) Accumulated translational forces
420  call unpack_dbuf(lbuf,fc(llocpar,:),pface)
421 ! 20) Accumulated torque forces
422  call unpack_dbuf(lbuf,tow(llocpar,:),pface)
423  IF(energy_eq) THEN
424 ! 21) Temperature
425  call unpack_dbuf(lbuf,des_t_s(llocpar),pface)
426 ! 22) Species composition
427  call unpack_dbuf(lbuf,des_x_s(llocpar,:),pface)
428  ENDIF
429 ! 23) User defined variable
430  IF(des_usr_var_size > 0) &
431  call unpack_dbuf(lbuf,des_usr_var(:,llocpar),pface)
432 ! 24) Particle orientation
433  IF(particle_orientation) &
434  call unpack_dbuf(lbuf,orientation(:,llocpar),pface)
435 
436 ! -- Higher order integration variables
437  IF (do_old) THEN
438 ! 25) Position (previous)
439  call unpack_dbuf(lbuf,des_pos_old(llocpar,:),pface)
440 ! 26) Translational velocity (previous)
441  call unpack_dbuf(lbuf,des_vel_old(llocpar,:),pface)
442 ! 27) Rotational velocity (previous)
443  call unpack_dbuf(lbuf,omega_old(llocpar,:),pface)
444 ! 28) Translational acceleration (previous)
445  call unpack_dbuf(lbuf,des_acc_old(llocpar,:),pface)
446 ! 29) Rotational acceleration (previous)
447  call unpack_dbuf(lbuf,rot_acc_old(llocpar,:),pface)
448  ENDIF
449 
450  IF(des_explicitly_coupled) THEN
451 ! 30) Explicit drag force
452  call unpack_dbuf(lbuf,drag_fc(llocpar,:),pface)
453 ! 31) Explicit convective heat transfer
454  IF(energy_eq) call unpack_dbuf(lbuf,conv_qs(llocpar),pface)
455 ! 32) Explicit heat of reaction
456  IF(any_species_eq) call unpack_dbuf(lbuf,rxns_qs(llocpar),pface)
457  ENDIF
458 
459 ! 33) Statistical weight
460  IF(mppic) call unpack_dbuf(lbuf,des_stat_wt(llocpar),pface)
461 
462  end do
463 
464 ! 34) Number of pair datasets
465  lbuf = lparcnt*iparticlepacketsize + ibufoffset
466  call unpack_dbuf(lbuf,num_neighborlists_sent,pface)
467 
468  do nn = 1, num_neighborlists_sent
469 ! 35) Global ID of packed particle.
470  call unpack_dbuf(lbuf,lparid,pface)
471 ! 36) DES grid IJK of cell receiving the particle.
472  call unpack_dbuf(lbuf,lparijk,pface)
473 
474 ! Locate the particle on the current process.
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 ! 37) Global ID of neighbor particle.
483  call unpack_dbuf(lbuf,lneighid,pface)
484 ! 38) DES grid IJK of cell containing the neighbor particle.
485  call unpack_dbuf(lbuf,lneighijk,pface)
486 
487 ! Locate the neighbor particle on the current process.
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 ! If the neighbor particle is a 'real' particle on this processor, then
500 ! the pair data may already exist. Check before adding it.
501 ! Create a new neighbor pair if it was not matched to an exiting pair.
502  cc = add_pair(llocpar,lneigh)
503 ! 39) Tangential collision history.
504  call unpack_dbuf(lbuf,pft_neighbor(:,cc),pface)
505  enddo
506 
507  END SUBROUTINE desmpi_unpack_parcross
508 
509 !----------------------------------------------------------------------!
510 ! Function: LOCATE_PAR !
511 ! Author: Pradeep Gopalakrishnan !
512 ! !
513 ! Purpose: Return the local index of the particle matching the passed !
514 ! global ID. The function returns TRUE if the particle is matched, !
515 ! otherwise it returns FALSE. !
516 !----------------------------------------------------------------------!
517  LOGICAL FUNCTION locate_par(pGLOBALID, pIJK, pLOCALNO)
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 ! Dummy arguments:
526 !---------------------------------------------------------------------//
527 ! Global ID of the particle
528  INTEGER, INTENT(IN) :: pGlobalID
529 ! IJK of DES grid cell containing the particle
530  INTEGER, INTENT(IN) :: pIJK
531 ! Local ID for the particle.
532  INTEGER, INTENT(OUT) :: pLocalNO
533 
534 ! Local variables
535 !---------------------------------------------------------------------//
536  INTEGER :: lpicloc, lcurpar
537 
538 ! Initialize the result.
539  locate_par = .false.
540 
541 ! Verify that the passied IJK value is within a valid range.
542  if(pijk < dg_ijkstart2 .or. pijk > dg_ijkend2) RETURN
543 
544 ! Loop the the particles in DES grid cell pIJK. Return to the calling
545 ! routine if the passed global ID matches the global ID of one of
546 ! the local particles.
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 ! Function: EXTEN_LOCATE_PAR !
561 ! Author: Pradeep Gopalakrishnan !
562 ! !
563 ! Purpose: Return the local index of the particle matching the passed !
564 ! global ID. The function returns TRUE if the particle is matched, !
565 ! otherwise it returns FALSE. !
566 !------------------------------------------------------------------------
567  LOGICAL FUNCTION exten_locate_par(pGlobalID, pIJK, pLocalNO)
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 ! Dummy variables:
580 !---------------------------------------------------------------------//
581 ! The global ID of the particle to be matched locally
582  INTEGER, INTENT(IN) :: pGlobalId
583 ! The DES grid cell index expected to contain the particle.
584  INTEGER, INTENT(IN) :: pIJK
585 ! The local ID of the matching particle.
586  INTEGER, INTENT(OUT) :: pLocalNo
587 
588 ! Local variables:
589 !---------------------------------------------------------------------//
590 ! Loop counters.
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 ! Unpack subroutine for single real variables !
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 ! Unpack subroutine for real arrays !
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 ! Unpack subroutine for single integer variables !
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 ! Unpack subroutine for integer arrays !
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 ! Unpack subroutine for logical variables !
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
subroutine unpack_i0(lbuf, idata, pface)
integer function dg_iof_lo(fijk)
Definition: desgrid_mod.f:307
subroutine des_mpi_stop(myid)
integer dg_ijkend2
Definition: desgrid_mod.f:57
double precision, dimension(:), allocatable des_t_s
subroutine, public desmpi_unpack_parcross(pface)
integer iparticlepacketsize
Definition: desmpi_mod.f:12
type(iap2), dimension(:), allocatable dg_pic
subroutine, public desmpi_unpack_ghostpar(pface)
subroutine unpack_i1(lbuf, idata, pface)
integer ispot
Definition: desmpi_mod.f:38
type(array), dimension(:), allocatable drecvbuf
Definition: desmpi_mod.f:28
integer function dg_kof_lo(fijk)
Definition: desgrid_mod.f:320
logical function exten_locate_par(pGlobalID, pIJK, pLocalNO)
double precision, dimension(:,:), allocatable des_x_s
Definition: des_rxns_mod.f:21
subroutine unpack_db0(lbuf, idata, pface)
logical any_species_eq
Definition: run_mod.f:118
Definition: run_mod.f:13
integer function dg_jof_lo(fijk)
Definition: desgrid_mod.f:295
integer dg_ijksize2
Definition: desgrid_mod.f:57
logical no_k
Definition: geometry_mod.f:28
integer dg_ijkstart2
Definition: desgrid_mod.f:57
double precision, dimension(:), allocatable conv_qs
integer, parameter ibufoffset
Definition: desmpi_mod.f:34
integer mype
Definition: compar_mod.f:24
logical function locate_par(pGLOBALID, pIJK, pLOCALNO)
logical energy_eq
Definition: run_mod.f:100
integer ighostpacketsize
Definition: desmpi_mod.f:11
integer, dimension(:), allocatable ineighproc
Definition: desmpi_mod.f:16
integer function dg_funijk(fi, fj, fk)
Definition: desgrid_mod.f:141
logical mppic
Definition: mfix_pic_mod.f:14
subroutine, public particle_grow(new_max_pip)
double precision, dimension(:), allocatable rxns_qs
subroutine unpack_db1(lbuf, idata, pface)
double precision, parameter pi
Definition: constant_mod.f:158
subroutine unpack_l0(lbuf, idata, pface)
double precision, dimension(:), allocatable des_stat_wt
Definition: mfix_pic_mod.f:54
double precision function, public add_pair(ii, jj)