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

Go to the source code of this file.

Functions/Subroutines

subroutine spmd_ds_isend (buf, size, itag, idest)
subroutine spmd_ds_rsend (buf, size, itag, idest)
subroutine spmd_ds_irecv (buf, size, itag, iprov)
subroutine spmd_ds_rrecv (buf, size, itag, iprov)
subroutine spmd_ds_iexch (sbuf, rbuf, size, msgoff)
subroutine spmd_ds_mexch (matr, n, iexch, rexch, msgoff, iads, iadr, nn)
subroutine spmd_ds_vexch (vect, n, iexch, rexch, msgoff, iads, iadr, nn, nv)
subroutine spmd_ds_vdesc (v, ndepl, vv, nddlc, nloc, msgoff)
subroutine spmd_iwlg (nddlt, nddlg, nddl, lsddl, iw, msgoff)

Function/Subroutine Documentation

◆ spmd_ds_iexch()

subroutine spmd_ds_iexch ( integer, dimension(size,*) sbuf,
integer, dimension(size,*) rbuf,
integer size,
integer msgoff )

Definition at line 183 of file spmd_dsreso.F.

184C-----------------------------------------------
185C I m p l i c i t T y p e s
186C-----------------------------------------------
187 USE spmd_comm_world_mod, ONLY : spmd_comm_world
188#include "implicit_f.inc"
189C-----------------------------------------------
190C M e s s a g e P a s s i n g
191C-----------------------------------------------
192#include "spmd.inc"
193C-----------------------------------------------
194C C o m m o n B l o c k s
195C-----------------------------------------------
196#include "com01_c.inc"
197#include "task_c.inc"
198C-----------------------------------------------
199C D u m m y A r g u m e n t s
200C-----------------------------------------------
201 INTEGER SIZE, SBUF(SIZE,*), RBUF(SIZE,*), MSGOFF
202C-----------------------------------------------
203C L o c a l V a r i a b l e s
204C-----------------------------------------------
205#ifdef MPI
206 INTEGER I, ITAG, REQ(2), IERR,
207 . TSTAT(MPI_STATUS_SIZE,2)
208C
209 DO i=1,nspmd
210 IF (ispmd==i-1) THEN
211 rbuf(1,i)=sbuf(1,i)
212 rbuf(2,i)=sbuf(2,i)
213 ELSE
214 itag=msgoff + nspmd*ispmd + i
215 CALL mpi_isend(sbuf(1,i), SIZE, mpi_integer, it_spmd(i),
216 . itag, spmd_comm_world, req(1), ierr)
217C
218 itag=msgoff + nspmd*(i-1) + ispmd+1
219 CALL mpi_irecv(rbuf(1,i), SIZE, mpi_integer, it_spmd(i),
220 . itag, spmd_comm_world, req(2), ierr)
221C
222 CALL mpi_waitall(2, req, tstat, ierr)
223 ENDIF
224 ENDDO
225C
226#endif
227 RETURN
subroutine mpi_isend(buf, cnt, datatype, dest, tag, comm, ireq, ierr)
Definition mpi.f:382
subroutine mpi_waitall(cnt, array_of_requests, status, ierr)
Definition mpi.f:536
subroutine mpi_irecv(buf, cnt, datatype, source, tag, comm, ireq, ierr)
Definition mpi.f:372

◆ spmd_ds_irecv()

subroutine spmd_ds_irecv ( integer, dimension(*) buf,
integer size,
integer itag,
integer iprov )

Definition at line 106 of file spmd_dsreso.F.

107C-----------------------------------------------
108C I m p l i c i t T y p e s
109C-----------------------------------------------
110 USE spmd_comm_world_mod, ONLY : spmd_comm_world
111#include "implicit_f.inc"
112C-----------------------------------------------
113C M e s s a g e P a s s i n g
114C-----------------------------------------------
115#include "spmd.inc"
116C-----------------------------------------------
117C C o m m o n B l o c k s
118C-----------------------------------------------
119#include "task_c.inc"
120C-----------------------------------------------
121C D u m m y A r g u m e n t s
122C-----------------------------------------------
123 INTEGER BUF(*), SIZE, ITAG, IPROV
124C-----------------------------------------------
125C L o c a l V a r i a b l e s
126C-----------------------------------------------
127#ifdef MPI
128 INTEGER IERR, ISTAT(MPI_STATUS_SIZE), LEN_STR, IERR_STR
129 CHARACTER STR_ERROR*(MPI_MAX_ERROR_STRING)
130C
131 CALL mpi_recv(buf, SIZE, mpi_integer, it_spmd(iprov), itag,
132 . spmd_comm_world, istat, ierr)
133C
134#endif
135 RETURN
subroutine mpi_recv(buf, cnt, datatype, source, tag, comm, status, ierr)
Definition mpi.f:461

◆ spmd_ds_isend()

subroutine spmd_ds_isend ( integer, dimension(*) buf,
integer size,
integer itag,
integer idest )

Definition at line 29 of file spmd_dsreso.F.

30C-----------------------------------------------
31C I m p l i c i t T y p e s
32C-----------------------------------------------
33 USE spmd_comm_world_mod, ONLY : spmd_comm_world
34#include "implicit_f.inc"
35C-----------------------------------------------
36C M e s s a g e P a s s i n g
37C-----------------------------------------------
38#include "spmd.inc"
39C-----------------------------------------------
40C C o m m o n B l o c k s
41C-----------------------------------------------
42#include "task_c.inc"
43C-----------------------------------------------
44C D u m m y A r g u m e n t s
45C-----------------------------------------------
46 INTEGER BUF(*), SIZE, ITAG, IDEST
47C-----------------------------------------------
48C L o c a l V a r i a b l e s
49C-----------------------------------------------
50#ifdef MPI
51 INTEGER IERR
52C
53 CALL mpi_send(buf, SIZE, mpi_integer, it_spmd(idest), itag,
54 . spmd_comm_world, ierr)
55* WRITE(*,*) 'Requete I envoyee - ITAG ',ITAG,' IERR ',IERR
56C
57#endif
58 RETURN
subroutine mpi_send(buf, cnt, datatype, dest, tag, comm, ierr)
Definition mpi.f:480

◆ spmd_ds_mexch()

subroutine spmd_ds_mexch ( matr,
integer n,
integer, dimension(nn,*) iexch,
rexch,
integer msgoff,
integer, dimension(*) iads,
integer, dimension(*) iadr,
integer nn )

Definition at line 235 of file spmd_dsreso.F.

237C-----------------------------------------------
238C I m p l i c i t T y p e s
239C-----------------------------------------------
240 USE spmd_comm_world_mod, ONLY : spmd_comm_world
241#include "implicit_f.inc"
242C-----------------------------------------------
243C M e s s a g e P a s s i n g
244C-----------------------------------------------
245#include "spmd.inc"
246C-----------------------------------------------
247C C o m m o n B l o c k s
248C-----------------------------------------------
249#include "com01_c.inc"
250#include "task_c.inc"
251C-----------------------------------------------
252C D u m m y A r g u m e n t s
253C-----------------------------------------------
254 INTEGER N, NN, IEXCH(NN,*), MSGOFF, IADS(*), IADR(*)
255 my_real
256 . matr(n,*), rexch(*)
257C-----------------------------------------------
258C L o c a l V a r i a b l e s
259C-----------------------------------------------
260#ifdef MPI
261 INTEGER I, IADC, LEN, ITAG(3), J, IR, IC, REQ(6),
262 . TSTAT(MPI_STATUS_SIZE,6), IERR, LENR
263 INTEGER, DIMENSION(:), ALLOCATABLE :: IROW, ICOL
264 my_real
265 . , DIMENSION(:), ALLOCATABLE :: val
266C
267 DO i=1,nspmd
268 IF (ispmd==i-1) THEN
269 iadc=iads(i)
270 len=iads(i+1)-iads(i)
271 DO j=1,len
272 ir=iexch(iadc+j-1,1)
273 ic=iexch(iadc+j-1,2)
274 matr(ir,ic)=matr(ir,ic)+rexch(iadc+j-1)
275 ENDDO
276 ELSE
277C Reception
278 len=iadr(i+1)-iadr(i)
279 lenr=len
280 ALLOCATE(irow(len), icol(len), val(len))
281 itag(1)=msgoff + nspmd*3*(i-1) + ispmd+1
282 itag(2)=msgoff + nspmd*3*(i-1) + nspmd+ispmd+1
283 itag(3)=msgoff + nspmd*3*(i-1) + 2*nspmd+ispmd+1
284 CALL mpi_irecv(irow, len, mpi_integer, it_spmd(i),
285 . itag(1), spmd_comm_world, req(1), ierr)
286 CALL mpi_irecv(icol, len, mpi_integer, it_spmd(i),
287 . itag(2), spmd_comm_world, req(2), ierr)
288 CALL mpi_irecv(val, len, real, it_spmd(i),
289 . itag(3), spmd_comm_world, req(3), ierr)
290C Envoi
291 iadc=iads(i)
292 len=iads(i+1)-iads(i)
293 itag(1)=msgoff + nspmd*3*ispmd + i
294 itag(2)=msgoff + nspmd*3*ispmd + nspmd+i
295 itag(3)=msgoff + nspmd*3*ispmd + 2*nspmd+i
296 CALL mpi_isend(iexch(iadc,1), len, mpi_integer, it_spmd(i),
297 . itag(1), spmd_comm_world, req(4), ierr)
298 CALL mpi_isend(iexch(iadc,2), len, mpi_integer, it_spmd(i),
299 . itag(2), spmd_comm_world, req(5), ierr)
300 CALL mpi_isend(rexch(iadc), len, real, it_spmd(i),
301 . itag(3), spmd_comm_world, req(6), ierr)
302C
303 CALL mpi_waitall(6, req, tstat, ierr)
304C
305 DO j=1,lenr
306 ir=irow(j)
307 ic=icol(j)
308 matr(ir,ic)=matr(ir,ic)+val(j)
309 ENDDO
310 DEALLOCATE(irow, icol, val)
311 ENDIF
312 ENDDO
313C
314#endif
315 RETURN
#define my_real
Definition cppsort.cpp:32

◆ spmd_ds_rrecv()

subroutine spmd_ds_rrecv ( buf,
integer size,
integer itag,
integer iprov )

Definition at line 145 of file spmd_dsreso.F.

146C-----------------------------------------------
147C I m p l i c i t T y p e s
148C-----------------------------------------------
149 USE spmd_comm_world_mod, ONLY : spmd_comm_world
150#include "implicit_f.inc"
151C-----------------------------------------------
152C M e s s a g e P a s s i n g
153C-----------------------------------------------
154#include "spmd.inc"
155C-----------------------------------------------
156C C o m m o n B l o c k s
157C-----------------------------------------------
158#include "task_c.inc"
159C-----------------------------------------------
160C D u m m y A r g u m e n t s
161C-----------------------------------------------
162 INTEGER SIZE, ITAG, IPROV
163 my_real
164 . buf(*)
165C-----------------------------------------------
166C L o c a l V a r i a b l e s
167C-----------------------------------------------
168#ifdef MPI
169 INTEGER IERR, ISTAT(MPI_STATUS_SIZE)
170C
171 CALL mpi_recv(buf, SIZE, real, it_spmd(iprov), itag,
172 . spmd_comm_world, istat, ierr)
173C
174#endif
175 RETURN

◆ spmd_ds_rsend()

subroutine spmd_ds_rsend ( buf,
integer size,
integer itag,
integer idest )

Definition at line 68 of file spmd_dsreso.F.

69C-----------------------------------------------
70C I m p l i c i t T y p e s
71C-----------------------------------------------
72 USE spmd_comm_world_mod, ONLY : spmd_comm_world
73#include "implicit_f.inc"
74C-----------------------------------------------
75C M e s s a g e P a s s i n g
76C-----------------------------------------------
77#include "spmd.inc"
78C-----------------------------------------------
79C C o m m o n B l o c k s
80C-----------------------------------------------
81#include "task_c.inc"
82C-----------------------------------------------
83C D u m m y A r g u m e n t s
84C-----------------------------------------------
85 INTEGER SIZE, ITAG, IDEST
87 . buf(*)
88C-----------------------------------------------
89C L o c a l V a r i a b l e s
90C-----------------------------------------------
91#ifdef MPI
92 INTEGER IERR
93C
94 CALL mpi_send(buf, SIZE, real, it_spmd(idest), itag,
95 . spmd_comm_world, ierr)
96C
97#endif
98 RETURN

◆ spmd_ds_vdesc()

subroutine spmd_ds_vdesc ( v,
integer ndepl,
vv,
integer nddlc,
integer nloc,
integer msgoff )

Definition at line 425 of file spmd_dsreso.F.

427C-----------------------------------------------
428C I m p l i c i t T y p e s
429C-----------------------------------------------
430 USE spmd_comm_world_mod, ONLY : spmd_comm_world
431#include "implicit_f.inc"
432C-----------------------------------------------
433C M e s s a g e P a s s i n g
434C-----------------------------------------------
435#include "spmd.inc"
436C-----------------------------------------------
437C C o m m o n B l o c k s
438C-----------------------------------------------
439#include "com01_c.inc"
440#include "task_c.inc"
441C-----------------------------------------------
442C D u m m y A r g u m e n t s
443C-----------------------------------------------
444 INTEGER NDEPL, NDDLC, NLOC, MSGOFF
445 my_real
446 . v(nloc,*), vv(nddlc,*)
447C-----------------------------------------------
448C L o c a l V a r i a b l e s
449C-----------------------------------------------
450#ifdef MPI
451 INTEGER ITAG, LEN, IERR, I, NPLOC, IPROV,
452 . ISTAT(MPI_STATUS_SIZE), J, KK, IAD1, K
453 my_real
454 . , DIMENSION(:,:), ALLOCATABLE :: vp
455C
456 IF (ispmd/=0.AND.mod(ispmd,dsncol)==0) THEN
457 itag=msgoff + ispmd
458 CALL mpi_send(nloc, 1, mpi_integer, it_spmd(1), itag,
459 . spmd_comm_world, ierr)
460C
461 itag=msgoff + nspmd + ispmd
462 len=nloc*ndepl
463 CALL mpi_send(v, len, real, it_spmd(1), itag,
464 . spmd_comm_world, ierr)
465 ENDIF
466C
467 IF (ispmd==0) THEN
468 DO i=1,dsnrow
469 IF (i==1) THEN
470 nploc=nloc
471 ALLOCATE(vp(nploc,ndepl))
472 DO j=1,ndepl
473 DO k=1,nploc
474 vp(k,j)=v(k,j)
475 ENDDO
476 ENDDO
477 ELSE
478 iprov=(i-1)*dsncol+1
479 itag=msgoff + iprov-1
480 CALL mpi_recv(nploc, 1, mpi_integer, it_spmd(iprov),
481 . itag, spmd_comm_world, istat, ierr)
482C
483 itag=msgoff + nspmd + iprov-1
484 ALLOCATE(vp(nploc,ndepl))
485 len=nploc*ndepl
486 CALL mpi_recv(vp, len, real, it_spmd(iprov), itag,
487 . spmd_comm_world, istat, ierr)
488 ENDIF
489C
490 DO j=1,ndepl
491 kk=0
492 iad1=(i-1)*dsnbloc
493 DO k=1,nploc
494 iad1=iad1+1
495 kk=kk+1
496 IF (kk>dsnbloc) THEN
497 iad1=iad1+dsnbloc*(dsnrow-1)
498 kk=1
499 ENDIF
500 vv(iad1,j)=vp(k,j)
501 ENDDO
502 ENDDO
503 DEALLOCATE(vp)
504 ENDDO
505C
506 DO i=2,nspmd
507 itag=msgoff + 2*nspmd + i-1
508 len=nddlc*ndepl
509 CALL mpi_send(vv, len, real, it_spmd(i), itag,
510 . spmd_comm_world, ierr)
511 ENDDO
512 ELSE
513 itag=msgoff + 2*nspmd + ispmd
514 len=nddlc*ndepl
515 CALL mpi_recv(vv, len, real, it_spmd(1), itag,
516 . spmd_comm_world, istat, ierr)
517 ENDIF
518C
519#endif
520 RETURN

◆ spmd_ds_vexch()

subroutine spmd_ds_vexch ( vect,
integer n,
integer, dimension(*) iexch,
rexch,
integer msgoff,
integer, dimension(*) iads,
integer, dimension(*) iadr,
integer nn,
integer nv )

Definition at line 323 of file spmd_dsreso.F.

325C-----------------------------------------------
326C I m p l i c i t T y p e s
327C-----------------------------------------------
328 USE spmd_comm_world_mod, ONLY : spmd_comm_world
329#include "implicit_f.inc"
330C-----------------------------------------------
331C M e s s a g e P a s s i n g
332C-----------------------------------------------
333#include "spmd.inc"
334C-----------------------------------------------
335C C o m m o n B l o c k s
336C-----------------------------------------------
337#include "com01_c.inc"
338#include "task_c.inc"
339C-----------------------------------------------
340C D u m m y A r g u m e n t s
341C-----------------------------------------------
342 INTEGER N, IEXCH(*), MSGOFF, IADS(*), IADR(*), NN, NV
343 my_real
344 . vect(n,*), rexch(nn,*)
345C-----------------------------------------------
346C L o c a l V a r i a b l e s
347C-----------------------------------------------
348#ifdef MPI
349 INTEGER I, IADC, LEN, ITAG, J, K, IR, REQ(4),
350 . TSTAT(MPI_STATUS_SIZE, 4), IERR, II, LEN2, LENR
351 INTEGER, DIMENSION(:), ALLOCATABLE :: IROW
352 my_real
353 . , DIMENSION(:,:), ALLOCATABLE :: val, vals
354C
355 DO i=1,nspmd
356 IF (ispmd==i-1) THEN
357 iadc=iads(i)
358 len=iads(i+1)-iads(i)
359 DO j=1,nv
360 DO k=1,len
361 ir=iexch(iadc+k-1)
362 vect(ir,j)=vect(ir,j)+rexch(iadc+k-1,j)
363 ENDDO
364 ENDDO
365 ELSE
366C Reception
367 len=iadr(i+1)-iadr(i)
368 lenr=len
369 ii=0
370 IF (len>0) THEN
371 ALLOCATE(irow(len), val(len,nv))
372 itag=msgoff + nspmd*2*(i-1) + ispmd+1
373 ii=ii+1
374 CALL mpi_irecv(irow, len, mpi_integer, it_spmd(i),
375 . itag, spmd_comm_world, req(ii), ierr)
376 itag=msgoff + nspmd*2*(i-1) + nspmd + ispmd+1
377 len2=len*nv
378 ii=ii+1
379 CALL mpi_irecv(val, len2, real, it_spmd(i),
380 . itag, spmd_comm_world, req(ii), ierr)
381 ENDIF
382C Envoi
383 iadc=iads(i)
384 len=iads(i+1)-iads(i)
385 IF (len>0) THEN
386 itag=msgoff + nspmd*2*ispmd + i
387 ii=ii+1
388 CALL mpi_isend(iexch(iadc), len, mpi_integer, it_spmd(i),
389 . itag, spmd_comm_world, req(ii), ierr)
390 ALLOCATE(vals(len,nv))
391 DO j=1,nv
392 DO k=1,len
393 vals(k,j)=rexch(iadc+k-1,j)
394 ENDDO
395 ENDDO
396 itag=msgoff + nspmd*2*ispmd + nspmd + i
397 len2=len*nv
398 ii=ii+1
399 CALL mpi_isend(vals, len2, real, it_spmd(i),
400 . itag, spmd_comm_world, req(ii), ierr)
401 ENDIF
402C
403 IF (ii>0) CALL mpi_waitall(ii, req, tstat, ierr)
404C
405 DO j=1,nv
406 DO k=1,lenr
407 ir=irow(k)
408 vect(ir,j)=vect(ir,j)+val(k,j)
409 ENDDO
410 ENDDO
411 IF (len>0) DEALLOCATE(vals)
412 IF (lenr>0) DEALLOCATE(irow, val)
413 ENDIF
414 ENDDO
415C
416#endif
417 RETURN

◆ spmd_iwlg()

subroutine spmd_iwlg ( integer nddlt,
integer nddlg,
integer nddl,
integer, dimension(*) lsddl,
integer, dimension(*) iw,
integer msgoff )

Definition at line 528 of file spmd_dsreso.F.

530C-----------------------------------------------
531C I m p l i c i t T y p e s
532C-----------------------------------------------
533 USE spmd_comm_world_mod, ONLY : spmd_comm_world
534#include "implicit_f.inc"
535C-----------------------------------------------
536C M e s s a g e P a s s i n g
537C-----------------------------------------------
538#include "spmd.inc"
539C-----------------------------------------------
540C C o m m o n B l o c k s
541C-----------------------------------------------
542#include "com01_c.inc"
543#include "task_c.inc"
544C-----------------------------------------------
545C D u m m y A r g u m e n t s
546C-----------------------------------------------
547 INTEGER NDDLT, NDDLG, NDDL, LSDDL(*), IW(*), MSGOFF
548C-----------------------------------------------
549C L o c a l V a r i a b l e s
550C-----------------------------------------------
551#ifdef MPI
552 INTEGER I, ITAG(2,NDDLT), II, NDDLPM, IRQTAG, NDDLP(NSPMD-1),
553 . ISTAT(MPI_STATUS_SIZE), IERR, J, JJ, N
554 INTEGER, DIMENSION(:,:), ALLOCATABLE :: LDDLP
555C
556 IF (ispmd==0) THEN
557 DO i=1,nddlt
558 itag(1,i)=0
559 ENDDO
560C
561 DO i=1,nddl
562 ii=lsddl(i)
563 itag(1,ii)=1
564 ENDDO
565C
566 nddlpm=0
567 DO i=1,nspmd-1
568 irqtag=msgoff + i
569 CALL mpi_recv(nddlp(i), 1, mpi_integer, it_spmd(i+1),
570 . irqtag, spmd_comm_world, istat, ierr)
571 nddlpm=max(nddlpm,nddlp(i))
572 ENDDO
573 ALLOCATE(lddlp(nddlpm,nspmd-1))
574 DO i=1,nspmd-1
575 irqtag=msgoff + nspmd + i
576 CALL mpi_recv(
577 . lddlp(1,i), nddlp(i), mpi_integer, it_spmd(i+1),
578 . irqtag, spmd_comm_world, istat, ierr)
579 DO j=1,nddlp(i)
580 jj=lddlp(j,i)
581 itag(1,jj)=1
582 ENDDO
583 ENDDO
584C
585 n=0
586 DO i=1,nddlt
587 IF (itag(1,i)==1) THEN
588 n=n+1
589 itag(2,i)=n
590 ENDIF
591 ENDDO
592 nddlg=n
593C
594 DO i=1,nddl
595 ii=lsddl(i)
596 iw(i)=itag(2,ii)
597 ENDDO
598C
599 DO i=1,nspmd-1
600 DO j=1,nddlp(i)
601 jj=lddlp(j,i)
602 lddlp(j,i)=itag(2,jj)
603 ENDDO
604 irqtag=msgoff + 2*nspmd + i
605 CALL mpi_send(nddlg, 1, mpi_integer, it_spmd(i+1),
606 . irqtag, spmd_comm_world, ierr)
607 irqtag=msgoff + 3*nspmd + i
608 CALL mpi_send(
609 . lddlp(1,i), nddlp(i), mpi_integer, it_spmd(i+1),
610 . irqtag, spmd_comm_world, ierr)
611 ENDDO
612 DEALLOCATE(lddlp)
613 ELSE
614 irqtag=msgoff + ispmd
615 CALL mpi_send(nddl, 1, mpi_integer, it_spmd(1),
616 . irqtag, spmd_comm_world, ierr)
617 irqtag=msgoff + nspmd + ispmd
618 CALL mpi_send(lsddl, nddl, mpi_integer, it_spmd(1),
619 . irqtag, spmd_comm_world, ierr)
620 irqtag=msgoff + 2*nspmd + ispmd
621 CALL mpi_recv(nddlg, 1, mpi_integer, it_spmd(1),
622 . irqtag, spmd_comm_world, istat, ierr)
623 irqtag=msgoff + 3*nspmd + ispmd
624 CALL mpi_recv(iw, nddl, mpi_integer, it_spmd(1),
625 . irqtag, spmd_comm_world, istat, ierr)
626 ENDIF
627C
628#endif
629 RETURN
#define max(a, b)
Definition macros.h:21