MFIX  2016-1
mpi_pack_des_mod.f
Go to the documentation of this file.
1 !----------------------------------------------------------------------!
2 ! Module: MPI_PACK_DES !
3 ! Author: Pradeep Gopalakrishnan !
4 ! !
5 ! Purpose: Contains routines for packing real and ghost particles !
6 ! into the MPI send buffers. !
7 !----------------------------------------------------------------------!
8  MODULE mpi_pack_des
9 
10  PRIVATE
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 ! Subroutine: DESMPI_PACK_GHOSTPAR !
25 ! Author: Pradeep Gopalakrishnan !
26 ! !
27 ! Purpose: Packs ghost particle in the send buffer. !
28 !----------------------------------------------------------------------!
29  SUBROUTINE desmpi_pack_ghostpar(PFACE)
30 
31 ! Global Variables:
32 !---------------------------------------------------------------------//
33 ! Size of ghost particle packets
34  use desmpi, only: ighostpacketsize
35 ! Flag indicating that the ghost particle was updated
36  use discretelement, only: ighost_updated
37 ! The MPI send buffer
38  use desmpi, only: dsendbuf
39 ! Buffer offset
40  use desmpi, only: ibufoffset
41 ! Runtime flag for solving the energy equations
42  use run, only: energy_eq
43 ! The neighbor processor's rank
44  use desmpi, only: ineighproc
45 ! DES grid cell containing each particle: current/previous
46  use discretelement, only: dg_pijkprv
47 ! The global ID for each particle
48  use discretelement, only: iglobal_id
49 ! Particle positions: current/previous
50  use discretelement, only: des_pos_new
51 ! Particle tangential velocities: current/previous
52  use discretelement, only: des_vel_new
53 ! Particle rotational velocities: current/previous
54  use discretelement, only: omega_new
55 ! Particle tempertures. current/previous
56  use des_thermo, only: des_t_s
57 ! Particle radius, volume
58  use discretelement, only: des_radius
59 ! Map to fluid grid cells and solids phase (I,J,K,IJK,M)
60  use discretelement, only: pijk
61 ! Number of particles on the process (max particle array size)
62  use discretelement, only: max_pip
63 ! User-defined variables for each particle.
64  use discretelement, only: des_usr_var, des_usr_var_size
65 ! Function to convert DES grid IJK to new proc value.
66  use desgrid, only: dg_ijkconv
67 ! Size of the send buffer
68  use desmpi, only: isendcnt
69 ! Offset for particles with cycle BCs (otherwise zero)
70  use desmpi, only: dcycl_offset
71 ! Map of particles to DES grid
72  use derived_types, only: dg_pic
73 ! Cell number of ghost particles
74  use desmpi, only: isendindices
75 
76 ! Global Constants:
77 !---------------------------------------------------------------------//
78  use constant, only: pi
79 ! Dimension of particle spatial arrays.
80  use discretelement, only: dimn
81  use functions, only: is_exiting
82  use functions, only: is_ghost, is_entering_ghost, is_exiting_ghost
83 
84  IMPLICIT NONE
85 
86 ! Dummy arguments:
87 !---------------------------------------------------------------------//
88 ! Processor boundary being packed (Top/Bottom/North/South/East/West)
89  INTEGER, INTENT(IN) :: PFACE
90 
91 ! Local variables
92 !---------------------------------------------------------------------//
93  integer :: lijk,lindx,ltot_ind,lpicloc,lpar_cnt,lcurpar
94  integer :: lbuf
95 !......................................................................!
96 
97  lpar_cnt = 0
98  ltot_ind = isendindices(1,pface)
99  do lindx = 2,ltot_ind+1
100  lijk = isendindices(lindx,pface)
101  do lpicloc =1,dg_pic(lijk)%isize
102  lbuf = lpar_cnt*ighostpacketsize + ibufoffset
103  lcurpar = dg_pic(lijk)%p(lpicloc)
104 
105 ! Do not send particle data for a ghost particle whose owner has not yet
106 ! updated the particle's data on this processor.
107  if((is_ghost(lcurpar) .or. &
108  is_entering_ghost(lcurpar) .or. &
109  is_exiting_ghost(lcurpar)) .and. &
110  .not.ighost_updated(lcurpar)) cycle
111 
112 ! 1) Global ID
113  call pack_dbuf(lbuf,iglobal_id(lcurpar),pface)
114 ! 2) DES grid IJK
115  call pack_dbuf(lbuf,dg_ijkconv(lijk,pface, &
116  ineighproc(pface)),pface)
117 ! 3) DES grid IJK - previous
118  call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lcurpar),pface, &
119  ineighproc(pface)),pface)
120 ! 4) Radius
121  call pack_dbuf(lbuf,des_radius(lcurpar),pface)
122 ! 5) Phase index
123  call pack_dbuf(lbuf,pijk(lcurpar,5),pface)
124 ! 6) Position
125  call pack_dbuf(lbuf,des_pos_new(lcurpar,:)+ &
126  dcycl_offset(pface,:),pface)
127 ! 7) Translational Velocity
128  call pack_dbuf(lbuf,des_vel_new(lcurpar,:),pface)
129 ! 8) Rotational Velocity
130  call pack_dbuf(lbuf,omega_new(lcurpar,:),pface)
131 ! 9) Exiting particle flag
132  call pack_dbuf(lbuf,merge(1,0,is_exiting(lcurpar).or.is_exiting_ghost(lcurpar)),pface)
133 ! 10) Temperature
134  IF(energy_eq) &
135  call pack_dbuf(lbuf,des_t_s(lcurpar),pface)
136 ! 11) User Variable
137  IF(des_usr_var_size > 0) &
138  call pack_dbuf(lbuf,des_usr_var(:,lcurpar),pface)
139 
140  lpar_cnt = lpar_cnt + 1
141  end do
142  end do
143  dsendbuf(1+mod(pface,2))%facebuf(1)=lpar_cnt
144  isendcnt(pface) = lpar_cnt*ighostpacketsize + ibufoffset
145 
146  end subroutine desmpi_pack_ghostpar
147 
148 !----------------------------------------------------------------------!
149 ! Subroutine: DESMPI_PACK_PARCROSS !
150 ! Author: Pradeep Gopalakrishnan !
151 ! !
152 ! Purpose: Packs real particle in the send buffer. !
153 !----------------------------------------------------------------------!
154  SUBROUTINE desmpi_pack_parcross(PFACE)
156 ! Global Variables:
157 !---------------------------------------------------------------------//
158 ! The MPI send buffer
159  use desmpi, only: dsendbuf
160 ! Buffer offset
161  use desmpi, only: ibufoffset
162 ! Runtime flag for solving the energy equations
163  use run, only: energy_eq
164 ! Runtime flag for solving species equations
165  use run, only: any_species_eq
166 ! Runtime flag for MPPIC solids
167  use mfix_pic, only: mppic
168 ! DES grid cell containing each particle: current/previous
169  use discretelement, only: dg_pijkprv
170 ! The neighbor processor's rank
171  use desmpi, only: ineighproc
172 ! The statistical weight of each particle.
173  use mfix_pic, only: des_stat_wt
174 ! The global ID for each particle
175  use discretelement, only: iglobal_id
176 ! Particle positions: current/previous
177  use discretelement, only: des_pos_new, des_pos_old
178 ! Particle tangential velocities: current/previous
179  use discretelement, only: des_vel_new, des_vel_old
180 ! Particle rotational velocities: current/previous
181  use discretelement, only: omega_new, omega_old
182 ! Particle orientation
183  use discretelement, only: particle_orientation,orientation
184 ! Particle radius, volume, density, mass
185  use discretelement, only: des_radius, pvol, ro_sol, pmass
186 ! Previous value for particle acceleration (tangential/rotational)
187  use discretelement, only: des_acc_old, rot_acc_old
188 ! Particle species composition
189  use des_rxns, only: des_x_s
190 ! Particle tempertures. current/previous
191  use des_thermo, only: des_t_s
192 ! Force arrays acting on the particle
193  use discretelement, only: fc, tow
194 ! One of the moment of inertia
195  use discretelement, only: omoi
196 ! Map to fluid grid cells and solids phase (I,J,K,IJK,M)
197  use discretelement, only: pijk
198 ! Flag to send/recv old (previous) values
199  use discretelement, only: do_old
200 ! Number of particles on the process (max particle array size)
201  use discretelement, only: pip, max_pip
202 ! Number of ghost particles on the current process
203  use discretelement, only: ighost_cnt
204 ! User-defined variables for each particle.
205  use discretelement, only: des_usr_var, des_usr_var_size
206 ! Particle pair (neighborhood) arrays:
207  use discretelement, only: neighbors, neighbor_index, neigh_num
208 ! Pair collision history information
209  use discretelement, only: pft_neighbor
210 ! Dimension of particle spatial arrays.
211  use discretelement, only: dimn
212 ! Flag indicating the the fluid-particle drag is explicitly coupled.
213  use discretelement, only: des_explicitly_coupled
214 ! Explicit particle drag force
215  use discretelement, only: drag_fc
216 ! Explict convection and HOR
217  use des_thermo, only: conv_qs, rxns_qs
218 
219  use desgrid, only: dg_ijkconv, icycoffset
220  use desmpi, only: dcycl_offset, isendcnt
221  use desmpi, only: irecvindices
222 
223  use desmpi, only: iparticlepacketsize
224  use desmpi, only: ipairpacketsize
225  use derived_types, only: dg_pic
226 
227  use functions
228 
229  implicit none
230 
231 ! Dummy arguments:
232 !---------------------------------------------------------------------//
233 ! Processor boundary being packed (Top/Bottom/North/South/East/West)
234  INTEGER, INTENT(IN) :: PFACE
235 
236 ! Local variables
237 !---------------------------------------------------------------------//
238  integer :: li, lj, lk
239  integer :: ltot_ind,lindx,cc
240  integer :: lneigh,lijk,&
241  lpicloc,lparcnt,lcurpar
242  integer :: lbuf,num_neighborlists_to_send
243 
244  logical, allocatable, dimension(:) :: going_to_send
245 
246 ! Location in the buffer where the number of pair data is specified.
247  integer :: num_neighborlists_send_buf_loc
248 !......................................................................!
249 
250 ! pack the particle crossing the boundary
251  ltot_ind = irecvindices(1,pface)
252  lparcnt = 0
253 
254  allocate(going_to_send(max_pip))
255  going_to_send(:) = .false.
256 
257  do lindx = 2,ltot_ind+1
258  lijk = irecvindices(lindx,pface)
259  do lpicloc = 1,dg_pic(lijk)%isize
260  lcurpar = dg_pic(lijk)%p(lpicloc)
261 
262 ! if ghost particle then cycle
263  if(is_ghost(lcurpar) .or. &
264  is_entering_ghost(lcurpar) .or. &
265  is_exiting_ghost(lcurpar)) cycle
266 
267  going_to_send(lcurpar) = .true.
268  lbuf = lparcnt*iparticlepacketsize + ibufoffset
269 
270 ! 1) Global ID
271  call pack_dbuf(lbuf,iglobal_id(lcurpar),pface)
272 ! 2) DES Grid IJK
273  call pack_dbuf(lbuf,dg_ijkconv(lijk,pface, &
274  ineighproc(pface)),pface)
275 ! 3) DES grid IJK - previous
276  call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lcurpar),pface, &
277  ineighproc(pface)),pface)
278 ! 4) Radius
279  call pack_dbuf(lbuf,des_radius(lcurpar),pface)
280 ! 5) Fluid cell I index with cycle offset
281  li = pijk(lcurpar,1) + icycoffset(pface,1)
282  call pack_dbuf(lbuf,li,pface)
283 ! 6) Fluid cell J index with cycle offset
284  lj = pijk(lcurpar,2) + icycoffset(pface,2)
285  call pack_dbuf(lbuf,lj,pface)
286 ! 7) Fluid cell K index with cycle offset
287  lk = pijk(lcurpar,3) + icycoffset(pface,3)
288  call pack_dbuf(lbuf,lk,pface)
289 ! 8) Fluid cell IJK on destination process
290  call pack_dbuf(lbuf,funijk_proc(li,lj,lk, &
291  ineighproc(pface)),pface)
292 ! 9) Particle solids phase index
293  call pack_dbuf(lbuf,pijk(lcurpar,5),pface)
294 ! 10) Entering particle flag.
295  call pack_dbuf(lbuf, is_entering(lcurpar).or.is_entering_ghost(lcurpar), pface)
296 ! 11) Exiting particle flag.
297  call pack_dbuf(lbuf, is_exiting(lcurpar).or.is_exiting_ghost(lcurpar), pface)
298 ! 12) Density
299  call pack_dbuf(lbuf,ro_sol(lcurpar),pface)
300 ! 13) Volume
301  call pack_dbuf(lbuf,pvol(lcurpar),pface)
302 ! 14) Mass
303  call pack_dbuf(lbuf,pmass(lcurpar),pface)
304 ! 15) 1/Moment of Inertia
305  call pack_dbuf(lbuf,omoi(lcurpar),pface)
306 ! 16) Position with cyclic shift
307  call pack_dbuf(lbuf,des_pos_new(lcurpar,:) + &
308  dcycl_offset(pface,:),pface)
309 ! 17) Translational velocity
310  call pack_dbuf(lbuf,des_vel_new(lcurpar,:),pface)
311 ! 18) Rotational velocity
312  call pack_dbuf(lbuf,omega_new(lcurpar,:),pface)
313 ! 19) Accumulated translational forces
314  call pack_dbuf(lbuf,fc(lcurpar,:),pface)
315 ! 20) Accumulated torque forces
316  call pack_dbuf(lbuf,tow(lcurpar,:),pface)
317  IF(energy_eq) THEN
318 ! 21) Temperature
319  call pack_dbuf(lbuf,des_t_s(lcurpar),pface)
320 ! 22) Species composition
321  call pack_dbuf(lbuf,des_x_s(lcurpar,:),pface)
322  ENDIF
323 ! 23) User defined variable
324  IF(des_usr_var_size > 0) &
325  call pack_dbuf(lbuf, des_usr_var(:,lcurpar),pface)
326 ! 24) Particle orientation
327  IF(particle_orientation) &
328  call pack_dbuf(lbuf,orientation(:,lcurpar),pface)
329 
330 ! -- Higher order integration variables
331  IF (do_old) THEN
332 ! 25) Position (previous)
333  call pack_dbuf(lbuf,des_pos_old(lcurpar,:) + &
334  dcycl_offset(pface,:),pface)
335 ! 26) Translational velocity (previous)
336  call pack_dbuf(lbuf,des_vel_old(lcurpar,:),pface)
337 ! 27) Rotational velocity (previous)
338  call pack_dbuf(lbuf,omega_old(lcurpar,:),pface)
339 ! 28) Translational acceleration (previous)
340  call pack_dbuf(lbuf,des_acc_old(lcurpar,:),pface)
341 ! 29) Rotational acceleration (previous)
342  call pack_dbuf(lbuf,rot_acc_old(lcurpar,:),pface)
343  ENDIF
344 
345  IF(des_explicitly_coupled) THEN
346 ! 30) Explicit drag force
347  call pack_dbuf(lbuf, drag_fc(lcurpar,:),pface)
348 ! 31) Explicit convective heat transfer
349  IF(energy_eq) call pack_dbuf(lbuf,conv_qs(lcurpar),pface)
350 ! 32) Explicit heat of reaction
351  IF(any_species_eq) call pack_dbuf(lbuf, rxns_qs(lcurpar),pface)
352  ENDIF
353 
354 ! PIC particles are removed and the number of particles on the processor
355 ! is decremented.
356  IF (mppic) THEN
357 ! 34) Statistical weight
358  call pack_dbuf(lbuf,des_stat_wt(lcurpar),pface)
359  call set_nonexistent(lcurpar)
360  pip = pip - 1
361 
362 ! DEM particles are converted to ghost particles. This action does not
363 ! change the number of particles.
364  ELSE
365  if (is_entering(lcurpar)) then
366  call set_entering_ghost(lcurpar)
367  elseif (is_exiting(lcurpar)) then
368  call set_exiting_ghost(lcurpar)
369  else
370  call set_ghost(lcurpar)
371  endif
372  ighost_cnt = ighost_cnt + 1
373  END IF
374 
375 ! Clear out the force array.
376  fc(lcurpar,:) = 0.
377  lparcnt = lparcnt + 1
378  end do
379  end do
380 
381 ! Calculate the location in buffer where the number of pair data is
382 ! stored and skip specifying the entry. After all the pair data is
383 ! packed, then this value is set.
384  lbuf = lparcnt*iparticlepacketsize + ibufoffset
385  num_neighborlists_send_buf_loc = lbuf
386  lbuf = lbuf+1
387 
388  num_neighborlists_to_send = 0
389  lcurpar = 1
390  do cc = 1, neigh_num
391  IF (0 .eq. neighbors(cc)) EXIT
392 
393  IF (cc.eq.neighbor_index(lcurpar)) THEN
394  lcurpar = lcurpar + 1
395  ENDIF
396 
397 ! Only packup pairing data for particles being transfered.
398  if (.not. going_to_send(lcurpar)) cycle
399 
400 ! Do not send pairing data if the pair no longer exists or if the
401 ! particle is exiting as it may be locatable during unpacking.
402  lneigh = neighbors(lcurpar)
403  if(is_nonexistent(lneigh)) cycle
404  if(is_exiting(lneigh)) cycle
405 
406 ! 35) Global ID of particle being packed.
407  call pack_dbuf(lbuf,iglobal_id(lcurpar),pface)
408 ! 36) DES grid IJK of cell receiving the particle.
409  call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lcurpar),pface, &
410  ineighproc(pface)),pface)
411 ! 37) Global ID of neighbor particle.
412  call pack_dbuf(lbuf,iglobal_id(lneigh),pface)
413 ! 38) DES grid IJK of cell containing the neighbor particle.
414  call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lneigh),pface, &
415  ineighproc(pface)),pface)
416 ! 39) Tangential collision history.
417  call pack_dbuf(lbuf,pft_neighbor(:,cc),pface)
418 ! Increment the number of pairs being sent.
419  num_neighborlists_to_send = num_neighborlists_to_send + 1
420  enddo
421 
422 ! Store the number of pair datasets being sent. This information is
423 ! stored before the pairing data so the receiving process knows the
424 ! amount of data to 'unpack.'
425  lbuf = num_neighborlists_send_buf_loc
426 ! 34) Number of pair datasets.
427  call pack_dbuf(lbuf,num_neighborlists_to_send,pface)
428 
429  dsendbuf(1+mod(pface,2))%facebuf(1) = lparcnt
430  isendcnt(pface) = lparcnt*iparticlepacketsize + &
431  num_neighborlists_to_send*ipairpacketsize + ibufoffset + 3
432 
433  deallocate(going_to_send)
434 
435  RETURN
436  END SUBROUTINE desmpi_pack_parcross
437 
438 !----------------------------------------------------------------------!
439 ! PACK SUBROUTINE FOR SINGLE REAL VARIABLES !
440 !----------------------------------------------------------------------!
441  subroutine pack_db0(lbuf,idata,pface)
442  use desmpi, only: dsendbuf
443  integer, intent(inout) :: lbuf
444  integer, intent(in) :: pface
445  double precision, intent(in) :: idata
446 
447  dsendbuf(1+mod(pface,2))%facebuf(lbuf) = idata
448  lbuf = lbuf + 1
449 
450  return
451  end subroutine pack_db0
452 
453 !----------------------------------------------------------------------!
454 ! Pack subroutine for real arrays !
455 !----------------------------------------------------------------------!
456  subroutine pack_db1(lbuf,idata,pface)
457  use desmpi, only: dsendbuf
458  integer, intent(inout) :: lbuf
459  integer, intent(in) :: pface
460  double precision, intent(in) :: idata(:)
461 
462  integer :: lsize
463 
464  lsize = size(idata)
465 
466  dsendbuf(1+mod(pface,2))%facebuf(lbuf:lbuf+lsize-1) = idata
467  lbuf = lbuf + lsize
468 
469  return
470  end subroutine pack_db1
471 
472 !----------------------------------------------------------------------!
473 ! Pack subroutine for single integer variables !
474 !----------------------------------------------------------------------!
475  subroutine pack_i0(lbuf,idata,pface)
476  use desmpi, only: dsendbuf
477  integer, intent(inout) :: lbuf
478  integer, intent(in) :: pface
479  integer, intent(in) :: idata
480 
481  dsendbuf(1+mod(pface,2))%facebuf(lbuf) = idata
482  lbuf = lbuf + 1
483 
484  return
485  end subroutine pack_i0
486 
487 !----------------------------------------------------------------------!
488 ! Pack subroutine for integer arrays !
489 !----------------------------------------------------------------------!
490  subroutine pack_i1(lbuf,idata,pface)
491  use desmpi, only: dsendbuf
492  integer, intent(inout) :: lbuf
493  integer, intent(in) :: pface
494  integer, intent(in) :: idata(:)
495 
496  integer :: lsize
497 
498  lsize = size(idata)
499 
500  dsendbuf(1+mod(pface,2))%facebuf(lbuf:lbuf+lsize-1) = idata
501  lbuf = lbuf + lsize
502 
503  return
504  end subroutine pack_i1
505 
506 !----------------------------------------------------------------------!
507 ! Pack subroutine for logical scalars !
508 !----------------------------------------------------------------------!
509  subroutine pack_l0(lbuf,ldata,pface)
510  use desmpi, only: dsendbuf
511 
512  integer, intent(inout) :: lbuf
513  integer, intent(in) :: pface
514  logical, intent(in) :: ldata
515 
516  dsendbuf(1+mod(pface,2))%facebuf(lbuf) = merge(1.0, 0.0, ldata)
517  lbuf = lbuf + 1
518 
519  return
520  end subroutine pack_l0
521 
522  end module mpi_pack_des
double precision, dimension(:), allocatable des_t_s
integer function dg_ijkconv(fijk, fface, fto_proc)
Definition: desgrid_mod.f:334
integer iparticlepacketsize
Definition: desmpi_mod.f:12
integer, dimension(:,:), allocatable icycoffset
Definition: desgrid_mod.f:67
subroutine, public desmpi_pack_parcross(PFACE)
subroutine pack_db1(lbuf, idata, pface)
type(iap2), dimension(:), allocatable dg_pic
subroutine pack_i1(lbuf, idata, pface)
subroutine, public desmpi_pack_ghostpar(PFACE)
type(array), dimension(:), allocatable dsendbuf
Definition: desmpi_mod.f:27
integer, dimension(:,:), allocatable isendindices
Definition: desmpi_mod.f:54
subroutine pack_l0(lbuf, ldata, pface)
subroutine pack_i0(lbuf, idata, pface)
double precision, dimension(:,:), allocatable des_x_s
Definition: des_rxns_mod.f:21
logical any_species_eq
Definition: run_mod.f:118
Definition: run_mod.f:13
double precision, dimension(:), allocatable conv_qs
integer ipairpacketsize
Definition: desmpi_mod.f:13
integer, parameter ibufoffset
Definition: desmpi_mod.f:34
double precision, dimension(:,:), allocatable dcycl_offset
Definition: desmpi_mod.f:20
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, dimension(:,:), allocatable irecvindices
Definition: desmpi_mod.f:55
logical mppic
Definition: mfix_pic_mod.f:14
integer, dimension(:), allocatable isendcnt
Definition: desmpi_mod.f:30
double precision, dimension(:), allocatable rxns_qs
double precision, parameter pi
Definition: constant_mod.f:158
double precision, dimension(:), allocatable des_stat_wt
Definition: mfix_pic_mod.f:54
subroutine pack_db0(lbuf, idata, pface)