OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
zmumps_sol_lr Module Reference

Functions/Subroutines

subroutine zmumps_sol_fwd_lr_su (inode, n, iwhdlr, npiv_global, nslaves, iw, ipos_init, liw, liell, wcb, lwcb, ld_wcbpiv, ld_wcbcb, ppiv_init, pcb_init, rhscomp, lrhscomp, nrhs, posinrhscomp_fwd, jbdeb, jbfin, mtype, keep, keep8, oocwrite_compatible_with_blr, iflag, ierror)
subroutine zmumps_sol_slave_lr_u (inode, iwhdlr, npiv_global, wcb, lwcb, ldx, ldy, ptrx_init, ptry_init, jbdeb, jbfin, mtype, keep, keep8, iflag, ierror)
subroutine zmumps_sol_fwd_blr_update (arraypiv, lpiv, lpivcol, ldpiv, pospiv, pospivcol, arraycb, lcb, ldcb, poscb, posdiag, nrhs_b, npiv, blr_panel, last_blr, current_blr, begs_blr_static, keep8, k34, k450, is_t2_slave, iflag, ierror)
subroutine zmumps_sol_bwd_lr_su (inode, iwhdlr, npiv_global, nslaves, liell, wcb, lwcb, nrhs_b, ptwcb, rhscomp, lrhscomp, nrhs, iposinrhscomp, jbdeb, mtype, keep, keep8, iflag, ierror)
subroutine zmumps_sol_bwd_blr_update (arraypiv, lpiv, lpivcol, ldpiv, pospiv, pospivcol, arraycb, lcb, ldcb, poscb, posdiag, nrhs_b, npiv, blr_panel, last_blr, current_blr, begs_blr_static, keep8, k34, k450, is_t2_slave, iflag, ierror)

Function/Subroutine Documentation

◆ zmumps_sol_bwd_blr_update()

subroutine zmumps_sol_lr::zmumps_sol_bwd_blr_update ( complex(kind=8), dimension(lpiv,lpivcol), intent(inout) arraypiv,
integer(8), intent(in) lpiv,
integer, intent(in) lpivcol,
integer, intent(in) ldpiv,
integer(8), intent(in) pospiv,
integer, intent(in) pospivcol,
complex(kind=8), dimension(lcb), intent(inout) arraycb,
integer(8), intent(in) lcb,
integer, intent(in) ldcb,
integer(8), intent(in) poscb,
integer(8), intent(in) posdiag,
integer, intent(in) nrhs_b,
integer, intent(in) npiv,
type(lrb_type), dimension(:), intent(in), target blr_panel,
integer, intent(in) last_blr,
integer, intent(in) current_blr,
integer, dimension(:) begs_blr_static,
integer(8), dimension(150), intent(in) keep8,
integer, intent(in) k34,
integer, intent(in) k450,
logical, intent(in) is_t2_slave,
integer, intent(inout) iflag,
integer, intent(inout) ierror )

Definition at line 528 of file zsol_lr.F.

537!$ USE OMP_LIB
538 INTEGER(8), INTENT(IN) :: LPIV, LCB, POSPIV, POSCB, POSDIAG
539 INTEGER,INTENT(IN) :: LPIVCOL, POSPIVCOL
540 COMPLEX(kind=8), INTENT(INOUT) :: ARRAYPIV(LPIV,LPIVCOL)
541 COMPLEX(kind=8), INTENT(INOUT) :: ARRAYCB(LCB)
542 INTEGER, INTENT(IN) :: LAST_BLR, NRHS_B, LDPIV, LDCB,
543 & CURRENT_BLR, NPIV, K34, K450
544 TYPE(LRB_TYPE), TARGET,INTENT(IN) ::
545 & BLR_PANEL(:)
546 INTEGER(8), INTENT(IN) :: KEEP8(150)
547 LOGICAL, INTENT(IN) :: IS_T2_SLAVE
548 INTEGER :: BEGS_BLR_STATIC(:)
549 INTEGER, INTENT(INOUT) :: IFLAG, IERROR
550 INTEGER :: I, K, M, N, IBEG_BLOCK, IEND_BLOCK
551 INTEGER :: KMAX
552 INTEGER(8) :: POSBLOCK
553 TYPE(LRB_TYPE), POINTER :: LRB
554 COMPLEX(kind=8), ALLOCATABLE, DIMENSION(:) :: TEMP_BLOCK
555 COMPLEX(kind=8), ALLOCATABLE, DIMENSION(:) :: DEST_ARRAY
556 INTEGER :: allocok
557 COMPLEX(kind=8) :: ONE, MONE, ZERO
558 parameter(one=(1.0d0,0.0d0), mone=(-1.0d0,0.0d0))
559 parameter(zero=(0.0d0,0.0d0))
560 kmax = -1
561 DO i = current_blr+1, last_blr
562 kmax = max(kmax, blr_panel(i-current_blr)%K)
563 ENDDO
564 IF (current_blr.LT.last_blr) THEN
565 n = blr_panel(1)%N
566 ELSE
567 RETURN
568 ENDIF
569 allocate(dest_array(n*nrhs_b),stat=allocok)
570 IF (allocok .GT. 0) THEN
571 iflag = -13
572 ierror = n * nrhs_b
573 GOTO 100
574 ENDIF
575 dest_array = zero
576#if defined(BLR_MT)
577!$OMP PARALLEL PRIVATE(TEMP_BLOCK,allocok)
578#endif
579 IF (kmax.GT.0) THEN
580 allocate(temp_block(kmax*nrhs_b), stat=allocok )
581 IF (allocok .GT. 0) THEN
582 iflag = -13
583 ierror = nrhs_b * kmax
584 write(*,*) 'Allocation problem in BLR routine
585 & ZMUMPS_SOL_BWD_BLR_UPDATE: ',
586 & 'not enough memory? memory requested = ', ierror
587 ENDIF
588 ENDIF
589#if defined(BLR_MT)
590!$OMP BARRIER
591#endif
592#if defined(BLR_MT)
593!$OMP DO SCHEDULE(DYNAMIC,1)
594!$OMP& PRIVATE(IBEG_BLOCK, IEND_BLOCK, LRB, K, M,
595!$OMP& POSBLOCK)
596!$OMP& REDUCTION(+:DEST_ARRAY)
597#endif
598 DO i = current_blr+1, last_blr
599 IF (iflag.LT.0) cycle
600 ibeg_block = begs_blr_static(i)
601 iend_block = begs_blr_static(i+1)-1
602 lrb => blr_panel(i-current_blr)
603 k = lrb%K
604 m = lrb%M
605 IF (lrb%ISLR) THEN
606 IF (k.GT.0) THEN
607 IF (is_t2_slave) THEN
608 posblock = poscb +int(ibeg_block-1,8)
609 CALL zgemm('T', 'N', k, nrhs_b, m, one,
610 & lrb%Q(1,1), m,
611 & arraycb(posblock), ldcb, zero,
612 & temp_block(1), k)
613 ELSE IF (ibeg_block.LE.npiv.AND.iend_block.GT.npiv) THEN
614 posblock = pospiv+int(ibeg_block-1,8)
615 CALL zgemm('T', 'N', k, nrhs_b, npiv-ibeg_block+1, one,
616 & lrb%Q(1,1), m,
617 & arraypiv(posblock,pospivcol), ldpiv,
618 & zero, temp_block(1), k)
619 CALL zgemm('T', 'N', k, nrhs_b, ibeg_block+m-npiv-1,
620 & one, lrb%Q(npiv-ibeg_block+2,1), m,
621 & arraycb(poscb), ldcb,
622 & one,
623 & temp_block(1), k)
624 ELSEIF (ibeg_block.LE.npiv) THEN
625 posblock = pospiv+int(ibeg_block-1,8)
626 CALL zgemm('T', 'N', k, nrhs_b, m, one,
627 & lrb%Q(1,1), m,
628 & arraypiv(posblock,pospivcol), ldpiv,
629 & zero, temp_block(1), k)
630 ELSE
631 posblock = poscb+int(ibeg_block-1-npiv,8)
632 CALL zgemm('T', 'N', k, nrhs_b, m, one,
633 & lrb%Q(1,1), m,
634 & arraycb(posblock), ldcb, zero,
635 & temp_block(1), k)
636 ENDIF
637 CALL zgemm('T', 'N', n, nrhs_b, k, mone,
638 & lrb%R(1,1), k,
639 & temp_block(1), k, one,
640 & dest_array(1), n)
641 ENDIF
642 ELSE
643 IF (is_t2_slave) THEN
644 posblock = poscb+int(ibeg_block-1,8)
645 CALL zgemm('T', 'N', n, nrhs_b, m, mone,
646 & lrb%Q(1,1), m, arraycb(posblock), ldcb,
647 & one, dest_array(1), n)
648 ELSE IF (ibeg_block.LE.npiv.AND.iend_block.GT.npiv) THEN
649 posblock = pospiv+int(ibeg_block-1,8)
650 CALL zgemm('T', 'N', n, nrhs_b, npiv-ibeg_block+1, mone,
651 & lrb%Q(1,1), m, arraypiv(posblock,pospivcol), ldpiv,
652 & one, dest_array(1), n)
653 CALL zgemm('T', 'N', n, nrhs_b, ibeg_block+m-npiv-1, mone,
654 & lrb%Q(npiv-ibeg_block+2,1), m, arraycb(poscb),
655 & ldcb, one, dest_array(1), n)
656 ELSEIF (ibeg_block.LE.npiv) THEN
657 posblock = pospiv+int(ibeg_block-1,8)
658 CALL zgemm('T', 'N', n, nrhs_b, m, mone,
659 & lrb%Q(1,1), m, arraypiv(posblock,pospivcol), ldpiv,
660 & one, dest_array(1), n)
661 ELSE
662 posblock = poscb+int(ibeg_block-1-npiv,8)
663 CALL zgemm('T', 'N', n, nrhs_b, m, mone,
664 & lrb%Q(1,1), m, arraycb(posblock), ldcb,
665 & one, dest_array(1), n)
666 ENDIF
667 ENDIF
668 ENDDO
669#if defined(BLR_MT)
670!$OMP END DO
671#endif
672 IF (kmax.GT.0) THEN
673 IF (allocated(temp_block)) deallocate(temp_block)
674 ENDIF
675#if defined(BLR_MT)
676!$OMP END PARALLEL
677#endif
678 IF (is_t2_slave) THEN
679 DO i=1,nrhs_b
680 call zaxpy(n, one, dest_array((i-1)*n+1), 1,
681 & arraypiv(posdiag+(i-1)*ldpiv,pospivcol), 1)
682 ENDDO
683 ELSE
684 DO i=1,nrhs_b
685 call zaxpy(n, one, dest_array((i-1)*n+1), 1,
686 & arraypiv(posdiag,pospivcol+i-1), 1)
687 ENDDO
688 ENDIF
689 100 CONTINUE
690 IF (allocated(dest_array)) DEALLOCATE(dest_array)
691 RETURN
subroutine zaxpy(n, za, zx, incx, zy, incy)
ZAXPY
Definition zaxpy.f:88
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:187
#define max(a, b)
Definition macros.h:21

◆ zmumps_sol_bwd_lr_su()

subroutine zmumps_sol_lr::zmumps_sol_bwd_lr_su ( integer, intent(in) inode,
integer, intent(in) iwhdlr,
integer, intent(in) npiv_global,
integer, intent(in) nslaves,
integer, intent(in) liell,
complex(kind=8), dimension(lwcb), intent(inout) wcb,
integer(8), intent(in) lwcb,
integer, intent(in) nrhs_b,
integer(8), intent(in) ptwcb,
complex(kind=8), dimension(lrhscomp,nrhs) rhscomp,
integer, intent(in) lrhscomp,
integer, intent(in) nrhs,
integer, intent(in) iposinrhscomp,
integer, intent(in) jbdeb,
integer, intent(in) mtype,
integer, dimension(500), intent(in) keep,
integer(8), dimension(150), intent(in) keep8,
integer, intent(inout) iflag,
integer, intent(inout) ierror )

Definition at line 379 of file zsol_lr.F.

386 INTEGER, INTENT(IN) :: INODE, IWHDLR, NPIV_GLOBAL, NSLAVES
387 INTEGER, INTENT(IN) :: MTYPE, LIELL, KEEP(500)
388 INTEGER(8), INTENT(IN) :: KEEP8(150)
389 INTEGER, INTENT(IN) :: IPOSINRHSCOMP, JBDEB, LRHSCOMP, NRHS
390 INTEGER(8), INTENT(IN) :: LWCB, PTWCB
391 INTEGER, INTENT(IN) :: NRHS_B
392 INTEGER, INTENT(INOUT) :: IFLAG, IERROR
393 COMPLEX(kind=8), INTENT(INOUT) :: WCB(LWCB)
394 COMPLEX(kind=8) RHSCOMP(LRHSCOMP,NRHS)
395 INTEGER :: I, NPARTSASS, NB_BLR, LAST_BLR,
396 & NELIM_PANEL, LD_WCB,
397 & DIAGSIZ_DYN, DIAGSIZ_STA, LDADIAG,
398 & IEND_BLR, IBEG_BLR, PCBINRHSCOMP
399 INTEGER(8) :: PCB_LAST, PWCB
400 INTEGER :: IPIV_PANEL
401 COMPLEX(kind=8), POINTER, DIMENSION(:) :: DIAG
402 TYPE(LRB_TYPE), POINTER, DIMENSION(:) :: BLR_PANEL
403 COMPLEX(kind=8) :: ONE, MONE, ZERO
404 parameter(one=(1.0d0,0.0d0), mone=(-1.0d0,0.0d0))
405 parameter(zero=(0.0d0,0.0d0))
406 IF ((mtype.EQ.1).AND.(keep(50).EQ.0)) THEN
407 IF (associated(blr_array(iwhdlr)%PANELS_U))
408 & THEN
409 npartsass=size(blr_array(iwhdlr)%PANELS_U)
410 nb_blr = size(blr_array(iwhdlr)%BEGS_BLR_STATIC) -1
411 ENDIF
412 ELSE
413 IF (associated(blr_array(iwhdlr)%PANELS_L))
414 & THEN
415 npartsass=size(blr_array(iwhdlr)%PANELS_L)
416 nb_blr = size(blr_array(iwhdlr)%BEGS_BLR_STATIC) -1
417 ELSE
418 WRITE(6,*) " Internal error in ZMUMPS_SOL_FWD_SU_MASTER"
419 ENDIF
420 ENDIF
421 pcbinrhscomp= iposinrhscomp + npiv_global
422 pcb_last = ptwcb + int(liell ,8)
423 pwcb = ptwcb + int(npiv_global,8)
424 ld_wcb = liell
425 DO i=npartsass,1,-1
426 ibeg_blr = blr_array(iwhdlr)%BEGS_BLR_DYNAMIC(i)
427 iend_blr = blr_array(iwhdlr)%BEGS_BLR_DYNAMIC(i+1) -1
428 diagsiz_dyn = blr_array(iwhdlr)%BEGS_BLR_DYNAMIC(i+1) -
429 & ibeg_blr
430 diagsiz_sta = blr_array(iwhdlr)%BEGS_BLR_STATIC(i+1) -
431 & ibeg_blr
432 IF (keep(50).NE.0) THEN
433 ldadiag = diagsiz_dyn
434 ELSE
435 ldadiag = diagsiz_sta
436 ENDIF
437 IF (diagsiz_dyn.EQ.0) cycle
438 nelim_panel = diagsiz_sta - diagsiz_dyn
439 ipiv_panel = iposinrhscomp + ibeg_blr -1
440 IF ( mtype .EQ. 1 .AND. keep(50).EQ.0) THEN
441 blr_panel => blr_array(iwhdlr)%PANELS_U(i)%LRB_PANEL
442 ELSE
443 blr_panel => blr_array(iwhdlr)%PANELS_L(i)%LRB_PANEL
444 END IF
445 IF (keep(50).EQ.0 .AND. nslaves.GT.0 .AND. mtype.NE.1) THEN
446 last_blr = npartsass
447 ELSE
448 last_blr = nb_blr
449 ENDIF
450 CALL zmumps_sol_bwd_blr_update (
451 & rhscomp, int(lrhscomp,8), nrhs, lrhscomp,
452 & int(iposinrhscomp,8), jbdeb,
453 & wcb, lwcb, ld_wcb, pwcb,
454 & int(ipiv_panel,8),
455 & nrhs_b, npiv_global,
456 & blr_panel, last_blr,
457 & i, blr_array(iwhdlr)%BEGS_BLR_STATIC,
458 & keep8, keep(34), keep(450), .false., iflag, ierror)
459 IF (iflag.LT.0) RETURN
460 diag => blr_array(iwhdlr)%DIAG_BLOCKS(i)%DIAG_BLOCK
461 IF (nelim_panel.GT.0) THEN
462 IF (mtype.EQ.1.AND.keep(50).EQ.0) THEN
463 IF (iend_blr.EQ.npiv_global) THEN
464 CALL zgemm('T', 'N', diagsiz_dyn, nrhs_b, nelim_panel,
465 & mone, diag(1+diagsiz_dyn), diagsiz_sta, wcb(pwcb),
466 & ld_wcb, one , rhscomp(ipiv_panel,jbdeb),lrhscomp)
467 ELSE
468 IF (iend_blr+1.LE.npiv_global .AND.
469 & iend_blr+nelim_panel.GT.npiv_global) THEN
470 CALL zgemm('T', 'N', diagsiz_dyn, nrhs_b,
471 & npiv_global-iend_blr,
472 & mone, diag(1+diagsiz_dyn), diagsiz_sta,
473 & rhscomp(ipiv_panel+diagsiz_dyn,jbdeb), lrhscomp,
474 & one, rhscomp(ipiv_panel,jbdeb), lrhscomp)
475 CALL zgemm('T', 'N', diagsiz_dyn, nrhs_b,
476 & iend_blr+nelim_panel-npiv_global,
477 & mone, diag(1+diagsiz_dyn+npiv_global-iend_blr),
478 & diagsiz_sta,
479 & wcb(pwcb), ld_wcb,
480 & one, rhscomp(ipiv_panel,jbdeb), lrhscomp)
481 ELSE
482 CALL zgemm('T', 'N', diagsiz_dyn, nrhs_b, nelim_panel,
483 & mone, diag(1+diagsiz_dyn), diagsiz_sta,
484 & rhscomp(ipiv_panel+diagsiz_dyn,jbdeb), lrhscomp,
485 & one, rhscomp(ipiv_panel,jbdeb), lrhscomp)
486 ENDIF
487 ENDIF
488 ELSE
489 IF (iend_blr.EQ.npiv_global) THEN
490 CALL zgemm('N', 'N', diagsiz_dyn, nrhs_b, nelim_panel,
491 & mone, diag(1+diagsiz_dyn*ldadiag), diagsiz_dyn,
492 & wcb(pwcb), ld_wcb, one,
493 & rhscomp(ipiv_panel,jbdeb), lrhscomp)
494 ELSE
495 IF (iend_blr+1.LE.npiv_global .AND.
496 & iend_blr+nelim_panel.GT.npiv_global) THEN
497 CALL zgemm('N', 'N', diagsiz_dyn, nrhs_b,
498 & npiv_global-iend_blr,
499 & mone, diag(1+diagsiz_dyn*ldadiag), diagsiz_dyn,
500 & rhscomp(ipiv_panel+diagsiz_dyn,jbdeb), lrhscomp,
501 & one, rhscomp(ipiv_panel,jbdeb), lrhscomp)
502 CALL zgemm('N', 'N', diagsiz_dyn, nrhs_b,
503 & iend_blr+nelim_panel-npiv_global,
504 & mone, diag(1+diagsiz_dyn*ldadiag +
505 & (npiv_global-iend_blr)*diagsiz_dyn),
506 & diagsiz_dyn,
507 & wcb(pwcb), ld_wcb,
508 & one, rhscomp(ipiv_panel,jbdeb), lrhscomp)
509 ELSE
510 CALL zgemm('N', 'N', diagsiz_dyn, nrhs_b, nelim_panel,
511 & mone, diag(1+diagsiz_dyn*ldadiag), diagsiz_dyn,
512 & rhscomp(ipiv_panel+diagsiz_dyn,jbdeb), lrhscomp,
513 & one, rhscomp(ipiv_panel,jbdeb), lrhscomp)
514 ENDIF
515 ENDIF
516 ENDIF
517 ENDIF
518 IF (iflag.LT.0) RETURN
520 & diag(1), size(diag), diagsiz_dyn, nelim_panel, liell,
521 & nrhs_b, wcb, lwcb,
522 & rhscomp, lrhscomp, nrhs,
523 & ipiv_panel, jbdeb,
524 & mtype, keep )
525 ENDDO
526 RETURN
subroutine zmumps_solve_bwd_lr_trsolve(diag, ldiag, npiv, nelim, liell, nrhs_b, w, lwc, rhscomp, lrhscomp, nrhs, ppivinrhscomp, jbdeb, mtype, keep)
Definition zsol_lr.F:700

◆ zmumps_sol_fwd_blr_update()

subroutine zmumps_sol_lr::zmumps_sol_fwd_blr_update ( complex(kind=8), dimension(lpiv,lpivcol), intent(inout) arraypiv,
integer(8), intent(in) lpiv,
integer, intent(in) lpivcol,
integer, intent(in) ldpiv,
integer(8), intent(in) pospiv,
integer, intent(in) pospivcol,
complex(kind=8), dimension(lcb), intent(inout) arraycb,
integer(8), intent(in) lcb,
integer, intent(in) ldcb,
integer(8), intent(in) poscb,
integer(8), intent(in) posdiag,
integer, intent(in) nrhs_b,
integer, intent(in) npiv,
type(lrb_type), dimension(:), intent(in), target blr_panel,
integer, intent(in) last_blr,
integer, intent(in) current_blr,
integer, dimension(:) begs_blr_static,
integer(8), dimension(150), intent(in) keep8,
integer, intent(in) k34,
integer, intent(in) k450,
logical, intent(in) is_t2_slave,
integer, intent(inout) iflag,
integer, intent(inout) ierror )

Definition at line 246 of file zsol_lr.F.

254!$ USE OMP_LIB
255 INTEGER(8), INTENT(IN) :: LPIV, LCB, POSPIV, POSCB, POSDIAG
256 INTEGER, INTENT(IN) :: LPIVCOL, POSPIVCOL
257 COMPLEX(kind=8), INTENT(INOUT) :: ARRAYPIV(LPIV,LPIVCOL)
258 COMPLEX(kind=8), INTENT(INOUT) :: ARRAYCB(LCB)
259 INTEGER, INTENT(IN) :: LAST_BLR, NRHS_B, LDPIV, LDCB,
260 & CURRENT_BLR, NPIV, K34, K450
261 TYPE(LRB_TYPE), TARGET,INTENT(IN) ::
262 & BLR_PANEL(:)
263 INTEGER :: BEGS_BLR_STATIC(:)
264 LOGICAL, INTENT(IN) :: IS_T2_SLAVE
265 INTEGER(8), INTENT(IN) :: KEEP8(150)
266 INTEGER, INTENT(INOUT) :: IFLAG, IERROR
267 INTEGER :: I, K, M, N, IBEG_BLOCK, IEND_BLOCK
268 INTEGER :: KMAX
269 INTEGER(8) :: POSBLOCK
270 INTEGER :: allocok
271 TYPE(LRB_TYPE), POINTER :: LRB
272 COMPLEX(kind=8), ALLOCATABLE,DIMENSION(:) :: TEMP_BLOCK
273 COMPLEX(kind=8) :: ONE, MONE, ZERO
274 parameter(one=(1.0d0,0.0d0), mone=(-1.0d0,0.0d0))
275 parameter(zero=(0.0d0,0.0d0))
276 kmax = -1
277 DO i = current_blr+1, last_blr
278 kmax = max(kmax, blr_panel(i-current_blr)%K)
279 ENDDO
280#if defined(BLR_MT)
281!$OMP PARALLEL PRIVATE(TEMP_BLOCK, allocok)
282#endif
283 IF (kmax.GT.0) THEN
284 allocate(temp_block(kmax*nrhs_b), stat=allocok )
285 IF (allocok .GT. 0) THEN
286 iflag = -13
287 ierror = nrhs_b * kmax
288 write(*,*) 'Allocation problem in BLR routine
289 & ZMUMPS_SOL_FWD_BLR_UPDATE: ',
290 & 'not enough memory? memory requested = ', ierror
291 ENDIF
292 ENDIF
293#if defined(BLR_MT)
294!$OMP BARRIER
295#endif
296#if defined(BLR_MT)
297!$OMP DO SCHEDULE(DYNAMIC,1)
298!$OMP& PRIVATE(IBEG_BLOCK, IEND_BLOCK, LRB, K, M, N,
299!$OMP& POSBLOCK)
300#endif
301 DO i = current_blr+1, last_blr
302 IF (iflag.LT.0) cycle
303 ibeg_block = begs_blr_static(i)
304 iend_block = begs_blr_static(i+1)-1
305 IF (ibeg_block .EQ. iend_block + 1) cycle
306 lrb => blr_panel(i-current_blr)
307 k = lrb%K
308 m = lrb%M
309 n = lrb%N
310 IF (lrb%ISLR) THEN
311 IF (k.GT.0) THEN
312 CALL zgemm('N', 'N', k, nrhs_b, n, one,
313 & lrb%R(1,1), k, arraypiv(posdiag,pospivcol), ldpiv,
314 & zero, temp_block(1), k)
315 IF (is_t2_slave) THEN
316 posblock = poscb+int(ibeg_block-1,8)
317 CALL zgemm('N', 'N', m, nrhs_b, k, mone,
318 & lrb%Q(1,1), m, temp_block(1), k,
319 & one, arraycb(posblock), ldcb)
320 ELSEIF (ibeg_block.LE.npiv.AND.iend_block.GT.npiv) THEN
321 posblock = pospiv+int(ibeg_block-1,8)
322 CALL zgemm('N', 'N', npiv-ibeg_block+1, nrhs_b, k,
323 & mone, lrb%Q(1,1), m, temp_block(1), k,
324 & one, arraypiv(posblock,pospivcol), ldpiv)
325 CALL zgemm('N', 'N', ibeg_block+m-npiv-1, nrhs_b, k,
326 & mone, lrb%Q(npiv-ibeg_block+2,1), m, temp_block(1),
327 & k, one, arraycb(poscb), ldcb)
328 ELSEIF (ibeg_block.LE.npiv) THEN
329 posblock = pospiv+int(ibeg_block-1,8)
330 CALL zgemm('N', 'N', m, nrhs_b, k, mone,
331 & lrb%Q(1,1), m, temp_block(1), k,
332 & one, arraypiv(posblock,pospivcol), ldpiv)
333 ELSE
334 posblock = poscb+int(ibeg_block-1-npiv,8)
335 CALL zgemm('N', 'N', m, nrhs_b, k, mone,
336 & lrb%Q(1,1), m, temp_block(1), k,
337 & one, arraycb(posblock), ldcb)
338 ENDIF
339 ENDIF
340 ELSE
341 IF (is_t2_slave) THEN
342 posblock = poscb + int(ibeg_block-1,8)
343 CALL zgemm('N', 'N', m, nrhs_b, n, mone,
344 & lrb%Q(1,1), m, arraypiv(posdiag,pospivcol), ldpiv,
345 & one, arraycb(posblock), ldcb)
346 ELSEIF (ibeg_block.LE.npiv.AND.iend_block.GT.npiv) THEN
347 posblock = pospiv+int(ibeg_block-1,8)
348 CALL zgemm('N', 'N', npiv-ibeg_block+1, nrhs_b, n, mone,
349 & lrb%Q(1,1), m, arraypiv(posdiag,pospivcol), ldpiv,
350 & one, arraypiv(posblock,pospivcol), ldpiv)
351 CALL zgemm('N', 'N', ibeg_block+m-npiv-1, nrhs_b, n, mone,
352 & lrb%Q(npiv-ibeg_block+2,1), m,
353 & arraypiv(posdiag,pospivcol),
354 & ldpiv, one, arraycb(poscb), ldcb)
355 ELSEIF (ibeg_block.LE.npiv) THEN
356 posblock = pospiv+int(ibeg_block-1,8)
357 CALL zgemm('N', 'N', m, nrhs_b, n, mone,
358 & lrb%Q(1,1), m, arraypiv(posdiag,pospivcol), ldpiv,
359 & one, arraypiv(posblock,pospivcol), ldpiv)
360 ELSE
361 posblock = poscb + int(ibeg_block-1-npiv,8)
362 CALL zgemm('N', 'N', m, nrhs_b, n, mone,
363 & lrb%Q(1,1), m, arraypiv(posdiag,pospivcol), ldpiv,
364 & one, arraycb(posblock), ldcb)
365 ENDIF
366 ENDIF
367 ENDDO
368#if defined(BLR_MT)
369!$OMP END DO
370#endif
371 IF (kmax.GT.0) THEN
372 IF (allocated(temp_block)) deallocate(temp_block)
373 ENDIF
374#if defined(BLR_MT)
375!$OMP END PARALLEL
376#endif
377 RETURN

◆ zmumps_sol_fwd_lr_su()

subroutine zmumps_sol_lr::zmumps_sol_fwd_lr_su ( integer, intent(in) inode,
integer, intent(in) n,
integer, intent(in) iwhdlr,
integer, intent(in) npiv_global,
integer, intent(in) nslaves,
integer, dimension(liw), intent(in) iw,
integer, intent(in) ipos_init,
integer, intent(in) liw,
integer, intent(in) liell,
complex(kind=8), dimension(lwcb), intent(inout) wcb,
integer(8), intent(in) lwcb,
integer, intent(in) ld_wcbpiv,
integer, intent(in) ld_wcbcb,
integer(8), intent(in) ppiv_init,
integer(8), intent(in) pcb_init,
complex(kind=8), 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,
integer(8), dimension(150), intent(in) keep8,
logical, intent(in) oocwrite_compatible_with_blr,
integer, intent(inout) iflag,
integer, intent(inout) ierror )

Definition at line 21 of file zsol_lr.F.

31 INTEGER, INTENT(IN) :: INODE, N, IWHDLR, NPIV_GLOBAL, NSLAVES
32 INTEGER, INTENT(IN) :: MTYPE, LIELL, KEEP(500)
33 INTEGER(8), INTENT(IN) :: KEEP8(150)
34 INTEGER, INTENT(IN) :: LIW, IPOS_INIT, LRHSCOMP
35 INTEGER, INTENT(IN) :: IW(LIW), POSINRHSCOMP_FWD(N)
36 INTEGER(8), INTENT(IN) :: LWCB, PPIV_INIT, PCB_INIT
37 INTEGER, INTENT(IN) :: LD_WCBPIV, LD_WCBCB, NRHS, JBDEB, JBFIN
38 COMPLEX(kind=8), INTENT(INOUT) :: WCB(LWCB)
39 INTEGER, INTENT(INOUT) :: IFLAG, IERROR
40 COMPLEX(kind=8), INTENT(INOUT) :: RHSCOMP(LRHSCOMP, NRHS)
41 LOGICAL, INTENT(IN) :: OOCWRITE_COMPATIBLE_WITH_BLR
42 INTEGER :: I, NPARTSASS, NB_BLR , NELIM, LDADIAG,
43 & DIAGSIZ_DYN, DIAGSIZ_STA, IBEG_BLR, IEND_BLR,
44 & LD_CB, NELIM_GLOBAL, NRHS_B, IPOS, KCB
45 INTEGER(8) :: PPIV, PCB
46 INTEGER :: LAST_BLR
47 COMPLEX(kind=8), POINTER, DIMENSION(:) :: DIAG
48 TYPE(LRB_TYPE), POINTER, DIMENSION(:) :: BLR_PANEL
49 COMPLEX(kind=8) :: ONE, MONE, ZERO
50 parameter(one=(1.0d0,0.0d0), mone=(-1.0d0,0.0d0))
51 parameter(zero=(0.0d0,0.0d0))
52 nrhs_b = jbfin-jbdeb+1
53 IF (mtype.EQ.1) THEN
54 IF (associated(blr_array(iwhdlr)%PANELS_L))
55 & THEN
56 npartsass=size(blr_array(iwhdlr)%PANELS_L)
57 nb_blr = size(blr_array(iwhdlr)%BEGS_BLR_STATIC) -1
58 ELSE
59 WRITE(6,*) " Internal error in ZMUMPS_SOL_FWD_SU_MASTER"
60 ENDIF
61 ELSE
62 IF (associated(blr_array(iwhdlr)%PANELS_U))
63 & THEN
64 npartsass=size(blr_array(iwhdlr)%PANELS_U)
65 nb_blr = size(blr_array(iwhdlr)%BEGS_BLR_STATIC) -1
66 ENDIF
67 ENDIF
68 IF (nslaves.EQ.0 .OR. (keep(50).eq.0 .and. mtype .NE.1)) THEN
69 last_blr = nb_blr
70 ELSE
71 last_blr = npartsass
72 ENDIF
73 ipos = ipos_init
74 ppiv = ppiv_init
75 nelim_global =
76 & blr_array(iwhdlr)%BEGS_BLR_STATIC(npartsass+1)
77 & - blr_array(iwhdlr)%BEGS_BLR_DYNAMIC(npartsass+1)
78 DO i=1, npartsass
79 ibeg_blr = blr_array(iwhdlr)%BEGS_BLR_DYNAMIC(i)
80 iend_blr = blr_array(iwhdlr)%BEGS_BLR_DYNAMIC(i+1) -1
81 diagsiz_dyn = blr_array(iwhdlr)%BEGS_BLR_DYNAMIC(i+1) -
82 & ibeg_blr
83 diagsiz_sta = blr_array(iwhdlr)%BEGS_BLR_STATIC(i+1) -
84 & ibeg_blr
85 IF (keep(50).NE.0) THEN
86 ldadiag = diagsiz_dyn
87 ELSE
88 ldadiag = diagsiz_sta
89 ENDIF
90 IF (iend_blr.EQ.npiv_global) THEN
91 pcb = pcb_init
92 ELSE
93 pcb = ppiv + int(diagsiz_dyn,8)
94 ENDIF
95 IF ( diagsiz_dyn.EQ.0) cycle
96 nelim = diagsiz_sta - diagsiz_dyn
97 IF ( mtype .EQ. 1 ) THEN
98 blr_panel => blr_array(iwhdlr)%PANELS_L(i)%LRB_PANEL
99 ELSE
100 blr_panel => blr_array(iwhdlr)%PANELS_U(i)%LRB_PANEL
101 END IF
102 diag => blr_array(iwhdlr)%DIAG_BLOCKS(i)%DIAG_BLOCK
103 CALL zmumps_solve_fwd_trsolve (diag(1), int(size(diag),8), 1_8,
104 & diagsiz_dyn , ldadiag, nrhs_b, wcb, lwcb, npiv_global,
105 & ppiv, mtype, keep)
106 IF (nelim.GT.0) THEN
107 kcb = int(pcb-ppiv_init+1)
108 IF (iend_blr.EQ.npiv_global) THEN
109 ld_cb = ld_wcbcb
110 ELSE
111 ld_cb = ld_wcbpiv
112 ENDIF
113 IF (mtype.EQ.1) THEN
114 IF (kcb.LE.npiv_global .AND.
115 & kcb+nelim-1.GT.npiv_global) THEN
116 CALL zgemm('T', 'N', npiv_global-kcb+1, nrhs_b,
117 & diagsiz_dyn, mone,
118 & diag(1+diagsiz_dyn*ldadiag), diagsiz_dyn,
119 & wcb(ppiv), ld_wcbpiv,
120 & one, wcb(pcb), ld_cb)
121 CALL zgemm('T', 'N', kcb+nelim-npiv_global-1,
122 & nrhs_b, diagsiz_dyn, mone,
123 & diag(1+diagsiz_dyn*ldadiag +
124 & (npiv_global-kcb+1)*diagsiz_dyn),
125 & diagsiz_dyn,
126 & wcb(ppiv), ld_wcbpiv,
127 & one, wcb(pcb_init), ld_wcbcb)
128 ELSE
129 CALL zgemm('T', 'N', nelim, nrhs_b, diagsiz_dyn, mone,
130 & diag(1+diagsiz_dyn*ldadiag), diagsiz_dyn,
131 & wcb(ppiv), ld_wcbpiv,
132 & one, wcb(pcb), ld_cb)
133 ENDIF
134 ELSE
135 IF (kcb.LE.npiv_global .AND.
136 & kcb+nelim-1.GT.npiv_global) THEN
137 CALL zgemm('N', 'N', npiv_global-kcb+1,
138 & nrhs_b, diagsiz_dyn, mone,
139 & diag(1+diagsiz_dyn), diagsiz_sta,
140 & wcb(ppiv), ld_wcbpiv,
141 & one, wcb(pcb), ld_cb)
142 CALL zgemm('N', 'N', kcb+nelim-npiv_global-1,
143 & nrhs_b, diagsiz_dyn, mone,
144 & diag(1+diagsiz_dyn+npiv_global-kcb+1),
145 & diagsiz_sta,
146 & wcb(ppiv), ld_wcbpiv,
147 & one, wcb(pcb_init), ld_wcbcb)
148 ELSE
149 CALL zgemm('N', 'N', nelim, nrhs_b, diagsiz_dyn, mone,
150 & diag(1+diagsiz_dyn), diagsiz_sta,
151 & wcb(ppiv), ld_wcbpiv,
152 & one, wcb(pcb), ld_cb)
153 ENDIF
154 ENDIF
155 ENDIF
156 CALL zmumps_sol_fwd_blr_update (
157 & wcb, lwcb, 1, ld_wcbpiv, ppiv_init, 1,
158 & wcb, lwcb, ld_wcbcb, pcb_init,
159 & ppiv,
160 & nrhs_b, npiv_global,
161 & blr_panel, last_blr, i,
162 & blr_array(iwhdlr)%BEGS_BLR_STATIC,
163 & keep8, keep(34), keep(450), .false.,
164 & iflag, ierror)
165 IF (iflag.LT.0) RETURN
167 & inode, n, diagsiz_dyn, liell, nelim, nslaves,
168 & ppiv,
169 & iw, ipos, liw,
170 & diag(1), int(size(diag),8), 1_8,
171 & wcb, lwcb, ld_wcbpiv,
172 & rhscomp, lrhscomp, nrhs,
173 & posinrhscomp_fwd, jbdeb, jbfin,
174 & mtype, keep, oocwrite_compatible_with_blr,
175 & .true.
176 & )
177 ppiv = ppiv + int(diagsiz_dyn,8)
178 ipos = ipos + diagsiz_dyn
179 ENDDO
180 RETURN
subroutine zmumps_solve_fwd_trsolve(a, la, apos, npiv, ldadiag, nrhs_b, wcb, lwcb, lda_wcb, ppiv_courant, mtype, keep)
Definition zsol_aux.F:1148
subroutine zmumps_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)
Definition zsol_aux.F:1381

◆ zmumps_sol_slave_lr_u()

subroutine zmumps_sol_lr::zmumps_sol_slave_lr_u ( integer, intent(in) inode,
integer, intent(in) iwhdlr,
integer, intent(in) npiv_global,
complex(kind=8), dimension(lwcb), intent(inout) wcb,
integer(8), intent(in) lwcb,
integer, intent(in) ldx,
integer, intent(in) ldy,
integer(8), intent(in) ptrx_init,
integer(8), intent(in) ptry_init,
integer, intent(in) jbdeb,
integer, intent(in) jbfin,
integer, intent(in) mtype,
integer, dimension(500), intent(in) keep,
integer(8), dimension(150), intent(in) keep8,
integer, intent(inout) iflag,
integer, intent(inout) ierror )

Definition at line 182 of file zsol_lr.F.

189 INTEGER, INTENT(IN) :: INODE, IWHDLR, NPIV_GLOBAL
190 INTEGER, INTENT(IN) :: MTYPE, KEEP(500)
191 INTEGER(8), INTENT(IN) :: KEEP8(150)
192 INTEGER(8), INTENT(IN) :: LWCB, PTRX_INIT, PTRY_INIT
193 INTEGER, INTENT(IN) :: LDX, LDY, JBDEB, JBFIN
194 COMPLEX(kind=8), INTENT(INOUT) :: WCB(LWCB)
195 INTEGER, INTENT(INOUT) :: IFLAG, IERROR
196 INTEGER :: I, NPARTSASS, NB_BLR , NRHS_B
197 INTEGER(8) :: PTRX, PTRY
198 TYPE(LRB_TYPE), POINTER, DIMENSION(:) :: BLR_PANEL
199 COMPLEX(kind=8) :: ONE, MONE, ZERO
200 parameter(one=(1.0d0,0.0d0), mone=(-1.0d0,0.0d0))
201 parameter(zero=(0.0d0,0.0d0))
202 nrhs_b = jbfin-jbdeb+1
203 IF (associated(blr_array(iwhdlr)%PANELS_L))
204 & THEN
205 npartsass=size(blr_array(iwhdlr)%PANELS_L)
206 nb_blr = size(blr_array(iwhdlr)%BEGS_BLR_STATIC)
207 nb_blr = nb_blr - 2
208 ELSE
209 WRITE(6,*) " Internal error 1 in ZMUMPS_SOL_SLAVE_LR_U"
210 CALL mumps_abort()
211 ENDIF
212 ptrx = ptrx_init
213 ptry = ptry_init
214 DO i = 1, npartsass
215 blr_panel => blr_array(iwhdlr)%PANELS_L(i)%LRB_PANEL
216 IF (associated(blr_panel)) THEN
217 IF (mtype.EQ.1) THEN
218 CALL zmumps_sol_fwd_blr_update (
219 & wcb, lwcb, 1, ldx, -99999_8, 1,
220 & wcb, lwcb, ldy, ptry,
221 & ptrx,
222 & nrhs_b, npiv_global,
223 & blr_panel, nb_blr, 0,
224 & blr_array(iwhdlr)%BEGS_BLR_STATIC(2:nb_blr+2),
225 & keep8, keep(34), keep(450), .true., iflag, ierror )
226 ELSE
227 CALL zmumps_sol_bwd_blr_update (
228 & wcb, lwcb, 1, ldy, -99999_8, 1,
229 & wcb, lwcb, ldx, ptrx,
230 & ptry,
231 & nrhs_b, npiv_global,
232 & blr_panel, nb_blr, 0,
233 & blr_array(iwhdlr)%BEGS_BLR_STATIC(2:nb_blr+2),
234 & keep8, keep(34), keep(450), .true., iflag, ierror )
235 ENDIF
236 IF (mtype .EQ. 1) THEN
237 ptrx = ptrx + blr_panel(1)%N
238 ELSE
239 ptry = ptry + blr_panel(1)%N
240 ENDIF
241 IF (iflag.LT.0) RETURN
242 ENDIF
243 ENDDO
244 RETURN
#define mumps_abort
Definition VE_Metis.h:25