Fix complex version.

psblas3-type-indexed
Salvatore Filippone 19 years ago
parent 8ca0d6fa93
commit 2c21b017b3

@ -85,7 +85,7 @@ SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
#ifdef Have_SLU_ #ifdef Have_SLU_
#include "zsp_defs.h" #include "zsp_defs.h"
#undef Have_SLU_
#define HANDLE_SIZE 8 #define HANDLE_SIZE 8
/* kind of integer to hold a pointer. Use int. /* kind of integer to hold a pointer. Use int.
This might need to be changed on 64-bit systems. */ This might need to be changed on 64-bit systems. */
@ -102,7 +102,6 @@ typedef struct {
int *perm_r; int *perm_r;
} factors_t; } factors_t;
#undef Have_SLU_
#else #else
@ -132,7 +131,7 @@ typedef struct {
void void
psb_zslu_factor_(int *n, int *nnz, psb_zslu_factor_(int *n, int *nnz,
double *values, int *rowind, int *colptr, doublecomplex *values, int *rowind, int *colptr,
#ifdef Have_SLU_ #ifdef Have_SLU_
fptr *f_factors, /* a handle containing the address fptr *f_factors, /* a handle containing the address
pointing to the factored matrices */ pointing to the factored matrices */
@ -181,8 +180,8 @@ psb_zslu_factor_(int *n, int *nnz,
for (i = 0; i < *nnz; ++i) --colptr[i]; for (i = 0; i < *nnz; ++i) --colptr[i];
for (i = 0; i <= *n; ++i) --rowind[i]; for (i = 0; i <= *n; ++i) --rowind[i];
dCreate_CompRow_Matrix(&A, *n, *n, *nnz, values, colptr, rowind, zCreate_CompRow_Matrix(&A, *n, *n, *nnz, values, colptr, rowind,
SLU_NR, SLU_D, SLU_GE); SLU_NC, SLU_Z, SLU_GE);
L = (SuperMatrix *) SUPERLU_MALLOC( sizeof(SuperMatrix) ); L = (SuperMatrix *) SUPERLU_MALLOC( sizeof(SuperMatrix) );
U = (SuperMatrix *) SUPERLU_MALLOC( sizeof(SuperMatrix) ); U = (SuperMatrix *) SUPERLU_MALLOC( sizeof(SuperMatrix) );
if ( !(perm_r = intMalloc(*n)) ) ABORT("Malloc fails for perm_r[]."); if ( !(perm_r = intMalloc(*n)) ) ABORT("Malloc fails for perm_r[].");
@ -205,13 +204,13 @@ psb_zslu_factor_(int *n, int *nnz,
panel_size = sp_ienv(1); panel_size = sp_ienv(1);
relax = sp_ienv(2); relax = sp_ienv(2);
dgstrf(&options, &AC, drop_tol, relax, panel_size, zgstrf(&options, &AC, drop_tol, relax, panel_size,
etree, NULL, 0, perm_c, perm_r, L, U, &stat, info); etree, NULL, 0, perm_c, perm_r, L, U, &stat, info);
if ( *info == 0 ) { if ( *info == 0 ) {
Lstore = (SCformat *) L->Store; Lstore = (SCformat *) L->Store;
Ustore = (NCformat *) U->Store; Ustore = (NCformat *) U->Store;
dQuerySpace(L, U, &mem_usage); zQuerySpace(L, U, &mem_usage);
#if 0 #if 0
printf("No of nonzeros in factor L = %d\n", Lstore->nnz); printf("No of nonzeros in factor L = %d\n", Lstore->nnz);
printf("No of nonzeros in factor U = %d\n", Ustore->nnz); printf("No of nonzeros in factor U = %d\n", Ustore->nnz);
@ -223,7 +222,7 @@ psb_zslu_factor_(int *n, int *nnz,
} else { } else {
printf("dgstrf() error returns INFO= %d\n", *info); printf("dgstrf() error returns INFO= %d\n", *info);
if ( *info <= *n ) { /* factorization completes */ if ( *info <= *n ) { /* factorization completes */
dQuerySpace(L, U, &mem_usage); zQuerySpace(L, U, &mem_usage);
printf("L\\U MB %.3f\ttotal MB needed %.3f\texpansions %d\n", printf("L\\U MB %.3f\ttotal MB needed %.3f\texpansions %d\n",
mem_usage.for_lu/1e6, mem_usage.total_needed/1e6, mem_usage.for_lu/1e6, mem_usage.total_needed/1e6,
mem_usage.expansions); mem_usage.expansions);
@ -256,7 +255,7 @@ psb_zslu_factor_(int *n, int *nnz,
void void
psb_zslu_solve_(int *itrans, int *n, int *nrhs, psb_zslu_solve_(int *itrans, int *n, int *nrhs,
double *b, int *ldb, doublecomplex *b, int *ldb,
#ifdef Have_SLU_ #ifdef Have_SLU_
fptr *f_factors, /* a handle containing the address fptr *f_factors, /* a handle containing the address
pointing to the factored matrices */ pointing to the factored matrices */
@ -272,7 +271,7 @@ psb_zslu_solve_(int *itrans, int *n, int *nrhs,
* *
*/ */
#ifdef Have_SLU_ #ifdef Have_SLU_
SuperMatrix A, AC, B; SuperMatrix B;
SuperMatrix *L, *U; SuperMatrix *L, *U;
int *perm_r; /* row permutations from partial pivoting */ int *perm_r; /* row permutations from partial pivoting */
int *perm_c; /* column permutation vector */ int *perm_c; /* column permutation vector */
@ -306,9 +305,14 @@ psb_zslu_solve_(int *itrans, int *n, int *nrhs,
perm_c = LUfactors->perm_c; perm_c = LUfactors->perm_c;
perm_r = LUfactors->perm_r; perm_r = LUfactors->perm_r;
dCreate_Dense_Matrix(&B, *n, *nrhs, b, *ldb, SLU_DN, SLU_D, SLU_GE); zCreate_Dense_Matrix(&B, *n, *nrhs, b, *ldb, SLU_DN, SLU_Z, SLU_GE);
/* Solve the system A*X=B, overwriting B with X. */ /* Solve the system A*X=B, overwriting B with X. */
dgstrs (trans, L, U, perm_c, perm_r, &B, &stat, info); zgstrs (trans, L, U, perm_c, perm_r, &B, &stat, info);
if (info != 0) {
if (B.Stype != SLU_DN) fprintf(stderr,"zgstrs error kind 1: SLU_DN\n");
if (B.Dtype != SLU_Z) fprintf(stderr,"zgstrs error kind 2: SLU_Z\n");
if (B.Mtype != SLU_GE) fprintf(stderr,"zgstrs error kind 3: SLU_GE\n");
}
Destroy_SuperMatrix_Store(&B); Destroy_SuperMatrix_Store(&B);
StatFree(&stat); StatFree(&stat);

@ -262,7 +262,7 @@ subroutine psb_zspsm(alpha,a,x,beta,y,desc_a,info,&
if(info.ne.0) then if(info.ne.0) then
info = 4010 info = 4010
ch_err='dcssm' ch_err='zcssm'
call psb_errpush(info,name,a_err=ch_err) call psb_errpush(info,name,a_err=ch_err)
goto 9999 goto 9999
end if end if

@ -42,7 +42,7 @@ C
PARAMETER (DEBUG=.FALSE.) PARAMETER (DEBUG=.FALSE.)
C .. Local Arrays .. C .. Local Arrays ..
CHARACTER*20 NAME CHARACTER*20 NAME
INTEGER INT_VAL(5) INTEGER INT_VAL(5), err_Act
NAME = 'DCSRSM\0' NAME = 'DCSRSM\0'
IERROR = 0 IERROR = 0

@ -1,3 +1,33 @@
C
C Parallel Sparse BLAS v2.0
C (C) Copyright 2006 Salvatore Filippone University of Rome Tor Vergata
C Alfredo Buttari University of Rome Tor Vergata
C
C Redistribution and use in source and binary forms, with or without
C modification, are permitted provided that the following conditions
C are met:
C 1. Redistributions of source code must retain the above copyright
C notice, this list of conditions and the following disclaimer.
C 2. Redistributions in binary form must reproduce the above copyright
C notice, this list of conditions, and the following disclaimer in the
C documentation and/or other materials provided with the distribution.
C 3. The name of the PSBLAS group or the names of its contributors may
C not be used to endorse or promote products derived from this
C software without specific written permission.
C
C THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
C ``AS IS'' AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED
C TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
C PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE PSBLAS GROUP OR ITS CONTRIBUTORS
C BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
C CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
C SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
C INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
C CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
C ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
C POSSIBILITY OF SUCH DAMAGE.
C
C
SUBROUTINE ZCSRSM(TRANST,M,N,UNITD,D,ALPHA,DESCRA,A,JA,IA, SUBROUTINE ZCSRSM(TRANST,M,N,UNITD,D,ALPHA,DESCRA,A,JA,IA,
* B,LDB,BETA,C,LDC,WORK,LWORK) * B,LDB,BETA,C,LDC,WORK,LWORK)
COMPLEX*16 ALPHA, BETA COMPLEX*16 ALPHA, BETA
@ -8,24 +38,31 @@
CHARACTER DESCRA*11 CHARACTER DESCRA*11
INTEGER I, K INTEGER I, K
CHARACTER DIAG, UPLO CHARACTER DIAG, UPLO
EXTERNAL XERBLA
IF ((ALPHA.NE.(1.D0, 0.D0)).OR. C .. Local Arrays ..
+ (BETA.NE.(0.D0, 0.D0))) THEN CHARACTER*20 NAME
CALL XERBLA('ZCSSM ',9) INTEGER INT_VAL(5),err_act
RETURN
ENDIF NAME = 'DCSRSM\0'
IERROR = 0
CALL FCPSB_ERRACTIONSAVE(ERR_ACT)
int_Val(1)=0
UPLO = '?' UPLO = '?'
IF (DESCRA(1:1).EQ.'T' .AND. DESCRA(2:2).EQ.'U') UPLO = 'U' IF (DESCRA(1:1).EQ.'T' .AND. DESCRA(2:2).EQ.'U') UPLO = 'U'
IF (DESCRA(1:1).EQ.'T' .AND. DESCRA(2:2).EQ.'L') UPLO = 'L' IF (DESCRA(1:1).EQ.'T' .AND. DESCRA(2:2).EQ.'L') UPLO = 'L'
IF (UPLO.EQ.'?') THEN IF (UPLO.EQ.'?') THEN
CALL XERBLA('ZCSSM ',9) int_val(1) = 7
RETURN IERROR=5
CALL FCPSB_ERRPUSH(IERROR,NAME,INT_VAL)
GOTO 9999
END IF END IF
IF (DESCRA(3:3).EQ.'N') DIAG = 'N' IF (DESCRA(3:3).EQ.'N') DIAG = 'N'
IF (DESCRA(3:3).EQ.'U') DIAG = 'U' IF (DESCRA(3:3).EQ.'U') DIAG = 'U'
IF(UNITD.EQ.'B') THEN IF(UNITD.EQ.'B') THEN
CALL XERBLA('ZCSSM ',9) IERROR=5
RETURN int_val(1) = 4
CALL FCPSB_ERRPUSH(IERROR,NAME,INT_VAL)
GOTO 9999
ENDIF ENDIF
IF (UNITD.EQ.'R') THEN IF (UNITD.EQ.'R') THEN
DO 40 I = 1, N DO 40 I = 1, N
@ -35,9 +72,24 @@
40 CONTINUE 40 CONTINUE
END IF END IF
DO 60 I = 1, N if ((alpha.ne.(1.d0,0.d0)) .or.(beta .ne.(0.0d0,0.d0))) then
CALL ZSRSV(UPLO,TRANST,DIAG,M,A,JA,IA,B(1,I),C(1,I)) if (lwork .lt. m) then
60 CONTINUE int_val(1) = 17
IERROR=5
CALL FCPSB_ERRPUSH(IERROR,NAME,INT_VAL)
GOTO 9999
END IF
DO I = 1, N
CALL ZSRSV(UPLO,TRANST,DIAG,M,A,JA,IA,B(1,I),work)
do k=1,m
c(k,i) = beta*c(k,i) + alpha*work(k)
enddo
enddo
else
DO I = 1, N
CALL ZSRSV(UPLO,TRANST,DIAG,M,A,JA,IA,B(1,I),c(1,i))
enddo
endif
IF (UNITD.EQ.'L') THEN IF (UNITD.EQ.'L') THEN
DO 45 I = 1, N DO 45 I = 1, N
DO 25 K = 1, M DO 25 K = 1, M
@ -45,7 +97,17 @@
25 CONTINUE 25 CONTINUE
45 CONTINUE 45 CONTINUE
END IF END IF
CALL FCPSB_ERRACTIONRESTORE(ERR_ACT)
RETURN RETURN
END
9999 CONTINUE
CALL FCPSB_ERRACTIONRESTORE(ERR_ACT)
IF ( ERR_ACT .NE. 0 ) THEN
CALL FCPSB_SERROR()
RETURN
ENDIF
RETURN
END

@ -231,6 +231,7 @@ C
C Error handling C Error handling
C C
IF(IERROR.NE.0) THEN IF(IERROR.NE.0) THEN
write(0,*) 'zcssm ierror',ierror
CALL FCPSB_ERRPUSH(IERROR,NAME,INT_VAL) CALL FCPSB_ERRPUSH(IERROR,NAME,INT_VAL)
GOTO 9999 GOTO 9999
@ -336,6 +337,7 @@ C
9999 CONTINUE 9999 CONTINUE
CALL FCPSB_ERRACTIONRESTORE(ERR_ACT) CALL FCPSB_ERRACTIONRESTORE(ERR_ACT)
IF ( ERR_ACT .NE. 0 ) THEN IF ( ERR_ACT .NE. 0 ) THEN
CALL FCPSB_SERROR() CALL FCPSB_SERROR()
ENDIF ENDIF

Loading…
Cancel
Save