File: /nfs/home/0/users/jenkins/mfix.git/model/des/mpi_pack_des_mod.f
1
2
3
4
5
6
7
8 MODULE MPI_PACK_DES
9
10
11 PRIVATE
12 PUBLIC :: DESMPI_PACK_PARCROSS, DESMPI_PACK_GHOSTPAR
13
14 interface pack_dbuf
15 module procedure pack_db0
16 module procedure pack_db1
17 module procedure pack_i0
18 module procedure pack_i1
19 module procedure pack_l0
20 end interface pack_dbuf
21
22 CONTAINS
23
24
25
26
27
28
29
30 SUBROUTINE DESMPI_PACK_GHOSTPAR(PFACE)
31
32
33
34
35 use desmpi, only: iGhostPacketSize
36
37 use discretelement, only: iGHOST_UPDATED
38
39 use desmpi, only: dSENDBUF
40
41 use desmpi, only: iBUFOFFSET
42
43 use run, only: ENERGY_EQ
44
45 use run, only: ANY_SPECIES_EQ
46
47 use mfix_pic, only: MPPIC
48
49 use desmpi, only: iNEIGHPROC
50
51 use desgrid, only: DG_IJKSIZE2
52
53 use discretelement, only: DG_PIJK, DG_PIJKPRV
54
55 use discretelement, only: iGLOBAL_ID
56
57 use discretelement, only: DES_POS_NEW, DES_POS_OLD
58
59 use discretelement, only: DES_VEL_NEW, DES_VEL_OLD
60
61 use discretelement, only: OMEGA_NEW, OMEGA_OLD
62
63 use des_rxns, only: DES_X_s
64
65 use des_thermo, only: DES_T_s_NEW, DES_T_s_OLD
66
67 use discretelement, only: DES_RADIUS, PVOL
68
69 use particle_filter, only: FILTER_SIZE
70
71 use particle_filter, only: FILTER_CELL, FILTER_WEIGHT
72
73 use discretelement, only: PEA
74
75 use discretelement, only: PIJK
76
77 use discretelement, only: DO_OLD
78
79 use discretelement, only: DO_NSEARCH
80
81 use discretelement, only: PIP, MAX_PIP
82
83 use discretelement, only: iGHOST_CNT
84
85 use discretelement, only: DES_USR_VAR, DES_USR_VAR_SIZE
86
87 use desgrid, only: dg_ijkconv
88
89 use desmpi, only: isendcnt
90
91 use desmpi, only: dcycl_offset
92
93 use discretelement, only: DG_PIC
94
95 use desmpi, only: iSENDINDICES
96
97
98
99
100 use constant, only: PI
101
102 use discretelement, only: DIMN
103
104
105 IMPLICIT NONE
106
107
108
109
110 INTEGER, INTENT(IN) :: PFACE
111
112
113
114 integer :: lijk,lindx,ltot_ind,lpicloc,lpar_cnt,lcurpar
115 integer :: lbuf
116
117
118 = 0
119 ltot_ind = isendindices(1,pface)
120 do lindx = 2,ltot_ind+1
121 lijk = isendindices(lindx,pface)
122 do lpicloc =1,dg_pic(lijk)%isize
123 lbuf = lpar_cnt*iGhostPacketSize + ibufoffset
124 lcurpar = dg_pic(lijk)%p(lpicloc)
125
126
127
128 if(pea(lcurpar,4) .and. .not.ighost_updated(lcurpar)) cycle
129
130
131 call pack_dbuf(lbuf,iglobal_id(lcurpar),pface)
132
133 call pack_dbuf(lbuf,dg_ijkconv(lijk,pface, &
134 ineighproc(pface)),pface)
135
136 call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lcurpar),pface, &
137 ineighproc(pface)),pface)
138
139 call pack_dbuf(lbuf,des_radius(lcurpar),pface)
140
141 call pack_dbuf(lbuf,pijk(lcurpar,5),pface)
142
143 call pack_dbuf(lbuf,des_pos_new(:,lcurpar)+ &
144 dcycl_offset(pface,:),pface)
145
146 call pack_dbuf(lbuf,des_vel_new(:,lcurpar),pface)
147
148 call pack_dbuf(lbuf,omega_new(:,lcurpar),pface)
149
150 call pack_dbuf(lbuf,merge(1,0,pea(lcurpar,3)),pface)
151
152 IF(ENERGY_EQ) &
153 call pack_dbuf(lbuf,des_t_s_new(lcurpar),pface)
154
155 IF(DES_USR_VAR_SIZE > 0) &
156 call pack_dbuf(lbuf,des_usr_var(:,lcurpar),pface)
157
158 IF(FILTER_SIZE > 0) THEN
159 call pack_dbuf(lbuf,filter_cell(:,lcurpar),pface)
160 call pack_dbuf(lbuf,filter_weight(:,lcurpar),pface)
161 ENDIF
162
163 lpar_cnt = lpar_cnt + 1
164 end do
165 end do
166 dsendbuf(1,pface)=lpar_cnt
167 isendcnt(pface) = lpar_cnt*iGhostPacketSize + ibufoffset
168
169 end subroutine desmpi_pack_ghostpar
170
171
172
173
174
175
176
177
178 SUBROUTINE DESMPI_PACK_PARCROSS(PFACE)
179
180
181
182
183 use desmpi, only: iSPOT
184
185 use discretelement, only: iGHOST_UPDATED
186
187 use desmpi, only: dSENDBUF
188
189 use desmpi, only: iBUFOFFSET
190
191 use run, only: ENERGY_EQ
192
193 use run, only: ANY_SPECIES_EQ
194
195 use mfix_pic, only: MPPIC
196
197 use desgrid, only: DG_IJKSIZE2
198
199 use discretelement, only: DG_PIJK, DG_PIJKPRV
200
201 use desmpi, only: iNEIGHPROC
202
203 use mfix_pic, only: DES_STAT_WT
204
205 use discretelement, only: iGLOBAL_ID
206
207 use discretelement, only: DES_POS_NEW, DES_POS_OLD
208
209 use discretelement, only: DES_VEL_NEW, DES_VEL_OLD
210
211 use discretelement, only: OMEGA_NEW, OMEGA_OLD
212
213 use discretelement, only: PARTICLE_ORIENTATION,ORIENTATION
214
215 use discretelement, only: DES_RADIUS, PVOL, RO_SOL, PMASS
216
217 use discretelement, only: DES_ACC_OLD, ROT_ACC_OLD
218
219 use des_rxns, only: DES_X_s
220
221 use des_thermo, only: DES_T_s_NEW, DES_T_s_OLD
222
223 use discretelement, only: FC, TOW
224
225 use discretelement, only: OMOI
226
227 use discretelement, only: PEA
228
229 use discretelement, only: PIJK
230
231 use discretelement, only: DO_OLD
232
233 use discretelement, only: DO_NSEARCH
234
235 use discretelement, only: PIP, MAX_PIP
236
237 use discretelement, only: iGHOST_CNT
238
239 use discretelement, only: DES_USR_VAR, DES_USR_VAR_SIZE
240
241 use discretelement, only: PAIR_NUM, PAIRS
242
243 use discretelement, only: PV_PAIR, PFN_PAIR, PFT_PAIR
244
245 use discretelement, only: DIMN
246
247 use compar, only: myPE
248
249 use discretelement, only: DES_EXPLICITLY_COUPLED
250
251 use discretelement, only: DRAG_FC
252
253 use particle_filter, only: FILTER_SIZE
254
255 use particle_filter, only: FILTER_CELL, FILTER_WEIGHT
256
257 use desgrid, only: dg_ijkconv, icycoffset
258 use desmpi, only: dcycl_offset, isendcnt
259 use discretelement, only: DG_PIC
260 use desmpi, only: iSENDINDICES
261 use desmpi, only: irecvindices
262
263 use desmpi, only: iParticlePacketSize
264 use desmpi, only: iPairPacketSize
265
266 use functions
267
268 implicit none
269
270
271
272
273 INTEGER, INTENT(IN) :: PFACE
274
275
276
277 integer :: li, lj, lk
278 integer :: ltot_ind,lindx,cc,ii
279 integer :: lneigh,lijk,&
280 lpicloc,lparcnt,lcurpar
281 integer :: lbuf,num_pairs_to_send
282
283 logical, allocatable, dimension(:) :: going_to_send
284
285
286 integer :: num_pairs_send_buf_loc
287
288
289
290 = irecvindices(1,pface)
291 lparcnt = 0
292
293 allocate(going_to_send(max_pip))
294 going_to_send(:) = .false.
295
296 do lindx = 2,ltot_ind+1
297 lijk = irecvindices(lindx,pface)
298 do lpicloc = 1,dg_pic(lijk)%isize
299 lcurpar = dg_pic(lijk)%p(lpicloc)
300
301 if (pea(lcurpar,4)) cycle
302
303 (lcurpar) = .true.
304 lbuf = lparcnt*iParticlePacketSize + ibufoffset
305
306
307 call pack_dbuf(lbuf,iglobal_id(lcurpar),pface)
308
309 call pack_dbuf(lbuf,dg_ijkconv(lijk,pface, &
310 ineighproc(pface)),pface)
311
312 call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lcurpar),pface, &
313 ineighproc(pface)),pface)
314
315 call pack_dbuf(lbuf,des_radius(lcurpar),pface)
316
317 = pijk(lcurpar,1) + icycoffset(pface,1)
318 call pack_dbuf(lbuf,li,pface)
319
320 = pijk(lcurpar,2) + icycoffset(pface,2)
321 call pack_dbuf(lbuf,lj,pface)
322
323 = pijk(lcurpar,3) + icycoffset(pface,3)
324 call pack_dbuf(lbuf,lk,pface)
325
326 call pack_dbuf(lbuf,funijk_proc(li,lj,lk, &
327 ineighproc(pface)),pface)
328
329 call pack_dbuf(lbuf,pijk(lcurpar,5),pface)
330
331 call pack_dbuf(lbuf, pea(lcurpar,2), pface)
332
333 call pack_dbuf(lbuf, pea(lcurpar,3), pface)
334
335 call pack_dbuf(lbuf,ro_sol(lcurpar),pface)
336
337 call pack_dbuf(lbuf,pvol(lcurpar),pface)
338
339 call pack_dbuf(lbuf,pmass(lcurpar),pface)
340
341 call pack_dbuf(lbuf,omoi(lcurpar),pface)
342
343 call pack_dbuf(lbuf,des_pos_new(:,lcurpar) + &
344 dcycl_offset(pface,:),pface)
345
346 call pack_dbuf(lbuf,des_vel_new(:,lcurpar),pface)
347
348 call pack_dbuf(lbuf,omega_new(:,lcurpar),pface)
349
350 call pack_dbuf(lbuf,fc(:,lcurpar),pface)
351
352 call pack_dbuf(lbuf,tow(:,lcurpar),pface)
353
354 IF(ENERGY_EQ) &
355 call pack_dbuf(lbuf,des_t_s_new(lcurpar),pface)
356
357 IF(ANY_SPECIES_EQ) &
358 call pack_dbuf(lbuf,des_x_s(lcurpar,:),pface)
359
360 IF(DES_EXPLICITLY_COUPLED) &
361 call pack_dbuf(lbuf, drag_fc(:,lcurpar),pface)
362
363 IF(DES_USR_VAR_SIZE > 0) &
364 call pack_dbuf(lbuf, des_usr_var(:,lcurpar),pface)
365
366 IF(PARTICLE_ORIENTATION) &
367 call pack_dbuf(lbuf,orientation(:,lcurpar),pface)
368
369
370 IF (DO_OLD) THEN
371
372 call pack_dbuf(lbuf,des_pos_old(:,lcurpar) + &
373 dcycl_offset(pface,:),pface)
374
375 call pack_dbuf(lbuf,des_vel_old(:,lcurpar),pface)
376
377 call pack_dbuf(lbuf,omega_old(:,lcurpar),pface)
378
379 call pack_dbuf(lbuf,des_acc_old(:,lcurpar),pface)
380
381 call pack_dbuf(lbuf,rot_acc_old(:,lcurpar),pface)
382
383 IF(ENERGY_EQ) &
384 call pack_dbuf(lbuf,des_t_s_old(lcurpar),pface)
385 ENDIF
386
387
388
389 IF (MPPIC) THEN
390
391 call pack_dbuf(lbuf,des_stat_wt(lcurpar),pface)
392 pea(lcurpar,1:4) = .false.
393 pip = pip - 1
394
395
396
397 ELSE
398 pea(lcurpar,4) = .true.
399 ighost_cnt = ighost_cnt + 1
400 END IF
401
402
403 (:,lcurpar) = 0.
404 lparcnt = lparcnt + 1
405 end do
406 end do
407
408
409
410
411 = lparcnt*iParticlePacketSize + ibufoffset
412 num_pairs_send_buf_loc = lbuf
413 lbuf = lbuf+1
414
415 num_pairs_to_send = 0
416 do cc = 1, pair_num
417 lcurpar = PAIRS(1,CC)
418
419 if (.not. going_to_send(lcurpar)) cycle
420
421
422
423 = PAIRS(2,CC)
424 if(.not.PEA(lneigh,1)) cycle
425 if(PEA(lneigh,3)) cycle
426
427 call pack_dbuf(lbuf,iglobal_id(lcurpar),pface)
428
429 call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lcurpar),pface, &
430 ineighproc(pface)),pface)
431
432 call pack_dbuf(lbuf,iglobal_id(lneigh),pface)
433
434 call pack_dbuf(lbuf,dg_ijkconv(dg_pijkprv(lneigh),pface, &
435 ineighproc(pface)),pface)
436
437 call pack_dbuf(lbuf,pv_PAIR(CC),pface)
438
439 call pack_dbuf(lbuf,PFN_PAIR(:,CC),pface)
440
441 call pack_dbuf(lbuf,PFT_PAIR(:,CC),pface)
442
443 = num_pairs_to_send + 1
444 enddo
445
446
447
448
449 = num_pairs_send_buf_loc
450
451 call pack_dbuf(lbuf,num_pairs_to_send,pface)
452
453
454 dsendbuf(1,pface) = lparcnt
455 isendcnt(pface) = lparcnt*iParticlePacketSize + &
456 num_pairs_to_send*iPairPacketSize + ibufoffset + 3
457
458 deallocate(going_to_send)
459
460 RETURN
461 END SUBROUTINE DESMPI_PACK_PARCROSS
462
463
464
465
466 subroutine pack_db0(lbuf,idata,pface)
467 use desmpi, only: dSENDBUF
468 integer, intent(inout) :: lbuf
469 integer, intent(in) :: pface
470 double precision, intent(in) :: idata
471
472 dsendbuf(lbuf,pface) = idata
473 lbuf = lbuf + 1
474
475 return
476 end subroutine pack_db0
477
478
479
480
481 subroutine pack_db1(lbuf,idata,pface)
482 use desmpi, only: dSENDBUF
483 integer, intent(inout) :: lbuf
484 integer, intent(in) :: pface
485 double precision, intent(in) :: idata(:)
486
487 integer :: lsize
488
489 lsize = size(idata)
490
491 dsendbuf(lbuf:lbuf+lsize-1,pface) = idata
492 lbuf = lbuf + lsize
493
494 return
495 end subroutine pack_db1
496
497
498
499
500
501 subroutine pack_i0(lbuf,idata,pface)
502 use desmpi, only: dSENDBUF
503 integer, intent(inout) :: lbuf
504 integer, intent(in) :: pface
505 integer, intent(in) :: idata
506
507 dsendbuf(lbuf,pface) = idata
508 lbuf = lbuf + 1
509
510 return
511 end subroutine pack_i0
512
513
514
515
516 subroutine pack_i1(lbuf,idata,pface)
517 use desmpi, only: dSENDBUF
518 integer, intent(inout) :: lbuf
519 integer, intent(in) :: pface
520 integer, intent(in) :: idata(:)
521
522 integer :: lsize
523
524 lsize = size(idata)
525
526 dsendbuf(lbuf:lbuf+lsize-1,pface) = idata
527 lbuf = lbuf + lsize
528
529 return
530 end subroutine pack_i1
531
532
533
534
535 subroutine pack_l0(lbuf,ldata,pface)
536 use desmpi, only: dSENDBUF
537
538 integer, intent(inout) :: lbuf
539 integer, intent(in) :: pface
540 logical, intent(in) :: ldata
541
542 dsendbuf(lbuf,pface) = merge(1.0, 0.0, ldata)
543 lbuf = lbuf + 1
544
545 return
546 end subroutine pack_l0
547
548
549 end module mpi_pack_des
550