1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 367e560aaSBarry Smith /* 467e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 567e560aaSBarry Smith */ 6289bc588SBarry Smith 770f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h" 8b4c862e0SSatish Balay #include "petscblaslapack.h" 9289bc588SBarry Smith 104a2ae208SSatish Balay #undef __FUNCT__ 114a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense" 12f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 131987afe7SBarry Smith { 141987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 15f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 1613f74950SBarry Smith PetscInt j; 170805154bSBarry Smith PetscBLASInt N,m,ldax,lday,one = 1; 18efee365bSSatish Balay PetscErrorCode ierr; 193a40ed3dSBarry Smith 203a40ed3dSBarry Smith PetscFunctionBegin; 210805154bSBarry Smith N = PetscBLASIntCast(X->rmap.n*X->cmap.n); 220805154bSBarry Smith m = PetscBLASIntCast(X->rmap.n); 230805154bSBarry Smith ldax = PetscBLASIntCast(x->lda); 240805154bSBarry Smith lday = PetscBLASIntCast(y->lda); 25a5ce6ee0Svictorle if (ldax>m || lday>m) { 26899cda47SBarry Smith for (j=0; j<X->cmap.n; j++) { 27f4df32b1SMatthew Knepley BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one); 28a5ce6ee0Svictorle } 29a5ce6ee0Svictorle } else { 30f4df32b1SMatthew Knepley BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one); 31a5ce6ee0Svictorle } 32efee365bSSatish Balay ierr = PetscLogFlops(2*N-1);CHKERRQ(ierr); 333a40ed3dSBarry Smith PetscFunctionReturn(0); 341987afe7SBarry Smith } 351987afe7SBarry Smith 364a2ae208SSatish Balay #undef __FUNCT__ 374a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 38dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 39289bc588SBarry Smith { 40899cda47SBarry Smith PetscInt N = A->rmap.n*A->cmap.n; 413a40ed3dSBarry Smith 423a40ed3dSBarry Smith PetscFunctionBegin; 43899cda47SBarry Smith info->rows_global = (double)A->rmap.n; 44899cda47SBarry Smith info->columns_global = (double)A->cmap.n; 45899cda47SBarry Smith info->rows_local = (double)A->rmap.n; 46899cda47SBarry Smith info->columns_local = (double)A->cmap.n; 474e220ebcSLois Curfman McInnes info->block_size = 1.0; 484e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 496de62eeeSBarry Smith info->nz_used = (double)N; 506de62eeeSBarry Smith info->nz_unneeded = (double)0; 514e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 524e220ebcSLois Curfman McInnes info->mallocs = 0; 537adad957SLisandro Dalcin info->memory = ((PetscObject)A)->mem; 544e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 554e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 564e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 573a40ed3dSBarry Smith PetscFunctionReturn(0); 58289bc588SBarry Smith } 59289bc588SBarry Smith 604a2ae208SSatish Balay #undef __FUNCT__ 614a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense" 62f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha) 6380cd9d93SLois Curfman McInnes { 64273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 65f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 66efee365bSSatish Balay PetscErrorCode ierr; 670805154bSBarry Smith PetscBLASInt one = 1,j,nz,lda = PetscBLASIntCast(a->lda); 6880cd9d93SLois Curfman McInnes 693a40ed3dSBarry Smith PetscFunctionBegin; 70899cda47SBarry Smith if (lda>A->rmap.n) { 710805154bSBarry Smith nz = PetscBLASIntCast(A->rmap.n); 72899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 73f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v+j*lda,&one); 74a5ce6ee0Svictorle } 75a5ce6ee0Svictorle } else { 760805154bSBarry Smith nz = PetscBLASIntCast(A->rmap.n*A->cmap.n); 77f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v,&one); 78a5ce6ee0Svictorle } 79efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 803a40ed3dSBarry Smith PetscFunctionReturn(0); 8180cd9d93SLois Curfman McInnes } 8280cd9d93SLois Curfman McInnes 831cbb95d3SBarry Smith #undef __FUNCT__ 841cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense" 851cbb95d3SBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscTruth *fl) 861cbb95d3SBarry Smith { 871cbb95d3SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 881cbb95d3SBarry Smith PetscInt i,j,m = A->rmap.n,N; 891cbb95d3SBarry Smith PetscScalar *v = a->v; 901cbb95d3SBarry Smith 911cbb95d3SBarry Smith PetscFunctionBegin; 921cbb95d3SBarry Smith *fl = PETSC_FALSE; 931cbb95d3SBarry Smith if (A->rmap.n != A->cmap.n) PetscFunctionReturn(0); 941cbb95d3SBarry Smith N = a->lda; 951cbb95d3SBarry Smith 961cbb95d3SBarry Smith for (i=0; i<m; i++) { 971cbb95d3SBarry Smith for (j=i+1; j<m; j++) { 981cbb95d3SBarry Smith if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0); 991cbb95d3SBarry Smith } 1001cbb95d3SBarry Smith } 1011cbb95d3SBarry Smith *fl = PETSC_TRUE; 1021cbb95d3SBarry Smith PetscFunctionReturn(0); 1031cbb95d3SBarry Smith } 1041cbb95d3SBarry Smith 105289bc588SBarry Smith /* ---------------------------------------------------------------*/ 1066831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 1076831982aSBarry Smith rather than put it in the Mat->row slot.*/ 1084a2ae208SSatish Balay #undef __FUNCT__ 1094a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense" 110dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo) 111289bc588SBarry Smith { 112a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 11381f479a6Svictorle PetscFunctionBegin; 114a07ab388SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 115a07ab388SBarry Smith #else 116c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1176849ba73SBarry Smith PetscErrorCode ierr; 1180805154bSBarry Smith PetscBLASInt n,m,info; 11955659b69SBarry Smith 1203a40ed3dSBarry Smith PetscFunctionBegin; 1210805154bSBarry Smith n = PetscBLASIntCast(A->cmap.n); 1220805154bSBarry Smith m = PetscBLASIntCast(A->rmap.n); 123289bc588SBarry Smith if (!mat->pivots) { 124899cda47SBarry Smith ierr = PetscMalloc((A->rmap.n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr); 125899cda47SBarry Smith ierr = PetscLogObjectMemory(A,A->rmap.n*sizeof(PetscBLASInt));CHKERRQ(ierr); 126289bc588SBarry Smith } 127f1af5d2fSBarry Smith A->factor = FACTOR_LU; 128899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 12971044d3cSBarry Smith LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info); 13029bbc08cSBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 13129bbc08cSBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 132899cda47SBarry Smith ierr = PetscLogFlops((2*A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr); 133a07ab388SBarry Smith #endif 1343a40ed3dSBarry Smith PetscFunctionReturn(0); 135289bc588SBarry Smith } 1366ee01492SSatish Balay 137*b24902e0SBarry Smith 138*b24902e0SBarry Smith #undef __FUNCT__ 139*b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense" 140*b24902e0SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 141*b24902e0SBarry Smith { 142*b24902e0SBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 143*b24902e0SBarry Smith PetscErrorCode ierr; 144*b24902e0SBarry Smith PetscInt lda = (PetscInt)mat->lda,j,m; 145*b24902e0SBarry Smith Mat newi = *newmat; 146*b24902e0SBarry Smith 147*b24902e0SBarry Smith PetscFunctionBegin; 148*b24902e0SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 149*b24902e0SBarry Smith l = (Mat_SeqDense*)newi->data; 150*b24902e0SBarry Smith if (lda>A->rmap.n) { 151*b24902e0SBarry Smith m = A->rmap.n; 152*b24902e0SBarry Smith for (j=0; j<A->cmap.n; j++) { 153*b24902e0SBarry Smith ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 154*b24902e0SBarry Smith } 155*b24902e0SBarry Smith } else { 156*b24902e0SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 157*b24902e0SBarry Smith } 158*b24902e0SBarry Smith } 159*b24902e0SBarry Smith newi->assembled = PETSC_TRUE; 160*b24902e0SBarry Smith PetscFunctionReturn(0); 161*b24902e0SBarry Smith } 162*b24902e0SBarry Smith 1634a2ae208SSatish Balay #undef __FUNCT__ 1644a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 165dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 16602cad45dSBarry Smith { 16702cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 1686849ba73SBarry Smith PetscErrorCode ierr; 16913f74950SBarry Smith PetscInt lda = (PetscInt)mat->lda,j,m; 17002cad45dSBarry Smith Mat newi; 17102cad45dSBarry Smith 1723a40ed3dSBarry Smith PetscFunctionBegin; 1737adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&newi);CHKERRQ(ierr); 174899cda47SBarry Smith ierr = MatSetSizes(newi,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);CHKERRQ(ierr); 1757adad957SLisandro Dalcin ierr = MatSetType(newi,((PetscObject)A)->type_name);CHKERRQ(ierr); 1765c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr); 17702cad45dSBarry Smith *newmat = newi; 178*b24902e0SBarry Smith ierr = MatDuplicateNoCreate_SeqDense(A,cpvalues,newmat);CHKERRQ(ierr); 179*b24902e0SBarry Smith PetscFunctionReturn(0); 180*b24902e0SBarry Smith } 181*b24902e0SBarry Smith 182*b24902e0SBarry Smith #undef __FUNCT__ 183*b24902e0SBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 184*b24902e0SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact) 185*b24902e0SBarry Smith { 186*b24902e0SBarry Smith PetscErrorCode ierr; 187*b24902e0SBarry Smith 188*b24902e0SBarry Smith PetscFunctionBegin; 189*b24902e0SBarry Smith ierr = MatDuplicateNoCreate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 1903a40ed3dSBarry Smith PetscFunctionReturn(0); 19102cad45dSBarry Smith } 19202cad45dSBarry Smith 1934a2ae208SSatish Balay #undef __FUNCT__ 1944a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 195dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact) 196289bc588SBarry Smith { 197dfbe8321SBarry Smith PetscErrorCode ierr; 1983a40ed3dSBarry Smith 1993a40ed3dSBarry Smith PetscFunctionBegin; 200*b24902e0SBarry Smith ierr = MatDuplicateNoCreate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 201*b24902e0SBarry Smith 202*b24902e0SBarry Smith PetscFunctionReturn(0); 203*b24902e0SBarry Smith } 204*b24902e0SBarry Smith 205*b24902e0SBarry Smith #undef __FUNCT__ 206*b24902e0SBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc" 207*b24902e0SBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,Mat *fact) 208*b24902e0SBarry Smith { 209*b24902e0SBarry Smith PetscErrorCode ierr; 210*b24902e0SBarry Smith 211*b24902e0SBarry Smith PetscFunctionBegin; 212*b24902e0SBarry Smith ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr); 213*b24902e0SBarry Smith ierr = MatSetSizes(*fact,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);CHKERRQ(ierr); 214*b24902e0SBarry Smith ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr); 2153a40ed3dSBarry Smith PetscFunctionReturn(0); 216289bc588SBarry Smith } 2176ee01492SSatish Balay 2184a2ae208SSatish Balay #undef __FUNCT__ 2194a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 220af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact) 221289bc588SBarry Smith { 22202cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data; 2236849ba73SBarry Smith PetscErrorCode ierr; 224899cda47SBarry Smith PetscInt lda1=mat->lda,lda2=l->lda, m=A->rmap.n,n=A->cmap.n, j; 2254482741eSBarry Smith MatFactorInfo info; 2263a40ed3dSBarry Smith 2273a40ed3dSBarry Smith PetscFunctionBegin; 22802cad45dSBarry Smith /* copy the numerical values */ 2291b807ce4Svictorle if (lda1>m || lda2>m ) { 2301b807ce4Svictorle for (j=0; j<n; j++) { 2311b807ce4Svictorle ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 2321b807ce4Svictorle } 2331b807ce4Svictorle } else { 234899cda47SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 2351b807ce4Svictorle } 23602cad45dSBarry Smith (*fact)->factor = 0; 2374482741eSBarry Smith ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr); 2383a40ed3dSBarry Smith PetscFunctionReturn(0); 239289bc588SBarry Smith } 2406ee01492SSatish Balay 2414a2ae208SSatish Balay #undef __FUNCT__ 2424a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense" 243dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo) 244289bc588SBarry Smith { 245a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 246a07ab388SBarry Smith PetscFunctionBegin; 247a07ab388SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 248a07ab388SBarry Smith #else 249c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2506849ba73SBarry Smith PetscErrorCode ierr; 2510805154bSBarry Smith PetscBLASInt info,n = PetscBLASIntCast(A->cmap.n); 25255659b69SBarry Smith 2533a40ed3dSBarry Smith PetscFunctionBegin; 254606d414cSSatish Balay ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 255752f5567SLois Curfman McInnes mat->pivots = 0; 25605b42c5fSBarry Smith 257899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 25871044d3cSBarry Smith LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info); 25977431f27SBarry Smith if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 260c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 261899cda47SBarry Smith ierr = PetscLogFlops((A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr); 262a07ab388SBarry Smith #endif 2633a40ed3dSBarry Smith PetscFunctionReturn(0); 264289bc588SBarry Smith } 265289bc588SBarry Smith 2664a2ae208SSatish Balay #undef __FUNCT__ 2670b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 268af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact) 2690b4b3355SBarry Smith { 270dfbe8321SBarry Smith PetscErrorCode ierr; 27115e8a5b3SHong Zhang MatFactorInfo info; 2720b4b3355SBarry Smith 2730b4b3355SBarry Smith PetscFunctionBegin; 27415e8a5b3SHong Zhang info.fill = 1.0; 27515e8a5b3SHong Zhang ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr); 2760b4b3355SBarry Smith PetscFunctionReturn(0); 2770b4b3355SBarry Smith } 2780b4b3355SBarry Smith 2790b4b3355SBarry Smith #undef __FUNCT__ 2804a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 281dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 282289bc588SBarry Smith { 283c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2846849ba73SBarry Smith PetscErrorCode ierr; 28587828ca2SBarry Smith PetscScalar *x,*y; 2860805154bSBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap.n); 28767e560aaSBarry Smith 2883a40ed3dSBarry Smith PetscFunctionBegin; 2891ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2901ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 291899cda47SBarry Smith ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 292c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 293ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 29480444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 295ae7cfcebSSatish Balay #else 29671044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 2974ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve"); 298ae7cfcebSSatish Balay #endif 2997a97a34bSBarry Smith } else if (A->factor == FACTOR_CHOLESKY){ 300ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 30180444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 302ae7cfcebSSatish Balay #else 30371044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 3044ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve"); 305ae7cfcebSSatish Balay #endif 306289bc588SBarry Smith } 30729bbc08cSBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 3081ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3091ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 310899cda47SBarry Smith ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr); 3113a40ed3dSBarry Smith PetscFunctionReturn(0); 312289bc588SBarry Smith } 3136ee01492SSatish Balay 3144a2ae208SSatish Balay #undef __FUNCT__ 3154a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 316dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 317da3a660dSBarry Smith { 318c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 319dfbe8321SBarry Smith PetscErrorCode ierr; 32087828ca2SBarry Smith PetscScalar *x,*y; 3210805154bSBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap.n); 32267e560aaSBarry Smith 3233a40ed3dSBarry Smith PetscFunctionBegin; 3241ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3251ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 326899cda47SBarry Smith ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 327752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 328da3a660dSBarry Smith if (mat->pivots) { 329ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 33080444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 331ae7cfcebSSatish Balay #else 33271044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 3334ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve"); 334ae7cfcebSSatish Balay #endif 3357a97a34bSBarry Smith } else { 336ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 33780444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 338ae7cfcebSSatish Balay #else 33971044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 3404ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve"); 341ae7cfcebSSatish Balay #endif 342da3a660dSBarry Smith } 3431ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3441ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 345899cda47SBarry Smith ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr); 3463a40ed3dSBarry Smith PetscFunctionReturn(0); 347da3a660dSBarry Smith } 3486ee01492SSatish Balay 3494a2ae208SSatish Balay #undef __FUNCT__ 3504a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 351dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 352da3a660dSBarry Smith { 353c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 354dfbe8321SBarry Smith PetscErrorCode ierr; 35587828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 356da3a660dSBarry Smith Vec tmp = 0; 3570805154bSBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap.n); 35867e560aaSBarry Smith 3593a40ed3dSBarry Smith PetscFunctionBegin; 3601ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3611ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 362899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 363da3a660dSBarry Smith if (yy == zz) { 36478b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 36552e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 36678b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 367da3a660dSBarry Smith } 368899cda47SBarry Smith ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 369752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 370da3a660dSBarry Smith if (mat->pivots) { 371ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 37280444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 373ae7cfcebSSatish Balay #else 37471044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 3752ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 376ae7cfcebSSatish Balay #endif 377a8c6a408SBarry Smith } else { 378ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 37980444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 380ae7cfcebSSatish Balay #else 38171044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 3822ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 383ae7cfcebSSatish Balay #endif 384da3a660dSBarry Smith } 3852dcb1b2aSMatthew Knepley if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 3862dcb1b2aSMatthew Knepley else {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);} 3871ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3881ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 389899cda47SBarry Smith ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr); 3903a40ed3dSBarry Smith PetscFunctionReturn(0); 391da3a660dSBarry Smith } 39267e560aaSBarry Smith 3934a2ae208SSatish Balay #undef __FUNCT__ 3944a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 395dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 396da3a660dSBarry Smith { 397c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3986849ba73SBarry Smith PetscErrorCode ierr; 39987828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 400da3a660dSBarry Smith Vec tmp; 4010805154bSBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap.n); 40267e560aaSBarry Smith 4033a40ed3dSBarry Smith PetscFunctionBegin; 404899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 4051ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4061ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 407da3a660dSBarry Smith if (yy == zz) { 40878b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 40952e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 41078b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 411da3a660dSBarry Smith } 412899cda47SBarry Smith ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 413752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 414da3a660dSBarry Smith if (mat->pivots) { 415ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 41680444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 417ae7cfcebSSatish Balay #else 41871044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 4192ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 420ae7cfcebSSatish Balay #endif 4213a40ed3dSBarry Smith } else { 422ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 42380444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 424ae7cfcebSSatish Balay #else 42571044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 4262ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 427ae7cfcebSSatish Balay #endif 428da3a660dSBarry Smith } 42990f02eecSBarry Smith if (tmp) { 4302dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 43190f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 4323a40ed3dSBarry Smith } else { 4332dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 43490f02eecSBarry Smith } 4351ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4361ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 437899cda47SBarry Smith ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr); 4383a40ed3dSBarry Smith PetscFunctionReturn(0); 439da3a660dSBarry Smith } 440289bc588SBarry Smith /* ------------------------------------------------------------------*/ 4414a2ae208SSatish Balay #undef __FUNCT__ 4424a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense" 44313f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 444289bc588SBarry Smith { 445c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 44687828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 447dfbe8321SBarry Smith PetscErrorCode ierr; 448899cda47SBarry Smith PetscInt m = A->rmap.n,i; 449aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 4500805154bSBarry Smith PetscBLASInt o = 1,bm = PetscBLASIntCast(m); 451bc1b551cSSatish Balay #endif 452289bc588SBarry Smith 4533a40ed3dSBarry Smith PetscFunctionBegin; 454289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 45571044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 4562dcb1b2aSMatthew Knepley ierr = VecSet(xx,zero);CHKERRQ(ierr); 457289bc588SBarry Smith } 4581ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4591ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 460b965ef7fSBarry Smith its = its*lits; 46177431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 462289bc588SBarry Smith while (its--) { 463fccaa45eSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 464289bc588SBarry Smith for (i=0; i<m; i++) { 465aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 466f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 467f1747703SBarry Smith not happy about returning a double complex */ 46813f74950SBarry Smith PetscInt _i; 46987828ca2SBarry Smith PetscScalar sum = b[i]; 470f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4713f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 472f1747703SBarry Smith } 473f1747703SBarry Smith xt = sum; 474f1747703SBarry Smith #else 47571044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 476f1747703SBarry Smith #endif 47755a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 478289bc588SBarry Smith } 479289bc588SBarry Smith } 480fccaa45eSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 481289bc588SBarry Smith for (i=m-1; i>=0; i--) { 482aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 483f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 484f1747703SBarry Smith not happy about returning a double complex */ 48513f74950SBarry Smith PetscInt _i; 48687828ca2SBarry Smith PetscScalar sum = b[i]; 487f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4883f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 489f1747703SBarry Smith } 490f1747703SBarry Smith xt = sum; 491f1747703SBarry Smith #else 49271044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 493f1747703SBarry Smith #endif 49455a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 495289bc588SBarry Smith } 496289bc588SBarry Smith } 497289bc588SBarry Smith } 4981ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 4991ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5003a40ed3dSBarry Smith PetscFunctionReturn(0); 501289bc588SBarry Smith } 502289bc588SBarry Smith 503289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5044a2ae208SSatish Balay #undef __FUNCT__ 5054a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 506dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 507289bc588SBarry Smith { 508c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 50987828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 510dfbe8321SBarry Smith PetscErrorCode ierr; 5110805154bSBarry Smith PetscBLASInt m, n,_One=1; 512ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 5133a40ed3dSBarry Smith 5143a40ed3dSBarry Smith PetscFunctionBegin; 5150805154bSBarry Smith m = PetscBLASIntCast(A->rmap.n); 5160805154bSBarry Smith n = PetscBLASIntCast(A->cmap.n); 517899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 5181ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5191ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 52071044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One); 5211ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5221ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 523899cda47SBarry Smith ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr); 5243a40ed3dSBarry Smith PetscFunctionReturn(0); 525289bc588SBarry Smith } 5266ee01492SSatish Balay 5274a2ae208SSatish Balay #undef __FUNCT__ 5284a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 529dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 530289bc588SBarry Smith { 531c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 53287828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 533dfbe8321SBarry Smith PetscErrorCode ierr; 5340805154bSBarry Smith PetscBLASInt m, n, _One=1; 5353a40ed3dSBarry Smith 5363a40ed3dSBarry Smith PetscFunctionBegin; 5370805154bSBarry Smith m = PetscBLASIntCast(A->rmap.n); 5380805154bSBarry Smith n = PetscBLASIntCast(A->cmap.n); 539899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 5401ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5411ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 54271044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 5431ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5441ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 545899cda47SBarry Smith ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->rmap.n);CHKERRQ(ierr); 5463a40ed3dSBarry Smith PetscFunctionReturn(0); 547289bc588SBarry Smith } 5486ee01492SSatish Balay 5494a2ae208SSatish Balay #undef __FUNCT__ 5504a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 551dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 552289bc588SBarry Smith { 553c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 55487828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 555dfbe8321SBarry Smith PetscErrorCode ierr; 5560805154bSBarry Smith PetscBLASInt m, n, _One=1; 5573a40ed3dSBarry Smith 5583a40ed3dSBarry Smith PetscFunctionBegin; 5590805154bSBarry Smith m = PetscBLASIntCast(A->rmap.n); 5600805154bSBarry Smith n = PetscBLASIntCast(A->cmap.n); 561899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 562600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5631ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5641ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 56571044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5661ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5671ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 568899cda47SBarry Smith ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr); 5693a40ed3dSBarry Smith PetscFunctionReturn(0); 570289bc588SBarry Smith } 5716ee01492SSatish Balay 5724a2ae208SSatish Balay #undef __FUNCT__ 5734a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 574dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 575289bc588SBarry Smith { 576c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 57787828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 578dfbe8321SBarry Smith PetscErrorCode ierr; 5790805154bSBarry Smith PetscBLASInt m, n, _One=1; 58087828ca2SBarry Smith PetscScalar _DOne=1.0; 5813a40ed3dSBarry Smith 5823a40ed3dSBarry Smith PetscFunctionBegin; 5830805154bSBarry Smith m = PetscBLASIntCast(A->rmap.n); 5840805154bSBarry Smith n = PetscBLASIntCast(A->cmap.n); 585899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 586600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5871ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5881ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 58971044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5901ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5911ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 592899cda47SBarry Smith ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr); 5933a40ed3dSBarry Smith PetscFunctionReturn(0); 594289bc588SBarry Smith } 595289bc588SBarry Smith 596289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5974a2ae208SSatish Balay #undef __FUNCT__ 5984a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 59913f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 600289bc588SBarry Smith { 601c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 60287828ca2SBarry Smith PetscScalar *v; 6036849ba73SBarry Smith PetscErrorCode ierr; 60413f74950SBarry Smith PetscInt i; 60567e560aaSBarry Smith 6063a40ed3dSBarry Smith PetscFunctionBegin; 607899cda47SBarry Smith *ncols = A->cmap.n; 608289bc588SBarry Smith if (cols) { 609899cda47SBarry Smith ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr); 610899cda47SBarry Smith for (i=0; i<A->cmap.n; i++) (*cols)[i] = i; 611289bc588SBarry Smith } 612289bc588SBarry Smith if (vals) { 613899cda47SBarry Smith ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 614289bc588SBarry Smith v = mat->v + row; 615899cda47SBarry Smith for (i=0; i<A->cmap.n; i++) {(*vals)[i] = *v; v += mat->lda;} 616289bc588SBarry Smith } 6173a40ed3dSBarry Smith PetscFunctionReturn(0); 618289bc588SBarry Smith } 6196ee01492SSatish Balay 6204a2ae208SSatish Balay #undef __FUNCT__ 6214a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 62213f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 623289bc588SBarry Smith { 624dfbe8321SBarry Smith PetscErrorCode ierr; 625606d414cSSatish Balay PetscFunctionBegin; 626606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 627606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 6283a40ed3dSBarry Smith PetscFunctionReturn(0); 629289bc588SBarry Smith } 630289bc588SBarry Smith /* ----------------------------------------------------------------*/ 6314a2ae208SSatish Balay #undef __FUNCT__ 6324a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 63313f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 634289bc588SBarry Smith { 635c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 63613f74950SBarry Smith PetscInt i,j,idx=0; 637d6dfbf8fSBarry Smith 6383a40ed3dSBarry Smith PetscFunctionBegin; 639289bc588SBarry Smith if (!mat->roworiented) { 640dbb450caSBarry Smith if (addv == INSERT_VALUES) { 641289bc588SBarry Smith for (j=0; j<n; j++) { 642cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6432515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 644899cda47SBarry Smith if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1); 64558804f6dSBarry Smith #endif 646289bc588SBarry Smith for (i=0; i<m; i++) { 647cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6482515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 649899cda47SBarry Smith if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1); 65058804f6dSBarry Smith #endif 651cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 652289bc588SBarry Smith } 653289bc588SBarry Smith } 6543a40ed3dSBarry Smith } else { 655289bc588SBarry Smith for (j=0; j<n; j++) { 656cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6572515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 658899cda47SBarry Smith if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1); 65958804f6dSBarry Smith #endif 660289bc588SBarry Smith for (i=0; i<m; i++) { 661cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6622515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 663899cda47SBarry Smith if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1); 66458804f6dSBarry Smith #endif 665cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 666289bc588SBarry Smith } 667289bc588SBarry Smith } 668289bc588SBarry Smith } 6693a40ed3dSBarry Smith } else { 670dbb450caSBarry Smith if (addv == INSERT_VALUES) { 671e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 672cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 6732515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 674899cda47SBarry Smith if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1); 67558804f6dSBarry Smith #endif 676e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 677cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 6782515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 679899cda47SBarry Smith if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1); 68058804f6dSBarry Smith #endif 681cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 682e8d4e0b9SBarry Smith } 683e8d4e0b9SBarry Smith } 6843a40ed3dSBarry Smith } else { 685289bc588SBarry Smith for (i=0; i<m; i++) { 686cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 6872515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 688899cda47SBarry Smith if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap.n-1); 68958804f6dSBarry Smith #endif 690289bc588SBarry Smith for (j=0; j<n; j++) { 691cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 6922515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 693899cda47SBarry Smith if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap.n-1); 69458804f6dSBarry Smith #endif 695cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 696289bc588SBarry Smith } 697289bc588SBarry Smith } 698289bc588SBarry Smith } 699e8d4e0b9SBarry Smith } 7003a40ed3dSBarry Smith PetscFunctionReturn(0); 701289bc588SBarry Smith } 702e8d4e0b9SBarry Smith 7034a2ae208SSatish Balay #undef __FUNCT__ 7044a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 70513f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 706ae80bb75SLois Curfman McInnes { 707ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 70813f74950SBarry Smith PetscInt i,j; 709ae80bb75SLois Curfman McInnes 7103a40ed3dSBarry Smith PetscFunctionBegin; 711ae80bb75SLois Curfman McInnes /* row-oriented output */ 712ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 71397e567efSBarry Smith if (indexm[i] < 0) {v += n;continue;} 714f12f8aa1SBarry Smith if (indexm[i] >= A->rmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap.n); 715ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 7166f31f424SBarry Smith if (indexn[j] < 0) {v++; continue;} 7176f31f424SBarry Smith if (indexn[j] >= A->cmap.n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap.n); 71897e567efSBarry Smith *v++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 719ae80bb75SLois Curfman McInnes } 720ae80bb75SLois Curfman McInnes } 7213a40ed3dSBarry Smith PetscFunctionReturn(0); 722ae80bb75SLois Curfman McInnes } 723ae80bb75SLois Curfman McInnes 724289bc588SBarry Smith /* -----------------------------------------------------------------*/ 725289bc588SBarry Smith 726e090d566SSatish Balay #include "petscsys.h" 72788e32edaSLois Curfman McInnes 7284a2ae208SSatish Balay #undef __FUNCT__ 7294a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense" 730f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, MatType type,Mat *A) 73188e32edaSLois Curfman McInnes { 73288e32edaSLois Curfman McInnes Mat_SeqDense *a; 73388e32edaSLois Curfman McInnes Mat B; 7346849ba73SBarry Smith PetscErrorCode ierr; 73513f74950SBarry Smith PetscInt *scols,i,j,nz,header[4]; 73613f74950SBarry Smith int fd; 73713f74950SBarry Smith PetscMPIInt size; 73813f74950SBarry Smith PetscInt *rowlengths = 0,M,N,*cols; 73987828ca2SBarry Smith PetscScalar *vals,*svals,*v,*w; 74019bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 74188e32edaSLois Curfman McInnes 7423a40ed3dSBarry Smith PetscFunctionBegin; 743d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 74429bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 745b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 7460752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 747552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 74888e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 74988e32edaSLois Curfman McInnes 75090ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 751f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 752f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr); 753be5d1d56SKris Buschelman ierr = MatSetType(*A,type);CHKERRQ(ierr); 7545c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr); 75590ace30eSBarry Smith B = *A; 75690ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 7574a41a4d3SSatish Balay v = a->v; 7584a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 7594a41a4d3SSatish Balay from row major to column major */ 760b72907f3SBarry Smith ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 76190ace30eSBarry Smith /* read in nonzero values */ 7624a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 7634a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 7644a41a4d3SSatish Balay for (j=0; j<N; j++) { 7654a41a4d3SSatish Balay for (i=0; i<M; i++) { 7664a41a4d3SSatish Balay *v++ =w[i*N+j]; 7674a41a4d3SSatish Balay } 7684a41a4d3SSatish Balay } 769606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 7706d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7716d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 77290ace30eSBarry Smith } else { 77388e32edaSLois Curfman McInnes /* read row lengths */ 77413f74950SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 7750752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 77688e32edaSLois Curfman McInnes 77788e32edaSLois Curfman McInnes /* create our matrix */ 778f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 779f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr); 780be5d1d56SKris Buschelman ierr = MatSetType(*A,type);CHKERRQ(ierr); 7815c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr); 78288e32edaSLois Curfman McInnes B = *A; 78388e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 78488e32edaSLois Curfman McInnes v = a->v; 78588e32edaSLois Curfman McInnes 78688e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 78713f74950SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr); 788b0a32e0cSBarry Smith cols = scols; 7890752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 79087828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 791b0a32e0cSBarry Smith vals = svals; 7920752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 79388e32edaSLois Curfman McInnes 79488e32edaSLois Curfman McInnes /* insert into matrix */ 79588e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 79688e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 79788e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 79888e32edaSLois Curfman McInnes } 799606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 800606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 801606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 80288e32edaSLois Curfman McInnes 8036d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8046d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 80590ace30eSBarry Smith } 8063a40ed3dSBarry Smith PetscFunctionReturn(0); 80788e32edaSLois Curfman McInnes } 80888e32edaSLois Curfman McInnes 809e090d566SSatish Balay #include "petscsys.h" 810932b0c3eSLois Curfman McInnes 8114a2ae208SSatish Balay #undef __FUNCT__ 8124a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 8136849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 814289bc588SBarry Smith { 815932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 816dfbe8321SBarry Smith PetscErrorCode ierr; 81713f74950SBarry Smith PetscInt i,j; 8182dcb1b2aSMatthew Knepley const char *name; 81987828ca2SBarry Smith PetscScalar *v; 820f3ef73ceSBarry Smith PetscViewerFormat format; 8215f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX) 8225f481a85SSatish Balay PetscTruth allreal = PETSC_TRUE; 8235f481a85SSatish Balay #endif 824932b0c3eSLois Curfman McInnes 8253a40ed3dSBarry Smith PetscFunctionBegin; 826435da068SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 827b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 828456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 8293a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 830fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 831b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 832899cda47SBarry Smith for (i=0; i<A->rmap.n; i++) { 83344cd7ae7SLois Curfman McInnes v = a->v + i; 83477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 835899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 836aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 837329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 838a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 839329f5518SBarry Smith } else if (PetscRealPart(*v)) { 840a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr); 8416831982aSBarry Smith } 84280cd9d93SLois Curfman McInnes #else 8436831982aSBarry Smith if (*v) { 844a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr); 8456831982aSBarry Smith } 84680cd9d93SLois Curfman McInnes #endif 8471b807ce4Svictorle v += a->lda; 84880cd9d93SLois Curfman McInnes } 849b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 85080cd9d93SLois Curfman McInnes } 851b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 8523a40ed3dSBarry Smith } else { 853b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 854aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 85547989497SBarry Smith /* determine if matrix has all real values */ 85647989497SBarry Smith v = a->v; 857899cda47SBarry Smith for (i=0; i<A->rmap.n*A->cmap.n; i++) { 858ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 85947989497SBarry Smith } 86047989497SBarry Smith #endif 861fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 8623a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 863899cda47SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap.n,A->cmap.n);CHKERRQ(ierr); 864899cda47SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap.n,A->cmap.n);CHKERRQ(ierr); 865fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 866ffac6cdbSBarry Smith } 867ffac6cdbSBarry Smith 868899cda47SBarry Smith for (i=0; i<A->rmap.n; i++) { 869932b0c3eSLois Curfman McInnes v = a->v + i; 870899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 871aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 87247989497SBarry Smith if (allreal) { 873f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr); 87447989497SBarry Smith } else { 875f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 87647989497SBarry Smith } 877289bc588SBarry Smith #else 878f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr); 879289bc588SBarry Smith #endif 8801b807ce4Svictorle v += a->lda; 881289bc588SBarry Smith } 882b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 883289bc588SBarry Smith } 884fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 885b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 886ffac6cdbSBarry Smith } 887b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 888da3a660dSBarry Smith } 889b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 8903a40ed3dSBarry Smith PetscFunctionReturn(0); 891289bc588SBarry Smith } 892289bc588SBarry Smith 8934a2ae208SSatish Balay #undef __FUNCT__ 8944a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 8956849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 896932b0c3eSLois Curfman McInnes { 897932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 8986849ba73SBarry Smith PetscErrorCode ierr; 89913f74950SBarry Smith int fd; 900899cda47SBarry Smith PetscInt ict,j,n = A->cmap.n,m = A->rmap.n,i,*col_lens,nz = m*n; 90187828ca2SBarry Smith PetscScalar *v,*anonz,*vals; 902f3ef73ceSBarry Smith PetscViewerFormat format; 903932b0c3eSLois Curfman McInnes 9043a40ed3dSBarry Smith PetscFunctionBegin; 905b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 90690ace30eSBarry Smith 907b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 908fb9695e5SSatish Balay if (format == PETSC_VIEWER_BINARY_NATIVE) { 90990ace30eSBarry Smith /* store the matrix as a dense matrix */ 91013f74950SBarry Smith ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 911552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 91290ace30eSBarry Smith col_lens[1] = m; 91390ace30eSBarry Smith col_lens[2] = n; 91490ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 9156f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 916606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 91790ace30eSBarry Smith 91890ace30eSBarry Smith /* write out matrix, by rows */ 91987828ca2SBarry Smith ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 92090ace30eSBarry Smith v = a->v; 92190ace30eSBarry Smith for (j=0; j<n; j++) { 922578230a0SSatish Balay for (i=0; i<m; i++) { 923578230a0SSatish Balay vals[j + i*n] = *v++; 92490ace30eSBarry Smith } 92590ace30eSBarry Smith } 9266f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 927606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 92890ace30eSBarry Smith } else { 92913f74950SBarry Smith ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 930552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 931932b0c3eSLois Curfman McInnes col_lens[1] = m; 932932b0c3eSLois Curfman McInnes col_lens[2] = n; 933932b0c3eSLois Curfman McInnes col_lens[3] = nz; 934932b0c3eSLois Curfman McInnes 935932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 936932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 9376f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 938932b0c3eSLois Curfman McInnes 939932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 940932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 941932b0c3eSLois Curfman McInnes ict = 0; 942932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 943932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 944932b0c3eSLois Curfman McInnes } 9456f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 946606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 947932b0c3eSLois Curfman McInnes 948932b0c3eSLois Curfman McInnes /* store nonzero values */ 94987828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 950932b0c3eSLois Curfman McInnes ict = 0; 951932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 952932b0c3eSLois Curfman McInnes v = a->v + i; 953932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 9541b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 955932b0c3eSLois Curfman McInnes } 956932b0c3eSLois Curfman McInnes } 9576f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 958606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 95990ace30eSBarry Smith } 9603a40ed3dSBarry Smith PetscFunctionReturn(0); 961932b0c3eSLois Curfman McInnes } 962932b0c3eSLois Curfman McInnes 9634a2ae208SSatish Balay #undef __FUNCT__ 9644a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 965dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 966f1af5d2fSBarry Smith { 967f1af5d2fSBarry Smith Mat A = (Mat) Aa; 968f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9696849ba73SBarry Smith PetscErrorCode ierr; 970899cda47SBarry Smith PetscInt m = A->rmap.n,n = A->cmap.n,color,i,j; 97187828ca2SBarry Smith PetscScalar *v = a->v; 972b0a32e0cSBarry Smith PetscViewer viewer; 973b0a32e0cSBarry Smith PetscDraw popup; 974329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 975f3ef73ceSBarry Smith PetscViewerFormat format; 976f1af5d2fSBarry Smith 977f1af5d2fSBarry Smith PetscFunctionBegin; 978f1af5d2fSBarry Smith 979f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 980b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 981b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 982f1af5d2fSBarry Smith 983f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 984fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 985f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 986b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 987f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 988f1af5d2fSBarry Smith x_l = j; 989f1af5d2fSBarry Smith x_r = x_l + 1.0; 990f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 991f1af5d2fSBarry Smith y_l = m - i - 1.0; 992f1af5d2fSBarry Smith y_r = y_l + 1.0; 993f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 994329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 995b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 996329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 997b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 998f1af5d2fSBarry Smith } else { 999f1af5d2fSBarry Smith continue; 1000f1af5d2fSBarry Smith } 1001f1af5d2fSBarry Smith #else 1002f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 1003b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1004f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 1005b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1006f1af5d2fSBarry Smith } else { 1007f1af5d2fSBarry Smith continue; 1008f1af5d2fSBarry Smith } 1009f1af5d2fSBarry Smith #endif 1010b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1011f1af5d2fSBarry Smith } 1012f1af5d2fSBarry Smith } 1013f1af5d2fSBarry Smith } else { 1014f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 1015f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 1016f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 1017f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 1018f1af5d2fSBarry Smith } 1019b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 1020b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 1021b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 1022f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1023f1af5d2fSBarry Smith x_l = j; 1024f1af5d2fSBarry Smith x_r = x_l + 1.0; 1025f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1026f1af5d2fSBarry Smith y_l = m - i - 1.0; 1027f1af5d2fSBarry Smith y_r = y_l + 1.0; 1028b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 1029b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1030f1af5d2fSBarry Smith } 1031f1af5d2fSBarry Smith } 1032f1af5d2fSBarry Smith } 1033f1af5d2fSBarry Smith PetscFunctionReturn(0); 1034f1af5d2fSBarry Smith } 1035f1af5d2fSBarry Smith 10364a2ae208SSatish Balay #undef __FUNCT__ 10374a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 1038dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 1039f1af5d2fSBarry Smith { 1040b0a32e0cSBarry Smith PetscDraw draw; 1041f1af5d2fSBarry Smith PetscTruth isnull; 1042329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 1043dfbe8321SBarry Smith PetscErrorCode ierr; 1044f1af5d2fSBarry Smith 1045f1af5d2fSBarry Smith PetscFunctionBegin; 1046b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1047b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1048abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 1049f1af5d2fSBarry Smith 1050f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1051899cda47SBarry Smith xr = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0; 1052f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 1053b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1054b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 1055f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 1056f1af5d2fSBarry Smith PetscFunctionReturn(0); 1057f1af5d2fSBarry Smith } 1058f1af5d2fSBarry Smith 10594a2ae208SSatish Balay #undef __FUNCT__ 10604a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 1061dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1062932b0c3eSLois Curfman McInnes { 1063dfbe8321SBarry Smith PetscErrorCode ierr; 10646805f65bSBarry Smith PetscTruth iascii,isbinary,isdraw; 1065932b0c3eSLois Curfman McInnes 10663a40ed3dSBarry Smith PetscFunctionBegin; 106732077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1068fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 1069fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 10700f5bd95cSBarry Smith 1071c45a1595SBarry Smith if (iascii) { 1072c45a1595SBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 10730f5bd95cSBarry Smith } else if (isbinary) { 10743a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1075f1af5d2fSBarry Smith } else if (isdraw) { 1076f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 10775cd90555SBarry Smith } else { 1078958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1079932b0c3eSLois Curfman McInnes } 10803a40ed3dSBarry Smith PetscFunctionReturn(0); 1081932b0c3eSLois Curfman McInnes } 1082289bc588SBarry Smith 10834a2ae208SSatish Balay #undef __FUNCT__ 10844a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1085dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1086289bc588SBarry Smith { 1087ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1088dfbe8321SBarry Smith PetscErrorCode ierr; 108990f02eecSBarry Smith 10903a40ed3dSBarry Smith PetscFunctionBegin; 1091aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1092899cda47SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap.n,mat->cmap.n); 1093a5a9c739SBarry Smith #endif 109405b42c5fSBarry Smith ierr = PetscFree(l->pivots);CHKERRQ(ierr); 10956857c123SSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1096606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 1097dbd8c25aSHong Zhang 1098dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 1099901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 11004ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11014ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11024ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11033a40ed3dSBarry Smith PetscFunctionReturn(0); 1104289bc588SBarry Smith } 1105289bc588SBarry Smith 11064a2ae208SSatish Balay #undef __FUNCT__ 11074a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1108fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout) 1109289bc588SBarry Smith { 1110c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11116849ba73SBarry Smith PetscErrorCode ierr; 111213f74950SBarry Smith PetscInt k,j,m,n,M; 111387828ca2SBarry Smith PetscScalar *v,tmp; 111448b35521SBarry Smith 11153a40ed3dSBarry Smith PetscFunctionBegin; 1116899cda47SBarry Smith v = mat->v; m = A->rmap.n; M = mat->lda; n = A->cmap.n; 1117e9695a30SBarry Smith if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */ 1118a5ce6ee0Svictorle if (m != n) { 1119634064b4SBarry Smith SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1120d64ed03dSBarry Smith } else { 1121d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1122289bc588SBarry Smith for (k=0; k<j; k++) { 11231b807ce4Svictorle tmp = v[j + k*M]; 11241b807ce4Svictorle v[j + k*M] = v[k + j*M]; 11251b807ce4Svictorle v[k + j*M] = tmp; 1126289bc588SBarry Smith } 1127289bc588SBarry Smith } 1128d64ed03dSBarry Smith } 11293a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1130d3e5ee88SLois Curfman McInnes Mat tmat; 1131ec8511deSBarry Smith Mat_SeqDense *tmatd; 113287828ca2SBarry Smith PetscScalar *v2; 1133ea709b57SSatish Balay 1134fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 11357adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr); 1136899cda47SBarry Smith ierr = MatSetSizes(tmat,A->cmap.n,A->rmap.n,A->cmap.n,A->rmap.n);CHKERRQ(ierr); 11377adad957SLisandro Dalcin ierr = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 11385c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1139fc4dec0aSBarry Smith } else { 1140fc4dec0aSBarry Smith tmat = *matout; 1141fc4dec0aSBarry Smith } 1142ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 11430de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1144d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 11451b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1146d3e5ee88SLois Curfman McInnes } 11476d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11486d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1149d3e5ee88SLois Curfman McInnes *matout = tmat; 115048b35521SBarry Smith } 11513a40ed3dSBarry Smith PetscFunctionReturn(0); 1152289bc588SBarry Smith } 1153289bc588SBarry Smith 11544a2ae208SSatish Balay #undef __FUNCT__ 11554a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1156dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1157289bc588SBarry Smith { 1158c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1159c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 116013f74950SBarry Smith PetscInt i,j; 116187828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 11629ea5d5aeSSatish Balay 11633a40ed3dSBarry Smith PetscFunctionBegin; 1164899cda47SBarry Smith if (A1->rmap.n != A2->rmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1165899cda47SBarry Smith if (A1->cmap.n != A2->cmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1166899cda47SBarry Smith for (i=0; i<A1->rmap.n; i++) { 11671b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1168899cda47SBarry Smith for (j=0; j<A1->cmap.n; j++) { 11693a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 11701b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 11711b807ce4Svictorle } 1172289bc588SBarry Smith } 117377c4ece6SBarry Smith *flg = PETSC_TRUE; 11743a40ed3dSBarry Smith PetscFunctionReturn(0); 1175289bc588SBarry Smith } 1176289bc588SBarry Smith 11774a2ae208SSatish Balay #undef __FUNCT__ 11784a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1179dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1180289bc588SBarry Smith { 1181c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1182dfbe8321SBarry Smith PetscErrorCode ierr; 118313f74950SBarry Smith PetscInt i,n,len; 118487828ca2SBarry Smith PetscScalar *x,zero = 0.0; 118544cd7ae7SLois Curfman McInnes 11863a40ed3dSBarry Smith PetscFunctionBegin; 11872dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 11887a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 11891ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1190899cda47SBarry Smith len = PetscMin(A->rmap.n,A->cmap.n); 1191899cda47SBarry Smith if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 119244cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 11931b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1194289bc588SBarry Smith } 11951ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 11963a40ed3dSBarry Smith PetscFunctionReturn(0); 1197289bc588SBarry Smith } 1198289bc588SBarry Smith 11994a2ae208SSatish Balay #undef __FUNCT__ 12004a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1201dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1202289bc588SBarry Smith { 1203c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 120487828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1205dfbe8321SBarry Smith PetscErrorCode ierr; 1206899cda47SBarry Smith PetscInt i,j,m = A->rmap.n,n = A->cmap.n; 120755659b69SBarry Smith 12083a40ed3dSBarry Smith PetscFunctionBegin; 120928988994SBarry Smith if (ll) { 12107a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 12111ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1212899cda47SBarry Smith if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1213da3a660dSBarry Smith for (i=0; i<m; i++) { 1214da3a660dSBarry Smith x = l[i]; 1215da3a660dSBarry Smith v = mat->v + i; 1216da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1217da3a660dSBarry Smith } 12181ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1219efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1220da3a660dSBarry Smith } 122128988994SBarry Smith if (rr) { 12227a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 12231ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1224899cda47SBarry Smith if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1225da3a660dSBarry Smith for (i=0; i<n; i++) { 1226da3a660dSBarry Smith x = r[i]; 1227da3a660dSBarry Smith v = mat->v + i*m; 1228da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1229da3a660dSBarry Smith } 12301ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1231efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1232da3a660dSBarry Smith } 12333a40ed3dSBarry Smith PetscFunctionReturn(0); 1234289bc588SBarry Smith } 1235289bc588SBarry Smith 12364a2ae208SSatish Balay #undef __FUNCT__ 12374a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1238dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1239289bc588SBarry Smith { 1240c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 124187828ca2SBarry Smith PetscScalar *v = mat->v; 1242329f5518SBarry Smith PetscReal sum = 0.0; 1243899cda47SBarry Smith PetscInt lda=mat->lda,m=A->rmap.n,i,j; 1244efee365bSSatish Balay PetscErrorCode ierr; 124555659b69SBarry Smith 12463a40ed3dSBarry Smith PetscFunctionBegin; 1247289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1248a5ce6ee0Svictorle if (lda>m) { 1249899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 1250a5ce6ee0Svictorle v = mat->v+j*lda; 1251a5ce6ee0Svictorle for (i=0; i<m; i++) { 1252a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1253a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1254a5ce6ee0Svictorle #else 1255a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1256a5ce6ee0Svictorle #endif 1257a5ce6ee0Svictorle } 1258a5ce6ee0Svictorle } 1259a5ce6ee0Svictorle } else { 1260899cda47SBarry Smith for (i=0; i<A->cmap.n*A->rmap.n; i++) { 1261aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1262329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1263289bc588SBarry Smith #else 1264289bc588SBarry Smith sum += (*v)*(*v); v++; 1265289bc588SBarry Smith #endif 1266289bc588SBarry Smith } 1267a5ce6ee0Svictorle } 1268064f8208SBarry Smith *nrm = sqrt(sum); 1269899cda47SBarry Smith ierr = PetscLogFlops(2*A->cmap.n*A->rmap.n);CHKERRQ(ierr); 12703a40ed3dSBarry Smith } else if (type == NORM_1) { 1271064f8208SBarry Smith *nrm = 0.0; 1272899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 12731b807ce4Svictorle v = mat->v + j*mat->lda; 1274289bc588SBarry Smith sum = 0.0; 1275899cda47SBarry Smith for (i=0; i<A->rmap.n; i++) { 127633a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1277289bc588SBarry Smith } 1278064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1279289bc588SBarry Smith } 1280899cda47SBarry Smith ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr); 12813a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1282064f8208SBarry Smith *nrm = 0.0; 1283899cda47SBarry Smith for (j=0; j<A->rmap.n; j++) { 1284289bc588SBarry Smith v = mat->v + j; 1285289bc588SBarry Smith sum = 0.0; 1286899cda47SBarry Smith for (i=0; i<A->cmap.n; i++) { 12871b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1288289bc588SBarry Smith } 1289064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1290289bc588SBarry Smith } 1291899cda47SBarry Smith ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr); 12923a40ed3dSBarry Smith } else { 129329bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1294289bc588SBarry Smith } 12953a40ed3dSBarry Smith PetscFunctionReturn(0); 1296289bc588SBarry Smith } 1297289bc588SBarry Smith 12984a2ae208SSatish Balay #undef __FUNCT__ 12994a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 13004e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscTruth flg) 1301289bc588SBarry Smith { 1302c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 130363ba0a88SBarry Smith PetscErrorCode ierr; 130467e560aaSBarry Smith 13053a40ed3dSBarry Smith PetscFunctionBegin; 1306b5a2b587SKris Buschelman switch (op) { 1307b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 13084e0d8c25SBarry Smith aij->roworiented = flg; 1309b5a2b587SKris Buschelman break; 1310512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 1311b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 13123971808eSMatthew Knepley case MAT_NEW_NONZERO_ALLOCATION_ERR: 13134e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 1314b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1315b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 131677e54ba9SKris Buschelman case MAT_SYMMETRIC: 131777e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 13189a4540c5SBarry Smith case MAT_HERMITIAN: 13199a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1320600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 1321290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 132277e54ba9SKris Buschelman break; 1323b5a2b587SKris Buschelman default: 1324600fe468SBarry Smith SETERRQ1(PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 13253a40ed3dSBarry Smith } 13263a40ed3dSBarry Smith PetscFunctionReturn(0); 1327289bc588SBarry Smith } 1328289bc588SBarry Smith 13294a2ae208SSatish Balay #undef __FUNCT__ 13304a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1331dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 13326f0a148fSBarry Smith { 1333ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 13346849ba73SBarry Smith PetscErrorCode ierr; 1335899cda47SBarry Smith PetscInt lda=l->lda,m=A->rmap.n,j; 13363a40ed3dSBarry Smith 13373a40ed3dSBarry Smith PetscFunctionBegin; 1338a5ce6ee0Svictorle if (lda>m) { 1339899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 1340a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1341a5ce6ee0Svictorle } 1342a5ce6ee0Svictorle } else { 1343899cda47SBarry Smith ierr = PetscMemzero(l->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 1344a5ce6ee0Svictorle } 13453a40ed3dSBarry Smith PetscFunctionReturn(0); 13466f0a148fSBarry Smith } 13476f0a148fSBarry Smith 13484a2ae208SSatish Balay #undef __FUNCT__ 13494a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 1350f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag) 13516f0a148fSBarry Smith { 1352ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1353899cda47SBarry Smith PetscInt n = A->cmap.n,i,j; 135487828ca2SBarry Smith PetscScalar *slot; 135555659b69SBarry Smith 13563a40ed3dSBarry Smith PetscFunctionBegin; 13576f0a148fSBarry Smith for (i=0; i<N; i++) { 13586f0a148fSBarry Smith slot = l->v + rows[i]; 13596f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 13606f0a148fSBarry Smith } 1361f4df32b1SMatthew Knepley if (diag != 0.0) { 13626f0a148fSBarry Smith for (i=0; i<N; i++) { 13636f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1364f4df32b1SMatthew Knepley *slot = diag; 13656f0a148fSBarry Smith } 13666f0a148fSBarry Smith } 13673a40ed3dSBarry Smith PetscFunctionReturn(0); 13686f0a148fSBarry Smith } 1369557bce09SLois Curfman McInnes 13704a2ae208SSatish Balay #undef __FUNCT__ 13714a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 1372dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 137364e87e97SBarry Smith { 1374c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 13753a40ed3dSBarry Smith 13763a40ed3dSBarry Smith PetscFunctionBegin; 1377899cda47SBarry Smith if (mat->lda != A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows"); 137864e87e97SBarry Smith *array = mat->v; 13793a40ed3dSBarry Smith PetscFunctionReturn(0); 138064e87e97SBarry Smith } 13810754003eSLois Curfman McInnes 13824a2ae208SSatish Balay #undef __FUNCT__ 13834a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 1384dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1385ff14e315SSatish Balay { 13863a40ed3dSBarry Smith PetscFunctionBegin; 138709b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 13883a40ed3dSBarry Smith PetscFunctionReturn(0); 1389ff14e315SSatish Balay } 13900754003eSLois Curfman McInnes 13914a2ae208SSatish Balay #undef __FUNCT__ 13924a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 139313f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 13940754003eSLois Curfman McInnes { 1395c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 13966849ba73SBarry Smith PetscErrorCode ierr; 139721a2c019SBarry Smith PetscInt i,j,*irow,*icol,nrows,ncols; 139887828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 13990754003eSLois Curfman McInnes Mat newmat; 14000754003eSLois Curfman McInnes 14013a40ed3dSBarry Smith PetscFunctionBegin; 140278b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 140378b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1404e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1405e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 14060754003eSLois Curfman McInnes 1407182d2002SSatish Balay /* Check submatrixcall */ 1408182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 140913f74950SBarry Smith PetscInt n_cols,n_rows; 1410182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 141121a2c019SBarry Smith if (n_rows != nrows || n_cols != ncols) { 141221a2c019SBarry Smith /* resize the result result matrix to match number of requested rows/columns */ 141321a2c019SBarry Smith ierr = MatSetSizes(*B,nrows,nrows,nrows,nrows);CHKERRQ(ierr); 141421a2c019SBarry Smith } 1415182d2002SSatish Balay newmat = *B; 1416182d2002SSatish Balay } else { 14170754003eSLois Curfman McInnes /* Create and fill new matrix */ 14187adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr); 1419f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 14207adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 14215c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1422182d2002SSatish Balay } 1423182d2002SSatish Balay 1424182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1425182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1426182d2002SSatish Balay 1427182d2002SSatish Balay for (i=0; i<ncols; i++) { 14286de62eeeSBarry Smith av = v + mat->lda*icol[i]; 1429182d2002SSatish Balay for (j=0; j<nrows; j++) { 1430182d2002SSatish Balay *bv++ = av[irow[j]]; 14310754003eSLois Curfman McInnes } 14320754003eSLois Curfman McInnes } 1433182d2002SSatish Balay 1434182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 14356d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14366d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14370754003eSLois Curfman McInnes 14380754003eSLois Curfman McInnes /* Free work space */ 143978b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 144078b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1441182d2002SSatish Balay *B = newmat; 14423a40ed3dSBarry Smith PetscFunctionReturn(0); 14430754003eSLois Curfman McInnes } 14440754003eSLois Curfman McInnes 14454a2ae208SSatish Balay #undef __FUNCT__ 14464a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 144713f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1448905e6a2fSBarry Smith { 14496849ba73SBarry Smith PetscErrorCode ierr; 145013f74950SBarry Smith PetscInt i; 1451905e6a2fSBarry Smith 14523a40ed3dSBarry Smith PetscFunctionBegin; 1453905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1454b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1455905e6a2fSBarry Smith } 1456905e6a2fSBarry Smith 1457905e6a2fSBarry Smith for (i=0; i<n; i++) { 14586a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1459905e6a2fSBarry Smith } 14603a40ed3dSBarry Smith PetscFunctionReturn(0); 1461905e6a2fSBarry Smith } 1462905e6a2fSBarry Smith 14634a2ae208SSatish Balay #undef __FUNCT__ 1464c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1465c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1466c0aa2d19SHong Zhang { 1467c0aa2d19SHong Zhang PetscFunctionBegin; 1468c0aa2d19SHong Zhang PetscFunctionReturn(0); 1469c0aa2d19SHong Zhang } 1470c0aa2d19SHong Zhang 1471c0aa2d19SHong Zhang #undef __FUNCT__ 1472c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1473c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1474c0aa2d19SHong Zhang { 1475c0aa2d19SHong Zhang PetscFunctionBegin; 1476c0aa2d19SHong Zhang PetscFunctionReturn(0); 1477c0aa2d19SHong Zhang } 1478c0aa2d19SHong Zhang 1479c0aa2d19SHong Zhang #undef __FUNCT__ 14804a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1481dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 14824b0e389bSBarry Smith { 14834b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 14846849ba73SBarry Smith PetscErrorCode ierr; 1485899cda47SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap.n,n=A->cmap.n, j; 14863a40ed3dSBarry Smith 14873a40ed3dSBarry Smith PetscFunctionBegin; 148833f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 148933f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1490cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 14913a40ed3dSBarry Smith PetscFunctionReturn(0); 14923a40ed3dSBarry Smith } 1493899cda47SBarry Smith if (m != B->rmap.n || n != B->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1494a5ce6ee0Svictorle if (lda1>m || lda2>m) { 14950dbb7854Svictorle for (j=0; j<n; j++) { 1496a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1497a5ce6ee0Svictorle } 1498a5ce6ee0Svictorle } else { 1499899cda47SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 1500a5ce6ee0Svictorle } 1501273d9f13SBarry Smith PetscFunctionReturn(0); 1502273d9f13SBarry Smith } 1503273d9f13SBarry Smith 15044a2ae208SSatish Balay #undef __FUNCT__ 15054a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1506dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A) 1507273d9f13SBarry Smith { 1508dfbe8321SBarry Smith PetscErrorCode ierr; 1509273d9f13SBarry Smith 1510273d9f13SBarry Smith PetscFunctionBegin; 1511273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 15123a40ed3dSBarry Smith PetscFunctionReturn(0); 15134b0e389bSBarry Smith } 15144b0e389bSBarry Smith 1515284134d9SBarry Smith #undef __FUNCT__ 1516284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense" 1517284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N) 1518284134d9SBarry Smith { 1519284134d9SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 152021a2c019SBarry Smith PetscErrorCode ierr; 1521284134d9SBarry Smith PetscFunctionBegin; 152221a2c019SBarry Smith /* this will not be called before lda, Mmax, and Nmax have been set */ 1523284134d9SBarry Smith m = PetscMax(m,M); 1524284134d9SBarry Smith n = PetscMax(n,N); 152521a2c019SBarry Smith if (m > a->Mmax) SETERRQ2(PETSC_ERR_SUP,"Cannot yet resize number rows of dense matrix larger then its initial size %d, requested %d",a->lda,(int)m); 1526284134d9SBarry Smith if (n > a->Nmax) SETERRQ2(PETSC_ERR_SUP,"Cannot yet resize number columns of dense matrix larger then its initial size %d, requested %d",a->Nmax,(int)n); 1527899cda47SBarry Smith A->rmap.n = A->rmap.n = m; 1528899cda47SBarry Smith A->cmap.n = A->cmap.N = n; 152921a2c019SBarry Smith if (a->changelda) a->lda = m; 153021a2c019SBarry Smith ierr = PetscMemzero(a->v,m*n*sizeof(PetscScalar));CHKERRQ(ierr); 1531284134d9SBarry Smith PetscFunctionReturn(0); 1532284134d9SBarry Smith } 1533170fe5c8SBarry Smith 1534284134d9SBarry Smith 1535a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/ 1536a9fe9ddaSSatish Balay #undef __FUNCT__ 1537a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1538a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1539a9fe9ddaSSatish Balay { 1540a9fe9ddaSSatish Balay PetscErrorCode ierr; 1541a9fe9ddaSSatish Balay 1542a9fe9ddaSSatish Balay PetscFunctionBegin; 1543a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1544a9fe9ddaSSatish Balay ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1545a9fe9ddaSSatish Balay } 1546a9fe9ddaSSatish Balay ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1547a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1548a9fe9ddaSSatish Balay } 1549a9fe9ddaSSatish Balay 1550a9fe9ddaSSatish Balay #undef __FUNCT__ 1551a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1552a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1553a9fe9ddaSSatish Balay { 1554ee16a9a1SHong Zhang PetscErrorCode ierr; 1555899cda47SBarry Smith PetscInt m=A->rmap.n,n=B->cmap.n; 1556ee16a9a1SHong Zhang Mat Cmat; 1557a9fe9ddaSSatish Balay 1558ee16a9a1SHong Zhang PetscFunctionBegin; 1559899cda47SBarry Smith if (A->cmap.n != B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->cmap.n %d != B->rmap.n %d\n",A->cmap.n,B->rmap.n); 1560ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1561ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1562ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1563ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1564ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1565ee16a9a1SHong Zhang *C = Cmat; 1566ee16a9a1SHong Zhang PetscFunctionReturn(0); 1567ee16a9a1SHong Zhang } 1568a9fe9ddaSSatish Balay 156998a3b096SSatish Balay #undef __FUNCT__ 1570a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1571a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1572a9fe9ddaSSatish Balay { 1573a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1574a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1575a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 15760805154bSBarry Smith PetscBLASInt m,n,k; 1577a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1578a9fe9ddaSSatish Balay 1579a9fe9ddaSSatish Balay PetscFunctionBegin; 15800805154bSBarry Smith m = PetscBLASIntCast(A->rmap.n); 15810805154bSBarry Smith n = PetscBLASIntCast(B->cmap.n); 15820805154bSBarry Smith k = PetscBLASIntCast(A->cmap.n); 1583a9fe9ddaSSatish Balay BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1584a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1585a9fe9ddaSSatish Balay } 1586a9fe9ddaSSatish Balay 1587a9fe9ddaSSatish Balay #undef __FUNCT__ 1588a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense" 1589a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1590a9fe9ddaSSatish Balay { 1591a9fe9ddaSSatish Balay PetscErrorCode ierr; 1592a9fe9ddaSSatish Balay 1593a9fe9ddaSSatish Balay PetscFunctionBegin; 1594a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1595a9fe9ddaSSatish Balay ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1596a9fe9ddaSSatish Balay } 1597a9fe9ddaSSatish Balay ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1598a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1599a9fe9ddaSSatish Balay } 1600a9fe9ddaSSatish Balay 1601a9fe9ddaSSatish Balay #undef __FUNCT__ 1602a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense" 1603a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1604a9fe9ddaSSatish Balay { 1605ee16a9a1SHong Zhang PetscErrorCode ierr; 1606899cda47SBarry Smith PetscInt m=A->cmap.n,n=B->cmap.n; 1607ee16a9a1SHong Zhang Mat Cmat; 1608a9fe9ddaSSatish Balay 1609ee16a9a1SHong Zhang PetscFunctionBegin; 1610899cda47SBarry Smith if (A->rmap.n != B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->rmap.n %d != B->rmap.n %d\n",A->rmap.n,B->rmap.n); 1611ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1612ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1613ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1614ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1615ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1616ee16a9a1SHong Zhang *C = Cmat; 1617ee16a9a1SHong Zhang PetscFunctionReturn(0); 1618ee16a9a1SHong Zhang } 1619a9fe9ddaSSatish Balay 1620a9fe9ddaSSatish Balay #undef __FUNCT__ 1621a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense" 1622a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1623a9fe9ddaSSatish Balay { 1624a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1625a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1626a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 16270805154bSBarry Smith PetscBLASInt m,n,k; 1628a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1629a9fe9ddaSSatish Balay 1630a9fe9ddaSSatish Balay PetscFunctionBegin; 16312fbe02b9SBarry Smith m = PetscBLASIntCast(A->cmap.n); 16320805154bSBarry Smith n = PetscBLASIntCast(B->cmap.n); 16332fbe02b9SBarry Smith k = PetscBLASIntCast(A->rmap.n); 16342fbe02b9SBarry Smith /* 16352fbe02b9SBarry Smith Note the m and n arguments below are the number rows and columns of A', not A! 16362fbe02b9SBarry Smith */ 1637a9fe9ddaSSatish Balay BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1638a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1639a9fe9ddaSSatish Balay } 1640985db425SBarry Smith 1641985db425SBarry Smith #undef __FUNCT__ 1642985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense" 1643985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 1644985db425SBarry Smith { 1645985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1646985db425SBarry Smith PetscErrorCode ierr; 1647985db425SBarry Smith PetscInt i,j,m = A->rmap.n,n = A->cmap.n,p; 1648985db425SBarry Smith PetscScalar *x; 1649985db425SBarry Smith MatScalar *aa = a->v; 1650985db425SBarry Smith 1651985db425SBarry Smith PetscFunctionBegin; 1652985db425SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1653985db425SBarry Smith 1654985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1655985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1656985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1657985db425SBarry Smith if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1658985db425SBarry Smith for (i=0; i<m; i++) { 1659985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1660985db425SBarry Smith for (j=1; j<n; j++){ 1661985db425SBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1662985db425SBarry Smith } 1663985db425SBarry Smith } 1664985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1665985db425SBarry Smith PetscFunctionReturn(0); 1666985db425SBarry Smith } 1667985db425SBarry Smith 1668985db425SBarry Smith #undef __FUNCT__ 1669985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense" 1670985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 1671985db425SBarry Smith { 1672985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1673985db425SBarry Smith PetscErrorCode ierr; 1674985db425SBarry Smith PetscInt i,j,m = A->rmap.n,n = A->cmap.n,p; 1675985db425SBarry Smith PetscScalar *x; 1676985db425SBarry Smith PetscReal atmp; 1677985db425SBarry Smith MatScalar *aa = a->v; 1678985db425SBarry Smith 1679985db425SBarry Smith PetscFunctionBegin; 1680985db425SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1681985db425SBarry Smith 1682985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1683985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1684985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1685985db425SBarry Smith if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1686985db425SBarry Smith for (i=0; i<m; i++) { 16879189402eSHong Zhang x[i] = PetscAbsScalar(aa[i]); 1688985db425SBarry Smith for (j=1; j<n; j++){ 1689985db425SBarry Smith atmp = PetscAbsScalar(aa[i+m*j]); 1690985db425SBarry Smith if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 1691985db425SBarry Smith } 1692985db425SBarry Smith } 1693985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1694985db425SBarry Smith PetscFunctionReturn(0); 1695985db425SBarry Smith } 1696985db425SBarry Smith 1697985db425SBarry Smith #undef __FUNCT__ 1698985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense" 1699985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 1700985db425SBarry Smith { 1701985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1702985db425SBarry Smith PetscErrorCode ierr; 1703985db425SBarry Smith PetscInt i,j,m = A->rmap.n,n = A->cmap.n,p; 1704985db425SBarry Smith PetscScalar *x; 1705985db425SBarry Smith MatScalar *aa = a->v; 1706985db425SBarry Smith 1707985db425SBarry Smith PetscFunctionBegin; 1708985db425SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1709985db425SBarry Smith 1710985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1711985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1712985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1713985db425SBarry Smith if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1714985db425SBarry Smith for (i=0; i<m; i++) { 1715985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1716985db425SBarry Smith for (j=1; j<n; j++){ 1717985db425SBarry Smith if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1718985db425SBarry Smith } 1719985db425SBarry Smith } 1720985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1721985db425SBarry Smith PetscFunctionReturn(0); 1722985db425SBarry Smith } 1723985db425SBarry Smith 17248d0534beSBarry Smith #undef __FUNCT__ 17258d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense" 17268d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 17278d0534beSBarry Smith { 17288d0534beSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 17298d0534beSBarry Smith PetscErrorCode ierr; 17308d0534beSBarry Smith PetscScalar *x; 17318d0534beSBarry Smith 17328d0534beSBarry Smith PetscFunctionBegin; 17338d0534beSBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 17348d0534beSBarry Smith 17358d0534beSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 17368d0534beSBarry Smith ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 17378d0534beSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 17388d0534beSBarry Smith PetscFunctionReturn(0); 17398d0534beSBarry Smith } 17408d0534beSBarry Smith 1741289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1742a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1743905e6a2fSBarry Smith MatGetRow_SeqDense, 1744905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1745905e6a2fSBarry Smith MatMult_SeqDense, 174697304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 17477c922b88SBarry Smith MatMultTranspose_SeqDense, 17487c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1749905e6a2fSBarry Smith MatSolve_SeqDense, 1750905e6a2fSBarry Smith MatSolveAdd_SeqDense, 17517c922b88SBarry Smith MatSolveTranspose_SeqDense, 175297304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense, 1753905e6a2fSBarry Smith MatLUFactor_SeqDense, 1754905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1755ec8511deSBarry Smith MatRelax_SeqDense, 1756ec8511deSBarry Smith MatTranspose_SeqDense, 175797304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1758905e6a2fSBarry Smith MatEqual_SeqDense, 1759905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1760905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1761905e6a2fSBarry Smith MatNorm_SeqDense, 1762c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense, 1763c0aa2d19SHong Zhang MatAssemblyEnd_SeqDense, 1764905e6a2fSBarry Smith 0, 1765905e6a2fSBarry Smith MatSetOption_SeqDense, 1766905e6a2fSBarry Smith MatZeroEntries_SeqDense, 176797304618SKris Buschelman /*25*/ MatZeroRows_SeqDense, 1768905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1769905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1770905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1771905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 177297304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense, 1773273d9f13SBarry Smith 0, 1774905e6a2fSBarry Smith 0, 1775905e6a2fSBarry Smith MatGetArray_SeqDense, 1776905e6a2fSBarry Smith MatRestoreArray_SeqDense, 177797304618SKris Buschelman /*35*/ MatDuplicate_SeqDense, 1778a5ae1ecdSBarry Smith 0, 1779a5ae1ecdSBarry Smith 0, 1780a5ae1ecdSBarry Smith 0, 1781a5ae1ecdSBarry Smith 0, 178297304618SKris Buschelman /*40*/ MatAXPY_SeqDense, 1783a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1784a5ae1ecdSBarry Smith 0, 17854b0e389bSBarry Smith MatGetValues_SeqDense, 1786a5ae1ecdSBarry Smith MatCopy_SeqDense, 1787985db425SBarry Smith /*45*/ MatGetRowMax_SeqDense, 1788a5ae1ecdSBarry Smith MatScale_SeqDense, 1789a5ae1ecdSBarry Smith 0, 1790a5ae1ecdSBarry Smith 0, 1791a5ae1ecdSBarry Smith 0, 1792521d7252SBarry Smith /*50*/ 0, 1793a5ae1ecdSBarry Smith 0, 1794a5ae1ecdSBarry Smith 0, 1795a5ae1ecdSBarry Smith 0, 1796a5ae1ecdSBarry Smith 0, 179797304618SKris Buschelman /*55*/ 0, 1798a5ae1ecdSBarry Smith 0, 1799a5ae1ecdSBarry Smith 0, 1800a5ae1ecdSBarry Smith 0, 1801a5ae1ecdSBarry Smith 0, 180297304618SKris Buschelman /*60*/ 0, 1803e03a110bSBarry Smith MatDestroy_SeqDense, 1804e03a110bSBarry Smith MatView_SeqDense, 1805357abbc8SBarry Smith 0, 180697304618SKris Buschelman 0, 180797304618SKris Buschelman /*65*/ 0, 180897304618SKris Buschelman 0, 180997304618SKris Buschelman 0, 181097304618SKris Buschelman 0, 181197304618SKris Buschelman 0, 1812985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqDense, 181397304618SKris Buschelman 0, 181497304618SKris Buschelman 0, 181597304618SKris Buschelman 0, 181697304618SKris Buschelman 0, 181797304618SKris Buschelman /*75*/ 0, 181897304618SKris Buschelman 0, 181997304618SKris Buschelman 0, 182097304618SKris Buschelman 0, 182197304618SKris Buschelman 0, 182297304618SKris Buschelman /*80*/ 0, 182397304618SKris Buschelman 0, 182497304618SKris Buschelman 0, 182597304618SKris Buschelman 0, 1826865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense, 1827865e5f61SKris Buschelman 0, 18281cbb95d3SBarry Smith MatIsHermitian_SeqDense, 1829865e5f61SKris Buschelman 0, 1830865e5f61SKris Buschelman 0, 1831865e5f61SKris Buschelman 0, 1832a9fe9ddaSSatish Balay /*90*/ MatMatMult_SeqDense_SeqDense, 1833a9fe9ddaSSatish Balay MatMatMultSymbolic_SeqDense_SeqDense, 1834a9fe9ddaSSatish Balay MatMatMultNumeric_SeqDense_SeqDense, 1835865e5f61SKris Buschelman 0, 1836865e5f61SKris Buschelman 0, 1837865e5f61SKris Buschelman /*95*/ 0, 1838a9fe9ddaSSatish Balay MatMatMultTranspose_SeqDense_SeqDense, 1839a9fe9ddaSSatish Balay MatMatMultTransposeSymbolic_SeqDense_SeqDense, 1840a9fe9ddaSSatish Balay MatMatMultTransposeNumeric_SeqDense_SeqDense, 1841284134d9SBarry Smith 0, 1842284134d9SBarry Smith /*100*/0, 1843284134d9SBarry Smith 0, 1844284134d9SBarry Smith 0, 1845284134d9SBarry Smith 0, 1846985db425SBarry Smith MatSetSizes_SeqDense, 1847985db425SBarry Smith 0, 1848985db425SBarry Smith 0, 1849985db425SBarry Smith 0, 1850985db425SBarry Smith 0, 1851985db425SBarry Smith 0, 1852985db425SBarry Smith /*110*/0, 1853985db425SBarry Smith 0, 18548d0534beSBarry Smith MatGetRowMin_SeqDense, 18558d0534beSBarry Smith MatGetColumnVector_SeqDense 1856985db425SBarry Smith }; 185790ace30eSBarry Smith 18584a2ae208SSatish Balay #undef __FUNCT__ 18594a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 18604b828684SBarry Smith /*@C 1861fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1862d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1863d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1864289bc588SBarry Smith 1865db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1866db81eaa0SLois Curfman McInnes 186720563c6bSBarry Smith Input Parameters: 1868db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 18690c775827SLois Curfman McInnes . m - number of rows 187018f449edSLois Curfman McInnes . n - number of columns 1871db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1872dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 187320563c6bSBarry Smith 187420563c6bSBarry Smith Output Parameter: 187544cd7ae7SLois Curfman McInnes . A - the matrix 187620563c6bSBarry Smith 1877b259b22eSLois Curfman McInnes Notes: 187818f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 187918f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1880b4fd4287SBarry Smith set data=PETSC_NULL. 188118f449edSLois Curfman McInnes 1882027ccd11SLois Curfman McInnes Level: intermediate 1883027ccd11SLois Curfman McInnes 1884dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1885d65003e9SLois Curfman McInnes 1886db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 188720563c6bSBarry Smith @*/ 1888be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 1889289bc588SBarry Smith { 1890dfbe8321SBarry Smith PetscErrorCode ierr; 18913b2fbd54SBarry Smith 18923a40ed3dSBarry Smith PetscFunctionBegin; 1893f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1894f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1895273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1896273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1897273d9f13SBarry Smith PetscFunctionReturn(0); 1898273d9f13SBarry Smith } 1899273d9f13SBarry Smith 19004a2ae208SSatish Balay #undef __FUNCT__ 1901afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation" 1902273d9f13SBarry Smith /*@C 1903273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1904273d9f13SBarry Smith 1905273d9f13SBarry Smith Collective on MPI_Comm 1906273d9f13SBarry Smith 1907273d9f13SBarry Smith Input Parameters: 1908273d9f13SBarry Smith + A - the matrix 1909273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1910273d9f13SBarry Smith 1911273d9f13SBarry Smith Notes: 1912273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1913273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1914284134d9SBarry Smith need not call this routine. 1915273d9f13SBarry Smith 1916273d9f13SBarry Smith Level: intermediate 1917273d9f13SBarry Smith 1918273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1919273d9f13SBarry Smith 1920273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1921273d9f13SBarry Smith @*/ 1922be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 1923273d9f13SBarry Smith { 1924dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar[]); 1925a23d5eceSKris Buschelman 1926a23d5eceSKris Buschelman PetscFunctionBegin; 1927a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1928a23d5eceSKris Buschelman if (f) { 1929a23d5eceSKris Buschelman ierr = (*f)(B,data);CHKERRQ(ierr); 1930a23d5eceSKris Buschelman } 1931a23d5eceSKris Buschelman PetscFunctionReturn(0); 1932a23d5eceSKris Buschelman } 1933a23d5eceSKris Buschelman 1934a23d5eceSKris Buschelman EXTERN_C_BEGIN 1935a23d5eceSKris Buschelman #undef __FUNCT__ 1936afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense" 1937be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 1938a23d5eceSKris Buschelman { 1939273d9f13SBarry Smith Mat_SeqDense *b; 1940dfbe8321SBarry Smith PetscErrorCode ierr; 1941273d9f13SBarry Smith 1942273d9f13SBarry Smith PetscFunctionBegin; 1943273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1944273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 19455afd5e0cSBarry Smith if (b->lda <= 0) b->lda = B->rmap.n; 19469e8f95c4SLisandro Dalcin if (!data) { /* petsc-allocated storage */ 19479e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 19485afd5e0cSBarry Smith ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 1949284134d9SBarry Smith ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 1950284134d9SBarry Smith ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 19519e8f95c4SLisandro Dalcin b->user_alloc = PETSC_FALSE; 1952273d9f13SBarry Smith } else { /* user-allocated storage */ 19539e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 1954273d9f13SBarry Smith b->v = data; 1955273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 1956273d9f13SBarry Smith } 1957273d9f13SBarry Smith PetscFunctionReturn(0); 1958273d9f13SBarry Smith } 1959a23d5eceSKris Buschelman EXTERN_C_END 1960273d9f13SBarry Smith 19611b807ce4Svictorle #undef __FUNCT__ 19621b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 19631b807ce4Svictorle /*@C 19641b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 19651b807ce4Svictorle 19661b807ce4Svictorle Input parameter: 19671b807ce4Svictorle + A - the matrix 19681b807ce4Svictorle - lda - the leading dimension 19691b807ce4Svictorle 19701b807ce4Svictorle Notes: 19711b807ce4Svictorle This routine is to be used in conjunction with MatSeqDenseSetPreallocation; 19721b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 19731b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 19741b807ce4Svictorle 19751b807ce4Svictorle Level: intermediate 19761b807ce4Svictorle 19771b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 19781b807ce4Svictorle 1979284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 19801b807ce4Svictorle @*/ 1981be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda) 19821b807ce4Svictorle { 19831b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 198421a2c019SBarry Smith 19851b807ce4Svictorle PetscFunctionBegin; 1986899cda47SBarry Smith if (lda < B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap.n); 19871b807ce4Svictorle b->lda = lda; 198821a2c019SBarry Smith b->changelda = PETSC_FALSE; 198921a2c019SBarry Smith b->Mmax = PetscMax(b->Mmax,lda); 19901b807ce4Svictorle PetscFunctionReturn(0); 19911b807ce4Svictorle } 19921b807ce4Svictorle 19930bad9183SKris Buschelman /*MC 1994fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 19950bad9183SKris Buschelman 19960bad9183SKris Buschelman Options Database Keys: 19970bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 19980bad9183SKris Buschelman 19990bad9183SKris Buschelman Level: beginner 20000bad9183SKris Buschelman 200189665df3SBarry Smith .seealso: MatCreateSeqDense() 200289665df3SBarry Smith 20030bad9183SKris Buschelman M*/ 20040bad9183SKris Buschelman 2005273d9f13SBarry Smith EXTERN_C_BEGIN 20064a2ae208SSatish Balay #undef __FUNCT__ 20074a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 2008be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B) 2009273d9f13SBarry Smith { 2010273d9f13SBarry Smith Mat_SeqDense *b; 2011dfbe8321SBarry Smith PetscErrorCode ierr; 20127c334f02SBarry Smith PetscMPIInt size; 2013273d9f13SBarry Smith 2014273d9f13SBarry Smith PetscFunctionBegin; 20157adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr); 201629bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 201755659b69SBarry Smith 2018899cda47SBarry Smith B->rmap.bs = B->cmap.bs = 1; 20196148ca0dSBarry Smith ierr = PetscMapSetUp(&B->rmap);CHKERRQ(ierr); 20206148ca0dSBarry Smith ierr = PetscMapSetUp(&B->cmap);CHKERRQ(ierr); 2021273d9f13SBarry Smith 202238f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr); 2023549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 202444cd7ae7SLois Curfman McInnes B->factor = 0; 202590f02eecSBarry Smith B->mapping = 0; 202644cd7ae7SLois Curfman McInnes B->data = (void*)b; 202718f449edSLois Curfman McInnes 2028a5ae1ecdSBarry Smith 202944cd7ae7SLois Curfman McInnes b->pivots = 0; 2030273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 2031273d9f13SBarry Smith b->v = 0; 2032899cda47SBarry Smith b->lda = B->rmap.n; 203321a2c019SBarry Smith b->changelda = PETSC_FALSE; 2034899cda47SBarry Smith b->Mmax = B->rmap.n; 2035899cda47SBarry Smith b->Nmax = B->cmap.n; 20364e220ebcSLois Curfman McInnes 2037*b24902e0SBarry Smith 2038*b24902e0SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_seqdense_petsc_C", 2039*b24902e0SBarry Smith "MatGetFactor_seqdense_petsc", 2040*b24902e0SBarry Smith MatGetFactor_seqdense_petsc);CHKERRQ(ierr); 2041a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 2042a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 2043a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 20444ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C", 20454ae313f4SHong Zhang "MatMatMult_SeqAIJ_SeqDense", 20464ae313f4SHong Zhang MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 20474ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C", 20484ae313f4SHong Zhang "MatMatMultSymbolic_SeqAIJ_SeqDense", 20494ae313f4SHong Zhang MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 20504ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C", 20514ae313f4SHong Zhang "MatMatMultNumeric_SeqAIJ_SeqDense", 20524ae313f4SHong Zhang MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 205317667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 20543a40ed3dSBarry Smith PetscFunctionReturn(0); 2055289bc588SBarry Smith } 2056c0aa2d19SHong Zhang 2057c0aa2d19SHong Zhang 2058273d9f13SBarry Smith EXTERN_C_END 2059