14 public :: hecmwst_contact_comm
24 type hecmwst_contact_comm
26 integer(kind=kint) :: n_neighbor_pe
27 integer(kind=kint),
pointer :: neighbor_pe(:)
28 integer(kind=kint) :: MPI_COMM
29 integer(kind=kint),
pointer :: ext_index(:)
30 integer(kind=kint),
pointer :: ext_item(:)
31 integer(kind=kint),
pointer :: int_index(:)
32 integer(kind=kint),
pointer :: int_item(:)
33 end type hecmwst_contact_comm
35 integer(kind=kint),
parameter :: op_overwrite = 46810
37 integer(kind=kint),
parameter :: DBG = 0
43 type (hecmwst_contact_comm),
intent(inout) :: concomm
45 integer(kind=kint),
intent(in) :: ndof, n_contact_dof
46 integer(kind=kint),
intent(in) :: contact_dofs(:)
47 integer(kind=kint),
pointer :: ext_index(:), ext_item(:), int_index(:), int_item(:)
48 integer(kind=kint),
allocatable :: n_ext_per_dom(:), n_int_per_dom(:), ext_item_remote(:)
49 integer(kind=kint),
allocatable :: statuses(:,:), requests(:)
50 integer(kind=kint) :: nn_int, np, ilag, icontact, irow, irank, idom, tag, idof, idx, irow_remote
51 integer(kind=kint) :: n_send, is, ie, len
52 if (hecmesh%n_neighbor_pe == 0)
then
53 concomm%n_neighbor_pe = 0
56 nn_int = hecmesh%nn_internal
59 allocate(n_ext_per_dom(hecmesh%n_neighbor_pe))
61 do ilag = 1, n_contact_dof
62 icontact = contact_dofs(ilag)
63 if (icontact <= nn_int*ndof) cycle
64 irow = (icontact+ndof-1) / ndof
65 irank = hecmesh%node_ID(2*irow)
66 call rank_to_idom(hecmesh, irank, idom)
67 n_ext_per_dom(idom) = n_ext_per_dom(idom) + 1
71 allocate(requests(hecmesh%n_neighbor_pe))
72 do idom = 1, hecmesh%n_neighbor_pe
73 irank = hecmesh%neighbor_pe(idom)
76 hecmesh%MPI_COMM, requests(idom))
78 allocate(n_int_per_dom(hecmesh%n_neighbor_pe))
79 do idom = 1, hecmesh%n_neighbor_pe
80 irank = hecmesh%neighbor_pe(idom)
83 hecmesh%MPI_COMM, statuses(:,1))
87 allocate(ext_index(0:hecmesh%n_neighbor_pe))
88 allocate(int_index(0:hecmesh%n_neighbor_pe))
91 do idom = 1, hecmesh%n_neighbor_pe
92 ext_index(idom) = ext_index(idom-1) + n_ext_per_dom(idom)
93 int_index(idom) = int_index(idom-1) + n_int_per_dom(idom)
96 allocate(ext_item(ext_index(hecmesh%n_neighbor_pe)))
97 allocate(ext_item_remote(ext_index(hecmesh%n_neighbor_pe)))
99 do ilag = 1, n_contact_dof
100 icontact = contact_dofs(ilag)
101 if (icontact <= nn_int*ndof) cycle
102 irow = (icontact+ndof-1) / ndof
103 idof = icontact - ndof*(irow-1)
104 irank = hecmesh%node_ID(2*irow)
105 call rank_to_idom(hecmesh, irank, idom)
106 n_ext_per_dom(idom) = n_ext_per_dom(idom) + 1
107 idx = ext_index(idom-1)+n_ext_per_dom(idom)
108 ext_item(idx) = icontact
109 irow_remote = hecmesh%node_ID(2*irow-1)
110 ext_item_remote(idx) = ndof*(irow_remote-1)+idof
112 deallocate(n_ext_per_dom)
113 deallocate(n_int_per_dom)
116 do idom = 1, hecmesh%n_neighbor_pe
117 irank = hecmesh%neighbor_pe(idom)
118 is = ext_index(idom-1)+1
125 hecmesh%MPI_COMM, requests(n_send))
127 allocate(int_item(int_index(hecmesh%n_neighbor_pe)))
128 do idom = 1, hecmesh%n_neighbor_pe
129 irank = hecmesh%neighbor_pe(idom)
130 is = int_index(idom-1)+1
136 hecmesh%MPI_COMM, statuses(:,1))
139 deallocate(statuses, requests)
141 write(0,*)
' DEBUG2: ext_index',ext_index(:)
142 write(0,*)
' DEBUG2: ext_item',ext_item(:)
143 write(0,*)
' DEBUG2: ext_item_remote',ext_item_remote(:)
144 write(0,*)
' DEBUG2: int_index',int_index(:)
145 write(0,*)
' DEBUG2: int_item',int_item(:)
147 deallocate(ext_item_remote)
149 concomm%n_neighbor_pe = hecmesh%n_neighbor_pe
150 allocate(concomm%neighbor_pe(concomm%n_neighbor_pe))
151 concomm%neighbor_pe(:) = hecmesh%neighbor_pe(:)
152 concomm%MPI_COMM = hecmesh%MPI_COMM
153 concomm%ext_index => ext_index
154 concomm%ext_item => ext_item
155 concomm%int_index => int_index
156 concomm%int_item => int_item
161 type (hecmwst_contact_comm),
intent(inout) :: concomm
162 if (concomm%n_neighbor_pe == 0)
return
163 if (
associated(concomm%neighbor_pe))
deallocate(concomm%neighbor_pe)
164 if (
associated(concomm%ext_index))
deallocate(concomm%ext_index)
165 if (
associated(concomm%ext_item))
deallocate(concomm%ext_item)
166 if (
associated(concomm%int_index))
deallocate(concomm%int_index)
167 if (
associated(concomm%int_item))
deallocate(concomm%int_item)
168 concomm%n_neighbor_pe = 0
174 type (hecmwst_contact_comm),
intent(in) :: concomm
175 real(kind=
kreal),
intent(inout) :: vec(:)
176 integer(kind=kint),
intent(in) :: op
177 if (concomm%n_neighbor_pe == 0)
return
178 call send_recv_contact_info_r(concomm%n_neighbor_pe, concomm%neighbor_pe, concomm%MPI_COMM, &
179 concomm%ext_index, concomm%ext_item, concomm%int_index, concomm%int_item, vec, op)
184 type (hecmwst_contact_comm),
intent(in) :: concomm
185 real(kind=
kreal),
intent(inout) :: vec(:)
186 integer(kind=kint) :: op
187 if (concomm%n_neighbor_pe == 0)
return
189 call send_recv_contact_info_r(concomm%n_neighbor_pe, concomm%neighbor_pe, concomm%MPI_COMM, &
190 concomm%int_index, concomm%int_item, concomm%ext_index, concomm%ext_item, vec, op)
195 type (hecmwst_contact_comm),
intent(in) :: concomm
196 integer(kind=kint),
intent(inout) :: vec(:)
197 integer(kind=kint),
intent(in) :: op
198 if (concomm%n_neighbor_pe == 0)
return
199 call send_recv_contact_info_i(concomm%n_neighbor_pe, concomm%neighbor_pe, concomm%MPI_COMM, &
200 concomm%ext_index, concomm%ext_item, concomm%int_index, concomm%int_item, vec, op)
205 type (hecmwst_contact_comm),
intent(in) :: concomm
206 integer(kind=kint),
intent(inout) :: vec(:)
207 integer(kind=kint) :: op
208 if (concomm%n_neighbor_pe == 0)
return
210 call send_recv_contact_info_i(concomm%n_neighbor_pe, concomm%neighbor_pe, concomm%MPI_COMM, &
211 concomm%int_index, concomm%int_item, concomm%ext_index, concomm%ext_item, vec, op)
216 type (hecmwst_contact_comm),
intent(in) :: concomm
217 real(kind=
kreal),
intent(inout) :: vec(:)
218 integer(kind=kint),
intent(in) :: op
225 type (hecmwst_contact_comm),
intent(in) :: concomm
226 integer(kind=kint),
intent(inout) :: vec(:)
227 integer(kind=kint),
intent(in) :: op
236 subroutine rank_to_idom(hecMESH, rank, idom)
239 integer(kind=kint),
intent(in) :: rank
240 integer(kind=kint),
intent(out) :: idom
241 integer(kind=kint) :: i
242 do i = 1, hecmesh%n_neighbor_pe
243 if (hecmesh%neighbor_pe(i) == rank)
then
248 stop
'ERROR: exp_rank not found in neighbor_pe'
249 end subroutine rank_to_idom
251 subroutine send_recv_contact_info_r(n_neighbor_pe, neighbor_pe, MPI_COMM, &
252 send_index, send_item, recv_index, recv_item, vec, op)
254 integer(kind=kint),
intent(in) :: n_neighbor_pe
255 integer(kind=kint),
intent(in) :: neighbor_pe(:)
256 integer(kind=kint),
intent(in) :: mpi_comm
257 integer(kind=kint),
pointer,
intent(in) :: send_index(:), send_item(:), recv_index(:), recv_item(:)
258 real(kind=
kreal),
intent(inout) :: vec(:)
259 integer(kind=kint),
intent(in) :: op
260 real(kind=
kreal),
allocatable :: send_buf(:), recv_buf(:)
261 integer(kind=kint) :: i, n_send, idom, irank, is, ie, len, tag
262 integer(kind=kint),
allocatable :: requests(:), statuses(:,:)
263 if (n_neighbor_pe == 0)
return
264 allocate(requests(n_neighbor_pe))
266 allocate(send_buf(send_index(n_neighbor_pe)))
267 allocate(recv_buf(recv_index(n_neighbor_pe)))
268 do i = 1, send_index(n_neighbor_pe)
269 send_buf(i) = vec(send_item(i))
272 do idom = 1, n_neighbor_pe
273 irank = neighbor_pe(idom)
274 is = send_index(idom-1)+1
275 ie = send_index(idom)
281 mpi_comm, requests(n_send))
283 do idom = 1, n_neighbor_pe
284 irank = neighbor_pe(idom)
285 is = recv_index(idom-1)+1
286 ie = recv_index(idom)
291 mpi_comm, statuses(:,1))
295 do i = 1, recv_index(n_neighbor_pe)
296 vec(recv_item(i)) = vec(recv_item(i)) + recv_buf(i)
299 do i = 1, recv_index(n_neighbor_pe)
300 vec(recv_item(i)) = vec(recv_item(i)) * recv_buf(i)
303 do i = 1, recv_index(n_neighbor_pe)
304 vec(recv_item(i)) = max(vec(recv_item(i)), recv_buf(i))
307 do i = 1, recv_index(n_neighbor_pe)
308 vec(recv_item(i)) = min(vec(recv_item(i)), recv_buf(i))
311 do i = 1, recv_index(n_neighbor_pe)
312 vec(recv_item(i)) = recv_buf(i)
318 write(0,*)
' DEBUG2: send_buf',send_buf(:)
319 write(0,*)
' DEBUG2: recv_buf',recv_buf(:)
323 end subroutine send_recv_contact_info_r
325 subroutine send_recv_contact_info_i(n_neighbor_pe, neighbor_pe, MPI_COMM, &
326 send_index, send_item, recv_index, recv_item, vec, op)
328 integer(kind=kint),
intent(in) :: n_neighbor_pe
329 integer(kind=kint),
intent(in) :: neighbor_pe(:)
330 integer(kind=kint),
intent(in) :: mpi_comm
331 integer(kind=kint),
pointer,
intent(in) :: send_index(:), send_item(:), recv_index(:), recv_item(:)
332 integer(kind=kint),
intent(inout) :: vec(:)
333 integer(kind=kint),
intent(in) :: op
334 integer(kind=kint),
allocatable :: send_buf(:), recv_buf(:)
335 integer(kind=kint) :: i, n_send, idom, irank, is, ie, len, tag
336 integer(kind=kint),
allocatable :: requests(:), statuses(:,:)
337 if (n_neighbor_pe == 0)
return
338 allocate(requests(n_neighbor_pe))
340 allocate(send_buf(send_index(n_neighbor_pe)))
341 allocate(recv_buf(recv_index(n_neighbor_pe)))
342 do i = 1, send_index(n_neighbor_pe)
343 send_buf(i) = vec(send_item(i))
346 do idom = 1, n_neighbor_pe
347 irank = neighbor_pe(idom)
348 is = send_index(idom-1)+1
349 ie = send_index(idom)
355 mpi_comm, requests(n_send))
357 do idom = 1, n_neighbor_pe
358 irank = neighbor_pe(idom)
359 is = recv_index(idom-1)+1
360 ie = recv_index(idom)
365 mpi_comm, statuses(:,1))
369 do i = 1, recv_index(n_neighbor_pe)
370 vec(recv_item(i)) = vec(recv_item(i)) + recv_buf(i)
373 do i = 1, recv_index(n_neighbor_pe)
374 vec(recv_item(i)) = vec(recv_item(i)) * recv_buf(i)
377 do i = 1, recv_index(n_neighbor_pe)
378 vec(recv_item(i)) = max(vec(recv_item(i)), recv_buf(i))
381 do i = 1, recv_index(n_neighbor_pe)
382 vec(recv_item(i)) = min(vec(recv_item(i)), recv_buf(i))
385 do i = 1, recv_index(n_neighbor_pe)
386 vec(recv_item(i)) = recv_buf(i)
392 write(0,*)
' DEBUG2: send_buf',send_buf(:)
393 write(0,*)
' DEBUG2: recv_buf',recv_buf(:)
397 end subroutine send_recv_contact_info_i