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

Go to the source code of this file.

Functions/Subroutines

subroutine spmd_x_section (nstrf, x, ms, weight, wa)
subroutine spmd_exch_sec (nstrf, x, ms, weight, xsec, fr_sec, iad_sec, lsend1, lrecv1, lsend2, lrecv2, weight_md)
subroutine spmd_wrt_cutd (nnod, nstrf, d, dr, rg_cut, iad_cut, nsize, nnodg, weight, iflg)
subroutine spmd_wrt_cutf (nnod, nstrf, secfcum, rg_cut, iad_cut, nsize, nnodg, weight, iflg)
subroutine spmd_exch_cut (nstrf, secfcum, iad_elem, fr_elem, size, lenr, nnod, weight)
subroutine spmd_sd_cut (secbufg, nnodg, secbuf1, secbuf2, nnod, fr_cut, iad_cut, iflg)

Function/Subroutine Documentation

◆ spmd_exch_cut()

subroutine spmd_exch_cut ( integer, dimension(*) nstrf,
secfcum,
integer, dimension(2,*) iad_elem,
integer, dimension(*) fr_elem,
integer size,
integer lenr,
integer nnod,
integer, dimension(*) weight )

Definition at line 1226 of file spmd_section.F.

1228C cumul de secfcum au noeuds frontieres puis mise a 0 sur procs additionnels
1229C-----------------------------------------------
1230C I m p l i c i t T y p e s
1231C-----------------------------------------------
1232 USE spmd_comm_world_mod, ONLY : spmd_comm_world
1233#include "implicit_f.inc"
1234C-----------------------------------------------------------------
1235C M e s s a g e P a s s i n g
1236C-----------------------------------------------
1237#include "spmd.inc"
1238C-----------------------------------------------
1239C C o m m o n B l o c k s
1240C-----------------------------------------------
1241#include "com01_c.inc"
1242#include "task_c.inc"
1243C-----------------------------------------------
1244C D u m m y A r g u m e n t s
1245C-----------------------------------------------
1246 INTEGER NSTRF(*), IAD_ELEM(2,*), FR_ELEM(*), WEIGHT(*),
1247 . SIZE, LENR, NNOD
1248 my_real
1249 . secfcum(7,*)
1250C-----------------------------------------------
1251C L o c a l V a r i a b l e s
1252C-----------------------------------------------
1253#ifdef MPI
1254 INTEGER MSGTYP,I,NOD,LOC_PROC,IERROR,
1255 . SIZ,J,L,NB_NOD,
1256 . STATUS(MPI_STATUS_SIZE),
1257 . IAD_SEND(NSPMD+1),IAD_RECV(NSPMD+1),
1258 . REQ_R(NSPMD),REQ_S(NSPMD),MSGOFF
1259 my_real
1260 . rbuf(size*lenr), sbuf(size*lenr)
1261 DATA msgoff/4005/
1262C-----------------------------------------------
1263C S o u r c e L i n e s
1264C-----------------------------------------------
1265 loc_proc = ispmd + 1
1266 l = 1
1267 iad_recv(1) = 1
1268 DO i=1,nspmd
1269 siz = size*(iad_elem(1,i+1)-iad_elem(1,i))
1270 IF(siz/=0)THEN
1271 msgtyp = msgoff
1272 CALL mpi_irecv(
1273 s rbuf(l),siz,real,it_spmd(i),msgtyp,
1274 g spmd_comm_world,req_r(i),ierror)
1275 l = l + siz
1276 ENDIF
1277 iad_recv(i+1) = l
1278 END DO
1279 l = 1
1280 iad_send(1) = 1
1281 DO i=1,nspmd
1282 IF(size==6) THEN
1283#include "vectorize.inc"
1284 DO j=iad_elem(1,i),iad_elem(1,i+1)-1
1285 nod = fr_elem(j)
1286 sbuf(l ) = secfcum(1,nod)
1287 sbuf(l+1) = secfcum(2,nod)
1288 sbuf(l+2) = secfcum(3,nod)
1289 sbuf(l+3) = secfcum(5,nod)
1290 sbuf(l+4) = secfcum(6,nod)
1291 sbuf(l+5) = secfcum(7,nod)
1292 l = l + SIZE
1293 END DO
1294 ELSE
1295#include "vectorize.inc"
1296 DO j=iad_elem(1,i),iad_elem(1,i+1)-1
1297 nod = fr_elem(j)
1298 sbuf(l ) = secfcum(1,nod)
1299 sbuf(l+1) = secfcum(2,nod)
1300 sbuf(l+2) = secfcum(3,nod)
1301 l = l + SIZE
1302 END DO
1303 END IF
1304 iad_send(i+1) = l
1305 ENDDO
1306C
1307C echange messages
1308C
1309 DO i=1,nspmd
1310C--------------------------------------------------------------------
1311 IF(iad_elem(1,i+1)-iad_elem(1,i)>0)THEN
1312 msgtyp = msgoff
1313 siz = iad_send(i+1)-iad_send(i)
1314 l = iad_send(i)
1315 CALL mpi_isend(
1316 s sbuf(l),siz,real,it_spmd(i),msgtyp,
1317 g spmd_comm_world,req_s(i),ierror)
1318 ENDIF
1319C--------------------------------------------------------------------
1320 ENDDO
1321C
1322 DO i = 1, nspmd
1323 nb_nod = iad_elem(1,i+1)-iad_elem(1,i)
1324 IF(nb_nod>0)THEN
1325 CALL mpi_wait(req_r(i),status,ierror)
1326 l = iad_recv(i)
1327 IF(size==6) THEN
1328#include "vectorize.inc"
1329 DO j=iad_elem(1,i),iad_elem(1,i+1)-1
1330 nod = fr_elem(j)
1331 secfcum(1,nod) = secfcum(1,nod) + rbuf(l)
1332 secfcum(2,nod) = secfcum(2,nod) + rbuf(l+1)
1333 secfcum(3,nod) = secfcum(3,nod) + rbuf(l+2)
1334 secfcum(5,nod) = secfcum(5,nod) + rbuf(l+3)
1335 secfcum(6,nod) = secfcum(6,nod) + rbuf(l+4)
1336 secfcum(7,nod) = secfcum(7,nod) + rbuf(l+5)
1337 l = l + SIZE
1338 END DO
1339 ELSE
1340#include "vectorize.inc"
1341 DO j=iad_elem(1,i),iad_elem(1,i+1)-1
1342 nod = fr_elem(j)
1343 secfcum(1,nod) = secfcum(1,nod) + rbuf(l)
1344 secfcum(2,nod) = secfcum(2,nod) + rbuf(l+1)
1345 secfcum(3,nod) = secfcum(3,nod) + rbuf(l+2)
1346 l = l + SIZE
1347 END DO
1348 END IF
1349 END IF
1350 END DO
1351C
1352C Remise a 0 de SECFCUM si P non main
1353C
1354 IF(size==6) THEN
1355#include "vectorize.inc"
1356 DO i = 1, nnod
1357 nod = nstrf(i)
1358 secfcum(1,nod) = secfcum(1,nod)*weight(nod)
1359 secfcum(2,nod) = secfcum(2,nod)*weight(nod)
1360 secfcum(3,nod) = secfcum(3,nod)*weight(nod)
1361 secfcum(5,nod) = secfcum(5,nod)*weight(nod)
1362 secfcum(6,nod) = secfcum(6,nod)*weight(nod)
1363 secfcum(7,nod) = secfcum(7,nod)*weight(nod)
1364 END DO
1365 ELSE
1366#include "vectorize.inc"
1367 DO i = 1, nnod
1368 nod = nstrf(i)
1369 secfcum(1,nod) = secfcum(1,nod)*weight(nod)
1370 secfcum(2,nod) = secfcum(2,nod)*weight(nod)
1371 secfcum(3,nod) = secfcum(3,nod)*weight(nod)
1372 END DO
1373 END IF
1374C
1375 DO i = 1, nspmd
1376 IF(iad_elem(1,i+1)-iad_elem(1,i)>0)THEN
1377 CALL mpi_wait(req_s(i),status,ierror)
1378 ENDIF
1379 ENDDO
1380C
1381#endif
1382 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_sec()

subroutine spmd_exch_sec ( integer, dimension(*) nstrf,
x,
ms,
integer, dimension(*) weight,
xsec,
integer, dimension(nspmd+1,*) fr_sec,
integer, dimension(4,*) iad_sec,
integer lsend1,
integer lrecv1,
integer lsend2,
integer lrecv2,
integer, dimension(*) weight_md )

Definition at line 412 of file spmd_section.F.

415C maj X (MS) N1, N2, N3 (+ NNODS) sur procs distants
416C-----------------------------------------------
417C I m p l i c i t T y p e s
418C-----------------------------------------------
419 USE spmd_comm_world_mod, ONLY : spmd_comm_world
420#include "implicit_f.inc"
421C-----------------------------------------------
422C M e s s a g e P a s s i n g
423C-----------------------------------------------
424#include "spmd.inc"
425C-----------------------------------------------
426C C o m m o n B l o c k s
427C-----------------------------------------------
428#include "com01_c.inc"
429#include "com04_c.inc"
430#include "task_c.inc"
431C-----------------------------------------------
432C D u m m y A r g u m e n t s
433C-----------------------------------------------
434 INTEGER NSTRF(*), WEIGHT(*), FR_SEC(NSPMD+1,*), IAD_SEC(4,*),
435 . LSEND1, LSEND2, LRECV1, LRECV2,WEIGHT_MD(*)
436 my_real
437 . x(3,*), ms(*), xsec(4,3,*)
438C-----------------------------------------------
439C L o c a l V a r i a b l e s
440C-----------------------------------------------
441#ifdef MPI
442 INTEGER LOC_PROC,A_AR,N,L,I,J,II,K,M,JJ, LEN,A_AR2,
443 . MSGTYP,MSGOFF,MSGOFF2,SIZ,IDEBR,IDEBS,ICC,IFRAM,
444 . IERROR, NBIRECV, NBISEND, INDEX, NBRBY, NBNOD,
445 . PMAIN, IDEB, LENS, LENR, K0, K2, N1, N2, N3, NNOD,
446 . NELSEG,NN,NB,
447 . IAD_SEND(NSPMD+1),IAD_RECV(NSPMD+1),
448 . REQ_R(NSPMD), REQ_S(NSPMD),
449 . IRINDEX(NSPMD), ISINDEX(NSPMD),
450 . STATUS(MPI_STATUS_SIZE),IAD_STMP(NSPMD)
451 DATA msgoff/4001/
452 DATA msgoff2/4002/
453 parameter(a_ar = 5)
454 parameter(a_ar2 = 13)
455 my_real
456 . mas, xxc, yyc ,zzc, dsec(nsect),
457 . sbuf(a_ar*lsend1),sbuf2(a_ar2*lsend2),
458 . rbuf(a_ar*lrecv1),rbuf2(a_ar2*lrecv2)
459C-----------------------------------------------
460C S o u r c e L i n e s
461C-----------------------------------------------
462 loc_proc = ispmd + 1
463C
464 nbirecv = 0
465 nbisend = 0
466 idebr = 1
467 idebs = 1
468 DO i = 1, nspmd
469 iad_recv(i) = idebr
470 IF(iad_sec(2,i)>0) THEN
471 msgtyp = msgoff
472 nbirecv = nbirecv + 1
473 irindex(nbirecv) = i
474 siz = iad_sec(2,i)*a_ar
475 CALL mpi_irecv(
476 s rbuf(idebr),siz,real,it_spmd(i),msgtyp,
477 g spmd_comm_world,req_r(nbirecv),ierror)
478 idebr = idebr + siz
479 ENDIF
480 iad_send(i) = idebs
481 IF(iad_sec(1,i)>0) THEN
482 nbisend = nbisend + 1
483 isindex(nbisend) = i
484 siz = iad_sec(1,i)*a_ar
485 idebs = idebs + siz
486 iad_stmp(i)=iad_send(i)
487 ENDIF
488 ENDDO
489 iad_recv(nspmd+1) = idebr
490C
491 ideb = 0
492C
493 k0=nstrf(25)
494 DO i=1,nsect
495 pmain = fr_sec(nspmd+1,i)
496 n1 = nstrf(k0+3)
497 n2 = nstrf(k0+4)
498 n3 = nstrf(k0+5)
499 nnod = nstrf(k0+6)
500 k2 = k0+30+nstrf(k0+14)
501 nelseg = nstrf(k0+7)+nstrf(k0+8)+nstrf(k0+9)+nstrf(k0+10)+
502 + nstrf(k0+11)+nstrf(k0+12)+nstrf(k0+13)
503 ifram = nstrf(k0+26)
504 IF(pmain>0.AND.loc_proc/=pmain) THEN
505 l = iad_stmp(pmain)
506 IF (ifram<=10.OR.n1/=0) THEN
507 IF(n1>0) THEN
508 IF(weight(n1)==1) THEN
509 sbuf(l) = 1
510 sbuf(l+1) = x(1,n1)
511 sbuf(l+2) = x(2,n1)
512 sbuf(l+3) = x(3,n1)
513 sbuf(l+4) = zero
514 l = l + a_ar
515 END IF
516 END IF
517 IF(n2>0) THEN
518 IF(weight(n2)==1) THEN
519 sbuf(l) = 2
520 sbuf(l+1) = x(1,n2)
521 sbuf(l+2) = x(2,n2)
522 sbuf(l+3) = x(3,n2)
523 sbuf(l+4) = zero
524 l = l + a_ar
525 END IF
526 END IF
527 IF(n3>0) THEN
528 IF(weight(n3)==1) THEN
529 sbuf(l) = 3
530 sbuf(l+1) = x(1,n3)
531 sbuf(l+2) = x(2,n3)
532 sbuf(l+3) = x(3,n3)
533 sbuf(l+4) = zero
534 l = l + a_ar
535 END IF
536 END IF
537 END IF
538 IF(mod(ifram,10)==1) THEN
539 xxc = zero
540 yyc = zero
541 zzc = zero
542 icc = 0
543 DO nn = 1, nnod
544 n = nstrf(k2+nn-1)
545 IF (weight_md(n)==1) THEN
546 xxc=xxc+x(1,n)
547 yyc=yyc+x(2,n)
548 zzc=zzc+x(3,n)
549 icc = icc + 1
550 END IF
551 END DO
552 IF(icc>0) THEN
553 sbuf(l) = 4
554 sbuf(l+1) = xxc
555 sbuf(l+2) = yyc
556 sbuf(l+3) = zzc
557 sbuf(l+4) = icc
558 l = l + a_ar
559 END IF
560 ELSEIF(mod(ifram,10)==2) THEN
561 xxc = zero
562 yyc = zero
563 zzc = zero
564 mas = zero
565 icc = 0
566 DO nn = 1, nnod
567 n = nstrf(k2+nn-1)
568 IF (weight_md(n)==1) THEN
569 xxc=xxc+x(1,n)*ms(n)
570 yyc=yyc+x(2,n)*ms(n)
571 zzc=zzc+x(3,n)*ms(n)
572 mas=mas+ms(n)
573 icc = icc + 1
574 END IF
575 END DO
576 IF(icc>0) THEN
577 sbuf(l) = 5
578 sbuf(l+1) = xxc
579 sbuf(l+2) = yyc
580 sbuf(l+3) = zzc
581 sbuf(l+4) = mas
582 l = l + a_ar
583 END IF
584 END IF
585 iad_stmp(pmain)=l
586 ELSE
587C pmain => stockage direct
588 IF (ifram<=10.OR.n1/=0) THEN
589 IF(n1>0) THEN
590 IF(weight(n1)==1) THEN
591 xsec(1,1,i) = x(1,n1)
592 xsec(1,2,i) = x(2,n1)
593 xsec(1,3,i) = x(3,n1)
594 END IF
595 END IF
596 IF(n2>0) THEN
597 IF(weight(n2)==1) THEN
598 xsec(2,1,i) = x(1,n2)
599 xsec(2,2,i) = x(2,n2)
600 xsec(2,3,i) = x(3,n2)
601 END IF
602 END IF
603 IF(n3>0) THEN
604 IF(weight(n3)==1) THEN
605 xsec(3,1,i) = x(1,n3)
606 xsec(3,2,i) = x(2,n3)
607 xsec(3,3,i) = x(3,n3)
608 END IF
609 END IF
610 END IF
611 IF(mod(ifram,10)==1) THEN
612 xxc = zero
613 yyc = zero
614 zzc = zero
615 icc = 0
616 DO nn = 1, nnod
617 n = nstrf(k2+nn-1)
618 IF (weight_md(n)==1) THEN
619 xxc=xxc+x(1,n)
620 yyc=yyc+x(2,n)
621 zzc=zzc+x(3,n)
622 icc = icc + 1
623 END IF
624 END DO
625C SBUF(L) = 4
626 xsec(4,1,i) = xxc
627 xsec(4,2,i) = yyc
628 xsec(4,3,i) = zzc
629 dsec(i) = icc
630 ELSEIF(mod(ifram,10)==2) THEN
631 xxc = zero
632 yyc = zero
633 zzc = zero
634 mas = zero
635 icc = 0
636 DO nn = 1, nnod
637 n = nstrf(k2+nn-1)
638 IF (weight_md(n)==1) THEN
639 xxc=xxc+x(1,n)*ms(n)
640 yyc=yyc+x(2,n)*ms(n)
641 zzc=zzc+x(3,n)*ms(n)
642 mas=mas+ms(n)
643 icc = icc + 1
644 END IF
645 END DO
646C SBUF(L) = 5
647 xsec(4,1,i) = xxc
648 xsec(4,2,i) = yyc
649 xsec(4,3,i) = zzc
650 dsec(i) = mas
651 END IF
652 END IF
653 k0=nstrf(k0+24)
654 END DO
655C
656 DO l = 1, nbisend
657 i = isindex(l)
658 siz = iad_stmp(i)-iad_send(i)
659 idebs = iad_send(i)
660 msgtyp = msgoff
661 CALL mpi_isend(
662 s sbuf(idebs),siz,real,it_spmd(i),msgtyp,
663 g spmd_comm_world,req_s(i),ierror)
664 ENDDO
665C
666 DO ii = 1, nbirecv
667 CALL mpi_waitany(nbirecv,req_r,index,status,ierror)
668 i = irindex(index)
669 ideb = iad_recv(i)
670 DO n = 1, nsect
671 pmain = fr_sec(nspmd+1,n)
672 IF(loc_proc==pmain) THEN
673 nb = fr_sec(i,n)
674 IF(nb>0) THEN
675 DO k = 1, nb
676 nn = nint(rbuf(ideb+(k-1)*a_ar ))
677 IF(nn==1) THEN
678 xsec(1,1,n) = rbuf(ideb+(k-1)*a_ar+1)
679 xsec(1,2,n) = rbuf(ideb+(k-1)*a_ar+2)
680 xsec(1,3,n) = rbuf(ideb+(k-1)*a_ar+3)
681 ELSEIF(nn==2) THEN
682 xsec(2,1,n) = rbuf(ideb+(k-1)*a_ar+1)
683 xsec(2,2,n) = rbuf(ideb+(k-1)*a_ar+2)
684 xsec(2,3,n) = rbuf(ideb+(k-1)*a_ar+3)
685 ELSEIF(nn==3) THEN
686 xsec(3,1,n) = rbuf(ideb+(k-1)*a_ar+1)
687 xsec(3,2,n) = rbuf(ideb+(k-1)*a_ar+2)
688 xsec(3,3,n) = rbuf(ideb+(k-1)*a_ar+3)
689 ELSEIF(nn==4.OR.nn==5) THEN
690 xsec(4,1,n) = xsec(4,1,n)+rbuf(ideb+(k-1)*a_ar+1)
691 xsec(4,2,n) = xsec(4,2,n)+rbuf(ideb+(k-1)*a_ar+2)
692 xsec(4,3,n) = xsec(4,3,n)+rbuf(ideb+(k-1)*a_ar+3)
693 dsec(n) = dsec(n) + rbuf(ideb+(k-1)*a_ar+4)
694 END IF
695 END DO
696 ideb = ideb + a_ar*nb
697 END IF
698 END IF
699 END DO
700 END DO
701 DO l = 1, nbisend
702 i = isindex(l)
703 CALL mpi_wait(req_s(i),status,ierror)
704 END DO
705C
706C Calcul des valeurs nodales sur proc main si besoin (test global ?)
707C
708 k0=nstrf(25)
709 DO n = 1, nsect
710 pmain = fr_sec(nspmd+1,n)
711 ifram = nstrf(k0+26)
712 IF(loc_proc==pmain) THEN
713 IF(mod(ifram,10)==1.OR.mod(ifram,10)==2)THEN
714 IF(dsec(n)/=0) THEN
715 xsec(4,1,n) = xsec(4,1,n)/dsec(n)
716 xsec(4,2,n) = xsec(4,2,n)/dsec(n)
717 xsec(4,3,n) = xsec(4,3,n)/dsec(n)
718 END IF
719 END IF
720 END IF
721 k0=nstrf(k0+24)
722 END DO
723C
724 nbirecv = 0
725 idebr = 1
726 DO i = 1, nspmd
727 iad_recv(i) = idebr
728 IF(iad_sec(4,i)>0) THEN
729 msgtyp = msgoff2
730 nbirecv = nbirecv + 1
731 irindex(nbirecv) = i
732 siz = iad_sec(4,i)*a_ar2
733 CALL mpi_irecv(
734 s rbuf2(idebr),siz,real,it_spmd(i),msgtyp,
735 g spmd_comm_world,req_r(nbirecv),ierror)
736 idebr = idebr + siz
737 ENDIF
738 ENDDO
739C reste a coder l'envoi des noeuds vers les procs concernes
740 nbisend = 0
741 IF(iad_sec(3,nspmd+1)>0) THEN
742 l = 0
743 k0=nstrf(25)
744 DO n = 1, nsect
745 pmain = fr_sec(nspmd+1,n)
746 n1 = nstrf(k0+3)
747 n2 = nstrf(k0+4)
748 n3 = nstrf(k0+5)
749 ifram = nstrf(k0+26)
750 IF(loc_proc==pmain) THEN
751 sbuf2(l+1) = n
752 IF (ifram<=10.OR.n1/=0) THEN
753 sbuf2(l+2) = xsec(1,1,n)
754 sbuf2(l+3) = xsec(1,2,n)
755 sbuf2(l+4) = xsec(1,3,n)
756 sbuf2(l+5) = xsec(2,1,n)
757 sbuf2(l+6) = xsec(2,2,n)
758 sbuf2(l+7) = xsec(2,3,n)
759 sbuf2(l+8) = xsec(3,1,n)
760 sbuf2(l+9) = xsec(3,2,n)
761 sbuf2(l+10)= xsec(3,3,n)
762 ELSE
763 sbuf2(l+2) = zero
764 sbuf2(l+3) = zero
765 sbuf2(l+4) = zero
766 sbuf2(l+5) = zero
767 sbuf2(l+6) = zero
768 sbuf2(l+7) = zero
769 sbuf2(l+8) = zero
770 sbuf2(l+9) = zero
771 sbuf2(l+10)= zero
772 END IF
773 IF(mod(ifram,10)==1.OR.mod(ifram,10)==2) THEN
774 sbuf2(l+11) = xsec(4,1,n)
775 sbuf2(l+12) = xsec(4,2,n)
776 sbuf2(l+13) = xsec(4,3,n)
777 ELSE
778 sbuf2(l+11) = zero
779 sbuf2(l+12) = zero
780 sbuf2(l+13) = zero
781 END IF
782 l = l + a_ar2
783 END IF
784 k0=nstrf(k0+24)
785 END DO
786C
787 DO i = 1, nspmd
788 IF(iad_sec(3,i)>0) THEN
789 msgtyp = msgoff2
790 nbisend = nbisend + 1
791 isindex(nbisend) = i
792 CALL mpi_isend(
793 s sbuf2,l,real,it_spmd(i),msgtyp,
794 g spmd_comm_world,req_s(i),ierror)
795 END IF
796 END DO
797 END IF
798C
799 DO ii = 1, nbirecv
800 CALL mpi_waitany(nbirecv,req_r,index,status,ierror)
801 i = irindex(index)
802 l = iad_recv(i)
803 nbnod = iad_sec(4,i)
804 DO j = 1, nbnod
805 n = nint(rbuf2(l))
806 xsec(1,1,n) = rbuf2(l+1)
807 xsec(1,2,n) = rbuf2(l+2)
808 xsec(1,3,n) = rbuf2(l+3)
809 xsec(2,1,n) = rbuf2(l+4)
810 xsec(2,2,n) = rbuf2(l+5)
811 xsec(2,3,n) = rbuf2(l+6)
812 xsec(3,1,n) = rbuf2(l+7)
813 xsec(3,2,n) = rbuf2(l+8)
814 xsec(3,3,n) = rbuf2(l+9)
815 xsec(4,1,n) = rbuf2(l+10)
816 xsec(4,2,n) = rbuf2(l+11)
817 xsec(4,3,n) = rbuf2(l+12)
818 l = l + a_ar2
819 END DO
820 END DO
821C
822 DO l = 1, nbisend
823 i = isindex(l)
824 CALL mpi_wait(req_s(i),status,ierror)
825 END DO
826C
827#endif
828 RETURN
subroutine mpi_waitany(cnt, array_of_requests, index, status, ierr)
Definition mpi.f:549

◆ spmd_sd_cut()

subroutine spmd_sd_cut ( secbufg,
integer nnodg,
secbuf1,
secbuf2,
integer nnod,
integer, dimension(*) fr_cut,
integer, dimension(*) iad_cut,
integer iflg )

Definition at line 1393 of file spmd_section.F.

1395C-----------------------------------------------
1396C I m p l i c i t T y p e s
1397C-----------------------------------------------
1398 USE spmd_comm_world_mod, ONLY : spmd_comm_world
1399#include "implicit_f.inc"
1400C-----------------------------------------------
1401C M e s s a g e P a s s i n g
1402C-----------------------------------------------
1403#include "spmd.inc"
1404C-----------------------------------------------
1405C C o m m o n B l o c k s
1406C-----------------------------------------------
1407#include "com01_c.inc"
1408#include "task_c.inc"
1409C-----------------------------------------------
1410C D u m m y A r g u m e n t s
1411C-----------------------------------------------
1412 INTEGER FR_CUT(*),IAD_CUT(*),
1413 . NNOD, NNODG, IFLG
1414 my_real
1415 . secbufg(*), secbuf1(*), secbuf2(*)
1416C-----------------------------------------------
1417C L o c a l V a r i a b l e s
1418C-----------------------------------------------
1419#ifdef MPI
1420 INTEGER LOC_PROC,N,L,I,K,LEN,P,NN,N0,OFFG,
1421 . MSGTYP,MSGOFF,SIZ,IDEBR,
1422 . IERROR, NBIRECV,IDEB,
1423 . IAD_RECV(NSPMD+1),
1424 . STATUS(MPI_STATUS_SIZE)
1425 DATA msgoff/4006/
1426 my_real
1427 . sbuf(6*iflg*nnodg),rbuf(6*iflg*nnod)
1428C-----------------------------------------------
1429C S o u r c e L i n e s
1430C-----------------------------------------------
1431 loc_proc = ispmd + 1
1432C
1433 IF(loc_proc==1) THEN
1434 n0 = 0
1435 offg = 6*nnodg
1436C Traitement P0
1437 IF(iad_cut(1)/=0) THEN
1438 nn = iad_cut(1)
1439 IF(iflg==2) THEN
1440 DO i = 1, nn
1441 n = fr_cut(i+n0)
1442 secbuf1(6*i-5) = secbufg((n-1)*6+1)
1443 secbuf1(6*i-4) = secbufg((n-1)*6+2)
1444 secbuf1(6*i-3) = secbufg((n-1)*6+3)
1445 secbuf1(6*i-2) = secbufg((n-1)*6+4)
1446 secbuf1(6*i-1) = secbufg((n-1)*6+5)
1447 secbuf1(6*i) = secbufg((n-1)*6+6)
1448 secbuf2(6*i-5) = secbufg(offg+(n-1)*6+1)
1449 secbuf2(6*i-4) = secbufg(offg+(n-1)*6+2)
1450 secbuf2(6*i-3) = secbufg(offg+(n-1)*6+3)
1451 secbuf2(6*i-2) = secbufg(offg+(n-1)*6+4)
1452 secbuf2(6*i-1) = secbufg(offg+(n-1)*6+5)
1453 secbuf2(6*i) = secbufg(offg+(n-1)*6+6)
1454 END DO
1455 ELSEIF(iflg==1) THEN
1456 DO i = 1, nn
1457 n = fr_cut(i+n0)
1458 secbuf1(6*i-5) = secbufg((n-1)*6+1)
1459 secbuf1(6*i-4) = secbufg((n-1)*6+2)
1460 secbuf1(6*i-3) = secbufg((n-1)*6+3)
1461 secbuf1(6*i-2) = secbufg((n-1)*6+4)
1462 secbuf1(6*i-1) = secbufg((n-1)*6+5)
1463 secbuf1(6*i) = secbufg((n-1)*6+6)
1464 END DO
1465 END IF
1466 n0 = n0 + nn
1467 END IF
1468C Traitement autres procs
1469 DO p = 2, nspmd
1470 IF(iad_cut(p)/=0) THEN
1471 l = 0
1472 nn = iad_cut(p)
1473 IF(iflg==2) THEN
1474 DO k = 1, nn
1475 n = fr_cut(k+n0)
1476 sbuf(l+1) = secbufg((n-1)*6+1)
1477 sbuf(l+2) = secbufg((n-1)*6+2)
1478 sbuf(l+3) = secbufg((n-1)*6+3)
1479 sbuf(l+4) = secbufg((n-1)*6+4)
1480 sbuf(l+5) = secbufg((n-1)*6+5)
1481 sbuf(l+6) = secbufg((n-1)*6+6)
1482 sbuf(l+7) = secbufg(offg+(n-1)*6+1)
1483 sbuf(l+8) = secbufg(offg+(n-1)*6+2)
1484 sbuf(l+9) = secbufg(offg+(n-1)*6+3)
1485 sbuf(l+10)= secbufg(offg+(n-1)*6+4)
1486 sbuf(l+11)= secbufg(offg+(n-1)*6+5)
1487 sbuf(l+12)= secbufg(offg+(n-1)*6+6)
1488 l = l + 12
1489 END DO
1490 ELSEIF(iflg==1) THEN
1491 DO k = 1, nn
1492 n = fr_cut(k+n0)
1493 sbuf(l+1) = secbufg((n-1)*6+1)
1494 sbuf(l+2) = secbufg((n-1)*6+2)
1495 sbuf(l+3) = secbufg((n-1)*6+3)
1496 sbuf(l+4) = secbufg((n-1)*6+4)
1497 sbuf(l+5) = secbufg((n-1)*6+5)
1498 sbuf(l+6) = secbufg((n-1)*6+6)
1499 l = l + 6
1500 END DO
1501 END IF
1502 n0 = n0 + nn
1503 msgtyp = msgoff
1504 CALL mpi_send(sbuf,l,real,it_spmd(p),msgtyp,
1505 . spmd_comm_world,ierror)
1506 END IF
1507 END DO
1508 ELSEIF(loc_proc/=1.AND.nnod>0) THEN
1509 msgtyp = msgoff
1510 siz = nnod*iflg*6
1511 CALL mpi_recv(
1512 s rbuf,siz,real,it_spmd(1),msgtyp,
1513 g spmd_comm_world,status,ierror)
1514C
1515 IF(iflg==2) THEN
1516 l = 0
1517 DO i = 1, nnod
1518 secbuf1(6*i-5) = rbuf(l+1)
1519 secbuf1(6*i-4) = rbuf(l+2)
1520 secbuf1(6*i-3) = rbuf(l+3)
1521 secbuf1(6*i-2) = rbuf(l+4)
1522 secbuf1(6*i-1) = rbuf(l+5)
1523 secbuf1(6*i) = rbuf(l+6)
1524 secbuf2(6*i-5) = rbuf(l+7)
1525 secbuf2(6*i-4) = rbuf(l+8)
1526 secbuf2(6*i-3) = rbuf(l+9)
1527 secbuf2(6*i-2) = rbuf(l+10)
1528 secbuf2(6*i-1) = rbuf(l+11)
1529 secbuf2(6*i) = rbuf(l+12)
1530 l = l + 12
1531 END DO
1532 ELSEIF(iflg==1) THEN
1533 l = 0
1534 DO i = 1, nnod
1535 secbuf1(6*i-5) = rbuf(l+1)
1536 secbuf1(6*i-4) = rbuf(l+2)
1537 secbuf1(6*i-3) = rbuf(l+3)
1538 secbuf1(6*i-2) = rbuf(l+4)
1539 secbuf1(6*i-1) = rbuf(l+5)
1540 secbuf1(6*i) = rbuf(l+6)
1541 l = l + 6
1542 END DO
1543 END IF
1544 END IF
1545C
1546#endif
1547 RETURN
subroutine mpi_recv(buf, cnt, datatype, source, tag, comm, status, ierr)
Definition mpi.f:461
subroutine mpi_send(buf, cnt, datatype, dest, tag, comm, ierr)
Definition mpi.f:480

◆ spmd_wrt_cutd()

subroutine spmd_wrt_cutd ( integer nnod,
integer, dimension(*) nstrf,
d,
dr,
integer, dimension(*) rg_cut,
integer, dimension(*) iad_cut,
integer nsize,
integer nnodg,
integer, dimension(*) weight,
integer iflg )

Definition at line 840 of file spmd_section.F.

842C-----------------------------------------------
843C I m p l i c i t T y p e s
844C-----------------------------------------------
845 USE spmd_comm_world_mod, ONLY : spmd_comm_world
846#include "implicit_f.inc"
847C-----------------------------------------------
848C M e s s a g e P a s s i n g
849C-----------------------------------------------
850#include "spmd.inc"
851C-----------------------------------------------
852C C o m m o n B l o c k s
853C-----------------------------------------------
854#include "com01_c.inc"
855#include "task_c.inc"
856C-----------------------------------------------
857C D u m m y A r g u m e n t s
858C-----------------------------------------------
859 INTEGER NSTRF(*), WEIGHT(*), RG_CUT(*), IAD_CUT(*),
860 . NNOD, NSIZE, NNODG, IFLG
861 my_real
862 . d(3,*), dr(3,*)
863C-----------------------------------------------
864C L o c a l V a r i a b l e s
865C-----------------------------------------------
866#ifdef MPI
867 INTEGER LOC_PROC,N,L,I,K,LEN,II,INDEX,NB,
868 . MSGTYP,MSGOFF,SIZ,IDEBR,
869 . IERROR, NBIRECV,IDEB,
870 . IAD_RECV(NSPMD+1),REQ_R(NSPMD),IRINDEX(NSPMD),
871 . STATUS(MPI_STATUS_SIZE)
872 DATA msgoff/4003/
873 my_real
874 . sbuf((3*iflg+1)*nnod),rbuf((3*iflg+1)*nsize),
875 . secbufg(3*iflg,nnodg)
876 real*4 r4
877C-----------------------------------------------
878C S o u r c e L i n e s
879C-----------------------------------------------
880 loc_proc = ispmd + 1
881C
882 IF(loc_proc/=1.AND.nnod>0) THEN
883 l = 0
884 IF(iflg==2) THEN
885 DO k = 1, nnod
886 n = nstrf(k)
887 IF(weight(n)==1) THEN
888 sbuf(l+1) = rg_cut(k)
889 sbuf(l+2) = d(1,n)
890 sbuf(l+3) = d(2,n)
891 sbuf(l+4) = d(3,n)
892 sbuf(l+5) = dr(1,n)
893 sbuf(l+6) = dr(2,n)
894 sbuf(l+7) = dr(3,n)
895 l = l + 7
896 END IF
897 END DO
898 ELSEIF(iflg==1) THEN
899 DO k = 1, nnod
900 n = nstrf(k)
901 IF(weight(n)==1) THEN
902 sbuf(l+1) = rg_cut(k)
903 sbuf(l+2) = d(1,n)
904 sbuf(l+3) = d(2,n)
905 sbuf(l+4) = d(3,n)
906 l = l + 4
907 END IF
908 END DO
909 END IF
910 msgtyp = msgoff
911 CALL mpi_send(sbuf,l,real,it_spmd(1),msgtyp,
912 . spmd_comm_world,ierror)
913 ELSEIF(loc_proc==1) THEN
914C P0
915 nbirecv = 0
916 idebr = 1
917 DO i = 2, nspmd
918 iad_recv(i) = idebr
919 IF(iad_cut(i)>0) THEN
920 msgtyp = msgoff
921 nbirecv = nbirecv + 1
922 irindex(nbirecv) = i
923 siz = iad_cut(i)*(1+iflg*3)
924 CALL mpi_irecv(
925 s rbuf(idebr),siz,real,it_spmd(i),msgtyp,
926 g spmd_comm_world,req_r(nbirecv),ierror)
927 idebr = idebr + siz
928 ENDIF
929 END DO
930 iad_recv(nspmd+1) = idebr
931C
932 IF(iflg==2) THEN
933 DO k = 1, nnod
934 n = nstrf(k)
935 IF(weight(n)==1) THEN
936 i = rg_cut(k)
937 secbufg(1,i) = d(1,n)
938 secbufg(2,i) = d(2,n)
939 secbufg(3,i) = d(3,n)
940 secbufg(4,i) = dr(1,n)
941 secbufg(5,i) = dr(2,n)
942 secbufg(6,i) = dr(3,n)
943 END IF
944 END DO
945 ELSEIF(iflg==1) THEN
946 DO k = 1, nnod
947 n = nstrf(k)
948 IF(weight(n)==1) THEN
949 i = rg_cut(k)
950 secbufg(1,i) = d(1,n)
951 secbufg(2,i) = d(2,n)
952 secbufg(3,i) = d(3,n)
953 END IF
954 END DO
955 END IF
956C
957 DO ii = 1, nbirecv
958 CALL mpi_waitany(nbirecv,req_r,index,status,ierror)
959 i = irindex(index)
960 ideb = iad_recv(i)
961 len = iflg*3+1
962 CALL mpi_get_count(status,real,siz,ierror)
963 nb = siz/len
964 IF(iflg==2) THEN
965 DO k = 1, nb
966 i = nint(rbuf(ideb+(k-1)*len))
967 secbufg(1,i) = rbuf(ideb+(k-1)*len+1)
968 secbufg(2,i) = rbuf(ideb+(k-1)*len+2)
969 secbufg(3,i) = rbuf(ideb+(k-1)*len+3)
970 secbufg(4,i) = rbuf(ideb+(k-1)*len+4)
971 secbufg(5,i) = rbuf(ideb+(k-1)*len+5)
972 secbufg(6,i) = rbuf(ideb+(k-1)*len+6)
973 END DO
974 ELSE
975 DO k = 1, nb
976 i = nint(rbuf(ideb+(k-1)*len))
977 secbufg(1,i) = rbuf(ideb+(k-1)*len+1)
978 secbufg(2,i) = rbuf(ideb+(k-1)*len+2)
979 secbufg(3,i) = rbuf(ideb+(k-1)*len+3)
980 END DO
981 END IF
982 END DO
983C ecriture p0
984 IF(iflg==2) THEN
985 DO i = 1,nnodg
986 r4 = secbufg(1,i)
987 CALL write_r_c(r4,1)
988 r4 = secbufg(2,i)
989 CALL write_r_c(r4,1)
990 r4 = secbufg(3,i)
991 CALL write_r_c(r4,1)
992 r4 = secbufg(4,i)
993 CALL write_r_c(r4,1)
994 r4 = secbufg(5,i)
995 CALL write_r_c(r4,1)
996 r4 = secbufg(6,i)
997 CALL write_r_c(r4,1)
998 END DO
999 ELSEIF(iflg==1) THEN
1000 DO i = 1,nnodg
1001 r4 = secbufg(1,i)
1002 CALL write_r_c(r4,1)
1003 r4 = secbufg(2,i)
1004 CALL write_r_c(r4,1)
1005 r4 = secbufg(3,i)
1006 CALL write_r_c(r4,1)
1007 r4 = zero
1008 CALL write_r_c(r4,1)
1009 CALL write_r_c(r4,1)
1010 CALL write_r_c(r4,1)
1011 END DO
1012 END IF
1013 END IF
1014C
1015#endif
1016 RETURN
subroutine mpi_get_count(status, datatype, cnt, ierr)
Definition mpi.f:296
void write_r_c(float *w, int *len)

◆ spmd_wrt_cutf()

subroutine spmd_wrt_cutf ( integer nnod,
integer, dimension(*) nstrf,
secfcum,
integer, dimension(*) rg_cut,
integer, dimension(*) iad_cut,
integer nsize,
integer nnodg,
integer, dimension(*) weight,
integer iflg )

Definition at line 1028 of file spmd_section.F.

1030C-----------------------------------------------
1031C I m p l i c i t T y p e s
1032C-----------------------------------------------
1033 USE spmd_comm_world_mod, ONLY : spmd_comm_world
1034#include "implicit_f.inc"
1035C-----------------------------------------------
1036C M e s s a g e P a s s i n g
1037C-----------------------------------------------
1038#include "spmd.inc"
1039C-----------------------------------------------
1040C C o m m o n B l o c k s
1041C-----------------------------------------------
1042#include "com01_c.inc"
1043#include "task_c.inc"
1044C-----------------------------------------------
1045C D u m m y A r g u m e n t s
1046C-----------------------------------------------
1047 INTEGER NSTRF(*), WEIGHT(*), RG_CUT(*), IAD_CUT(*),
1048 . NNOD, NSIZE, NNODG, IFLG
1049 my_real
1050 . secfcum(7,*)
1051C-----------------------------------------------
1052C L o c a l V a r i a b l e s
1053C-----------------------------------------------
1054#ifdef MPI
1055 INTEGER LOC_PROC,N,L,I,K,LEN,II,INDEX,NB,
1056 . MSGTYP,MSGOFF,SIZ,IDEBR,
1057 . IERROR, NBIRECV,IDEB,
1058 . IAD_RECV(NSPMD+1),REQ_R(NSPMD),IRINDEX(NSPMD),
1059 . STATUS(MPI_STATUS_SIZE)
1060 DATA msgoff/4004/
1061 my_real
1062 . sbuf((3*iflg+1)*nnod),rbuf((3*iflg+1)*nsize),
1063 . secbufg(3*iflg,nnodg)
1064 real*4 r4
1065C-----------------------------------------------
1066C S o u r c e L i n e s
1067C-----------------------------------------------
1068 loc_proc = ispmd + 1
1069C
1070 IF(loc_proc/=1.AND.nnod>0) THEN
1071 l = 0
1072 IF(iflg==2) THEN
1073 DO k = 1, nnod
1074 n = nstrf(k)
1075 IF(weight(n)==1) THEN
1076 sbuf(l+1) = rg_cut(k)
1077 sbuf(l+2) = secfcum(1,n)
1078 sbuf(l+3) = secfcum(2,n)
1079 sbuf(l+4) = secfcum(3,n)
1080 sbuf(l+5) = secfcum(5,n)
1081 sbuf(l+6) = secfcum(6,n)
1082 sbuf(l+7) = secfcum(7,n)
1083 l = l + 7
1084 END IF
1085 END DO
1086 ELSEIF(iflg==1) THEN
1087 DO k = 1, nnod
1088 n = nstrf(k)
1089 IF(weight(n)==1) THEN
1090 sbuf(l+1) = rg_cut(k)
1091 sbuf(l+2) = secfcum(1,n)
1092 sbuf(l+3) = secfcum(2,n)
1093 sbuf(l+4) = secfcum(3,n)
1094 l = l + 4
1095 END IF
1096 END DO
1097 END IF
1098 msgtyp = msgoff
1099 CALL mpi_send(sbuf,l,real,it_spmd(1),msgtyp,
1100 . spmd_comm_world,ierror)
1101 ELSEIF(loc_proc==1) THEN
1102C P0
1103 nbirecv = 0
1104 idebr = 1
1105 DO i = 2, nspmd
1106 iad_recv(i) = idebr
1107 IF(iad_cut(i)>0) THEN
1108 msgtyp = msgoff
1109 nbirecv = nbirecv + 1
1110 irindex(nbirecv) = i
1111 siz = iad_cut(i)*(1+iflg*3)
1112 CALL mpi_irecv(
1113 s rbuf(idebr),siz,real,it_spmd(i),msgtyp,
1114 g spmd_comm_world,req_r(nbirecv),ierror)
1115 idebr = idebr + siz
1116 ENDIF
1117 END DO
1118 iad_recv(nspmd+1) = idebr
1119C
1120 IF(iflg==2) THEN
1121 DO i = 1, nnodg
1122 secbufg(1,i) = zero
1123 secbufg(2,i) = zero
1124 secbufg(3,i) = zero
1125 secbufg(4,i) = zero
1126 secbufg(5,i) = zero
1127 secbufg(6,i) = zero
1128 END DO
1129C
1130 DO k = 1, nnod
1131 n = nstrf(k)
1132 i = rg_cut(k)
1133 secbufg(1,i) = secfcum(1,n)
1134 secbufg(2,i) = secfcum(2,n)
1135 secbufg(3,i) = secfcum(3,n)
1136 secbufg(4,i) = secfcum(5,n)
1137 secbufg(5,i) = secfcum(6,n)
1138 secbufg(6,i) = secfcum(7,n)
1139 END DO
1140 ELSEIF(iflg==1) THEN
1141 DO i = 1, nnodg
1142 secbufg(1,i) = zero
1143 secbufg(2,i) = zero
1144 secbufg(3,i) = zero
1145 END DO
1146C
1147 DO k = 1, nnod
1148 n = nstrf(k)
1149 i = rg_cut(k)
1150 secbufg(1,i) = secfcum(1,n)
1151 secbufg(2,i) = secfcum(2,n)
1152 secbufg(3,i) = secfcum(3,n)
1153 END DO
1154 END IF
1155C
1156 DO ii = 1, nbirecv
1157 CALL mpi_waitany(nbirecv,req_r,index,status,ierror)
1158 CALL mpi_get_count(status,real,nb,ierror)
1159 i = irindex(index)
1160 ideb = iad_recv(i)
1161 len = iflg*3+1
1162 nb = nb / len
1163 IF(iflg==2) THEN
1164 DO k = 1, nb
1165 i = nint(rbuf(ideb+(k-1)*len))
1166 secbufg(1,i) = secbufg(1,i) + rbuf(ideb+(k-1)*len+1)
1167 secbufg(2,i) = secbufg(2,i) + rbuf(ideb+(k-1)*len+2)
1168 secbufg(3,i) = secbufg(3,i) + rbuf(ideb+(k-1)*len+3)
1169 secbufg(4,i) = secbufg(4,i) + rbuf(ideb+(k-1)*len+4)
1170 secbufg(5,i) = secbufg(5,i) + rbuf(ideb+(k-1)*len+5)
1171 secbufg(6,i) = secbufg(6,i) + rbuf(ideb+(k-1)*len+6)
1172 END DO
1173 ELSE
1174 DO k = 1, nb
1175 i = nint(rbuf(ideb+(k-1)*len))
1176 secbufg(1,i) = secbufg(1,i) + rbuf(ideb+(k-1)*len+1)
1177 secbufg(2,i) = secbufg(2,i) + rbuf(ideb+(k-1)*len+2)
1178 secbufg(3,i) = secbufg(3,i) + rbuf(ideb+(k-1)*len+3)
1179 END DO
1180 END IF
1181 END DO
1182C ecriture p0
1183 IF(iflg==2) THEN
1184 DO i = 1,nnodg
1185 r4 = secbufg(1,i)
1186 CALL write_r_c(r4,1)
1187 r4 = secbufg(2,i)
1188 CALL write_r_c(r4,1)
1189 r4 = secbufg(3,i)
1190 CALL write_r_c(r4,1)
1191 r4 = secbufg(4,i)
1192 CALL write_r_c(r4,1)
1193 r4 = secbufg(5,i)
1194 CALL write_r_c(r4,1)
1195 r4 = secbufg(6,i)
1196 CALL write_r_c(r4,1)
1197 END DO
1198 ELSEIF(iflg==1) THEN
1199 DO i = 1,nnodg
1200 r4 = secbufg(1,i)
1201 CALL write_r_c(r4,1)
1202 r4 = secbufg(2,i)
1203 CALL write_r_c(r4,1)
1204 r4 = secbufg(3,i)
1205 CALL write_r_c(r4,1)
1206 r4 = zero
1207 CALL write_r_c(r4,1)
1208 CALL write_r_c(r4,1)
1209 CALL write_r_c(r4,1)
1210 END DO
1211 END IF
1212 END IF
1213C
1214#endif
1215 RETURN

◆ spmd_x_section()

subroutine spmd_x_section ( integer, dimension(*) nstrf,
x,
ms,
integer, dimension(*) weight,
wa )

Definition at line 29 of file spmd_section.F.

30C maj X (MS) N1, N2, N3 (+ NNODS) sur procs distants
31C-----------------------------------------------
32C I m p l i c i t T y p e s
33C-----------------------------------------------
34 USE spmd_comm_world_mod, ONLY : spmd_comm_world
35#include "implicit_f.inc"
36C-----------------------------------------------
37C M e s s a g e P a s s i n g
38C-----------------------------------------------
39#include "spmd.inc"
40C-----------------------------------------------
41C C o m m o n B l o c k s
42C-----------------------------------------------
43#include "com01_c.inc"
44#include "com04_c.inc"
45#include "task_c.inc"
46#include "scr03_c.inc"
47C-----------------------------------------------
48C D u m m y A r g u m e n t s
49C-----------------------------------------------
50 INTEGER NSTRF(*), WEIGHT(*)
52 . x(3,*), ms(*), wa(*)
53C-----------------------------------------------
54C L o c a l V a r i a b l e s
55C-----------------------------------------------
56#ifdef MPI
57 INTEGER MSGOFF,INFO,I,J,K,L,NELSEG,NNOD,NELC,NELTG,
58 . NN,P,ATID,ATAG,ALEN,IDEB,SIZ,K0,K1,K2,A_AR,IFRAM,
59 . REM_PROC,MSGTYP,REM_PROC2,MSGTYP2,N1,N2,N3,NOD,
60 . SENDTO(PARASIZ),RECVFR(PARASIZ),LOC_PROC,NB_NOD,
61 . REQ_R(PARASIZ),BUFSIZ,IALL
62 INTEGER STATUS(MPI_STATUS_SIZE),IERROR
63 SAVE sendto,recvfr,bufsiz
64 DATA msgoff/4000/
65C-----------------------------------------------
66C S o u r c e L i n e s
67C-----------------------------------------------
68 loc_proc = ispmd + 1
69 a_ar = 4
70C
71 IF(codvers<42) THEN
72 IF (ncycle==1) THEN
73C
74 bufsiz = nsect * 3 * a_ar
75C
76 DO i = 1,2*nsect
77 wa(i) = zero
78 ENDDO
79C
80 k1 = 1
81 DO i=1,nsect
82 nelc = nstrf(k1)
83 n1 = nstrf(k1+1)
84 n2 = nstrf(k1+2)
85 n3 = nstrf(k1+3)
86 neltg = nstrf(k1+4+3*nelc)
87 IF (weight(n1)==1.OR.weight(n2)==1.OR.weight(n3)==1)
88 + THEN
89 wa(i) = one
90 ENDIF
91 IF (nelc+neltg/=0.AND.(weight(n1)==0.OR.weight(n2)==0
92 + .OR.weight(n3)==0.)) THEN
93 wa(i+nsect) = one
94 ENDIF
95 k1=k1+4+3*nstrf(k1)
96 k1=k1+1+2*nstrf(k1)
97 END DO
98C
99 DO i = 1, nspmd
100 sendto(i) = 0
101 recvfr(i) = 0
102 ENDDO
103 siz = 2*nsect
104 ideb = siz + 1
105 DO i = 1, nspmd
106 rem_proc = mod(loc_proc + i-1,nspmd)+1
107 msgtyp = msgoff
108 rem_proc2 = mod(loc_proc+nspmd-i-1,nspmd)+1
109 msgtyp2 = msgoff
110
111 IF(rem_proc/=loc_proc.OR.rem_proc2/=loc_proc) THEN
112 CALL mpi_sendrecv(
113 s wa,siz,real,it_spmd(rem_proc),msgtyp,
114 r wa(ideb),siz,real,it_spmd(rem_proc2),msgtyp2,
115 g spmd_comm_world,status,ierror)
116 ENDIF
117 IF (rem_proc2/=loc_proc) THEN
118 DO j = 1, nsect
119 sendto(rem_proc2) = sendto(rem_proc2) +
120 + nint(wa(j)*wa(ideb+nsect+j-1))
121 recvfr(rem_proc2) = recvfr(rem_proc2) +
122 + nint(wa(nsect+j)*wa(ideb+j-1))
123 ENDDO
124 ENDIF
125 ENDDO
126 ENDIF
127C
128 DO i = 1, nspmd
129 IF(recvfr(i)>0) THEN
130 msgtyp = msgoff
131 CALL mpi_irecv(
132 s wa(1+i*bufsiz),bufsiz,real,it_spmd(i),msgtyp,
133 g spmd_comm_world,req_r(i),ierror)
134 ENDIF
135 ENDDO
136C
137 l = 0
138 k1 = 1
139 DO i=1,nsect
140 n1 = nstrf(k1+1)
141 n2 = nstrf(k1+2)
142 n3 = nstrf(k1+3)
143 IF (weight(n1)==1)THEN
144 wa(l+1) = n1
145 wa(l+2) = x(1,n1)
146 wa(l+3) = x(2,n1)
147 wa(l+4) = x(3,n1)
148 l = l + 4
149 ENDIF
150 IF (weight(n2)==1)THEN
151 wa(l+1) = n2
152 wa(l+2) = x(1,n2)
153 wa(l+3) = x(2,n2)
154 wa(l+4) = x(3,n2)
155 l = l + 4
156 ENDIF
157 IF (weight(n3)==1)THEN
158 wa(l+1) = n3
159 wa(l+2) = x(1,n3)
160 wa(l+3) = x(2,n3)
161 wa(l+4) = x(3,n3)
162 l = l + 4
163 ENDIF
164 k1=k1+4+3*nstrf(k1)
165 k1=k1+1+2*nstrf(k1)
166 END DO
167C
168 siz = l
169 DO i=1,nspmd
170 IF(sendto(i)>0) THEN
171 msgtyp = msgoff
172 CALL mpi_send(
173 s wa,siz,real,it_spmd(i),msgtyp,
174 g spmd_comm_world,ierror)
175 ENDIF
176 ENDDO
177C
178 DO i = 1, nspmd
179 IF(recvfr(i)>0) THEN
180 ideb = 1+bufsiz*i
181
182 CALL mpi_wait(req_r(i),status,ierror)
183 CALL mpi_get_count(status,real,siz,ierror)
184 nb_nod=siz/a_ar
185
186 l = ideb
187 DO j = 1, nb_nod
188 nod = nint(wa(l))
189 x(1,nod) = wa(l+1)
190 x(2,nod) = wa(l+2)
191 x(3,nod) = wa(l+3)
192 l = l + a_ar
193 ENDDO
194 ENDIF
195 ENDDO
196 ELSE
197 IF (ncycle==1) THEN
198C
199 bufsiz = 0
200 k0=nstrf(25)
201 DO i = 1, nsect
202 n1 = nstrf(k0+3)
203 nnod = nstrf(k0+6)
204 ifram = nstrf(k0+26)
205 IF (ifram<=10.OR.n1/=0) THEN
206 bufsiz = bufsiz + 3*a_ar
207 ENDIF
208 IF(mod(ifram,10)==1) THEN
209 bufsiz = bufsiz + nnod*a_ar
210 ELSEIF( mod(ifram,10)==2) THEN
211 bufsiz = bufsiz + 2*nnod*a_ar
212 ENDIF
213 k0=nstrf(k0+24)
214 ENDDO
215C
216 DO i = 1,2*nsect
217 wa(i) = zero
218 ENDDO
219C
220 k0=nstrf(25)
221 DO i=1,nsect
222 n1 = nstrf(k0+3)
223 n2 = nstrf(k0+4)
224 n3 = nstrf(k0+5)
225 nnod = nstrf(k0+6)
226 k2 = k0+30+nstrf(k0+14)
227 nelseg = nstrf(k0+7)+nstrf(k0+8)+nstrf(k0+9)+nstrf(k0+10)+
228 + nstrf(k0+11)+nstrf(k0+12)+nstrf(k0+13)
229 ifram = nstrf(k0+26)
230 IF (ifram<=10.OR.n1/=0) THEN
231 IF (weight(n1)==1.OR.weight(n2)==1.OR.weight(n3)==1)
232 + THEN
233 wa(i) = one
234 ENDIF
235 IF (nelseg/=0.AND.(weight(n1)==0.OR.weight(n2)==0
236 + .OR.weight(n3)==zero)) THEN
237 wa(i+nsect) = one
238 ENDIF
239 ENDIF
240 IF(mod(ifram,10)==1.OR.mod(ifram,10)==2) THEN
241 iall = 1
242 DO nn = 1, nnod
243 IF (weight(nstrf(k2+nn-1))==1) THEN
244 wa(i) = one
245 ELSE
246 iall = zero
247 ENDIF
248 ENDDO
249 IF (nelseg/=0.AND.iall==0) THEN
250 wa(i+nsect) = one
251 ENDIF
252 ENDIF
253 k0=nstrf(k0+24)
254 END DO
255C
256 DO i = 1, nspmd
257 sendto(i) = 0
258 recvfr(i) = 0
259 ENDDO
260 siz = 2*nsect
261 ideb = siz + 1
262 DO i = 1, nspmd
263 rem_proc = mod(loc_proc + i-1,nspmd)+1
264 msgtyp = msgoff
265 rem_proc2 = mod(loc_proc+nspmd-i-1,nspmd)+1
266 msgtyp2 = msgoff
267
268 IF(rem_proc/=loc_proc.OR.rem_proc2/=loc_proc) THEN
269 CALL mpi_sendrecv(
270 s wa,siz,real,it_spmd(rem_proc),msgtyp,
271 r wa(ideb),siz,real,it_spmd(rem_proc2),msgtyp2,
272 g spmd_comm_world,status,ierror)
273 ENDIF
274
275 IF (rem_proc2/=loc_proc) THEN
276 DO j = 1, nsect
277 sendto(rem_proc2) = sendto(rem_proc2) +
278 + nint(wa(j)*wa(ideb+nsect+j-1))
279 recvfr(rem_proc2) = recvfr(rem_proc2) +
280 + nint(wa(nsect+j)*wa(ideb+j-1))
281 ENDDO
282 ENDIF
283 ENDDO
284 ENDIF
285C
286 DO i = 1, nspmd
287 IF(recvfr(i)>0) THEN
288 msgtyp = msgoff
289 CALL mpi_irecv(
290 s wa(1+i*bufsiz),bufsiz,real,it_spmd(i),msgtyp,
291 g spmd_comm_world,req_r(i),ierror)
292 ENDIF
293 ENDDO
294C
295 l = 0
296 k0=nstrf(25)
297 DO i = 1, nsect
298 n1 = nstrf(k0+3)
299 n2 = nstrf(k0+4)
300 n3 = nstrf(k0+5)
301 nnod = nstrf(k0+6)
302 k2 = k0+30+nstrf(k0+14)
303 nelseg = nstrf(k0+7)+nstrf(k0+8)+nstrf(k0+9)+nstrf(k0+10)+
304 + nstrf(k0+11)+nstrf(k0+12)+nstrf(k0+13)
305 ifram = nstrf(k0+26)
306C
307 IF (ifram<=10.OR.n1/=0) THEN
308 IF (weight(n1)==1) THEN
309 wa(l+1) = n1
310 wa(l+2) = x(1,n1)
311 wa(l+3) = x(2,n1)
312 wa(l+4) = x(3,n1)
313 l = l + 4
314 ENDIF
315 IF (weight(n2)==1) THEN
316 wa(l+1) = n2
317 wa(l+2) = x(1,n2)
318 wa(l+3) = x(2,n2)
319 wa(l+4) = x(3,n2)
320 l = l + 4
321 ENDIF
322 IF (weight(n3)==1) THEN
323 wa(l+1) = n3
324 wa(l+2) = x(1,n3)
325 wa(l+3) = x(2,n3)
326 wa(l+4) = x(3,n3)
327 l = l + 4
328 ENDIF
329 ENDIF
330C
331 IF(mod(ifram,10)==1) THEN
332 DO nn = 1, nnod
333 n3 = nstrf(k2+nn-1)
334 IF (weight(n3)==1) THEN
335 wa(l+1) = n3
336 wa(l+2) = x(1,n3)
337 wa(l+3) = x(2,n3)
338 wa(l+4) = x(3,n3)
339 l = l + 4
340 ENDIF
341 ENDDO
342C
343 ELSEIF( mod(ifram,10)==2) THEN
344 DO nn = 1, nnod
345 n3 = nstrf(k2+nn-1)
346 IF (weight(n3)==1) THEN
347 wa(l+1) = n3
348 wa(l+2) = x(1,n3)
349 wa(l+3) = x(2,n3)
350 wa(l+4) = x(3,n3)
351 l = l + 4
352C
353 wa(l+1) = -n3
354 wa(l+2) = ms(n3)
355 wa(l+3) = ms(n3)
356 wa(l+4) = ms(n3)
357 l = l + 4
358 ENDIF
359 ENDDO
360 ENDIF
361 k0=nstrf(k0+24)
362 ENDDO
363C
364 siz = l
365 DO i=1,nspmd
366 IF(sendto(i)>0) THEN
367 msgtyp = msgoff
368
369 CALL mpi_send(
370 s wa,siz,real,it_spmd(i),msgtyp,
371 g spmd_comm_world,ierror)
372
373 ENDIF
374 ENDDO
375C
376 DO i = 1, nspmd
377 IF(recvfr(i)>0) THen
378 ideb = 1+bufsiz*i
379
380 CALL mpi_wait(req_r(i),status,ierror)
381 CALL mpi_get_count(status,real,siz,ierror)
382 nb_nod=siz/a_ar
383
384 l = ideb
385 DO j = 1, nb_nod
386 nod = nint(wa(l))
387 IF (nod>0) THEN
388 x(1,nod) = wa(l+1)
389 x(2,nod) = wa(l+2)
390 x(3,nod) = wa(l+3)
391 ELSE
392 ms(-nod) = wa(l+1)
393 ENDIF
394 l = l + a_ar
395 ENDDO
396 ENDIF
397 ENDDO
398C
399 ENDIF
400C
401#endif
402 RETURN