MFIX  2016-1
mpi_funs_des_mod.f
Go to the documentation of this file.
1 !----------------------------------------------------------------------!
2 ! Module: MPI_FUNS_DES !
3 ! Author: Pradeep Gopalakrishnan !
4 ! !
5 ! Purpose: This module contains the subroutines and functions for MPI !
6 ! communications in DES simulations. !
7 !----------------------------------------------------------------------!
8  module mpi_funs_des
9 
10  contains
11 
12 !----------------------------------------------------------------------!
13 ! Subroutine: DES_PAR_EXCHANGE !
14 ! Author: Pradeep Gopalakrishnan !
15 ! !
16 ! Purpose: This subroutine controls entire exchange of particles !
17 ! between processors. !
18 ! !
19 ! Steps: !
20 ! 1) Bin the particles to the DES grid. !
21 ! 2) Check if the send and recv buffer size is large enough !
22 ! 3) Pack and send active particles located in ghost cells to the !
23 ! processors that own the ghost cells. The exchange occurs in !
24 ! the following order to take care of particles crossing at corners !
25 ! (e.g., crossing directly into the northwest block): !
26 ! a.) top-bottom interface !
27 ! b.) north-south interface !
28 ! c.) east-west interface !
29 ! 4) Bin the particles (if required) !
30 ! 5) Pack and send particles adjacent to neighboring processes. The !
31 ! exchange occurs in the following order: !
32 ! a.) east-west interface !
33 ! b.) north-south interface !
34 ! c.) top-bottom interface !
35 ! !
36 ! Comments: The DO_NSEARCH flag should be set before calling !
37 ! DES_PAR_EXCHANGE; When DO_NSEARCH is true, ghost particles are !
38 ! updated and later used to generate the PAIR lists. !
39 !----------------------------------------------------------------------!
40  subroutine des_par_exchange()
41 
42  use discretelement, only: do_nsearch, des_periodic_walls
43 
44  use mfix_pic, only: mppic
45  use desmpi, only: iexchflag
46  use desmpi, only: dsendbuf, drecvbuf
47  use discretelement, only: ighost_updated
48  use desmpi, only: ispot
49 
50  use geometry, only: no_k
51 
52 ! Module procedures
53 !---------------------------------------------------------------------//
56 
59 
62 
63  use desgrid, only: desgrid_pic
65  use compar, only: numpes
66 
67  implicit none
68 
69 ! Local variables:
70 !---------------------------------------------------------------------//
71 ! Loop counters.
72  integer :: linter, lface, ii
73 ! Number of calls since the buffer was last checked.
74  integer, save :: lcheckbuf = 0
75 
76 !......................................................................!
77 
78  IF (.not.((numpes>1) .OR. des_periodic_walls)) RETURN
79 
80 ! Check that the send/recv buffer is sufficient every 100 calls to
81 ! avoid the related global communications.
82  if (mod(lcheckbuf,100) == 0) then
83  call desmpi_check_sendrecvbuf(check_global=.true.)
84  lcheckbuf = 0
85  elseif (mod(lcheckbuf,5) == 0) then
86  call desmpi_check_sendrecvbuf(check_global=.false.)
87  end if
88  lcheckbuf = lcheckbuf + 1
89 
90 ! call particle crossing the boundary exchange in T-B,N-S,E-W order
91  do ii=1, size(dsendbuf)
92  dsendbuf(ii)%facebuf(1) = 0
93  drecvbuf(ii)%facebuf(1) = 0
94  end do
95 
96  ispot = 1
97  do linter = merge(2,3,no_k),1,-1
98  do lface = linter*2-1,linter*2
99  if(.not.iexchflag(lface))cycle
100  call desmpi_pack_parcross(lface)
101  call desmpi_sendrecv_init(lface)
102  end do
103  do lface = linter*2-1,linter*2
104  if(.not.iexchflag(lface)) cycle
105  call desmpi_sendrecv_wait(lface)
106  call desmpi_unpack_parcross(lface)
107  end do
108 ! update pic this is required for particles crossing corner cells
109  do lface = linter*2-1,linter*2
110  if(dsendbuf(1+mod(lface,2))%facebuf(1).gt.0 .or. &
111  drecvbuf(1+mod(lface,2))%facebuf(1).gt.0) then
112  call desgrid_pic(plocate=.false.)
113  exit
114  end if
115  end do
116  end do
117  call des_mpi_barrier
118 
119 ! call des_dbgmpi(5)
120 
121 
122  IF(.NOT.mppic) THEN
123 ! call ghost particle exchange in E-W, N-S, T-B order
124 
125  do ii=1, size(dsendbuf)
126  dsendbuf(ii)%facebuf(1) = 0
127  drecvbuf(ii)%facebuf(1) = 0
128  end do
129 
130  ighost_updated(:) = .false.
131  ispot = 1
132  do linter = 1,merge(2,3,no_k)
133  do lface = linter*2-1,linter*2
134  if(.not.iexchflag(lface))cycle
135  call desmpi_pack_ghostpar(lface)
136  call desmpi_sendrecv_init(lface)
137  end do
138  do lface = linter*2-1,linter*2
139  if(.not.iexchflag(lface)) cycle
140  call desmpi_sendrecv_wait(lface)
141  call desmpi_unpack_ghostpar(lface)
142  end do
143 
144 ! Rebin particles to the DES grid as ghost particles may be moved.
145  do lface = linter*2-1,linter*2
146  if(dsendbuf(1+mod(lface,2))%facebuf(1) .gt.0.or.&
147  drecvbuf(1+mod(lface,2))%facebuf(1).gt.0) then
148  call desgrid_pic(plocate=.false.)
149  exit
150  end if
151  end do
152  end do
153  call desmpi_cleanup
154  call des_mpi_barrier
155  ENDIF ! end if(.not.mppic)
156 
157 ! call des_dbgmpi(2)
158 ! call des_dbgmpi(3)
159 ! call des_dbgmpi(4)
160 ! call des_dbgmpi(6)
161 ! call des_dbgmpi(7)
162 
163  END SUBROUTINE des_par_exchange
164 
165 
166 !----------------------------------------------------------------------!
167 ! Subroutine: DESMPI_CHECK_SENDRECVBUF !
168 ! Author: Pradeep Gopalakrishnan !
169 ! !
170 ! Purpose: Checks if the sendrecvbuf size is large enough. If the !
171 ! buffers are not sufficent, they are resized. !
172 !----------------------------------------------------------------------!
173  SUBROUTINE desmpi_check_sendrecvbuf(check_global)
175  use derived_types, only: dg_pic
176  use discretelement, only: dimn
177  use desmpi, only: imaxbuf
178  use desmpi, only: ibufoffset
179  use desmpi, only: dsendbuf, drecvbuf
180  use desmpi, only: isendindices
181  use desmpi, only: ighostpacketsize
182 
183  use mpi_utility, only: global_all_max
184  use error_manager
185  implicit none
186 
187  logical, intent(in) :: check_global
188 
189 ! Local variables:
190 !---------------------------------------------------------------------//
191 ! Loop counters
192  INTEGER :: lface, lindx, lijk
193 ! Particle count in send/recv region on current face
194  INTEGER :: lparcnt
195 ! Max particle count in send/recv region over all proc faces.
196  INTEGER :: lmaxcnt
197 ! Total number of DES grid cells on lface in send/recv
198  INTEGER :: ltot_ind
199 ! Previous Buffer
200  INTEGER :: pBUF
201 ! Growth factor when resizing send/recv buffers.
202  REAL :: lfactor = 0.5
203  DOUBLE PRECISION, PARAMETER :: ONEMBo8 = 131072.0
204 !......................................................................!
205 
206  lmaxcnt = 0
207  do lface = 1,6
208  ltot_ind = isendindices(1,lface)
209  lparcnt = 0
210  do lindx = 2,ltot_ind+1
211  lijk = isendindices(lindx,lface)
212  lparcnt = lparcnt + dg_pic(lijk)%isize
213  enddo
214  if(lparcnt > lmaxcnt) lmaxcnt = lparcnt
215  enddo
216 
217  if(check_global) call global_all_max(lmaxcnt)
218 
219  if (imaxbuf < (1.0+0.5*lfactor)*lmaxcnt*ighostpacketsize) then
220  pbuf = imaxbuf
221  imaxbuf = (1.0+lfactor)*lmaxcnt*ighostpacketsize
222  do lface = 1,2*dimn
223  if(allocated(dsendbuf(1+mod(lface,2))%facebuf)) then
224  deallocate(dsendbuf(1+mod(lface,2))%&
225  facebuf,drecvbuf(1+mod(lface,2))%facebuf)
226  endif
227  allocate(dsendbuf(1+mod(lface,2))%facebuf(imaxbuf),&
228  drecvbuf(1+mod(lface,2))%facebuf(imaxbuf))
229  end do
230 
231  WRITE(err_msg, 1000) imaxbuf/onembo8, &
232  100.0d0+100.0d0*dble(imaxbuf-pbuf)/dble(pbuf)
233  CALL flush_err_msg(header=.false., footer=.false.)
234 
235  1000 FORMAT(/'Resizeing DES MPI buffers: ',f7.1,' MB (+',f5.1'%)')
236 
237  endif
238 
239  END SUBROUTINE desmpi_check_sendrecvbuf
240 
241 !----------------------------------------------------------------------!
242 ! Subroutine: DESMPI_CLEANUP !
243 ! Author: Pradeep Gopalakrishnan !
244 ! !
245 ! Purpose: Cleans the ghost particle array positions. !
246 !----------------------------------------------------------------------!
247  SUBROUTINE desmpi_cleanup
249  use discretelement, only: dimn
250  use discretelement, only: des_pos_new, des_pos_old
251  use discretelement, only: des_vel_new, des_vel_old
252  use discretelement, only: omega_new
253  use discretelement, only: particle_orientation
254  use discretelement, only: orientation,init_orientation
255  use discretelement, only: fc
256  use discretelement, only: do_old
257  use discretelement, only: pip
258  use discretelement, only: ighost_cnt
259  use discretelement, only: des_usr_var_size, des_usr_var
260  use derived_types, only: dg_pic
261  use discretelement, only: pijk
262 
263  use run, only: energy_eq,any_species_eq
264  use des_thermo, only: des_t_s
265 
266  use des_rxns, only: des_x_s
267 
268  use discretelement, only: ighost_updated
269  use functions, only: set_nonexistent
270  use desmpi, only: irecvindices
271  use desmpi, only: iexchflag
272 
273  use param, only: dimension_n_s
274  use param1, only: zero
275 
276  implicit none
277 
278 ! Local variables:
279 !---------------------------------------------------------------------//
280  integer ltot_ind,lface,lindx,lijk,lcurpar,lpicloc
281 
282  do lface = 1,dimn*2
283  if(.not.iexchflag(lface))cycle
284  ltot_ind = irecvindices(1,lface)
285  do lindx = 2,ltot_ind+1
286  lijk = irecvindices(lindx,lface)
287  do lpicloc =1,dg_pic(lijk)%isize
288  lcurpar = dg_pic(lijk)%p(lpicloc)
289  if(ighost_updated(lcurpar)) cycle
290  pip = pip - 1
291  ighost_cnt = ighost_cnt-1
292  call set_nonexistent(lcurpar)
293  fc(lcurpar,:) = 0.0
294  des_pos_new(lcurpar,:)=0
295  pijk(lcurpar,:) = -10
296  IF (do_old) THEN
297  des_pos_old(lcurpar,:)=0
298  des_vel_old(lcurpar,:)=0
299  ENDIF
300  des_vel_new(lcurpar,:)=0
301  omega_new(lcurpar,:)=0
302 
303  IF(particle_orientation) &
304  orientation(1:3,lcurpar) = init_orientation
305 
306  if(energy_eq) des_t_s(lcurpar)=0
307 
308  if(any_species_eq)then
309  des_x_s(lcurpar,1:dimension_n_s)= 0
310  endif
311 
312  IF(des_usr_var_size > 0)&
313  des_usr_var(:,lcurpar)= 0
314 
315  end do
316  end do
317  end do
318  END SUBROUTINE desmpi_cleanup
319 
320 
321 !----------------------------------------------------------------------!
322 ! Subroutine: DESMPI_SEND_RECV_FIELD_VARS !
323 ! Author: J.Musser !
324 ! !
325 ! Purpose: Collect data from ghost cells for interpolated cases and !
326 ! preform halo exchange go sync ghost cell data. !
327 !----------------------------------------------------------------------!
328  SUBROUTINE desmpi_send_recv_field_vars
330 ! Flags for solving energy, species and solids density
332 ! Flag for explicitly coupled simulations
333  use discretelement, only: des_explicitly_coupled
334 ! Flag for use interpolation
335  use particle_filter, only: des_interp_on
336 ! Gas phae volume fraction and bluk density
337  use fldvar, only: ep_g, rop_g
338 ! Solids bulk density and material density
339  use fldvar, only: rop_s, ro_s
340 ! Gas phase mass and species mass source terms
341  use des_rxns, only: des_r_gp, des_r_gc
342  use des_rxns, only: des_sum_r_g, des_r_phase
343 ! Gas phase energy equation source terms
344  use des_thermo, only: conv_sc
345  use des_rxns, only: des_hor_g
346 ! MPI function for collecting interpolated data from ghost cells.
347  use sendrecvnode, only: des_collect_gdata
348 ! MPI wrapper for halo exchange.
349  use sendrecv, only: send_recv
350 
351  implicit none
352 
353 ! Local variables:
354 !---------------------------------------------------------------------//
355 
356  CALL send_recv(ep_g,2)
357  CALL send_recv(rop_g,2)
358 
359  CALL send_recv(rop_s,2)
360  IF(any(solve_ros)) CALL send_recv(ro_s,2)
361 
362 ! Add in interpolated data from ghost cells and halo exchange
363  IF(.NOT.des_explicitly_coupled) THEN
364 
365  IF(des_interp_on) THEN
366  IF(energy_eq) CALL des_collect_gdata(conv_sc)
367  IF(any_species_eq) THEN
368  CALL des_collect_gdata(des_r_gp)
369  CALL des_collect_gdata(des_r_gc)
370  CALL des_collect_gdata(des_r_phase)
371  CALL des_collect_gdata(des_hor_g)
372  CALL des_collect_gdata(des_sum_r_g)
373  ENDIF
374  ENDIF
375 
376  IF(energy_eq) CALL send_recv(conv_sc, 2)
377  IF(any_species_eq) THEN
378  CALL send_recv(des_r_gp, 2)
379  CALL send_recv(des_r_gc, 2)
380  CALL send_recv(des_r_phase, 2)
381  CALL send_recv(des_hor_g, 2)
382  CALL send_recv(des_sum_r_g, 2)
383  ENDIF
384  ENDIF
385 
386  RETURN
387  END SUBROUTINE desmpi_send_recv_field_vars
388 
389  END MODULE mpi_funs_des
double precision, dimension(:), allocatable conv_sc
double precision, dimension(:), allocatable ep_g
Definition: fldvar_mod.f:15
double precision, dimension(:), allocatable des_t_s
subroutine desmpi_sendrecv_init(pface, pdebug)
subroutine, public desmpi_unpack_parcross(pface)
double precision, dimension(:), allocatable des_hor_g
Definition: des_rxns_mod.f:46
subroutine desgrid_pic(plocate)
Definition: desgrid_mod.f:711
subroutine desmpi_sendrecv_wait(pface, pdebug)
subroutine, public desmpi_pack_parcross(PFACE)
subroutine des_par_exchange()
logical, dimension(dim_m) solve_ros
Definition: run_mod.f:250
logical, dimension(:), allocatable iexchflag
Definition: desmpi_mod.f:17
double precision, dimension(:), allocatable des_sum_r_g
Definition: des_rxns_mod.f:42
subroutine desmpi_cleanup
type(iap2), dimension(:), allocatable dg_pic
subroutine, public desmpi_unpack_ghostpar(pface)
subroutine, public desmpi_pack_ghostpar(PFACE)
subroutine des_mpi_barrier
type(array), dimension(:), allocatable dsendbuf
Definition: desmpi_mod.f:27
integer ispot
Definition: desmpi_mod.f:38
integer numpes
Definition: compar_mod.f:24
type(array), dimension(:), allocatable drecvbuf
Definition: desmpi_mod.f:28
integer, dimension(:,:), allocatable isendindices
Definition: desmpi_mod.f:54
double precision, dimension(:,:), allocatable des_r_gc
Definition: des_rxns_mod.f:40
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
Definition: param_mod.f:2
double precision, dimension(:,:), allocatable ro_s
Definition: fldvar_mod.f:45
double precision, dimension(:,:), allocatable des_r_phase
Definition: des_rxns_mod.f:44
logical no_k
Definition: geometry_mod.f:28
integer imaxbuf
Definition: desmpi_mod.f:37
subroutine desmpi_check_sendrecvbuf(check_global)
integer, parameter ibufoffset
Definition: desmpi_mod.f:34
logical energy_eq
Definition: run_mod.f:100
integer ighostpacketsize
Definition: desmpi_mod.f:11
character(len=line_length), dimension(line_count) err_msg
subroutine desmpi_send_recv_field_vars
integer, dimension(:,:), allocatable irecvindices
Definition: desmpi_mod.f:55
double precision, dimension(:,:), allocatable rop_s
Definition: fldvar_mod.f:51
logical mppic
Definition: mfix_pic_mod.f:14
double precision, dimension(:,:), allocatable des_r_gp
Definition: des_rxns_mod.f:38
integer dimension_n_s
Definition: param_mod.f:21
double precision, dimension(:), allocatable rop_g
Definition: fldvar_mod.f:38
double precision, parameter zero
Definition: param1_mod.f:27
subroutine flush_err_msg(DEBUG, HEADER, FOOTER, ABORT, LOG, CALL_TREE)