File: RELATIVE:/../../../mfix.git/model/des/sendrecvnode_mod.f
1
2
3
4
5
6
7
8
9
10
11
12
13
14 module sendrecvnode
15
16
17
18
19 use parallel_mpi
20 use mpi_utility
21 use discretelement
22 use compar
23 use physprop
24 use sendrecv
25 use desmpi_wrapper
26 use desgrid
27 use functions
28
29
30 integer :: itotalneigh,itotalindx
31 integer,dimension(:),allocatable :: itoproc,iprocsumindx,istartindx
32
33
34 integer,dimension(:),allocatable :: isendnodes
35 double precision,dimension(:),allocatable :: dsendnodebuf,drecvnodebuf
36 integer,dimension(:),allocatable :: irecvreqnode,isendreqnode
37
38
39 contains
40
41
42
43
44
45
46
47 subroutine des_setnodeindices
48
49 implicit none
50
51
52
53 integer lijkproc,liproc,ljproc,lkproc
54 integer li,lj,lk
55 integer li2,lj2,lk2
56 integer liproc_start,liproc_end,ljproc_start,ljproc_end,lkproc_start,lkproc_end
57 integer lci,lcj,lck,lproc,lcount
58 integer linode,ljnode,lknode
59 integer linode_start,linode_end,ljnode_start,ljnode_end,lknode_start,lknode_end
60 logical lpresent
61
62
63
64 = iofproc(mype)
65 ljproc = jofproc(mype)
66 lkproc = kofproc(mype)
67
68
69 if(des_periodic_walls_x.and.nodesi.gt.1) then
70 liproc_start=liproc-1
71 liproc_end=liproc+1
72 else
73 liproc_start =max(liproc-1,0)
74 liproc_end=min(liproc+1,nodesi-1)
75 end if
76 if(des_periodic_walls_y.and.nodesj.gt.1) then
77 ljproc_start=ljproc-1
78 ljproc_end=ljproc+1
79 else
80 ljproc_start =max(ljproc-1,0)
81 ljproc_end=min(ljproc+1,nodesj-1)
82 end if
83 if(des_periodic_walls_z.and.nodesk.gt.1) then
84 lkproc_start=lkproc-1
85 lkproc_end=lkproc+1
86 else
87 lkproc_start =max(lkproc-1,0)
88 lkproc_end=min(lkproc+1,nodesk-1)
89 end if
90 itotalneigh = (liproc_end-liproc_start+1)*(ljproc_end-ljproc_start+1)* &
91 (lkproc_end-lkproc_start+1)-1
92
93 allocate (itoproc(itotalneigh),iprocsumindx(itotalneigh),istartindx(itotalneigh+1), &
94 irecvreqnode(itotalneigh),isendreqnode(itotalneigh))
95
96
97
98 = 0
99 itoproc(:)=-1
100 iprocsumindx(:) =0
101 do lk = lkproc_start,lkproc_end
102 do lj = ljproc_start,ljproc_end
103 do li = liproc_start,liproc_end
104 li2 = mod(li,nodesi);if(li2.lt.0)li2=nodesi-1
105 lj2 = mod(lj,nodesj);if(lj2.lt.0)lj2=nodesj-1
106 lk2 = mod(lk,nodesk);if(lk2.lt.0)lk2=nodesk-1
107 lijkproc = procijk(li2,lj2,lk2)
108 if (lijkproc.eq.mype) cycle
109
110 = .false.
111 do lproc = 1,itotalneigh
112 if (lijkproc .eq.itoproc(lproc)) then
113 lpresent = .true.
114 exit
115 end if
116 end do
117 if(.not.lpresent) then
118 itotalneigh = itotalneigh + 1
119 lproc = itotalneigh
120 end if
121 itoproc(lproc) = lijkproc
122 lci=(liproc-li);lcj=(ljproc-lj);lck=(lkproc-lk)
123 linode_start = istart2; linode_end=iend1
124 ljnode_start = jstart2; ljnode_end=jend1
125 lknode_start = kstart2; lknode_end=kend1
126 if(lci.eq.1) linode_end = istart2
127 if(lci.eq.-1) linode_start = iend1
128 if(lcj.eq.1) ljnode_end = jstart2
129 if(lcj.eq.-1) ljnode_start = jend1
130 if(lck.eq.1) lknode_end = kstart2
131 if(lck.eq.-1) lknode_start = kend1
132 do lknode = lknode_start,lknode_end
133 do linode = linode_start,linode_end
134 do ljnode = ljnode_start,ljnode_end
135 IF(DEAD_CELL_AT(linode,ljnode,lknode)) CYCLE
136 iprocsumindx(lproc) = iprocsumindx(lproc) + 1
137 end do
138 end do
139 end do
140 end do
141 end do
142 end do
143
144 do lproc =1,itotalneigh+1
145 istartindx(lproc)=sum(iprocsumindx(1:lproc-1))+1
146 end do
147 itotalindx=istartindx(itotalneigh+1)-1
148
149
150 allocate(isendnodes(itotalindx),dsendnodebuf(itotalindx),drecvnodebuf(itotalindx))
151
152
153 (:)=0
154 do lk = lkproc_start,lkproc_end
155 do lj = ljproc_start,ljproc_end
156 do li = liproc_start,liproc_end
157 li2 = mod(li,nodesi);if(li2.lt.0)li2=nodesi-1
158 lj2 = mod(lj,nodesj);if(lj2.lt.0)lj2=nodesj-1
159 lk2 = mod(lk,nodesk);if(lk2.lt.0)lk2=nodesk-1
160 lijkproc = procijk(li2,lj2,lk2)
161 if (lijkproc.eq.mype) cycle
162
163 do lproc =1,itotalneigh
164 if(lijkproc.eq.itoproc(lproc)) then
165 exit
166 end if
167 end do
168 lci=(liproc-li);lcj=(ljproc-lj);lck=(lkproc-lk)
169 linode_start = istart2; linode_end=iend1
170 ljnode_start = jstart2; ljnode_end=jend1
171 lknode_start = kstart2; lknode_end=kend1
172 if(lci.eq.1) linode_end = istart2
173 if(lci.eq.-1) linode_start = iend1
174 if(lcj.eq.1) ljnode_end = jstart2
175 if(lcj.eq.-1) ljnode_start = jend1
176 if(lck.eq.1) lknode_end = kstart2
177 if(lck.eq.-1) lknode_start = kend1
178 lcount = istartindx(lproc)+iprocsumindx(lproc)
179 do lknode = lknode_start,lknode_end
180 do linode = linode_start,linode_end
181 do ljnode = ljnode_start,ljnode_end
182 IF(DEAD_CELL_AT(linode,ljnode,lknode)) CYCLE
183 isendnodes(lcount)=funijk(linode,ljnode,lknode)
184 iprocsumindx(lproc)=iprocsumindx(lproc)+1
185 lcount = lcount+1
186 end do
187 end do
188 end do
189 end do
190 end do
191 end do
192
193
194 end subroutine des_setnodeindices
195
196
197
198
199
200
201
202
203
204
205 subroutine des_exchangenode(pvar,padd)
206
207 implicit none
208
209
210
211 double precision,dimension(:),intent(inout) ::pvar
212 logical:: padd
213
214
215
216 character(len=80), parameter :: name = 'des_exchangenode'
217 integer :: lindx,lcount,lcount2,lneigh,ltag,lerr
218 integer :: lstart,lend,ltotal
219
220
221
222 do lcount = 1,itotalneigh
223 lneigh = itoproc(lcount)
224 lstart = istartindx(lcount);lend=istartindx(lcount+1)-1
225 do lcount2 = lstart,lend
226 dsendnodebuf(lcount2) = pvar(isendnodes(lcount2))
227 end do
228 ltag = message_tag(lneigh,mype)
229 ltotal = lend-lstart+1
230 call des_mpi_irecv(drecvnodebuf(lstart:lend),ltotal, &
231 lneigh,ltag,irecvreqnode(lcount),lerr)
232 call mpi_check( name //':mpi_irecv ', lerr )
233 ltag = message_tag(mype,lneigh)
234 call des_mpi_isend(dsendnodebuf(lstart:lend),ltotal, &
235 lneigh,ltag,isendreqnode(lcount),lerr)
236 call mpi_check( name //':mpi_irecv ', lerr )
237 end do
238
239 do lcount = 1,itotalneigh
240 call des_mpi_wait(isendreqnode(lcount),lerr)
241 call mpi_check( name //':mpi_wait-send', lerr )
242 call des_mpi_wait(irecvreqnode(lcount),lerr)
243 call mpi_check( name //':mpi_wait-recv', lerr )
244 end do
245
246
247 if (padd) then
248 do lcount = 1,itotalindx
249 lindx = isendnodes(lcount)
250 pvar(lindx) = pvar(lindx) + drecvnodebuf(lcount)
251 end do
252 else
253 do lcount = 1,itotalindx
254 lindx = isendnodes(lcount)
255 pvar(lindx) = drecvnodebuf(lcount)
256 end do
257 end if
258 return
259
260 contains
261
262 integer function message_tag(lsource,ldest)
263 implicit none
264 integer, intent(in) :: lsource,ldest
265 message_tag = lsource+numpes*ldest+200
266 end function message_tag
267
268 end subroutine des_exchangenode
269
270
271
272
273
274 subroutine des_dbgnodesr()
275
276
277 use functions
278
279 implicit none
280
281
282
283 character (255) filename
284 integer ijk
285 integer lcount,lcount2,lstart,lend
286
287
288
289 write(filename,'("dbg_nodesr",I4.4,".dat")') mype
290 open(44,file=filename,convert='big_endian')
291 do lcount = 1,itotalneigh
292 lstart = istartindx(lcount);lend=istartindx(lcount+1)-1
293 write(44,*) "--------------------------------------------"
294 write(44,*) "to proc start end",itoproc(lcount),lstart,lend
295 write(44,*) "--------------------------------------------"
296 do lcount2 = lstart,lend
297 ijk = isendnodes(lcount2)
298 write(44,*)ijk,i_of(ijk),j_of(ijk),k_of(ijk)
299 end do
300 write(44,*) "--------------------------------------------"
301 end do
302 close (44)
303 end subroutine des_dbgnodesr
304
305 end module
306
307
308