OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
spmd_exch_sub.F File Reference
#include "implicit_f.inc"
#include "com01_c.inc"
#include "task_c.inc"
#include "parit_c.inc"
#include "scr02_c.inc"
#include "spmd.inc"

Go to the source code of this file.

Functions/Subroutines

subroutine spmd_sub_boundaries (nloc_dmg, iad_elem, fr_elem)
subroutine spmd_exch_sub_poff (nloc_dmg)
subroutine spmd_exch_sub_pon (nloc_dmg)

Function/Subroutine Documentation

◆ spmd_exch_sub_poff()

subroutine spmd_exch_sub_poff ( type (nlocal_str_), target nloc_dmg)

Definition at line 230 of file spmd_exch_sub.F.

231C-----------------------------------------------
232C M o d u l e s
233C-----------------------------------------------
235C-----------------------------------------------
236C I m p l i c i t T y p e s
237C-----------------------------------------------
238 USE spmd_comm_world_mod, ONLY : spmd_comm_world
239#include "implicit_f.inc"
240C-----------------------------------------------
241C M e s s a g e P a s s i n g
242C-----------------------------------------------
243#include "spmd.inc"
244C-----------------------------------------------
245C C o m m o n B l o c k s
246C-----------------------------------------------
247#include "com01_c.inc"
248#include "task_c.inc"
249#include "scr02_c.inc"
250C-----------------------------------------------
251C D u m m y A r g u m e n t s
252C-----------------------------------------------
253 TYPE (NLOCAL_STR_) ,TARGET :: NLOC_DMG
254C-----------------------------------------------
255C L o c a l V a r i a b l e s
256C-----------------------------------------------
257#ifdef MPI
258 INTEGER I,J,K,L,M,ND
259 INTEGER NDDL,NSIZESUB,NN,NPOS,SIZEA
260
261C MPI VARIABLES
262
263 INTEGER MSGTYP,NOD,LOC_PROC, SIZ,NB_NOD,NB,MAXLEV,VALUE
264
265 INTEGER STATUS(MPI_STATUS_SIZE),IERROR
266
267 INTEGER MSGOFF1
268 INTEGER IAD_SEND(NSPMD+1),IAD_RECV(NSPMD+1)
269 INTEGER SEND_SIZ(NSPMD),RECV_SIZ(NSPMD)
270 INTEGER SIZ_SEND,SIZ_RECV
271
272 INTEGER REQ_S(NSPMD)
273 INTEGER REQ_R(NSPMD)
274
275 my_real, DIMENSION(:), ALLOCATABLE :: send_buf,recv_buf
276
277 DATA msgoff1/1282/
278C=======================================================================
279 nsizesub=nloc_dmg%NNOD
280 sizea = nloc_dmg%L_NLOC
281C-------------------------------
282C 1. SEND & RECEIVE Buffer size
283C-------------------------------
284 siz = nloc_dmg%IAD_SIZE(nspmd+1)-nloc_dmg%IAD_SIZE(1)
285 ALLOCATE(recv_buf(siz))
286 ALLOCATE(send_buf(siz))
287c--------------------------------------------------------------------
288c MPI receive
289c--------------------------------------------------------------------
290 l = 1
291 iad_recv(1) = 1
292 DO i=1,nspmd
293 siz = nloc_dmg%IAD_SIZE(i+1)-nloc_dmg%IAD_SIZE(i)
294 IF (siz > 0)THEN
295 msgtyp = msgoff1
296 CALL mpi_irecv( recv_buf(l),siz,real,it_spmd(i),msgtyp,
297 g spmd_comm_world,req_r(i),ierror)
298 l = l + siz
299 ENDIF
300 iad_recv(i+1) = l
301 END DO
302C --------------------------------------------------------------------
303C Prepare send buffers
304C --------------------------------------------------------------------
305 k=1
306 iad_send(1)=k
307 DO i=1,nspmd
308 DO j=nloc_dmg%IAD_ELEM(i),nloc_dmg%IAD_ELEM(i+1)-1
309 nn = nloc_dmg%FR_ELEM(j)
310 npos = nloc_dmg%POSI(nn)
311 nddl = nloc_dmg%POSI(nn+1)-nloc_dmg%POSI(nn)
312
313 DO l=1,nddl
314 send_buf(k) = nloc_dmg%FNL(npos+l-1,1)
315 k = k + 1
316 IF (nodadt > 0) THEN
317 send_buf(k) = nloc_dmg%STIFNL(npos+l-1,1)
318 k = k + 1
319 ENDIF
320 ENDDO
321 ENDDO
322 iad_send(i+1)=k
323 ENDDO
324c--------------------------------------------------------------------
325c MPI send
326c--------------------------------------------------------------------
327 DO i=1,nspmd
328 siz = iad_send(i+1)-iad_send(i)
329 IF (siz > 0)THEN
330 l = iad_send(i)
331 msgtyp = msgoff1
332 CALL mpi_isend(
333 s send_buf(l),siz,real,it_spmd(i),msgtyp,
334 g spmd_comm_world,req_s(i),ierror)
335 ENDIF
336 ENDDO
337c--------------------------------------------------------------------
338c Receive buffer & compute receive_size
339c-------------------------------------------------------------------
340 k=1
341 DO i = 1, nspmd
342 siz = nloc_dmg%IAD_SIZE(i+1)-nloc_dmg%IAD_SIZE(i)
343 IF (siz > 0)THEN
344 CALL mpi_wait(req_r(i),status,ierror)
345
346 k = iad_recv(i)
347
348 DO j=nloc_dmg%IAD_ELEM(i),nloc_dmg%IAD_ELEM(i+1)-1
349 nn = nloc_dmg%FR_ELEM(j)
350 npos = nloc_dmg%POSI(nn)
351 nddl = nloc_dmg%POSI(nn+1) - nloc_dmg%POSI(nn)
352
353 DO l=1,nddl
354 nloc_dmg%FNL(npos+l-1,1)=nloc_dmg%FNL(npos+l-1,1)+recv_buf(k)
355 k=k+1
356 IF (nodadt > 0) THEN
357 nloc_dmg%STIFNL(npos+l-1,1)=nloc_dmg%STIFNL(npos+l-1,1)+recv_buf(k)
358 k = k + 1
359 ENDIF
360 END DO
361
362 END DO
363 ENDIF
364 END DO
365
366
367c--------------------------------------------------------------------
368c Terminate send
369c-------------------------------------------------------------------
370 DO i = 1, nspmd
371 siz = iad_send(i+1)-iad_send(i)
372 IF (siz > 0)THEN
373 CALL mpi_wait(req_s(i),status,ierror)
374 ENDIF
375 ENDDO
376
377 DEALLOCATE(send_buf)
378 DEALLOCATE(recv_buf)
379#endif
380
381 RETURN
#define my_real
Definition cppsort.cpp:32
subroutine mpi_isend(buf, cnt, datatype, dest, tag, comm, ireq, ierr)
Definition mpi.f:382
subroutine mpi_wait(ireq, status, ierr)
Definition mpi.f:525
subroutine mpi_irecv(buf, cnt, datatype, source, tag, comm, ireq, ierr)
Definition mpi.f:372

◆ spmd_exch_sub_pon()

subroutine spmd_exch_sub_pon ( type (nlocal_str_), target nloc_dmg)

Definition at line 394 of file spmd_exch_sub.F.

395!$COMMENT
396! SPMD_EXCH_SUB_PON
397! communication of neighbouring values for
398! NLOCAL option and parith/ON
399!
400! SPMD_EXCH_SUB_PON organization :
401! - prepare the non-blocking received
402! - fill the sending buffer
403! - send the sending buffer
404!$ENDCOMMENT
405C-----------------------------------------------
406C M o d u l e s
407C-----------------------------------------------
409C-----------------------------------------------
410C I m p l i c i t T y p e s
411C-----------------------------------------------
412 USE spmd_comm_world_mod, ONLY : spmd_comm_world
413#include "implicit_f.inc"
414C-----------------------------------------------------------------
415C M e s s a g e P a s s i n g
416C-----------------------------------------------
417#include "spmd.inc"
418C-----------------------------------------------
419C C o m m o n B l o c k s
420C-----------------------------------------------
421#include "com01_c.inc"
422#include "task_c.inc"
423#include "scr02_c.inc"
424C-----------------------------------------------
425C D u m m y A r g u m e n t s
426C-----------------------------------------------
427 TYPE (NLOCAL_STR_) ,TARGET :: NLOC_DMG
428C-----------------------------------------------
429C L o c a l V a r i a b l e s
430C-----------------------------------------------
431#ifdef MPI
432 INTEGER I,J,K,L,M,ND,IJK
433 INTEGER NDDL,NSIZESUB,NN,NPOS,SIZEA
434
435C MPI VARIABLES
436
437 INTEGER MSGTYP,NOD,LOC_PROC, SIZ,NB_NOD,NB,MAXLEV,VALUE
438
439 INTEGER STATUS(MPI_STATUS_SIZE),IERROR
440
441 INTEGER MSGOFF1
442 INTEGER IAD_SEND(NSPMD+1),IAD_RECV(NSPMD+1)
443 INTEGER SEND_SIZ(NSPMD),RECV_SIZ(NSPMD)
444 INTEGER SIZ_SEND,SIZ_RECV
445
446 INTEGER :: SIZ_S,SIZ_R
447 INTEGER REQ_S(NSPMD)
448 INTEGER REQ_R(NSPMD)
449
450 my_real, DIMENSION(:), ALLOCATABLE :: send_buf,recv_buf
451
452 DATA msgoff1/1283/
453C=======================================================================
454 nsizesub=nloc_dmg%NNOD
455 sizea = nloc_dmg%L_NLOC
456C-------------------------------
457C 1. SEND & RECEIVE Buffer size
458C-------------------------------
459 siz_s = nloc_dmg%FR_NBCC(1,nspmd+1)
460 siz_r = nloc_dmg%FR_NBCC(2,nspmd+1)
461 ALLOCATE(recv_buf(siz_r))
462 ALLOCATE(send_buf(siz_s))
463
464 recv_buf(1:siz_r) = -123456.
465 send_buf(1:siz_s) = 123456.
466c--------------------------------------------------------------------
467c MPI receive
468c--------------------------------------------------------------------
469 l = 1
470 iad_recv(1) = 1
471 DO i=1,nspmd
472 siz = nloc_dmg%FR_NBCC(2,i)
473 IF (siz > 0)THEN
474 msgtyp = msgoff1
475 CALL mpi_irecv( recv_buf(l),siz,real,it_spmd(i),msgtyp,
476 g spmd_comm_world,req_r(i),ierror)
477 l = l + siz
478 ENDIF
479 iad_recv(i+1) = l
480 END DO
481
482C --------------------------------------------------------------------
483C Prepare send buffers
484C --------------------------------------------------------------------
485 k=1
486 iad_send(1)=k
487 DO i=1,nspmd
488 DO j = nloc_dmg%IADSDP(i), nloc_dmg%IADSDP(i+1)-1
489 nn = nloc_dmg%FR_ELEM_S(j)
490 npos = nloc_dmg%POSI(nn)
491 nddl = nloc_dmg%POSI(nn+1)-nloc_dmg%POSI(nn)
492 ijk = nloc_dmg%ISENDSP(j)
493 DO l=1,nddl
494 send_buf(k) = nloc_dmg%FSKY(ijk,l)
495 k = k + 1
496 IF (nodadt > 0) THEN
497 send_buf(k) = nloc_dmg%STSKY(ijk,l)
498 k = k + 1
499 ENDIF
500 ENDDO
501
502 ENDDO
503 iad_send(i+1)=k
504 ENDDO
505c--------------------------------------------------------------------
506c MPI send
507c--------------------------------------------------------------------
508 DO i=1,nspmd
509 siz = nloc_dmg%FR_NBCC(1,i)
510 IF (siz > 0)THEN
511 l = iad_send(i)
512 msgtyp = msgoff1
513 CALL mpi_isend(
514 s send_buf(l),siz,real,it_spmd(i),msgtyp,
515 g spmd_comm_world,req_s(i),ierror)
516 ENDIF
517 ENDDO
518c--------------------------------------------------------------------
519c Receive buffer & compute receive_size
520c-------------------------------------------------------------------
521 k=1
522 DO i = 1, nspmd
523 siz = nloc_dmg%FR_NBCC(2,i)
524 IF (siz > 0)THEN
525 CALL mpi_wait(req_r(i),status,ierror)
526
527 k = iad_recv(i)
528 DO j=nloc_dmg%IADRCP(i),nloc_dmg%IADRCP(i+1)-1
529 nn = nloc_dmg%FR_ELEM_R(j)
530 npos = nloc_dmg%POSI(nn)
531 nddl = nloc_dmg%POSI(nn+1) - nloc_dmg%POSI(nn)
532
533 ijk = nloc_dmg%IRECSP(j)
534 DO l=1,nddl
535 nloc_dmg%FSKY(ijk,l) = recv_buf(k)
536 k = k + 1
537 IF (nodadt > 0) THEN
538 nloc_dmg%STSKY(ijk,l) = recv_buf(k)
539 k = k + 1
540 ENDIF
541 ENDDO
542 END DO
543 ENDIF
544 END DO
545
546c--------------------------------------------------------------------
547c Terminate send
548c-------------------------------------------------------------------
549 DO i = 1, nspmd
550 siz = nloc_dmg%FR_NBCC(1,i)
551 IF (siz > 0)THEN
552 CALL mpi_wait(req_s(i),status,ierror)
553 ENDIF
554 ENDDO
555
556 DEALLOCATE(send_buf)
557 DEALLOCATE(recv_buf)
558#endif
559
560 RETURN

◆ spmd_sub_boundaries()

subroutine spmd_sub_boundaries ( type (nlocal_str_), target nloc_dmg,
integer, dimension(2,*) iad_elem,
integer, dimension(*) fr_elem )

Definition at line 31 of file spmd_exch_sub.F.

32C-----------------------------------------------
33C M o d u l e s
34C-----------------------------------------------
36C-----------------------------------------------
37C I m p l i c i t T y p e s
38C-----------------------------------------------
39 USE spmd_comm_world_mod, ONLY : spmd_comm_world
40#include "implicit_f.inc"
41C-----------------------------------------------
42C C o m m o n B l o c k s
43C-----------------------------------------------
44#include "com01_c.inc"
45#include "task_c.inc"
46#include "parit_c.inc"
47#include "scr02_c.inc"
48C-----------------------------------------------
49C D u m m y A r g u m e n t s
50C-----------------------------------------------
51 INTEGER IAD_ELEM(2,*), FR_ELEM(*)
52 TYPE (NLOCAL_STR_) ,TARGET :: NLOC_DMG
53C-----------------------------------------------
54C L o c a l V a r i a b l e s
55C-----------------------------------------------
56 INTEGER, DIMENSION(:), ALLOCATABLE :: COUNT,COUNT2
57 INTEGER P,J,TOTAL_NODES,NCOUNT,NOD,NN,NPOS,NDDL
58
59 INTEGER :: LSD,LRD,CC,LOC_PROC,KK
60C-----------------------------------------------
61C
62
63C ---------------------------------
64C SPMD Non Local nodes boundaries
65C ---------------------------------
66
67
68 ! If already allocated, deallocate non-local IAD_ELEM, IAD_SIZE and FR_ELEM tables
69 IF (ALLOCATED(nloc_dmg%IAD_ELEM)) DEALLOCATE (nloc_dmg%IAD_ELEM)
70 IF (ALLOCATED(nloc_dmg%IAD_SIZE)) DEALLOCATE (nloc_dmg%IAD_SIZE)
71 IF (ALLOCATED(nloc_dmg%FR_ELEM )) DEALLOCATE (nloc_dmg%FR_ELEM)
72C
73 ! Allocation of the counter tables
74 ALLOCATE(count(nspmd)) ! For nodes
75 ALLOCATE(count2(nspmd)) ! For d.o.fs
76c
77 ! Initialization of the counter tables and the total boundaries non-local nodes
78 count(1:nspmd) = 0
79 count2(1:nspmd) = 0
80 total_nodes = 0
81c
82 kk = 1
83 IF (nodadt > 0) kk = 2
84c
85 ! Loop over domains
86 DO p=1,nspmd
87 DO j=iad_elem(1,p),iad_elem(1,p+1)-1 ! Loop over nodes at the boundaries
88 nod = fr_elem(j) ! Global number of the nodes
89 nn = nloc_dmg%IDXI(nod) ! Corresponding non-local node number
90 IF (nloc_dmg%IDXI(nod) > 0)THEN ! If the node belongs to the subset
91 npos = nloc_dmg%POSI(nn) ! Position of its first d.o.f in A
92 nddl = nloc_dmg%POSI(nn+1) - npos ! Number of additional d.o.fs
93 count(p) = count(p) +1 ! Updating the counter tables
94 count2(p)= count2(p)+kk*nddl
95 ENDIF
96 ENDDO
97 total_nodes = total_nodes + count(p) ! Total number of non-local d.o.fs at the boundaries
98 ENDDO
99c
100 ! Allocation of the non-local IAD_ELEM and FR_ELEM tables
101 ALLOCATE(nloc_dmg%IAD_ELEM(nspmd+1))
102 ALLOCATE(nloc_dmg%IAD_SIZE(nspmd+1))
103 ALLOCATE(nloc_dmg%FR_ELEM(total_nodes))
104c
105 ! Filling the non-local IAD_ELEM table
106 nloc_dmg%IAD_ELEM(1) = 1
107 nloc_dmg%IAD_SIZE(1) = 1
108c
109 ! Filling the IAD_ELEM/IAD_SIZE tables of the subset
110 DO p=2,nspmd+1
111 nloc_dmg%IAD_ELEM(p) = nloc_dmg%IAD_ELEM(p-1) + count(p-1)
112 nloc_dmg%IAD_SIZE(p) = nloc_dmg%IAD_SIZE(p-1) + count2(p-1)
113 ENDDO
114c
115 ! Filling the non-local FR_ELEM table
116 ncount=0
117 ! Loop over domains
118 DO p=1,nspmd
119 DO j=iad_elem(1,p),iad_elem(1,p+1)-1 ! Loop nodes at the boundaries
120 nod = fr_elem(j) ! Global number of the nodes
121 nn = nloc_dmg%IDXI(nod) ! Corresponding non-local node number
122 IF (nn > 0)THEN ! If the nodes belongs to the subset
123 ncount = ncount+1 ! Updating the counter
124 nloc_dmg%FR_ELEM(ncount) = nn ! Filling the FR_ELEM table of the subset
125 ENDIF
126 ENDDO
127 ENDDO
128
129 !WRITE(*,*) 'IAD_ELEM = ', NLOC_DMG%IAD_ELEM
130 !WRITE(*,*) 'IAD_SIZE = ', NLOC_DMG%IAD_SIZE
131 !WRITE(*,*) 'FR_ELEM = ', NLOC_DMG%FR_ELEM
132
133 IF (ALLOCATED(count)) DEALLOCATE(count)
134 IF (ALLOCATED(count2)) DEALLOCATE(count2)
135
136C -----------------------------
137C SPMD PON SKYLINE boundaries
138C -----------------------------
139 IF(iparit==1)THEN
140! ------------------------
141! compute the number of send / rcv
142 lsd = 0
143 lrd = 0
144 loc_proc = ispmd + 1
145 DO p=1,nspmd
146 DO j=nloc_dmg%IAD_ELEM(p),nloc_dmg%IAD_ELEM(p+1)-1
147 nn = nloc_dmg%FR_ELEM(j)
148 nddl = nloc_dmg%POSI(nn+1) - nloc_dmg%POSI(nn)
149 DO cc = nloc_dmg%ADDCNE(nn),nloc_dmg%ADDCNE(nn+1)-1
150 IF( nloc_dmg%PROCNE(cc)==loc_proc)THEN
151 lsd=lsd+1
152 ELSEIF(nloc_dmg%PROCNE(cc)==p) THEN
153 lrd=lrd+1
154 ENDIF
155
156 ENDDO
157 ENDDO
158 ENDDO
159! ------------------------
160! allocation of ISENDSP and IRECSP
161 ALLOCATE( nloc_dmg%ISENDSP(lsd) )
162 ALLOCATE( nloc_dmg%IRECSP(lrd) )
163 nloc_dmg%ISENDSP(1:lsd) = 0
164 nloc_dmg%IRECSP(1:lrd) = 0
165
166! ------------------------
167! fill ISENDSP and IRECSP
168 IF( .NOT.ALLOCATED( nloc_dmg%IADSDP ) ) ALLOCATE( nloc_dmg%IADSDP(nspmd+1) )
169 IF( .NOT.ALLOCATED( nloc_dmg%IADRCP ) ) ALLOCATE( nloc_dmg%IADRCP(nspmd+1) )
170 IF( .NOT.ALLOCATED( nloc_dmg%FR_NBCC ) ) ALLOCATE( nloc_dmg%FR_NBCC(2,nspmd+1) )
171
172 IF( .NOT.ALLOCATED( nloc_dmg%FR_ELEM_S ) ) ALLOCATE( nloc_dmg%FR_ELEM_S(lsd) )
173 IF( .NOT.ALLOCATED( nloc_dmg%FR_ELEM_R ) ) ALLOCATE( nloc_dmg%FR_ELEM_R(lrd) )
174 nloc_dmg%FR_ELEM_S(1:lsd) = 0
175 nloc_dmg%FR_ELEM_R(1:lrd) = 0
176
177 nloc_dmg%IADSDP(1:nspmd+1) = 0
178 nloc_dmg%IADRCP(1:nspmd+1) = 0
179 nloc_dmg%FR_NBCC(1:2,1:nspmd+1) = 0
180
181 lsd = 1
182 lrd = 1
183 loc_proc = ispmd + 1
184 DO p=1,nspmd
185 nloc_dmg%IADSDP(p)=lsd
186 nloc_dmg%IADRCP(p)=lrd
187
188 DO j=nloc_dmg%IAD_ELEM(p),nloc_dmg%IAD_ELEM(p+1)-1
189 nn = nloc_dmg%FR_ELEM(j)
190 nddl = nloc_dmg%POSI(nn+1) - nloc_dmg%POSI(nn)
191
192 DO cc = nloc_dmg%ADDCNE(nn),nloc_dmg%ADDCNE(nn+1)-1
193 IF( nloc_dmg%PROCNE(cc)==loc_proc)THEN
194 nloc_dmg%FR_NBCC(1,p) = nloc_dmg%FR_NBCC(1,p)+kk*nddl
195 nloc_dmg%ISENDSP(lsd)=cc
196 nloc_dmg%FR_ELEM_S(lsd)=nn
197 lsd=lsd+1
198 ELSEIF(nloc_dmg%PROCNE(cc)==p) THEN
199 nloc_dmg%FR_NBCC(2,p) = nloc_dmg%FR_NBCC(2,p)+kk*nddl
200 nloc_dmg%IRECSP(lrd)=cc
201 nloc_dmg%FR_ELEM_R(lrd)=nn
202 lrd=lrd+1
203 ENDIF
204
205 ENDDO
206 ENDDO
207 ENDDO
208! ------------------------
209 nloc_dmg%FR_NBCC(1:2,nspmd+1) = 0
210 DO j=1,nspmd
211 nloc_dmg%FR_NBCC(1,nspmd+1) = nloc_dmg%FR_NBCC(1,nspmd+1) + nloc_dmg%FR_NBCC(1,j)
212 nloc_dmg%FR_NBCC(2,nspmd+1) = nloc_dmg%FR_NBCC(2,nspmd+1) + nloc_dmg%FR_NBCC(2,j)
213 ENDDO
214
215 nloc_dmg%IADSDP(nspmd+1)=lsd
216 nloc_dmg%IADRCP(nspmd+1)=lrd
217 ENDIF
end diagonal values have been computed in the(sparse) matrix id.SOL