File: /nfs/home/0/users/jenkins/mfix.git/model/des/mpi_funs_des_mod.f
1
2
3
4
5
6
7
8 module mpi_funs_des
9
10 contains
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40 subroutine des_par_exchange()
41
42 use discretelement, only: DO_NSEARCH
43
44
45 use mfix_pic, only: MPPIC
46 use desmpi, only: iEXCHFLAG
47 use desmpi, only: dSENDBUF, dRECVBUF
48 use discretelement, only: iGHOST_UPDATED
49 use desmpi, only: iSPOT
50
51 use geometry, only: NO_K
52
53
54
55 use mpi_pack_des, only: desmpi_pack_parcross
56 use mpi_unpack_des, only: desmpi_unpack_parcross
57
58 use mpi_pack_des, only: desmpi_pack_ghostpar
59 use mpi_unpack_des, only: desmpi_unpack_ghostpar
60
61 use mpi_comm_des, only: desmpi_sendrecv_init
62 use mpi_comm_des, only: desmpi_sendrecv_wait
63
64 use desgrid, only: desgrid_pic
65 use desmpi_wrapper, only: des_mpi_barrier
66
67 implicit none
68
69
70
71
72 integer :: linter, lface
73
74 integer, save :: lcheckbuf = 0
75
76
77
78
79 call desgrid_pic(plocate=.true.)
80
81
82
83 if (mod(lcheckbuf,100) .eq. 0) then
84 call desmpi_check_sendrecvbuf
85 lcheckbuf = 0
86 end if
87 lcheckbuf = lcheckbuf + 1
88
89
90 (1,:) = 0; drecvbuf(1,:) =0
91 ispot = 1
92 do linter = merge(2,3,NO_K),1,-1
93 do lface = linter*2-1,linter*2
94 if(.not.iexchflag(lface))cycle
95 call desmpi_pack_parcross(lface)
96 call desmpi_sendrecv_init(lface)
97 end do
98 do lface = linter*2-1,linter*2
99 if(.not.iexchflag(lface)) cycle
100 call desmpi_sendrecv_wait(lface)
101 call desmpi_unpack_parcross(lface)
102 end do
103
104 do lface = linter*2-1,linter*2
105 if(dsendbuf(1,lface).gt.0.or.drecvbuf(1,lface).gt.0) then
106 call desgrid_pic(plocate=.false.)
107 exit
108 end if
109 end do
110 end do
111 call des_mpi_barrier
112
113
114
115
116 IF(.NOT.MPPIC) THEN
117
118 (1,:) = 0; drecvbuf(1,:) =0
119 ighost_updated(:) = .false.
120 ispot = 1
121 do linter = 1,merge(2,3,NO_K)
122 do lface = linter*2-1,linter*2
123 if(.not.iexchflag(lface))cycle
124 call desmpi_pack_ghostpar(lface)
125 call desmpi_sendrecv_init(lface)
126 end do
127 do lface = linter*2-1,linter*2
128 if(.not.iexchflag(lface)) cycle
129 call desmpi_sendrecv_wait(lface)
130 call desmpi_unpack_ghostpar(lface)
131 end do
132
133
134 do lface = linter*2-1,linter*2
135 if(dsendbuf(1,lface).gt.0.or.drecvbuf(1,lface).gt.0) then
136 call desgrid_pic(plocate=.false.)
137 exit
138 end if
139 end do
140 end do
141 if(do_nsearch) call desmpi_cleanup
142 call des_mpi_barrier
143 ENDIF
144
145
146
147
148
149
150
151 END SUBROUTINE DES_PAR_EXCHANGE
152
153
154
155
156
157
158
159
160
161 SUBROUTINE DESMPI_CHECK_SENDRECVBUF
162
163 use discretelement, only: DIMN, dg_pic
164 use desmpi, only: iMAXBUF
165 use desmpi, only: iBUFOFFSET
166 use desmpi, only: dSENDBUF, dRECVBUF
167 use desmpi, only: iSENDINDICES
168
169 use mpi_utility, only: global_all_max
170
171 implicit none
172
173
174
175
176 INTEGER :: lface, lindx, lijk
177
178 INTEGER :: lpacketsize
179
180 INTEGER :: lparcnt
181
182 INTEGER :: lmaxcnt
183
184 INTEGER :: ltot_ind
185
186 REAL :: lfactor = 1.5
187
188
189
190 = 0
191 lpacketsize = 2*dimn + 3 + 5
192 do lface = 1,2*dimn
193 ltot_ind = isendindices(1,lface)
194 lparcnt = 0
195 do lindx = 2,ltot_ind+1
196 lijk = isendindices(lindx,lface)
197 lparcnt = lparcnt + dg_pic(lijk)%isize
198 enddo
199 if(lparcnt.gt.lmaxcnt) lmaxcnt = lparcnt
200 enddo
201
202 call global_all_max(lmaxcnt)
203 if (imaxbuf .lt. lmaxcnt*lpacketsize+ibufoffset) then
204 imaxbuf = lmaxcnt*lpacketsize*lfactor
205 if(allocated(dsendbuf)) deallocate(dsendbuf,drecvbuf)
206 allocate(dsendbuf(imaxbuf,2*dimn),drecvbuf(imaxbuf,2*dimn))
207 endif
208
209 END SUBROUTINE DESMPI_CHECK_SENDRECVBUF
210
211
212
213
214
215
216
217 SUBROUTINE DESMPI_CLEANUP
218
219 use discretelement, only: DIMN
220 use discretelement, only: PEA
221 use discretelement, only: DES_POS_NEW, DES_POS_OLD
222 use discretelement, only: DES_VEL_NEW, DES_VEL_OLD
223 use discretelement, only: OMEGA_NEW, OMEGA_OLD
224 use discretelement, only: PARTICLE_ORIENTATION,ORIENTATION,INIT_ORIENTATION
225 use discretelement, only: FC
226 use discretelement, only: DO_OLD
227 use discretelement, only: PIP
228 use discretelement, only: iGHOST_CNT
229 use discretelement, only: DES_USR_VAR_SIZE, DES_USR_VAR
230 use discretelement, only: dg_pic, pijk
231
232 use run, only: ENERGY_EQ,ANY_SPECIES_EQ
233 use des_thermo, only: DES_T_s_NEW, DES_T_s_OLD
234
235 use des_rxns, only: DES_X_s
236
237 use discretelement, only: iGHOST_UPDATED
238 use desmpi, only: iRECVINDICES
239 use desmpi, only: iEXCHFLAG
240
241 use param, only: DIMENSION_N_s
242 use param1, only: ZERO
243
244 implicit none
245
246
247
248 integer ltot_ind,lface,lindx,lijk,lcurpar,lpicloc
249
250 do lface = 1,dimn*2
251 if(.not.iexchflag(lface))cycle
252 ltot_ind = irecvindices(1,lface)
253 do lindx = 2,ltot_ind+1
254 lijk = irecvindices(lindx,lface)
255 do lpicloc =1,dg_pic(lijk)%isize
256 lcurpar = dg_pic(lijk)%p(lpicloc)
257 if(ighost_updated(lcurpar)) cycle
258 pip = pip - 1
259 ighost_cnt = ighost_cnt-1
260 pea(lcurpar,1:4) = .false.
261 fc(:,lcurpar) = 0.0
262 des_pos_new(:,lcurpar)=0
263 pijk(lcurpar,:) = -10
264 IF (DO_OLD) THEN
265 des_pos_old(:,lcurpar)=0
266 des_vel_old(:,lcurpar)=0
267 ENDIF
268 des_vel_new(:,lcurpar)=0
269 omega_new(:,lcurpar)=0
270
271 IF(PARTICLE_ORIENTATION) ORIENTATION(1:3,lcurpar) = INIT_ORIENTATION
272
273 if(ENERGY_EQ) then
274 des_t_s_new(lcurpar)=0
275 des_t_s_old(lcurpar)=0
276 endif
277
278 if(ANY_SPECIES_EQ)then
279 des_x_s(lcurpar,1:dimension_n_s)= 0
280 endif
281
282 IF(DES_USR_VAR_SIZE > 0)&
283 des_usr_var(:,lcurpar)= 0
284
285 end do
286 end do
287 end do
288 END SUBROUTINE DESMPI_CLEANUP
289
290
291 END MODULE MPI_FUNS_DES
292