OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
ssol_aux.F File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine smumps_freetopso (n, keep28, iwcb, liww, w, lwc, poswcb, iwposcb, ptricb, ptracb)
subroutine smumps_compso (n, keep28, iwcb, liww, w, lwc, poswcb, iwposcb, ptricb, ptracb)
subroutine smumps_sol_x (a, nz8, n, irn, icn, z, keep, keep8, eff_size_schur, sym_perm)
subroutine smumps_scal_x (a, nz8, n, irn, icn, z, keep, keep8, colsca, eff_size_schur, sym_perm)
subroutine smumps_sol_y (a, nz8, n, irn, icn, rhs, x, r, w, keep, keep8)
subroutine smumps_sol_mulr (n, r, w)
subroutine smumps_sol_b (n, kase, x, est, w, iw, grain)
subroutine smumps_qd2 (mtype, n, nz8, aspk, irn, icn, lhs, wrhs, w, rhs, keep, keep8)
subroutine smumps_eltqd2 (mtype, n, nelt, eltptr, leltvar, eltvar, na_elt8, a_elt, lhs, wrhs, w, rhs, keep, keep8)
subroutine smumps_sol_x_elt (mtype, n, nelt, eltptr, leltvar, eltvar, na_elt8, a_elt, w, keep, keep8)
subroutine smumps_sol_scalx_elt (mtype, n, nelt, eltptr, leltvar, eltvar, na_elt8, a_elt, w, keep, keep8, colsca)
subroutine smumps_eltyd (mtype, n, nelt, eltptr, leltvar, eltvar, na_elt8, a_elt, saverhs, x, y, w, k50)
subroutine smumps_solve_get_ooc_node (inode, ptrfac, keep, a, la, step, keep8, n, must_be_permuted, ierr)
subroutine smumps_build_mapping_info (id)
subroutine smumps_sol_omega (n, rhs, x, y, r_w, c_w, iw, iflag, omega, noiter, testconv, lp, arret, grain)
subroutine smumps_sol_lcond (n, rhs, x, y, d, r_w, c_w, iw, kase, omega, erx, cond, lp, keep, keep8)
subroutine smumps_sol_cpy_fs2rhscomp (jbdeb, jbfin, nbrows, keep, rhscomp, nrhs, lrhscomp, first_row_rhscomp, w, ld_w, first_row_w)
subroutine smumps_sol_bwd_gthr (jbdeb, jbfin, j1, j2, rhscomp, nrhs, lrhscomp, w, ld_w, first_row_w, iw, liw, keep, n, posinrhscomp_bwd)
subroutine smumps_sol_q (mtype, iflag, n, lhs, wrhs, w, res, givnorm, anorm, xnorm, sclnrm, mprint, icntl, keep, keep8)
subroutine smumps_solve_fwd_trsolve (a, la, apos, npiv, ldadiag, nrhs_b, wcb, lwcb, lda_wcb, ppiv_courant, mtype, keep)
subroutine smumps_solve_bwd_trsolve (a, la, apos, npiv, ldadiag, nrhs_b, wcb, lwcb, lda_wcb, ppiv_courant, mtype, keep)
subroutine smumps_solve_fwd_panels (a, la, apos, npiv, iw, nrhs_b, wcb, lwcb, lda_wcb, ppiv_courant, mtype, keep)
subroutine smumps_solve_bwd_panels (a, la, apos, npiv, iw, nrhs_b, wcb, lwcb, lda_wcb, ppiv_courant, mtype, keep)
subroutine smumps_solve_gemm_update (a, la, apos1, nx, lda, ny, nrhs_b, wcb, lwcb, ptrx, ldx, ptry, ldy, mtype, keep, coef_y)
subroutine smumps_sol_ld_and_reload_panel (inode, n, npiv, liell, nelim, nslaves, ppiv_courant, iw, ipos, liw, a, la, apos, wcb, lwcb, ld_wcbpiv, rhscomp, lrhscomp, nrhs, posinrhscomp_fwd, jbdeb, jbfin, mtype, keep, oocwrite_compatible_with_blr, ignore_k459)
subroutine smumps_sol_ld_and_reload (inode, n, npiv, liell, nelim, nslaves, ppiv_courant, iw, ipos, liw, a, la, apos, wcb, lwcb, ld_wcbpiv, rhscomp, lrhscomp, nrhs, posinrhscomp_fwd, jbdeb, jbfin, mtype, keep, oocwrite_compatible_with_blr, ignore_k459)
subroutine smumps_set_scaling_loc (scaling_data, n, iloc, liloc, comm, myid, i_am_slave, master, nb_bytes, nb_bytes_max, k16_8, lp, lpok, icntl, info)

Function/Subroutine Documentation

◆ smumps_build_mapping_info()

subroutine smumps_build_mapping_info ( type(smumps_struc), target id)

Definition at line 774 of file ssol_aux.F.

776 IMPLICIT NONE
777 include 'mpif.h'
778 TYPE(SMUMPS_STRUC), TARGET :: id
779 INTEGER, ALLOCATABLE, DIMENSION(:) :: LOCAL_LIST
780 INTEGER :: I,IERR,TMP,NSTEPS,N_LOCAL_LIST
781 INTEGER :: MASTER,TAG_SIZE,TAG_LIST
782 INTEGER :: STATUS(MPI_STATUS_SIZE)
783 LOGICAL :: I_AM_SLAVE
784 parameter(master=0, tag_size=85,tag_list=86)
785 i_am_slave = (id%MYID .NE. master
786 & .OR. ((id%MYID.EQ.master).AND.(id%KEEP(46).EQ.1)))
787 nsteps = id%KEEP(28)
788 ALLOCATE(local_list(nsteps),stat=ierr)
789 IF(ierr.GT.0) THEN
790 WRITE(*,*)'Problem in solve: error allocating LOCAL_LIST'
791 CALL mumps_abort()
792 END IF
793 n_local_list = 0
794 IF(i_am_slave) THEN
795 DO i=1,nsteps
796 IF(id%PTLUST_S(i).NE.0) THEN
797 n_local_list = n_local_list + 1
798 local_list(n_local_list) = i
799 END IF
800 END DO
801 IF(id%MYID.NE.master) THEN
802 CALL mpi_send(n_local_list, 1,
803 & mpi_integer, master, tag_size, id%COMM,ierr)
804 CALL mpi_send(local_list, n_local_list,
805 & mpi_integer, master, tag_list, id%COMM,ierr)
806 DEALLOCATE(local_list)
807 ALLOCATE(id%IPTR_WORKING(1),
808 & id%WORKING(1),
809 & stat=ierr)
810 IF(ierr.GT.0) THEN
811 WRITE(*,*)'Problem in solve: error allocating ',
812 & 'IPTR_WORKING and WORKING'
813 CALL mumps_abort()
814 END IF
815 END IF
816 END IF
817 IF(id%MYID.EQ.master) THEN
818 ALLOCATE(id%IPTR_WORKING(id%NPROCS+1), stat=ierr)
819 IF(ierr.GT.0) THEN
820 WRITE(*,*)'Problem in solve: error allocating IPTR_WORKING'
821 CALL mumps_abort()
822 END IF
823 id%IPTR_WORKING = 0
824 id%IPTR_WORKING(1) = 1
825 id%IPTR_WORKING(master+2) = n_local_list
826 DO i=1, id%NPROCS-1
827 CALL mpi_recv(tmp, 1, mpi_integer, mpi_any_source,
828 & tag_size, id%COMM, status, ierr)
829 id%IPTR_WORKING(status(mpi_source)+2) = tmp
830 END DO
831 DO i=2, id%NPROCS+1
832 id%IPTR_WORKING(i) = id%IPTR_WORKING(i)
833 & + id%IPTR_WORKING(i-1)
834 END DO
835 ALLOCATE(id%WORKING(id%IPTR_WORKING(id%NPROCS+1)-1),stat=ierr)
836 IF(ierr.GT.0) THEN
837 WRITE(*,*)'Problem in solve: error allocating LOCAL_LIST'
838 CALL mumps_abort()
839 END IF
840 tmp = master + 1
841 IF (i_am_slave) THEN
842 id%WORKING(id%IPTR_WORKING(tmp):id%IPTR_WORKING(tmp+1)-1)
843 & = local_list(1:id%IPTR_WORKING(tmp+1)
844 & -id%IPTR_WORKING(tmp))
845 ENDIF
846 DO i=1,id%NPROCS-1
847 CALL mpi_recv(local_list, nsteps, mpi_integer,
848 & mpi_any_source, tag_list, id%COMM, status, ierr)
849 tmp = status(mpi_source)+1
850 id%WORKING(id%IPTR_WORKING(tmp):id%IPTR_WORKING(tmp+1)-1)
851 & = local_list(1:id%IPTR_WORKING(tmp+1)-
852 & id%IPTR_WORKING(tmp))
853 END DO
854 DEALLOCATE(local_list)
855 END IF
#define mumps_abort
Definition VE_Metis.h:25
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
initmumps id

◆ smumps_compso()

subroutine smumps_compso ( integer n,
integer keep28,
integer, dimension(liww) iwcb,
integer liww,
real, dimension(lwc) w,
integer(8), intent(in) lwc,
integer(8), intent(inout) poswcb,
integer iwposcb,
integer, dimension(keep28) ptricb,
integer(8), dimension(keep28) ptracb )

Definition at line 35 of file ssol_aux.F.

37 IMPLICIT NONE
38 INTEGER(8), INTENT(IN) :: LWC
39 INTEGER(8), INTENT(INOUT) :: POSWCB
40 INTEGER N,LIWW,IWPOSCB,KEEP28
41 INTEGER IWCB(LIWW),PTRICB(KEEP28)
42 INTEGER(8) :: PTRACB(KEEP28)
43 REAL W(LWC)
44 INTEGER IPTIW,SIZFI,LONGI
45 INTEGER(8) :: IPTA, LONGR, SIZFR, I8
46 INTEGER :: I
47 iptiw = iwposcb
48 ipta = poswcb
49 longi = 0
50 longr = 0_8
51 IF ( iptiw .EQ. liww ) RETURN
5210 CONTINUE
53 IF (iwcb(iptiw+2).EQ.0) THEN
54 sizfr = int(iwcb(iptiw+1),8)
55 sizfi = 2
56 IF (longi.NE.0) THEN
57 DO 20 i=0,longi-1
58 iwcb(iptiw + sizfi - i) = iwcb(iptiw - i)
59 20 CONTINUE
60 DO 30 i8=0,longr-1
61 w(ipta + sizfr - i8) = w(ipta - i8)
62 30 CONTINUE
63 ENDIF
64 DO 40 i=1,keep28
65 IF ((ptricb(i).LE.(iptiw+1)).AND.
66 & (ptricb(i).GT.iwposcb) ) THEN
67 ptricb(i) = ptricb(i) + sizfi
68 ptracb(i) = ptracb(i) + sizfr
69 ENDIF
7040 CONTINUE
71 iwposcb = iwposcb + sizfi
72 iptiw = iptiw + sizfi
73 poswcb = poswcb + sizfr
74 ipta = ipta + sizfr
75 ELSE
76 sizfr = int(iwcb(iptiw+1),8)
77 sizfi = 2
78 iptiw = iptiw + sizfi
79 longi = longi + sizfi
80 ipta = ipta + sizfr
81 longr = longr + sizfr
82 ENDIF
83 IF (iptiw.NE.liww) GOTO 10
84 RETURN

◆ smumps_eltqd2()

subroutine smumps_eltqd2 ( integer mtype,
integer n,
integer nelt,
integer, dimension(nelt+1) eltptr,
integer leltvar,
integer, dimension(leltvar) eltvar,
integer(8), intent(in) na_elt8,
real, dimension(na_elt8) a_elt,
real, dimension( n ) lhs,
real, dimension( n ) wrhs,
real, dimension(n) w,
real, dimension( n ) rhs,
integer, dimension(500) keep,
integer(8), dimension(150) keep8 )

Definition at line 506 of file ssol_aux.F.

509 IMPLICIT NONE
510 INTEGER MTYPE, N, NELT, LELTVAR
511 INTEGER(8), INTENT(IN) :: NA_ELT8
512 INTEGER ELTPTR(NELT+1), ELTVAR(LELTVAR)
513 INTEGER KEEP(500)
514 INTEGER(8) KEEP8(150)
515 REAL A_ELT(NA_ELT8)
516 REAL LHS( N ), WRHS( N ), RHS( N )
517 REAL W(N)
518 CALL smumps_mv_elt(n, nelt, eltptr, eltvar, a_elt,
519 & lhs, rhs, keep(50), mtype )
520 rhs = wrhs - rhs
521 CALL smumps_sol_x_elt( mtype, n,
522 & nelt, eltptr, leltvar, eltvar, na_elt8, a_elt,
523 & w, keep,keep8 )
524 RETURN
subroutine smumps_sol_x_elt(mtype, n, nelt, eltptr, leltvar, eltvar, na_elt8, a_elt, w, keep, keep8)
Definition ssol_aux.F:529
subroutine smumps_mv_elt(n, nelt, eltptr, eltvar, a_elt, x, y, k50, mtype)
Definition ssol_matvec.F:16

◆ smumps_eltyd()

subroutine smumps_eltyd ( integer mtype,
integer n,
integer nelt,
integer, dimension( nelt + 1 ) eltptr,
integer leltvar,
integer, dimension( leltvar ) eltvar,
integer(8) na_elt8,
real, dimension( na_elt8 ) a_elt,
real, dimension(n) saverhs,
real, dimension( n ) x,
real, dimension( n ) y,
real, dimension(n) w,
integer k50 )

Definition at line 650 of file ssol_aux.F.

653 IMPLICIT NONE
654 INTEGER N, NELT, K50, MTYPE, LELTVAR
655 INTEGER(8) :: NA_ELT8
656 INTEGER ELTPTR( NELT + 1 ), ELTVAR( LELTVAR )
657 REAL A_ELT( NA_ELT8 ), X( N ), Y( N ),
658 & SAVERHS(N)
659 REAL W(N)
660 INTEGER IEL, I , J, K, SIZEI, IELPTR
661 REAL ZERO
662 REAL TEMP
663 REAL TEMP2
664 parameter( zero = 0.0e0 )
665 y = saverhs
666 w = zero
667 k = 1
668 DO iel = 1, nelt
669 sizei = eltptr( iel + 1 ) - eltptr( iel )
670 ielptr = eltptr( iel ) - 1
671 IF ( k50 .eq. 0 ) THEN
672 IF ( mtype .eq. 1 ) THEN
673 DO j = 1, sizei
674 temp = x( eltvar( ielptr + j ) )
675 DO i = 1, sizei
676 y( eltvar( ielptr + i ) ) =
677 & y( eltvar( ielptr + i ) ) -
678 & a_elt( k ) * temp
679 w( eltvar( ielptr + i ) ) =
680 & w( eltvar( ielptr + i ) ) +
681 & abs( a_elt( k ) * temp )
682 k = k + 1
683 END DO
684 END DO
685 ELSE
686 DO j = 1, sizei
687 temp = y( eltvar( ielptr + j ) )
688 temp2 = w( eltvar( ielptr + j ) )
689 DO i = 1, sizei
690 temp = temp -
691 & a_elt( k ) * x( eltvar( ielptr + i ) )
692 temp2 = temp2 + abs(
693 & a_elt( k ) * x( eltvar( ielptr + i ) ) )
694 k = k + 1
695 END DO
696 y( eltvar( ielptr + j ) ) = temp
697 w( eltvar( ielptr + j ) ) = temp2
698 END DO
699 END IF
700 ELSE
701 DO j = 1, sizei
702 y( eltvar( ielptr + j ) ) =
703 & y( eltvar( ielptr + j ) ) -
704 & a_elt( k ) * x( eltvar( ielptr + j ) )
705 w( eltvar( ielptr + j ) ) =
706 & w( eltvar( ielptr + j ) ) + abs(
707 & a_elt( k ) * x( eltvar( ielptr + j ) ) )
708 k = k + 1
709 DO i = j+1, sizei
710 y( eltvar( ielptr + i ) ) =
711 & y( eltvar( ielptr + i ) ) -
712 & a_elt( k ) * x( eltvar( ielptr + j ) )
713 y( eltvar( ielptr + j ) ) =
714 & y( eltvar( ielptr + j ) ) -
715 & a_elt( k ) * x( eltvar( ielptr + i ) )
716 w( eltvar( ielptr + i ) ) =
717 & w( eltvar( ielptr + i ) ) + abs(
718 & a_elt( k ) * x( eltvar( ielptr + j ) ) )
719 w( eltvar( ielptr + j ) ) =
720 & w( eltvar( ielptr + j ) ) + abs(
721 & a_elt( k ) * x( eltvar( ielptr + i ) ) )
722 k = k + 1
723 END DO
724 END DO
725 END IF
726 END DO
727 RETURN

◆ smumps_freetopso()

subroutine smumps_freetopso ( integer n,
integer keep28,
integer, dimension(liww) iwcb,
integer liww,
real, dimension(lwc) w,
integer(8), intent(in) lwc,
integer(8), intent(inout) poswcb,
integer iwposcb,
integer, dimension(keep28) ptricb,
integer(8), dimension(keep28) ptracb )

Definition at line 14 of file ssol_aux.F.

17 IMPLICIT NONE
18 INTEGER(8), INTENT(IN) :: LWC
19 INTEGER(8), INTENT(INOUT) :: POSWCB
20 INTEGER N,LIWW,IWPOSCB, KEEP28
21 INTEGER IWCB(LIWW),PTRICB(KEEP28)
22 INTEGER(8) :: PTRACB(KEEP28)
23 REAL W(LWC)
24 INTEGER SIZFI, SIZFR
25 IF ( iwposcb .eq. liww ) RETURN
26 DO WHILE ( iwcb( iwposcb + 2 ) .eq. 0 )
27 sizfr = iwcb( iwposcb + 1 )
28 sizfi = 2
29 iwposcb = iwposcb + sizfi
30 poswcb = poswcb + sizfr
31 IF ( iwposcb .eq. liww ) RETURN
32 END DO
33 RETURN

◆ smumps_qd2()

subroutine smumps_qd2 ( integer mtype,
integer n,
integer(8), intent(in) nz8,
real, dimension( nz8 ), intent(in) aspk,
integer, dimension( nz8 ), intent(in) irn,
integer, dimension( nz8 ), intent(in) icn,
real, dimension( n ), intent(in) lhs,
real, dimension( n ), intent(in) wrhs,
real, dimension( n ), intent(out) w,
real, dimension( n ), intent(out) rhs,
integer, dimension(500) keep,
integer(8), dimension(150) keep8 )

Definition at line 420 of file ssol_aux.F.

422 IMPLICIT NONE
423 INTEGER MTYPE, N
424 INTEGER(8), INTENT(IN) :: NZ8
425 INTEGER, INTENT(IN) :: IRN( NZ8 ), ICN( NZ8 )
426 INTEGER KEEP(500)
427 INTEGER(8) KEEP8(150)
428 REAL, INTENT(IN) :: ASPK( NZ8 )
429 REAL, INTENT(IN) :: LHS( N ), WRHS( N )
430 REAL, INTENT(OUT):: RHS( N )
431 REAL, INTENT(OUT):: W( N )
432 INTEGER I, J
433 INTEGER(8) :: K8
434 REAL, PARAMETER :: DZERO = 0.0e0
435 DO i = 1, n
436 w(i) = dzero
437 rhs(i) = wrhs(i)
438 ENDDO
439 IF ( keep(50) .EQ. 0 ) THEN
440 IF (mtype .EQ. 1) THEN
441 IF (keep(264).EQ.0) THEN
442 DO k8 = 1_8, nz8
443 i = irn(k8)
444 j = icn(k8)
445 IF ((i .LE. 0) .OR. (i .GT. n) .OR. (j .LE. 0) .OR.
446 & (j .GT. n)) cycle
447 rhs(i) = rhs(i) - aspk(k8) * lhs(j)
448 w(i) = w(i) + abs(aspk(k8))
449 ENDDO
450 ELSE
451 DO k8 = 1_8, nz8
452 i = irn(k8)
453 j = icn(k8)
454 rhs(i) = rhs(i) - aspk(k8) * lhs(j)
455 w(i) = w(i) + abs(aspk(k8))
456 ENDDO
457 ENDIF
458 ELSE
459 IF (keep(264).EQ.0) THEN
460 DO k8 = 1_8, nz8
461 i = irn(k8)
462 j = icn(k8)
463 IF ((i .LE. 0) .OR. (i .GT. n) .OR. (j .LE. 0) .OR.
464 & (j .GT. n)) cycle
465 rhs(j) = rhs(j) - aspk(k8) * lhs(i)
466 w(j) = w(j) + abs(aspk(k8))
467 ENDDO
468 ELSE
469 DO k8 = 1_8, nz8
470 i = irn(k8)
471 j = icn(k8)
472 rhs(j) = rhs(j) - aspk(k8) * lhs(i)
473 w(j) = w(j) + abs(aspk(k8))
474 ENDDO
475 ENDIF
476 ENDIF
477 ELSE
478 IF (keep(264).EQ.0) THEN
479 DO k8 = 1_8, nz8
480 i = irn(k8)
481 j = icn(k8)
482 IF ((i .LE. 0) .OR. (i .GT. n) .OR. (j .LE. 0) .OR.
483 & (j .GT. n)) cycle
484 rhs(i) = rhs(i) - aspk(k8) * lhs(j)
485 w(i) = w(i) + abs(aspk(k8))
486 IF (j.NE.i) THEN
487 rhs(j) = rhs(j) - aspk(k8) * lhs(i)
488 w(j) = w(j) + abs(aspk(k8))
489 ENDIF
490 ENDDO
491 ELSE
492 DO k8 = 1_8, nz8
493 i = irn(k8)
494 j = icn(k8)
495 rhs(i) = rhs(i) - aspk(k8) * lhs(j)
496 w(i) = w(i) + abs(aspk(k8))
497 IF (j.NE.i) THEN
498 rhs(j) = rhs(j) - aspk(k8) * lhs(i)
499 w(j) = w(j) + abs(aspk(k8))
500 ENDIF
501 ENDDO
502 ENDIF
503 ENDIF
504 RETURN

◆ smumps_scal_x()

subroutine smumps_scal_x ( real, dimension(nz8), intent(in) a,
integer(8), intent(in) nz8,
integer, intent(in) n,
integer, dimension(nz8), intent(in) irn,
integer, dimension(nz8), intent(in) icn,
real, dimension(n), intent(out) z,
integer, dimension(500), intent(in) keep,
integer(8), dimension(150), intent(in) keep8,
real, dimension(n), intent(in) colsca,
integer, intent(in) eff_size_schur,
integer, dimension(n), intent(in) sym_perm )

Definition at line 171 of file ssol_aux.F.

174 INTEGER, INTENT(IN) :: N, KEEP(500)
175 INTEGER(8), INTENT(IN) :: NZ8
176 INTEGER(8), INTENT(IN) :: KEEP8(150)
177 INTEGER, INTENT(IN) :: IRN(NZ8), ICN(NZ8)
178 REAL, INTENT(IN) :: A(NZ8)
179 REAL, INTENT(IN) :: COLSCA(N)
180 REAL, INTENT(OUT) :: Z(N)
181 INTEGER, INTENT(IN) :: EFF_SIZE_SCHUR, SYM_PERM(N)
182 REAL, PARAMETER :: ZERO = 0.0e0
183 INTEGER :: I, J
184 INTEGER(8) :: K
185 LOGICAL :: SKIP_COLinSchur
186 DO 10 i = 1, n
187 z(i) = zero
188 10 CONTINUE
189 skip_colinschur = (eff_size_schur.GT.0)
190 IF (keep(50) .EQ.0) THEN
191 DO k = 1_8, nz8
192 i = irn(k)
193 j = icn(k)
194 IF ((i .LT. 1) .OR. (i .GT. n)) cycle
195 IF ((j .LT. 1) .OR. (j .GT. n)) cycle
196 IF ( skip_colinschur.AND.
197 & (sym_perm(j).GT.n-eff_size_schur)) cycle
198 IF ( skip_colinschur.AND.
199 & (sym_perm(i).GT.n-eff_size_schur)) cycle
200 z(i) = z(i) + abs(a(k)*colsca(j))
201 ENDDO
202 ELSE
203 DO k = 1, nz8
204 i = irn(k)
205 j = icn(k)
206 IF ((i .LT. 1) .OR. (i .GT. n)) cycle
207 IF ((j .LT. 1) .OR. (j .GT. n)) cycle
208 IF ( skip_colinschur.AND.
209 & ( (sym_perm(i).GT.n-eff_size_schur)
210 & .OR.
211 & (sym_perm(j).GT.n-eff_size_schur)
212 & )
213 & ) cycle
214 z(i) = z(i) + abs(a(k)*colsca(j))
215 IF (j.NE.i) THEN
216 z(j) = z(j) + abs(a(k)*colsca(i))
217 ENDIF
218 ENDDO
219 ENDIF
220 RETURN

◆ smumps_set_scaling_loc()

subroutine smumps_set_scaling_loc ( type (scaling_data_t), intent(inout) scaling_data,
integer, intent(in) n,
integer, dimension(liloc), intent(in) iloc,
integer, intent(in) liloc,
integer, intent(in) comm,
integer, intent(in) myid,
logical, intent(in) i_am_slave,
integer, intent(in) master,
integer(8), intent(inout) nb_bytes,
integer(8), intent(inout) nb_bytes_max,
integer(8), intent(in) k16_8,
integer, intent(in) lp,
logical, intent(in) lpok,
integer, dimension(60), intent(in) icntl,
integer, dimension(80), intent(inout) info )

Definition at line 1657 of file ssol_aux.F.

1660 IMPLICIT NONE
1661 type scaling_data_t
1662 sequence
1663 REAL, dimension(:), pointer :: SCALING
1664 REAL, dimension(:), pointer :: SCALING_LOC
1665 end type scaling_data_t
1666 type (scaling_data_t), INTENT(INOUT) :: scaling_data
1667 INTEGER, INTENT(IN) :: N, LILOC, COMM, MYID, MASTER, LP
1668 INTEGER, INTENT(IN) :: ILOC(LILOC)
1669 INTEGER(8), INTENT(INOUT) :: NB_BYTES, NB_BYTES_MAX
1670 INTEGER(8), INTENT(IN) :: K16_8
1671 LOGICAL, INTENT(IN) :: I_AM_SLAVE, LPOK
1672 INTEGER, INTENT(INOUT) :: INFO(80)
1673 INTEGER, INTENT(IN) :: ICNTL(60)
1674 REAL, POINTER, DIMENSION(:) :: SCALING
1675 INTEGER :: I, IERR_MPI, allocok
1676 include 'mpif.h'
1677 NULLIFY(scaling_data%SCALING_LOC)
1678 IF (i_am_slave) THEN
1679 ALLOCATE(scaling_data%SCALING_LOC(max(1,liloc)),
1680 & stat=allocok)
1681 IF (allocok > 0) THEN
1682 info(1)=-13
1683 info(2)=max(1,liloc)
1684 GOTO 35
1685 ENDIF
1686 nb_bytes = nb_bytes + int(max(1,liloc),8)*k16_8
1687 nb_bytes_max = max(nb_bytes_max,nb_bytes)
1688 ENDIF
1689 IF (myid .NE. master) THEN
1690 ALLOCATE(scaling(n), stat=allocok)
1691 IF (allocok > 0) THEN
1692 IF (lpok) THEN
1693 WRITE(lp,*) 'Error allocating temporary scaling array'
1694 ENDIF
1695 info(1)=-13
1696 info(2)=n
1697 GOTO 35
1698 ENDIF
1699 nb_bytes = nb_bytes + int(n,8)*k16_8
1700 nb_bytes_max = max(nb_bytes_max,nb_bytes)
1701 ELSE
1702 scaling => scaling_data%SCALING
1703 ENDIF
1704 35 CONTINUE
1705 CALL mumps_propinfo( icntl(1), info(1),
1706 & comm, myid )
1707 IF (info(1) .LT. 0) GOTO 90
1708 CALL mpi_bcast( scaling(1), n, mpi_real,
1709 & master, comm, ierr_mpi)
1710 IF ( i_am_slave ) THEN
1711 DO i = 1, liloc
1712 IF (iloc(i) .GE. 1 .AND. iloc(i) .LE. n) THEN
1713 scaling_data%SCALING_LOC(i) = scaling(iloc(i))
1714 ENDIF
1715 ENDDO
1716 ENDIF
1717 90 CONTINUE
1718 IF (myid.NE. master) THEN
1719 IF (associated(scaling)) THEN
1720 DEALLOCATE(scaling)
1721 nb_bytes = nb_bytes - int(n,8)*k16_8
1722 ENDIF
1723 ENDIF
1724 NULLIFY(scaling)
1725 IF (info(1) .LT. 0) THEN
1726 IF (associated(scaling_data%SCALING_LOC)) THEN
1727 DEALLOCATE(scaling_data%SCALING_LOC)
1728 NULLIFY(scaling_data%SCALING_LOC)
1729 ENDIF
1730 ENDIF
1731 RETURN
subroutine mumps_propinfo(icntl, info, comm, id)
#define max(a, b)
Definition macros.h:21
subroutine mpi_bcast(buffer, cnt, datatype, root, comm, ierr)
Definition mpi.f:205

◆ smumps_sol_b()

subroutine smumps_sol_b ( integer, intent(in) n,
integer, intent(inout) kase,
real, dimension(n) x,
real, intent(inout) est,
real, dimension(n) w,
integer, dimension(n) iw,
integer, intent(in) grain )

Definition at line 303 of file ssol_aux.F.

304 INTEGER, intent(in) :: N
305 INTEGER, intent(inout) :: KASE
306 INTEGER IW(N)
307 REAL W(N), X(N)
308 REAL, intent(inout) :: EST
309 INTEGER, intent(in) :: GRAIN
310 INTRINSIC abs, nint, real, sign
311 INTEGER SMUMPS_IXAMAX
312 EXTERNAL smumps_ixamax
313 INTEGER ITMAX
314 parameter(itmax = 5)
315 INTEGER I, ITER, J, JLAST, JUMP
316 REAL ALTSGN
317 REAL TEMP
318 SAVE iter, j, jlast, jump
319 REAL ZERO, ONE
320 parameter( zero = 0.0e0 )
321 parameter( one = 1.0e0 )
322 REAL, PARAMETER :: RZERO = 0.0e0
323 REAL, PARAMETER :: RONE = 1.0e0
324 IF (kase .EQ. 0) THEN
325 DO 10 i = 1, n
326 x(i) = one / real(n)
327 10 CONTINUE
328 kase = 1
329 jump = 1
330 RETURN
331 ENDIF
332 SELECT CASE (jump)
333 CASE (1)
334 GOTO 20
335 CASE(2)
336 GOTO 40
337 CASE(3)
338 GOTO 70
339 CASE(4)
340 GOTO 120
341 CASE(5)
342 GOTO 160
343 CASE DEFAULT
344 END SELECT
345 20 CONTINUE
346 IF (n .EQ. 1) THEN
347 w(1) = x(1)
348 est = abs(w(1))
349 GOTO 190
350 ENDIF
351 DO 30 i = 1, n
352 x(i) = sign( rone,real(x(i)) )
353 iw(i) = nint(real(x(i)))
354 30 CONTINUE
355 kase = 2
356 jump = 2
357 RETURN
358 40 CONTINUE
359 j = smumps_ixamax(n, x, 1, grain)
360 iter = 2
361 50 CONTINUE
362 DO 60 i = 1, n
363 x(i) = zero
364 60 CONTINUE
365 x(j) = one
366 kase = 1
367 jump = 3
368 RETURN
369 70 CONTINUE
370 DO 80 i = 1, n
371 w(i) = x(i)
372 80 CONTINUE
373 DO 90 i = 1, n
374 IF (nint(sign(rone, real(x(i)))) .NE. iw(i)) GOTO 100
375 90 CONTINUE
376 GOTO 130
377 100 CONTINUE
378 DO 110 i = 1, n
379 x(i) = sign(rone, real(x(i)))
380 iw(i) = nint(real(x(i)))
381 110 CONTINUE
382 kase = 2
383 jump = 4
384 RETURN
385 120 CONTINUE
386 jlast = j
387 j = smumps_ixamax(n, x, 1, grain)
388 IF ((abs(x(jlast)) .NE. abs(x(j))) .AND. (iter .LT. itmax)) THEN
389 iter = iter + 1
390 GOTO 50
391 ENDIF
392 130 CONTINUE
393 est = rzero
394 DO 140 i = 1, n
395 est = est + abs(w(i))
396 140 CONTINUE
397 altsgn = rone
398 DO 150 i = 1, n
399 x(i) = altsgn * (rone + real(i - 1) / real(n - 1))
400 altsgn = -altsgn
401 150 CONTINUE
402 kase = 1
403 jump = 5
404 RETURN
405 160 CONTINUE
406 temp = rzero
407 DO 170 i = 1, n
408 temp = temp + abs(x(i))
409 170 CONTINUE
410 temp = 2.0e0 * temp / real(3 * n)
411 IF (temp .GT. est) THEN
412 DO 180 i = 1, n
413 w(i) = x(i)
414 180 CONTINUE
415 est = temp
416 ENDIF
417 190 kase = 0
418 RETURN
integer function smumps_ixamax(n, x, incx, grain)

◆ smumps_sol_bwd_gthr()

subroutine smumps_sol_bwd_gthr ( integer, intent(in) jbdeb,
integer, intent(in) jbfin,
integer, intent(in) j1,
integer, intent(in) j2,
real, dimension(lrhscomp,nrhs), intent(inout) rhscomp,
integer, intent(in) nrhs,
integer, intent(in) lrhscomp,
real, dimension(ld_w*(jbfin-jbdeb+1)) w,
integer, intent(in) ld_w,
integer, intent(in) first_row_w,
integer, dimension(liw), intent(in) iw,
integer, intent(in) liw,
integer, dimension(500), intent(in) keep,
integer, intent(in) n,
integer, dimension(n), intent(in) posinrhscomp_bwd )

Definition at line 1060 of file ssol_aux.F.

1063 INTEGER, INTENT(IN) :: JBDEB, JBFIN, J1, J2
1064 INTEGER, INTENT(IN) :: NRHS, LRHSCOMP
1065 INTEGER, INTENT(IN) :: FIRST_ROW_W, LD_W, LIW
1066 INTEGER, INTENT(IN) :: IW(LIW)
1067 INTEGER, INTENT(IN) :: KEEP(500)
1068 REAL, INTENT(INOUT) :: RHSCOMP(LRHSCOMP,NRHS)
1069 REAL :: W(LD_W*(JBFIN-JBDEB+1))
1070 INTEGER, INTENT(IN) :: N
1071 INTEGER, INTENT(IN) :: POSINRHSCOMP_BWD(N)
1072 INTEGER :: ISHIFT, JJ, K, IPOSINRHSCOMP
1073!$OMP PARALLEL DO PRIVATE(JJ,ISHIFT,IPOSINRHSCOMP), IF
1074!$OMP& ((JBFIN-JBDEB+1 > 2*KEEP(362) .AND.
1075!$OMP& (JBFIN-JBDEB+1)*(J2-KEEP(253)-J1+1)>2*KEEP(363)))
1076 DO k=jbdeb, jbfin
1077 ishift = first_row_w+(k-jbdeb)*ld_w
1078 DO jj = j1, j2-keep(253)
1079 iposinrhscomp = abs(posinrhscomp_bwd(iw(jj)))
1080 w(ishift+jj-j1)= rhscomp(iposinrhscomp,k)
1081 ENDDO
1082 ENDDO
1083!$OMP END PARALLEL DO
1084 RETURN

◆ smumps_sol_cpy_fs2rhscomp()

subroutine smumps_sol_cpy_fs2rhscomp ( integer jbdeb,
integer jbfin,
integer nbrows,
integer, dimension(500), intent(in) keep,
real, dimension(lrhscomp,nrhs), intent(inout) rhscomp,
integer nrhs,
integer lrhscomp,
integer first_row_rhscomp,
real, dimension(ld_w*(jbfin-jbdeb+1)) w,
integer ld_w,
integer first_row_w )

Definition at line 1037 of file ssol_aux.F.

1040 INTEGER :: JBDEB, JBFIN, NBROWS
1041 INTEGER :: NRHS, LRHSCOMP
1042 INTEGER :: FIRST_ROW_RHSCOMP
1043 INTEGER, INTENT(IN) :: KEEP(500)
1044 REAL, INTENT(INOUT) :: RHSCOMP(LRHSCOMP,NRHS)
1045 INTEGER :: LD_W, FIRST_ROW_W
1046 REAL :: W(LD_W*(JBFIN-JBDEB+1))
1047 INTEGER :: JJ, K, ISHIFT
1048!$OMP PARALLEL DO PRIVATE(ISHIFT, JJ), IF
1049!$OMP& (JBFIN-JBDEB+1 > 2*KEEP(362) .AND.
1050!$OMP& NBROWS * (JBFIN-JBDEB+1) > 2*KEEP(363))
1051 DO k = jbdeb, jbfin
1052 ishift = first_row_w + ld_w * (k-jbdeb)
1053 DO jj = 0, nbrows-1
1054 rhscomp(first_row_rhscomp+jj,k) = w(ishift+jj)
1055 END DO
1056 END DO
1057!$OMP END PARALLEL DO
1058 RETURN

◆ smumps_sol_lcond()

subroutine smumps_sol_lcond ( integer n,
real, dimension(n) rhs,
real, dimension(n) x,
real, dimension(n) y,
real, dimension(n) d,
real, dimension(n,2) r_w,
real, dimension(n) c_w,
integer, dimension(n,2) iw,
integer kase,
real, dimension(2) omega,
real erx,
real, dimension(2) cond,
integer lp,
integer, dimension(500) keep,
integer(8), dimension(150) keep8 )

Definition at line 935 of file ssol_aux.F.

939 IMPLICIT NONE
940 INTEGER N, KASE, KEEP(500)
941 INTEGER(8) KEEP8(150)
942 INTEGER IW(N,2)
943 REAL RHS(N)
944 REAL X(N), Y(N)
945 REAL D(N)
946 REAL R_W(N,2)
947 REAL C_W(N)
948 INTEGER LP
949 REAL COND(2),OMEGA(2)
950 LOGICAL LCOND1, LCOND2
951 INTEGER JUMP, I, IMAX
952 REAL ERX, DXMAX
953 REAL DXIMAX
954 REAL, PARAMETER :: ZERO = 0.0e0
955 REAL, PARAMETER :: ONE = 1.0e0
956 INTEGER SMUMPS_IXAMAX
957 INTRINSIC abs, max
958 SAVE lcond1, lcond2, jump, dximax, dxmax
959 IF (kase .EQ. 0) THEN
960 lcond1 = .false.
961 lcond2 = .false.
962 cond(1) = one
963 cond(2) = one
964 erx = zero
965 jump = 1
966 ENDIF
967 SELECT CASE (jump)
968 CASE (1)
969 GOTO 30
970 CASE(2)
971 GOTO 10
972 CASE(3)
973 GOTO 110
974 CASE(4)
975 GOTO 150
976 CASE(5)
977 GOTO 35
978 CASE DEFAULT
979 END SELECT
980 10 CONTINUE
981 30 CONTINUE
982 35 CONTINUE
983 imax = smumps_ixamax(n, x, 1, keep(361))
984 dxmax = abs(x(imax))
985 DO i = 1, n
986 IF (iw(i, 1) .EQ. 1) THEN
987 r_w(i, 1) = r_w(i, 1) + abs(rhs(i))
988 r_w(i, 2) = zero
989 lcond1 = .true.
990 ELSE
991 r_w(i, 2) = r_w(i, 2) * dxmax + r_w(i, 1)
992 r_w(i, 1) = zero
993 lcond2 = .true.
994 ENDIF
995 ENDDO
996 DO i = 1, n
997 c_w(i) = x(i) * d(i)
998 ENDDO
999 imax = smumps_ixamax(n, c_w(1), 1, keep(361))
1000 dximax = abs(c_w(imax))
1001 IF (.NOT.lcond1) GOTO 130
1002 100 CONTINUE
1003 CALL smumps_sol_b(n, kase, y, cond(1), c_w, iw(1, 2), keep(361))
1004 IF (kase .EQ. 0) GOTO 120
1005 IF (kase .EQ. 1) CALL smumps_sol_mulr(n, y, d)
1006 IF (kase .EQ. 2) CALL smumps_sol_mulr(n, y, r_w)
1007 jump = 3
1008 RETURN
1009 110 CONTINUE
1010 IF (kase .EQ. 1) CALL smumps_sol_mulr(n, y, r_w)
1011 IF (kase .EQ. 2) CALL smumps_sol_mulr(n, y, d)
1012 GOTO 100
1013 120 CONTINUE
1014 IF (dximax .GT. zero) cond(1) = cond(1) / dximax
1015 erx = omega(1) * cond(1)
1016 130 CONTINUE
1017 IF (.NOT.lcond2) GOTO 170
1018 kase = 0
1019 140 CONTINUE
1020 CALL smumps_sol_b(n, kase, y, cond(2), c_w, iw(1, 2), keep(361))
1021 IF (kase .EQ. 0) GOTO 160
1022 IF (kase .EQ. 1) CALL smumps_sol_mulr(n, y, d)
1023 IF (kase .EQ. 2) CALL smumps_sol_mulr(n, y, r_w(1, 2))
1024 jump = 4
1025 RETURN
1026 150 CONTINUE
1027 IF (kase .EQ. 1) CALL smumps_sol_mulr(n, y, r_w(1, 2))
1028 IF (kase .EQ. 2) CALL smumps_sol_mulr(n, y, d)
1029 GOTO 140
1030 160 IF (dximax .GT. zero) THEN
1031 cond(2) = cond(2) / dximax
1032 ENDIF
1033 erx = erx + omega(2) * cond(2)
1034 170 CONTINUE
1035 RETURN
subroutine smumps_sol_mulr(n, r, w)
Definition ssol_aux.F:294
subroutine smumps_sol_b(n, kase, x, est, w, iw, grain)
Definition ssol_aux.F:304

◆ smumps_sol_ld_and_reload()

subroutine smumps_sol_ld_and_reload ( integer, intent(in) inode,
integer, intent(in) n,
integer, intent(in) npiv,
integer, intent(in) liell,
integer, intent(in) nelim,
integer, intent(in) nslaves,
integer(8), intent(in) ppiv_courant,
integer, dimension(liw), intent(in) iw,
integer, intent(in) ipos,
integer, intent(in) liw,
real, dimension( la ), intent(in) a,
integer(8), intent(in) la,
integer(8), intent(in) apos,
real, dimension( lwcb ), intent(in) wcb,
integer(8), intent(in) lwcb,
integer, intent(in) ld_wcbpiv,
real, dimension(lrhscomp, nrhs), intent(inout) rhscomp,
integer, intent(in) lrhscomp,
integer, intent(in) nrhs,
integer, dimension(n), intent(in) posinrhscomp_fwd,
integer, intent(in) jbdeb,
integer, intent(in) jbfin,
integer, intent(in) mtype,
integer, dimension(500), intent(in) keep,
logical, intent(in) oocwrite_compatible_with_blr,
logical, intent(in) ignore_k459 )

Definition at line 1498 of file ssol_aux.F.

1509 USE smumps_ooc
1510 INTEGER, INTENT(IN) :: MTYPE, INODE, N, NPIV, LIELL,
1511 & NELIM, NSLAVES
1512 INTEGER, INTENT(IN) :: LRHSCOMP, NRHS, LIW, JBDEB, JBFIN
1513 INTEGER, INTENT(IN) :: IW(LIW), IPOS, POSINRHSCOMP_FWD(N)
1514 INTEGER(8), INTENT(IN) :: LWCB, APOS, LA, PPIV_COURANT
1515 INTEGER, INTENT(IN) :: LD_WCBPIV
1516 INTEGER, INTENT(IN) :: KEEP(500)
1517 REAL, INTENT(IN) :: WCB( LWCB )
1518 REAL, INTENT(IN) :: A( LA )
1519 REAL, INTENT(INOUT) :: RHSCOMP(LRHSCOMP, NRHS)
1520 LOGICAL, INTENT(IN) :: OOCWRITE_COMPATIBLE_WITH_BLR
1521 LOGICAL, INTENT(IN) :: IGNORE_K459
1522 INTEGER :: TempNROW, J1, J3, PANEL_SIZE, TYPEF
1523 INTEGER :: IPOSINRHSCOMP, JJ, K, NBK, LDAJ,
1524 & LDAJ_ini, NBK_ini, LDAJ_FIRST_PANEL, NRHS_B
1525 INTEGER(8) :: IFR8 , APOS1, APOS2, APOSOFF, IFR_ini8,
1526 & POSWCB1, POSWCB2
1527 REAL :: VALPIV, A11, A22, A12, DETPIV
1528!$ LOGICAL :: OMP_FLAG
1529 REAL ONE
1530 parameter(one = 1.0e0)
1531 nrhs_b = jbfin-jbdeb+1
1532 IF ( mtype .EQ. 1 .OR. keep(50) .NE. 0 ) THEN
1533 j1 = ipos + 1
1534 j3 = ipos + npiv
1535 ELSE
1536 j1 = ipos + liell + 1
1537 j3 = ipos + liell + npiv
1538 END IF
1539 iposinrhscomp = posinrhscomp_fwd(iw(j1))
1540 IF ( keep(50) .eq. 0 ) THEN
1541!$ OMP_FLAG=(NRHS_B.GE.KEEP(362).AND.NRHS_B*NPIV.GE.KEEP(363))
1542!$OMP PARALLEL DO PRIVATE(IFR8) IF (OMP_FLAG)
1543 DO k=jbdeb,jbfin
1544 ifr8 = ppiv_courant + (k-jbdeb)*ld_wcbpiv
1545 rhscomp(iposinrhscomp:iposinrhscomp+npiv-1, k) =
1546 & wcb(ifr8:ifr8+int(npiv-1,8))
1547 ENDDO
1548!$OMP END PARALLEL DO
1549 ELSE
1550 ifr8 = ppiv_courant - 1_8
1551 IF (keep(201).EQ.1.AND.oocwrite_compatible_with_blr) THEN
1552 IF (mtype.EQ.1) THEN
1553 IF ((mtype.EQ.1).AND.nslaves.NE.0) THEN
1554 tempnrow= npiv+nelim
1555 ldaj_first_panel=tempnrow
1556 ELSE
1557 tempnrow= liell
1558 ldaj_first_panel=tempnrow
1559 ENDIF
1560 typef=typef_l
1561 ELSE
1562 tempnrow= npiv
1563 ldaj_first_panel=liell
1564 typef= typef_u
1565 ENDIF
1566 panel_size = smumps_ooc_panel_size( ldaj_first_panel )
1567 ldaj = tempnrow
1568 ELSE
1569 IF ( keep(459) .GT. 1 .AND. keep(50) .NE. 0
1570 & .AND. .NOT. ignore_k459 ) THEN
1571 CALL mumps_ldltpanel_nbtarget( npiv, panel_size, keep )
1572 ldaj = panel_size
1573 ELSE
1574 panel_size = -1
1575 ldaj = npiv
1576 ENDIF
1577 ENDIF
1578 IF (keep(201).EQ.1.AND.oocwrite_compatible_with_blr) THEN
1579 nbk = 0
1580 ENDIF
1581 ifr_ini8 = ppiv_courant - 1_8
1582 ldaj_ini = ldaj
1583 IF (keep(201).EQ.1.AND.oocwrite_compatible_with_blr)
1584 & nbk_ini = nbk
1585!$ OMP_FLAG = ( JBFIN-JBDEB+1.GE.KEEP(362) .AND.
1586!$ & ((J3-J1+1)*(JBFIN-JBDEB+1) .GE. KEEP(363)))
1587!$OMP PARALLEL DO PRIVATE(JJ,IFR8,NBK,APOS1,APOS2,APOSOFF,VALPIV,
1588!$OMP& POSWCB1, POSWCB2,A11,A22,A12,DETPIV,LDAJ) IF(OMP_FLAG)
1589 DO k = jbdeb, jbfin
1590 ifr8 = ifr_ini8 + int(k-jbdeb,8)*int(ld_wcbpiv,8)
1591 nbk = nbk_ini
1592 apos1 = apos
1593 ldaj = ldaj_ini
1594 jj = j1
1595 DO
1596 IF (jj .GT. j3) EXIT
1597 ifr8 = ifr8 + 1_8
1598 IF (iw(jj+liell) .GT. 0) THEN
1599 valpiv = one/a( apos1 )
1600 rhscomp(iposinrhscomp+jj-j1 , k ) =
1601 & wcb( ifr8 ) * valpiv
1602 IF (keep(201).EQ.1.AND.oocwrite_compatible_with_blr)
1603 & THEN
1604 nbk = nbk+1
1605 IF (nbk.EQ.panel_size) THEN
1606 nbk = 0
1607 ldaj = ldaj - panel_size
1608 ENDIF
1609 ENDIF
1610 apos1 = apos1 + int(ldaj + 1,8)
1611 jj = jj+1
1612 ELSE
1613 IF (keep(201).EQ.1.AND.oocwrite_compatible_with_blr)
1614 & THEN
1615 nbk = nbk+1
1616 ENDIF
1617 apos2 = apos1+int(ldaj+1,8)
1618 IF (keep(201).EQ.1.AND.oocwrite_compatible_with_blr)
1619 & THEN
1620 aposoff = apos1+int(ldaj,8)
1621 ELSE
1622 aposoff=apos1+1_8
1623 ENDIF
1624 a11 = a(apos1)
1625 a22 = a(apos2)
1626 a12 = a(aposoff)
1627 detpiv = a11*a22 - a12**2
1628 a22 = a11/detpiv
1629 a11 = a(apos2)/detpiv
1630 a12 = -a12/detpiv
1631 poswcb1 = ifr8
1632 poswcb2 = poswcb1+1_8
1633 rhscomp(iposinrhscomp+jj-j1,k) =
1634 & wcb(poswcb1)*a11
1635 & + wcb(poswcb2)*a12
1636 rhscomp(iposinrhscomp+jj-j1+1,k) =
1637 & wcb(poswcb1)*a12
1638 & + wcb(poswcb2)*a22
1639 IF (keep(201).EQ.1.AND.oocwrite_compatible_with_blr)
1640 & THEN
1641 nbk = nbk+1
1642 IF (nbk.GE.panel_size) THEN
1643 ldaj = ldaj - nbk
1644 nbk = 0
1645 ENDIF
1646 ENDIF
1647 apos1 = apos2 + int(ldaj + 1,8)
1648 jj = jj+2
1649 ifr8 = ifr8+1_8
1650 ENDIF
1651 ENDDO
1652 ENDDO
1653!$OMP END PARALLEL DO
1654 END IF
1655 RETURN
integer function, public smumps_ooc_panel_size(nnmax)
subroutine mumps_ldltpanel_nbtarget(npiv, nb_target, keep)

◆ smumps_sol_ld_and_reload_panel()

subroutine smumps_sol_ld_and_reload_panel ( integer, intent(in) inode,
integer, intent(in) n,
integer, intent(in) npiv,
integer, intent(in) liell,
integer, intent(in) nelim,
integer, intent(in) nslaves,
integer(8), intent(in) ppiv_courant,
integer, dimension(liw), intent(in) iw,
integer, intent(in) ipos,
integer, intent(in) liw,
real, dimension( la ), intent(in) a,
integer(8), intent(in) la,
integer(8), intent(in) apos,
real, dimension( lwcb ), intent(in) wcb,
integer(8), intent(in) lwcb,
integer, intent(in) ld_wcbpiv,
real, dimension(lrhscomp, nrhs), intent(inout) rhscomp,
integer, intent(in) lrhscomp,
integer, intent(in) nrhs,
integer, dimension(n), intent(in) posinrhscomp_fwd,
integer, intent(in) jbdeb,
integer, intent(in) jbfin,
integer, intent(in) mtype,
integer, dimension(500), intent(in) keep,
logical, intent(in) oocwrite_compatible_with_blr,
logical, intent(in) ignore_k459 )

Definition at line 1368 of file ssol_aux.F.

1379 USE smumps_ooc
1380 IMPLICIT NONE
1381 INTEGER, INTENT(IN) :: MTYPE, INODE, N, NPIV, LIELL,
1382 & NELIM, NSLAVES
1383 INTEGER, INTENT(IN) :: LRHSCOMP, NRHS, LIW, JBDEB, JBFIN
1384 INTEGER, INTENT(IN) :: IW(LIW), IPOS, POSINRHSCOMP_FWD(N)
1385 INTEGER(8), INTENT(IN) :: LWCB, APOS, LA, PPIV_COURANT
1386 INTEGER, INTENT(IN) :: LD_WCBPIV
1387 INTEGER, INTENT(IN) :: KEEP(500)
1388 REAL, INTENT(IN) :: WCB( LWCB )
1389 REAL, INTENT(IN) :: A( LA )
1390 REAL, INTENT(INOUT) :: RHSCOMP(LRHSCOMP, NRHS)
1391 LOGICAL, INTENT(IN) :: OOCWRITE_COMPATIBLE_WITH_BLR
1392 LOGICAL, INTENT(IN) :: IGNORE_K459
1393 INTEGER :: J1, J3
1394 INTEGER :: IPOSINRHSCOMP, JJ, K, NBK,
1395 & LDAJ, NRHS_B
1396 INTEGER(8) :: IFR8 , APOS1, APOS2, APOSOFF, IFR_ini8,
1397 & POSWCB1, POSWCB2
1398 REAL :: VALPIV, A11, A22, A12, DETPIV
1399 INTEGER, PARAMETER :: PANEL_TABSIZE = 20
1400 INTEGER(8) :: PANEL_POS(PANEL_TABSIZE)
1401 INTEGER :: PANEL_COL(PANEL_TABSIZE)
1402 INTEGER :: IPANEL, ICOL, NBPANELS, NB_TARGET
1403 LOGICAL :: SKIP_IT
1404 LOGICAL :: OMP_FLAG
1405 REAL ONE
1406 parameter(one = 1.0e0)
1407 IF ( npiv.EQ. 0 ) RETURN
1408 nrhs_b = jbfin-jbdeb+1
1409 IF ( mtype .EQ. 1 .OR. keep(50) .NE. 0 ) THEN
1410 j1 = ipos + 1
1411 j3 = ipos + npiv
1412 ELSE
1413 j1 = ipos + liell + 1
1414 j3 = ipos + liell + npiv
1415 END IF
1416 iposinrhscomp = posinrhscomp_fwd(iw(j1))
1417 IF ( keep(50) .eq. 0 ) THEN
1418 omp_flag = .false.
1419!$ OMP_FLAG=(int(NRHS_B,8)*int(NPIV,8).GE.int(KEEP(363),8))
1420 IF (omp_flag) THEN
1421!$OMP PARALLEL DO PRIVATE(IFR8) COLLAPSE(2)
1422 DO k = jbdeb, jbfin
1423 DO ifr8 = 0_8, int(npiv-1,8)
1424 rhscomp(iposinrhscomp+ifr8, k) =
1425 & wcb(ppiv_courant+(k-jbdeb)*ld_wcbpiv+ifr8)
1426 ENDDO
1427 ENDDO
1428!$OMP END PARALLEL DO
1429 ELSE
1430 DO k = jbdeb, jbfin
1431 DO ifr8 = 0_8, int(npiv-1,8)
1432 rhscomp(iposinrhscomp+ifr8, k) =
1433 & wcb(ppiv_courant+(k-jbdeb)*ld_wcbpiv+ifr8)
1434 ENDDO
1435 ENDDO
1436 ENDIF
1437 ELSE
1438 CALL mumps_ldltpanel_panelinfos( npiv, keep, iw(ipos+liell+1),
1439 & nb_target, nbpanels, panel_col, panel_pos, panel_tabsize,
1440 & ignore_k459 )
1441 ifr_ini8 = ppiv_courant
1442!$ OMP_FLAG = ( JBFIN-JBDEB+1.GE.KEEP(362) .AND.
1443!$ & ((J3-J1+1)*(JBFIN-JBDEB+1) .GE. KEEP(363)))
1444!$OMP PARALLEL DO PRIVATE(JJ,IFR8,NBK,APOS1,APOS2,APOSOFF,VALPIV,
1445!$OMP& IPANEL,ICOL,
1446!$OMP& POSWCB1,POSWCB2,A11,A22,A12,DETPIV,LDAJ,SKIP_IT)
1447!$OMP& IF(OMP_FLAG)
1448 DO k = jbdeb, jbfin
1449 DO jj = j1, j3
1450 ipanel = (jj-j1)/nb_target + 1
1451 IF ( jj-j1+1 .LT. panel_col(ipanel) ) ipanel = ipanel -1
1452 icol = jj-j1+1 - panel_col(ipanel) + 1
1453 ldaj = panel_col(ipanel+1) - panel_col(ipanel)
1454 apos1 = apos-1_8+panel_pos( ipanel ) + int(icol-1,8) *
1455 & int(ldaj+1,8)
1456 ifr8 = ifr_ini8 + int(k-jbdeb,8)*int(ld_wcbpiv,8) +
1457 & int(jj-j1,8)
1458 IF ( jj .NE. j1 ) THEN
1459 IF ( iw(liell+jj-1) .LT. 0 ) THEN
1460 skip_it = .true.
1461 ELSE
1462 skip_it = .false.
1463 ENDIF
1464 ELSE
1465 skip_it = .false.
1466 ENDIF
1467 IF (skip_it) THEN
1468 ELSE IF ( iw(jj+liell) .GT. 0 ) THEN
1469 valpiv = one/a( apos1 )
1470 rhscomp(iposinrhscomp+jj-j1 , k ) =
1471 & wcb( ifr8 ) * valpiv
1472 apos1 = apos1 + int(ldaj + 1,8)
1473 ELSE
1474 apos2 = apos1+int(ldaj+1,8)
1475 aposoff=apos1+1_8
1476 a11 = a(apos1)
1477 a22 = a(apos2)
1478 a12 = a(aposoff)
1479 detpiv = a11*a22 - a12**2
1480 a22 = a11/detpiv
1481 a11 = a(apos2)/detpiv
1482 a12 = -a12/detpiv
1483 poswcb1 = ifr8
1484 poswcb2 = poswcb1+1_8
1485 rhscomp(iposinrhscomp+jj-j1,k) =
1486 & wcb(poswcb1)*a11
1487 & + wcb(poswcb2)*a12
1488 rhscomp(iposinrhscomp+jj-j1+1,k) =
1489 & wcb(poswcb1)*a12
1490 & + wcb(poswcb2)*a22
1491 ENDIF
1492 ENDDO
1493 ENDDO
1494!$OMP END PARALLEL DO
1495 END IF
1496 RETURN
subroutine mumps_ldltpanel_panelinfos(npiv, keep, iw, nb_target, nbpanels, panel_col, panel_pos, panel_tabsize, ignore_k459)

◆ smumps_sol_mulr()

subroutine smumps_sol_mulr ( integer, intent(in) n,
real, dimension(n), intent(inout) r,
real, dimension(n), intent(in) w )

Definition at line 293 of file ssol_aux.F.

294 INTEGER, intent(in) :: N
295 REAL, intent(in) :: W(N)
296 REAL, intent(inout) :: R(N)
297 INTEGER I
298 DO 10 i = 1, n
299 r(i) = r(i) * w(i)
300 10 CONTINUE
301 RETURN

◆ smumps_sol_omega()

subroutine smumps_sol_omega ( integer n,
real, dimension(n) rhs,
real, dimension(n) x,
real, dimension(n) y,
real, dimension(n,2) r_w,
real, dimension(n) c_w,
integer, dimension(n,2) iw,
integer iflag,
real, dimension(2) omega,
integer noiter,
logical testconv,
integer lp,
real arret,
integer, intent(in) grain )

Definition at line 857 of file ssol_aux.F.

861 IMPLICIT NONE
862 INTEGER N, IFLAG
863 INTEGER IW(N,2)
864 REAL RHS(N)
865 REAL X(N), Y(N)
866 REAL R_W(N,2)
867 REAL C_W(N)
868 INTEGER LP, NOITER
869 LOGICAL TESTConv
870 REAL OMEGA(2)
871 REAL ARRET
872 INTEGER, intent(in) :: GRAIN
873 REAL, PARAMETER :: CGCE=0.2e0
874 REAL, PARAMETER :: CTAU=1.0e3
875 INTEGER I, IMAX
876 REAL OM1, OM2, DXMAX
877 REAL TAU, DD
878 REAL OLDOMG(2)
879 REAL, PARAMETER :: ZERO=0.0e0
880 REAL, PARAMETER :: ONE=1.0e0
881 INTEGER SMUMPS_IXAMAX
882 INTRINSIC abs, max
883 SAVE om1, oldomg
884 imax = smumps_ixamax(n, x, 1, grain)
885 dxmax = abs(x(imax))
886 omega(1) = zero
887 omega(2) = zero
888 DO i = 1, n
889 tau = (r_w(i, 2) * dxmax + abs(rhs(i))) * real(n) * ctau
890 dd = r_w(i, 1) + abs(rhs(i))
891 IF (dd .GT. tau * epsilon(ctau)) THEN
892 omega(1) = max(omega(1), abs(y(i)) / dd)
893 iw(i, 1) = 1
894 ELSE
895 IF (tau .GT. zero) THEN
896 omega(2) = max(omega(2),
897 & abs(y(i)) / (dd + r_w(i, 2) * dxmax))
898 ENDIF
899 iw(i, 1) = 2
900 ENDIF
901 ENDDO
902 IF (testconv) THEN
903 om2 = omega(1) + omega(2)
904 IF (om2 .LT. arret ) THEN
905 iflag = 1
906 GOTO 70
907 ENDIF
908 IF (noiter .GE. 1) THEN
909 IF (om2 .GT. om1 * cgce) THEN
910 IF (om2 .GT. om1) THEN
911 omega(1) = oldomg(1)
912 omega(2) = oldomg(2)
913 DO i = 1, n
914 x(i) = c_w(i)
915 ENDDO
916 iflag = 2
917 GOTO 70
918 ENDIF
919 iflag = 3
920 GOTO 70
921 ENDIF
922 ENDIF
923 DO i = 1, n
924 c_w(i) = x(i)
925 ENDDO
926 oldomg(1) = omega(1)
927 oldomg(2) = omega(2)
928 om1 = om2
929 ENDIF
930 iflag = 0
931 RETURN
932 70 CONTINUE
933 RETURN
subroutine arret(nn)
Definition arret.F:87

◆ smumps_sol_q()

subroutine smumps_sol_q ( integer mtype,
integer iflag,
integer n,
real, dimension(n) lhs,
real, dimension(n) wrhs,
real, dimension(n) w,
real, dimension(n) res,
logical givnorm,
real anorm,
real xnorm,
real sclnrm,
integer mprint,
integer, dimension(60) icntl,
integer, dimension(500) keep,
integer(8), dimension(150) keep8 )

Definition at line 1086 of file ssol_aux.F.

1089 INTEGER MTYPE,N,IFLAG,ICNTL(60), KEEP(500)
1090 INTEGER(8) KEEP8(150)
1091 REAL RES(N),LHS(N)
1092 REAL WRHS(N)
1093 REAL W(N)
1094 REAL RESMAX,RESL2,XNORM, SCLNRM
1095 REAL ANORM,DZERO
1096 LOGICAL GIVNORM,PROK
1097 INTEGER MPRINT, MP
1098 INTEGER K
1099 INTRINSIC abs, max, sqrt
1100 mp = icntl(2)
1101 prok = (mprint .GT. 0)
1102 dzero = 0.0e0
1103 IF (.NOT.givnorm) anorm = dzero
1104 resmax = dzero
1105 resl2 = dzero
1106 DO 40 k = 1, n
1107 resmax = max(resmax, abs(res(k)))
1108 resl2 = resl2 + abs(res(k)) * abs(res(k))
1109 IF (.NOT.givnorm) anorm = max(anorm, w(k))
1110 40 CONTINUE
1111 xnorm = dzero
1112 DO 50 k = 1, n
1113 xnorm = max(xnorm, abs(lhs(k)))
1114 50 CONTINUE
1115 IF ( xnorm .EQ. dzero .OR. (exponent(xnorm) .LT.
1116 & minexponent(xnorm) + keep(122) )
1117 & .OR.
1118 & ( exponent(anorm)+exponent(xnorm) .LT.
1119 & minexponent(xnorm) + keep(122) )
1120 & .OR.
1121 & ( exponent(anorm) + exponent(xnorm) -exponent(resmax)
1122 & .LT. minexponent(xnorm) + keep(122) )
1123 & ) THEN
1124 IF (mod(iflag/2,2) .EQ. 0) THEN
1125 iflag = iflag + 2
1126 ENDIF
1127 IF ((mp .GT. 0) .AND. (icntl(4) .GE. 2)) WRITE( mp, * )
1128 & ' max-NORM of computed solut. is zero or close to zero. '
1129 ENDIF
1130 IF (resmax .EQ. dzero) THEN
1131 sclnrm = dzero
1132 ELSE
1133 sclnrm = resmax / (anorm * xnorm)
1134 ENDIF
1135 resl2 = sqrt(resl2)
1136 IF (prok) WRITE( mprint, 90 ) resmax, resl2, anorm, xnorm,
1137 & sclnrm
1138 90 FORMAT (/' RESIDUAL IS ............ (MAX-NORM) =',1pd9.2/
1139 & ' .. (2-NORM) =',1pd9.2/
1140 & ' RINFOG(4):NORM OF input Matrix (MAX-NORM)=',1pd9.2/
1141 & ' RINFOG(5):NORM OF Computed SOLUT (MAX-NORM)=',1pd9.2/
1142 & ' RINFOG(6):SCALED RESIDUAL ...... (MAX-NORM)=',1pd9.2)
1143 RETURN

◆ smumps_sol_scalx_elt()

subroutine smumps_sol_scalx_elt ( integer mtype,
integer n,
integer nelt,
integer, dimension(nelt+1) eltptr,
integer leltvar,
integer, dimension(leltvar) eltvar,
integer(8), intent(in) na_elt8,
real, dimension(na_elt8) a_elt,
real, dimension(n) w,
integer, dimension(500) keep,
integer(8), dimension(150) keep8,
real, dimension(n) colsca )

Definition at line 585 of file ssol_aux.F.

588 IMPLICIT NONE
589 INTEGER MTYPE, N, NELT, LELTVAR
590 INTEGER(8), INTENT(IN) :: NA_ELT8
591 INTEGER ELTPTR(NELT+1), ELTVAR(LELTVAR)
592 INTEGER KEEP(500)
593 INTEGER(8) KEEP8(150)
594 REAL COLSCA(N)
595 REAL A_ELT(NA_ELT8)
596 REAL W(N)
597 REAL TEMP, TEMP2
598 INTEGER I, J, IEL, SIZEI, IELPTR
599 INTEGER(8) :: K8
600 REAL DZERO
601 parameter(dzero = 0.0e0)
602 w = dzero
603 k8 = 1_8
604 DO iel = 1, nelt
605 sizei = eltptr( iel + 1 ) - eltptr( iel )
606 ielptr = eltptr( iel ) - 1
607 IF ( keep(50).EQ.0 ) THEN
608 IF (mtype.EQ.1) THEN
609 DO j = 1, sizei
610 temp2 = abs(colsca(eltvar( ielptr + j) ))
611 DO i = 1, sizei
612 w( eltvar( ielptr + i) ) =
613 & w( eltvar( ielptr + i) )
614 & + abs(a_elt( k8 )) * temp2
615 k8 = k8 + 1_8
616 END DO
617 END DO
618 ELSE
619 DO j = 1, sizei
620 temp = w( eltvar( ielptr + j ) )
621 temp2= abs(colsca(eltvar( ielptr + j) ))
622 DO i = 1, sizei
623 temp = temp + abs(a_elt( k8 )) * temp2
624 k8 = k8 + 1_8
625 END DO
626 w(eltvar( ielptr + j )) =
627 & w(eltvar( ielptr + j )) + temp
628 END DO
629 ENDIF
630 ELSE
631 DO j = 1, sizei
632 w(eltvar( ielptr + j )) =
633 & w(eltvar( ielptr + j )) +
634 & abs( a_elt( k8 )*colsca(eltvar( ielptr + j)) )
635 k8 = k8 + 1_8
636 DO i = j+1, sizei
637 w(eltvar( ielptr + j )) =
638 & w(eltvar( ielptr + j )) +
639 & abs(a_elt( k8 )*colsca(eltvar( ielptr + j)))
640 w(eltvar( ielptr + i ) ) =
641 & w(eltvar( ielptr + i )) +
642 & abs(a_elt( k8 )*colsca(eltvar( ielptr + i)))
643 k8 = k8 + 1_8
644 END DO
645 ENDDO
646 ENDIF
647 ENDDO
648 RETURN

◆ smumps_sol_x()

subroutine smumps_sol_x ( real, dimension(nz8), intent(in) a,
integer(8), intent(in) nz8,
integer, intent(in) n,
integer, dimension(nz8), intent(in) irn,
integer, dimension(nz8), intent(in) icn,
real, dimension(n), intent(out) z,
integer, dimension(500), intent(in) keep,
integer(8), dimension(150), intent(in) keep8,
integer, intent(in) eff_size_schur,
integer, dimension(n), intent(in) sym_perm )

Definition at line 86 of file ssol_aux.F.

88 INTEGER, INTENT(IN) :: N, KEEP(500)
89 INTEGER(8), INTENT(IN) :: NZ8
90 INTEGER(8), INTENT(IN) :: KEEP8(150)
91 INTEGER, INTENT(IN) :: IRN(NZ8), ICN(NZ8)
92 REAL, INTENT(IN) :: A(NZ8)
93 REAL, INTENT(OUT) :: Z(N)
94 INTEGER, INTENT(IN) :: EFF_SIZE_SCHUR, SYM_PERM(N)
95 INTEGER :: I, J
96 LOGICAL :: SKIP_COLinSchur
97 REAL, PARAMETER :: ZERO = 0.0e0
98 INTEGER(8) :: K
99 INTRINSIC abs
100 DO 10 i = 1, n
101 z(i) = zero
102 10 CONTINUE
103 skip_colinschur = (eff_size_schur.GT.0)
104 IF (keep(264).EQ.0) THEN
105 IF (keep(50) .EQ.0) THEN
106 DO k = 1_8, nz8
107 i = irn(k)
108 j = icn(k)
109 IF ((i .LT. 1) .OR. (i .GT. n)) cycle
110 IF ((j .LT. 1) .OR. (j .GT. n)) cycle
111 IF ( skip_colinschur.AND.
112 & (sym_perm(j).GT.n-eff_size_schur)) cycle
113 IF ( skip_colinschur.AND.
114 & (sym_perm(i).GT.n-eff_size_schur)) cycle
115 z(i) = z(i) + abs(a(k))
116 ENDDO
117 ELSE
118 DO k = 1_8, nz8
119 i = irn(k)
120 j = icn(k)
121 IF ((i .LT. 1) .OR. (i .GT. n)) cycle
122 IF ((j .LT. 1) .OR. (j .GT. n)) cycle
123 IF ( skip_colinschur.AND.
124 & ( (sym_perm(i).GT.n-eff_size_schur)
125 & .OR.
126 & (sym_perm(j).GT.n-eff_size_schur)
127 & )
128 & ) cycle
129 z(i) = z(i) + abs(a(k))
130 IF (j.NE.i) THEN
131 z(j) = z(j) + abs(a(k))
132 ENDIF
133 ENDDO
134 ENDIF
135 ELSE
136 IF (keep(50) .EQ.0) THEN
137 IF (skip_colinschur) THEN
138 DO k = 1_8, nz8
139 j = icn(k)
140 IF ( sym_perm(j).GT.n-eff_size_schur ) cycle
141 i = irn(k)
142 IF ( sym_perm(i).GT.n-eff_size_schur ) cycle
143 z(i) = z(i) + abs(a(k))
144 ENDDO
145 ELSE
146 DO k = 1_8, nz8
147 i = irn(k)
148 j = icn(k)
149 z(i) = z(i) + abs(a(k))
150 ENDDO
151 ENDIF
152 ELSE
153 DO k = 1_8, nz8
154 i = irn(k)
155 j = icn(k)
156 IF ( skip_colinschur.AND.
157 & ( (sym_perm(i).GT.n-eff_size_schur)
158 & .OR.
159 & (sym_perm(j).GT.n-eff_size_schur)
160 & )
161 & ) cycle
162 z(i) = z(i) + abs(a(k))
163 IF (j.NE.i) THEN
164 z(j) = z(j) + abs(a(k))
165 ENDIF
166 ENDDO
167 ENDIF
168 ENDIF
169 RETURN

◆ smumps_sol_x_elt()

subroutine smumps_sol_x_elt ( integer mtype,
integer n,
integer nelt,
integer, dimension(nelt+1) eltptr,
integer leltvar,
integer, dimension(leltvar) eltvar,
integer(8), intent(in) na_elt8,
real, dimension(na_elt8) a_elt,
real, dimension(n) w,
integer, dimension(500) keep,
integer(8), dimension(150) keep8 )

Definition at line 526 of file ssol_aux.F.

529 IMPLICIT NONE
530 INTEGER MTYPE, N, NELT, LELTVAR
531 INTEGER(8), INTENT(IN) :: NA_ELT8
532 INTEGER ELTPTR(NELT+1), ELTVAR(LELTVAR)
533 INTEGER KEEP(500)
534 INTEGER(8) KEEP8(150)
535 REAL A_ELT(NA_ELT8)
536 REAL TEMP
537 REAL W(N)
538 INTEGER I, J, IEL, SIZEI, IELPTR
539 INTEGER(8) :: K8
540 REAL DZERO
541 parameter(dzero = 0.0e0)
542 w = dzero
543 k8 = 1_8
544 DO iel = 1, nelt
545 sizei = eltptr( iel + 1 ) - eltptr( iel )
546 ielptr = eltptr( iel ) - 1
547 IF ( keep(50).EQ.0 ) THEN
548 IF (mtype.EQ.1) THEN
549 DO j = 1, sizei
550 DO i = 1, sizei
551 w( eltvar( ielptr + i) ) =
552 & w( eltvar( ielptr + i) )
553 & + abs(a_elt( k8 ))
554 k8 = k8 + 1_8
555 END DO
556 END DO
557 ELSE
558 DO j = 1, sizei
559 temp = w( eltvar( ielptr + j ) )
560 DO i = 1, sizei
561 temp = temp + abs( a_elt(k8))
562 k8 = k8 + 1_8
563 END DO
564 w(eltvar( ielptr + j )) =
565 & w(eltvar( ielptr + j )) + temp
566 END DO
567 ENDIF
568 ELSE
569 DO j = 1, sizei
570 w(eltvar( ielptr + j )) =
571 & w(eltvar( ielptr + j )) + abs(a_elt( k8 ))
572 k8 = k8 + 1_8
573 DO i = j+1, sizei
574 w(eltvar( ielptr + j )) =
575 & w(eltvar( ielptr + j )) + abs(a_elt( k8 ))
576 w(eltvar( ielptr + i ) ) =
577 & w(eltvar( ielptr + i )) + abs(a_elt( k8 ))
578 k8 = k8 + 1_8
579 END DO
580 ENDDO
581 ENDIF
582 ENDDO
583 RETURN

◆ smumps_sol_y()

subroutine smumps_sol_y ( real, dimension(nz8), intent(in) a,
integer(8), intent(in) nz8,
integer, intent(in) n,
integer, dimension(nz8), intent(in) irn,
integer, dimension(nz8), intent(in) icn,
real, dimension(n), intent(in) rhs,
real, dimension(n), intent(in) x,
real, dimension(n), intent(out) r,
real, dimension(n), intent(out) w,
integer, dimension(500), intent(in) keep,
integer(8), dimension(150), intent(in) keep8 )

Definition at line 222 of file ssol_aux.F.

224 IMPLICIT NONE
225 INTEGER, INTENT(IN) :: N, KEEP(500)
226 INTEGER(8), INTENT(IN) :: NZ8
227 INTEGER(8), INTENT(IN) :: KEEP8(150)
228 INTEGER, INTENT(IN) :: IRN(NZ8), ICN(NZ8)
229 REAL, INTENT(IN) :: A(NZ8), RHS(N), X(N)
230 REAL, INTENT(OUT) :: W(N)
231 REAL, INTENT(OUT) :: R(N)
232 INTEGER I, J
233 INTEGER(8) :: K8
234 REAL, PARAMETER :: ZERO = 0.0e0
235 REAL D
236 DO i = 1, n
237 r(i) = rhs(i)
238 w(i) = zero
239 ENDDO
240 IF (keep(264).EQ.0) THEN
241 IF (keep(50) .EQ.0) THEN
242 DO k8 = 1_8, nz8
243 i = irn(k8)
244 j = icn(k8)
245 IF ((i .GT. n) .OR. (j .GT. n) .OR. (i .LT. 1) .OR.
246 & (j .LT. 1)) cycle
247 d = a(k8) * x(j)
248 r(i) = r(i) - d
249 w(i) = w(i) + abs(d)
250 ENDDO
251 ELSE
252 DO k8 = 1_8, nz8
253 i = irn(k8)
254 j = icn(k8)
255 IF ((i .GT. n) .OR. (j .GT. n) .OR. (i .LT. 1) .OR.
256 & (j .LT. 1)) cycle
257 d = a(k8) * x(j)
258 r(i) = r(i) - d
259 w(i) = w(i) + abs(d)
260 IF (i.NE.j) THEN
261 d = a(k8) * x(i)
262 r(j) = r(j) - d
263 w(j) = w(j) + abs(d)
264 ENDIF
265 ENDDO
266 ENDIF
267 ELSE
268 IF (keep(50) .EQ.0) THEN
269 DO k8 = 1_8, nz8
270 i = irn(k8)
271 j = icn(k8)
272 d = a(k8) * x(j)
273 r(i) = r(i) - d
274 w(i) = w(i) + abs(d)
275 ENDDO
276 ELSE
277 DO k8 = 1_8, nz8
278 i = irn(k8)
279 j = icn(k8)
280 d = a(k8) * x(j)
281 r(i) = r(i) - d
282 w(i) = w(i) + abs(d)
283 IF (i.NE.j) THEN
284 d = a(k8) * x(i)
285 r(j) = r(j) - d
286 w(j) = w(j) + abs(d)
287 ENDIF
288 ENDDO
289 ENDIF
290 ENDIF
291 RETURN

◆ smumps_solve_bwd_panels()

subroutine smumps_solve_bwd_panels ( real, dimension(la), intent(in) a,
integer(8), intent(in) la,
integer(8), intent(in) apos,
integer, intent(in) npiv,
integer, dimension(npiv), intent(in) iw,
integer, intent(in) nrhs_b,
real, dimension(lwcb), intent(inout) wcb,
integer(8), intent(in) lwcb,
integer, intent(in) lda_wcb,
integer(8), intent(in) ppiv_courant,
integer, intent(in) mtype,
integer, dimension(500), intent(in) keep )

Definition at line 1270 of file ssol_aux.F.

1274 INTEGER, INTENT(IN) :: MTYPE, NPIV, KEEP(500)
1275 INTEGER, INTENT(IN) :: IW(NPIV)
1276 INTEGER, INTENT(IN) :: NRHS_B, LDA_WCB
1277 INTEGER(8), INTENT(IN) :: LA, APOS, LWCB, PPIV_COURANT
1278 REAL, INTENT(IN) :: A(LA)
1279 REAL, INTENT(INOUT) :: WCB(LWCB)
1280 INTEGER, PARAMETER :: PANEL_TABSIZE = 20
1281 INTEGER(8) :: PANEL_POS(PANEL_TABSIZE)
1282 INTEGER :: PANEL_COL(PANEL_TABSIZE)
1283 INTEGER :: IPANEL, NBPANELS, NB_TARGET
1284 INTEGER :: NBROWS_PANEL, NBCOLS_PANEL
1285 INTEGER(8) :: PPIV_PANEL
1286 INTEGER :: MTYPE_TEMP
1287 REAL, PARAMETER :: ONE=1.0e0
1288 IF (keep(459) .LE. 1) THEN
1289 WRITE(*,*) " Internal error 1 in SMUMPS_SOLVE_BWD_PANELS"
1290 CALL mumps_abort()
1291 ENDIF
1292 IF ( keep(459)+1 .GT. panel_tabsize ) THEN
1293 WRITE(*,*) " Internal error 2 in SMUMPS_SOLVE_BWD_PANELS"
1294 CALL mumps_abort()
1295 ENDIF
1296 CALL mumps_ldltpanel_panelinfos( npiv, keep, iw,
1297 &nb_target, nbpanels, panel_col, panel_pos, panel_tabsize,
1298 &.false. )
1299 DO ipanel = nbpanels, 1, -1
1300 nbcols_panel = panel_col( ipanel+1 ) - panel_col( ipanel )
1301 nbrows_panel = npiv - panel_col( ipanel ) + 1
1302 ppiv_panel = ppiv_courant + panel_col( ipanel ) - 1
1303 IF ( nbrows_panel .GT. nbcols_panel ) THEN
1304 mtype_temp = 0
1305 CALL smumps_solve_gemm_update( a, la,
1306 & apos-1_8+panel_pos(ipanel)+
1307 & int(nbcols_panel,8)*int(nbcols_panel,8),
1308 & nbrows_panel-nbcols_panel, nbcols_panel,
1309 & nbcols_panel,
1310 & nrhs_b, wcb, lwcb, ppiv_panel+nbcols_panel, lda_wcb,
1311 & ppiv_panel, lda_wcb,
1312 & mtype_temp, keep, one )
1313 ENDIF
1314 CALL smumps_solve_bwd_trsolve (a, la,
1315 & apos+panel_pos(ipanel)-1_8,
1316 & nbcols_panel, nbcols_panel,
1317 & nrhs_b, wcb, lwcb, lda_wcb, ppiv_panel, mtype, keep)
1318 ENDDO
1319 RETURN
subroutine smumps_solve_bwd_trsolve(a, la, apos, npiv, ldadiag, nrhs_b, wcb, lwcb, lda_wcb, ppiv_courant, mtype, keep)
Definition ssol_aux.F:1185
subroutine smumps_solve_gemm_update(a, la, apos1, nx, lda, ny, nrhs_b, wcb, lwcb, ptrx, ldx, ptry, ldy, mtype, keep, coef_y)
Definition ssol_aux.F:1326

◆ smumps_solve_bwd_trsolve()

subroutine smumps_solve_bwd_trsolve ( real, dimension(la), intent(in) a,
integer(8), intent(in) la,
integer(8), intent(in) apos,
integer, intent(in) npiv,
integer, intent(in) ldadiag,
integer, intent(in) nrhs_b,
real, dimension(lwcb), intent(inout) wcb,
integer(8), intent(in) lwcb,
integer, intent(in) lda_wcb,
integer(8), intent(in) ppiv_courant,
integer, intent(in) mtype,
integer, dimension(500), intent(in) keep )

Definition at line 1183 of file ssol_aux.F.

1185 INTEGER, INTENT(IN) :: MTYPE, LDADIAG, NPIV, KEEP(500)
1186 INTEGER, INTENT(IN) :: NRHS_B, LDA_WCB
1187 INTEGER(8), INTENT(IN) :: LA, APOS, LWCB, PPIV_COURANT
1188 REAL, INTENT(IN) :: A(LA)
1189 REAL, INTENT(INOUT) :: WCB(LWCB)
1190 REAL ONE
1191 parameter(one = 1.0e0)
1192 IF (mtype .eq. 1 ) THEN
1193#if defined(MUMPS_USE_BLAS2)
1194 IF ( nrhs_b == 1 ) THEN
1195 CALL strsv( 'L', 't', 'n', NPIV, A(APOS), LDADIAG,
1196 & WCB(PPIV_COURANT), 1 )
1197 ELSE
1198#endif
1199 CALL strsm( 'l','l','t','n', npiv, nrhs_b, one,
1200 & a(apos), ldadiag, wcb(ppiv_courant),
1201 & lda_wcb )
1202#if defined(MUMPS_USE_BLAS2)
1203 ENDIF
1204#endif
1205 ELSE
1206#if defined(MUMPS_USE_BLAS2)
1207 IF ( nrhs_b == 1 ) THEN
1208 CALL strsv( 'U', 'N', 'U', npiv, a(apos), ldadiag,
1209 & wcb(ppiv_courant), 1 )
1210 ELSE
1211#endif
1212 CALL strsm( 'L','U','N','U', npiv, nrhs_b, one,
1213 & a(apos), ldadiag, wcb(ppiv_courant),
1214 & lda_wcb )
1215#if defined(MUMPS_USE_BLAS2)
1216 ENDIF
1217#endif
1218 ENDIF
1219 RETURN
subroutine strsv(uplo, trans, diag, n, a, lda, x, incx)
STRSV
Definition strsv.f:149
subroutine strsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRSM
Definition strsm.f:181

◆ smumps_solve_fwd_panels()

subroutine smumps_solve_fwd_panels ( real, dimension(la), intent(in) a,
integer(8), intent(in) la,
integer(8), intent(in) apos,
integer, intent(in) npiv,
integer, dimension(npiv), intent(in) iw,
integer, intent(in) nrhs_b,
real, dimension(lwcb), intent(inout) wcb,
integer(8), intent(in) lwcb,
integer, intent(in) lda_wcb,
integer(8), intent(in) ppiv_courant,
integer, intent(in) mtype,
integer, dimension(500), intent(in) keep )

Definition at line 1221 of file ssol_aux.F.

1225 INTEGER, INTENT(IN) :: MTYPE, NPIV, KEEP(500)
1226 INTEGER, INTENT(IN) :: IW(NPIV)
1227 INTEGER, INTENT(IN) :: NRHS_B, LDA_WCB
1228 INTEGER(8), INTENT(IN) :: LA, APOS, LWCB, PPIV_COURANT
1229 REAL, INTENT(IN) :: A(LA)
1230 REAL, INTENT(INOUT) :: WCB(LWCB)
1231 INTEGER :: NB_TARGET
1232 INTEGER :: NBPANELS
1233 INTEGER :: NBROWS_PANEL, NBCOLS_PANEL, ICOL_BEG, ICOL_END
1234 INTEGER(8) :: PANEL_APOS, PPIV_PANEL
1235 REAL, PARAMETER :: ONE=1.0e0
1236 IF (keep(459) .LE. 1) THEN
1237 WRITE(*,*) " Internal error in SMUMPS_SOLVE_FWD_PANELS"
1238 CALL mumps_abort()
1239 ENDIF
1240 CALL mumps_ldltpanel_nbtarget( npiv, nb_target, keep )
1241 panel_apos = apos
1242 nbpanels = 0
1243 icol_beg = 1
1244 nbrows_panel = npiv
1245 ppiv_panel = ppiv_courant
1246 DO WHILE ( icol_beg .LE. npiv )
1247 nbpanels = nbpanels + 1
1248 icol_end = min(nb_target * nbpanels, npiv)
1249 IF ( iw(icol_end) .LT. 0 ) icol_end=icol_end+1
1250 nbcols_panel = icol_end - icol_beg + 1
1251 CALL smumps_solve_fwd_trsolve (a, la, panel_apos,
1252 & nbcols_panel, nbcols_panel,
1253 & nrhs_b, wcb, lwcb, lda_wcb, ppiv_panel, mtype, keep)
1254 IF ( nbrows_panel .GT. nbcols_panel ) THEN
1255 CALL smumps_solve_gemm_update( a, la,
1256 & panel_apos + int(nbcols_panel,8) * int(nbcols_panel,8),
1257 & nbcols_panel, nbcols_panel, nbrows_panel-nbcols_panel,
1258 & nrhs_b, wcb, lwcb, ppiv_panel, lda_wcb,
1259 & ppiv_panel+nbcols_panel, lda_wcb,
1260 & mtype, keep, one )
1261 ENDIF
1262 icol_beg = icol_end + 1
1263 panel_apos = panel_apos + int(nbcols_panel,8) *
1264 & int(nbrows_panel,8)
1265 nbrows_panel = nbrows_panel - nbcols_panel
1266 ppiv_panel = ppiv_panel + nbcols_panel
1267 ENDDO
1268 RETURN
#define min(a, b)
Definition macros.h:20
subroutine smumps_solve_fwd_trsolve(a, la, apos, npiv, ldadiag, nrhs_b, wcb, lwcb, lda_wcb, ppiv_courant, mtype, keep)
Definition ssol_aux.F:1147

◆ smumps_solve_fwd_trsolve()

subroutine smumps_solve_fwd_trsolve ( real, dimension(la), intent(in) a,
integer(8), intent(in) la,
integer(8), intent(in) apos,
integer, intent(in) npiv,
integer, intent(in) ldadiag,
integer, intent(in) nrhs_b,
real, dimension(lwcb), intent(inout) wcb,
integer(8), intent(in) lwcb,
integer, intent(in) lda_wcb,
integer(8), intent(in) ppiv_courant,
integer, intent(in) mtype,
integer, dimension(500), intent(in) keep )

Definition at line 1145 of file ssol_aux.F.

1147 INTEGER, INTENT(IN) :: MTYPE, LDADIAG, NPIV, KEEP(500)
1148 INTEGER, INTENT(IN) :: NRHS_B, LDA_WCB
1149 INTEGER(8), INTENT(IN) :: LA, APOS, LWCB, PPIV_COURANT
1150 REAL, INTENT(IN) :: A(LA)
1151 REAL, INTENT(INOUT) :: WCB(LWCB)
1152 REAL ONE
1153 parameter(one = 1.0e0)
1154 IF (keep(50).NE.0 .OR. mtype .eq. 1 ) THEN
1155#if defined(MUMPS_USE_BLAS2)
1156 IF ( nrhs_b == 1 ) THEN
1157 CALL strsv( 'U', 'T', 'U', npiv, a(apos), ldadiag,
1158 & wcb(ppiv_courant), 1 )
1159 ELSE
1160#endif
1161 CALL strsm( 'L','U','T','U', npiv, nrhs_b, one,
1162 & a(apos), ldadiag, wcb(ppiv_courant),
1163 & lda_wcb )
1164#if defined(MUMPS_USE_BLAS2)
1165 ENDIF
1166#endif
1167 ELSE
1168#if defined(MUMPS_USE_BLAS2)
1169 IF ( nrhs_b == 1 ) THEN
1170 CALL strsv( 'L', 'N', 'N', npiv, a(apos), ldadiag,
1171 & wcb(ppiv_courant), 1 )
1172 ELSE
1173#endif
1174 CALL strsm( 'L','l','n','n', NPIV, NRHS_B, ONE,
1175 & A(APOS), LDADIAG, WCB(PPIV_COURANT),
1176 & LDA_WCB )
1177#if defined(MUMPS_USE_BLAS2)
1178 ENDIF
1179#endif
1180 ENDIF
1181 RETURN

◆ smumps_solve_gemm_update()

subroutine smumps_solve_gemm_update ( real, dimension(la), intent(in) a,
integer(8), intent(in) la,
integer(8), intent(in) apos1,
integer, intent(in) nx,
integer, intent(in) lda,
integer, intent(in) ny,
integer, intent(in) nrhs_b,
real, dimension(lwcb), intent(inout) wcb,
integer(8), intent(in) lwcb,
integer(8), intent(in) ptrx,
integer, intent(in) ldx,
integer(8), intent(in) ptry,
integer, intent(in) ldy,
integer, intent(in) mtype,
integer, dimension(500), intent(in) keep,
real, intent(in) coef_y )

Definition at line 1321 of file ssol_aux.F.

1326 INTEGER, INTENT(IN) :: MTYPE, NY, NX, KEEP(500)
1327 INTEGER, INTENT(IN) :: NRHS_B, LDY, LDA, LDX
1328 INTEGER(8), INTENT(IN) :: LA, APOS1, LWCB, PTRX,
1329 & PTRY
1330 REAL, INTENT(IN) :: A(LA)
1331 REAL, INTENT(INOUT) :: WCB(LWCB)
1332 REAL, INTENT(IN) :: COEF_Y
1333 REAL ALPHA, ZERO, ONE
1334 parameter(zero = 0.0e0, one = 1.0e0, alpha=-1.0e0)
1335 IF ( nx .NE. 0 .AND. ny.NE.0 ) THEN
1336 IF ( mtype .eq. 1 ) THEN
1337#if defined(MUMPS_USE_BLAS2)
1338 IF ( nrhs_b == 1 ) THEN
1339 CALL sgemv('T', nx, ny, alpha, a(apos1),
1340 & lda, wcb(ptrx), 1, coef_y,
1341 & wcb(ptry), 1)
1342 ELSE
1343#endif
1344 CALL sgemm('T', 'N', ny, nrhs_b, nx, alpha,
1345 & a(apos1), lda, wcb(ptrx), ldx, coef_y,
1346 & wcb(ptry), ldy)
1347#if defined(MUMPS_USE_BLAS2)
1348 END IF
1349#endif
1350 ELSE
1351#if defined(MUMPS_USE_BLAS2)
1352 IF ( nrhs_b == 1 ) THEN
1353 CALL sgemv('N',ny, nx, alpha, a(apos1),
1354 & lda, wcb(ptrx), 1,
1355 & coef_y, wcb(ptry), 1 )
1356 ELSE
1357#endif
1358 CALL sgemm('N', 'N', ny, nrhs_b, nx, alpha,
1359 & a(apos1), lda, wcb(ptrx), ldx,
1360 & coef_y, wcb(ptry), ldy)
1361#if defined(MUMPS_USE_BLAS2)
1362 END IF
1363#endif
1364 END IF
1365 END IF
1366 RETURN
#define alpha
Definition eval.h:35
subroutine sgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
SGEMV
Definition sgemv.f:156
subroutine sgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
SGEMM
Definition sgemm.f:187

◆ smumps_solve_get_ooc_node()

subroutine smumps_solve_get_ooc_node ( integer inode,
integer(8), dimension(keep(28)) ptrfac,
integer, dimension(500) keep,
real, dimension(la) a,
integer(8) la,
integer, dimension(n) step,
integer(8), dimension(150) keep8,
integer n,
logical must_be_permuted,
integer ierr )

Definition at line 729 of file ssol_aux.F.

732 USE smumps_ooc
733 IMPLICIT NONE
734 INTEGER INODE,KEEP(500),N
735 INTEGER(8) KEEP8(150)
736 INTEGER(8) :: LA
737 INTEGER(8) :: PTRFAC(KEEP(28))
738 INTEGER STEP(N)
739 INTEGER IERR
740 REAL A(LA)
741 INTEGER RETURN_VALUE
742 LOGICAL MUST_BE_PERMUTED
743 return_value=smumps_solve_is_inode_in_mem(inode,ptrfac,
744 & keep(28),a,la,ierr)
745 IF(return_value.EQ.ooc_node_not_in_mem)THEN
746 IF(ierr.LT.0)THEN
747 RETURN
748 ENDIF
749 CALL smumps_solve_alloc_factor_space(inode,ptrfac,
750 & keep,keep8,a,ierr)
751 IF(ierr.LT.0)THEN
752 RETURN
753 ENDIF
754 CALL smumps_read_ooc(
755 & a(ptrfac(step(inode))),
756 & inode,ierr
757 & )
758 IF(ierr.LT.0)THEN
759 RETURN
760 ENDIF
761 ELSE
762 IF(ierr.LT.0)THEN
763 RETURN
764 ENDIF
765 ENDIF
766 IF(return_value.NE.ooc_node_permuted)THEN
767 must_be_permuted=.true.
769 ELSE
770 must_be_permuted=.false.
771 ENDIF
772 RETURN
integer ooc_node_not_in_mem
Definition smumps_ooc.F:24
integer function smumps_solve_is_inode_in_mem(inode, ptrfac, nsteps, a, la, ierr)
integer ooc_node_permuted
Definition smumps_ooc.F:24
subroutine smumps_solve_modify_state_node(inode)
subroutine, public smumps_solve_alloc_factor_space(inode, ptrfac, keep, keep8, a, ierr)
subroutine, public smumps_read_ooc(dest, inode, ierr)
Definition smumps_ooc.F:394