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; 17899cda47SBarry Smith PetscBLASInt N = (PetscBLASInt)X->rmap.n*X->cmap.n,m=(PetscBLASInt)X->rmap.n,ldax = x->lda,lday=y->lda,one = 1; 18efee365bSSatish Balay PetscErrorCode ierr; 193a40ed3dSBarry Smith 203a40ed3dSBarry Smith PetscFunctionBegin; 21a5ce6ee0Svictorle if (ldax>m || lday>m) { 22899cda47SBarry Smith for (j=0; j<X->cmap.n; j++) { 23f4df32b1SMatthew Knepley BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one); 24a5ce6ee0Svictorle } 25a5ce6ee0Svictorle } else { 26f4df32b1SMatthew Knepley BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one); 27a5ce6ee0Svictorle } 28efee365bSSatish Balay ierr = PetscLogFlops(2*N-1);CHKERRQ(ierr); 293a40ed3dSBarry Smith PetscFunctionReturn(0); 301987afe7SBarry Smith } 311987afe7SBarry Smith 324a2ae208SSatish Balay #undef __FUNCT__ 334a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 34dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 35289bc588SBarry Smith { 36899cda47SBarry Smith PetscInt N = A->rmap.n*A->cmap.n; 373a40ed3dSBarry Smith 383a40ed3dSBarry Smith PetscFunctionBegin; 39899cda47SBarry Smith info->rows_global = (double)A->rmap.n; 40899cda47SBarry Smith info->columns_global = (double)A->cmap.n; 41899cda47SBarry Smith info->rows_local = (double)A->rmap.n; 42899cda47SBarry Smith info->columns_local = (double)A->cmap.n; 434e220ebcSLois Curfman McInnes info->block_size = 1.0; 444e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 456de62eeeSBarry Smith info->nz_used = (double)N; 466de62eeeSBarry Smith info->nz_unneeded = (double)0; 474e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 484e220ebcSLois Curfman McInnes info->mallocs = 0; 494e220ebcSLois Curfman McInnes info->memory = A->mem; 504e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 514e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 524e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 533a40ed3dSBarry Smith PetscFunctionReturn(0); 54289bc588SBarry Smith } 55289bc588SBarry Smith 564a2ae208SSatish Balay #undef __FUNCT__ 574a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense" 58f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha) 5980cd9d93SLois Curfman McInnes { 60273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 61f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 624ce68768SBarry Smith PetscBLASInt one = 1,lda = a->lda,j,nz; 63efee365bSSatish Balay PetscErrorCode ierr; 6480cd9d93SLois Curfman McInnes 653a40ed3dSBarry Smith PetscFunctionBegin; 66899cda47SBarry Smith if (lda>A->rmap.n) { 67899cda47SBarry Smith nz = (PetscBLASInt)A->rmap.n; 68899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 69f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v+j*lda,&one); 70a5ce6ee0Svictorle } 71a5ce6ee0Svictorle } else { 72899cda47SBarry Smith nz = (PetscBLASInt)A->rmap.n*A->cmap.n; 73f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v,&one); 74a5ce6ee0Svictorle } 75efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 763a40ed3dSBarry Smith PetscFunctionReturn(0); 7780cd9d93SLois Curfman McInnes } 7880cd9d93SLois Curfman McInnes 791cbb95d3SBarry Smith #undef __FUNCT__ 801cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense" 811cbb95d3SBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscTruth *fl) 821cbb95d3SBarry Smith { 831cbb95d3SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 841cbb95d3SBarry Smith PetscInt i,j,m = A->rmap.n,N; 851cbb95d3SBarry Smith PetscScalar *v = a->v; 861cbb95d3SBarry Smith 871cbb95d3SBarry Smith PetscFunctionBegin; 881cbb95d3SBarry Smith *fl = PETSC_FALSE; 891cbb95d3SBarry Smith if (A->rmap.n != A->cmap.n) PetscFunctionReturn(0); 901cbb95d3SBarry Smith N = a->lda; 911cbb95d3SBarry Smith 921cbb95d3SBarry Smith for (i=0; i<m; i++) { 931cbb95d3SBarry Smith for (j=i+1; j<m; j++) { 941cbb95d3SBarry Smith if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0); 951cbb95d3SBarry Smith } 961cbb95d3SBarry Smith } 971cbb95d3SBarry Smith *fl = PETSC_TRUE; 981cbb95d3SBarry Smith PetscFunctionReturn(0); 991cbb95d3SBarry Smith } 1001cbb95d3SBarry Smith 101289bc588SBarry Smith /* ---------------------------------------------------------------*/ 1026831982aSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 1036831982aSBarry Smith rather than put it in the Mat->row slot.*/ 1044a2ae208SSatish Balay #undef __FUNCT__ 1054a2ae208SSatish Balay #define __FUNCT__ "MatLUFactor_SeqDense" 106dfbe8321SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo) 107289bc588SBarry Smith { 108a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 10981f479a6Svictorle PetscFunctionBegin; 110a07ab388SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 111a07ab388SBarry Smith #else 112c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1136849ba73SBarry Smith PetscErrorCode ierr; 114899cda47SBarry Smith PetscBLASInt n = (PetscBLASInt)A->cmap.n,m = (PetscBLASInt)A->rmap.n,info; 11555659b69SBarry Smith 1163a40ed3dSBarry Smith PetscFunctionBegin; 117289bc588SBarry Smith if (!mat->pivots) { 118899cda47SBarry Smith ierr = PetscMalloc((A->rmap.n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr); 119899cda47SBarry Smith ierr = PetscLogObjectMemory(A,A->rmap.n*sizeof(PetscBLASInt));CHKERRQ(ierr); 120289bc588SBarry Smith } 121f1af5d2fSBarry Smith A->factor = FACTOR_LU; 122899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 12371044d3cSBarry Smith LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info); 12429bbc08cSBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 12529bbc08cSBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 126899cda47SBarry Smith ierr = PetscLogFlops((2*A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr); 127a07ab388SBarry Smith #endif 1283a40ed3dSBarry Smith PetscFunctionReturn(0); 129289bc588SBarry Smith } 1306ee01492SSatish Balay 1314a2ae208SSatish Balay #undef __FUNCT__ 1324a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 133dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 13402cad45dSBarry Smith { 13502cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 1366849ba73SBarry Smith PetscErrorCode ierr; 13713f74950SBarry Smith PetscInt lda = (PetscInt)mat->lda,j,m; 13802cad45dSBarry Smith Mat newi; 13902cad45dSBarry Smith 1403a40ed3dSBarry Smith PetscFunctionBegin; 141f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&newi);CHKERRQ(ierr); 142899cda47SBarry Smith ierr = MatSetSizes(newi,A->rmap.n,A->cmap.n,A->rmap.n,A->cmap.n);CHKERRQ(ierr); 1435c5985e7SKris Buschelman ierr = MatSetType(newi,A->type_name);CHKERRQ(ierr); 1445c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr); 1455609ef8eSBarry Smith if (cpvalues == MAT_COPY_VALUES) { 146a5ce6ee0Svictorle l = (Mat_SeqDense*)newi->data; 147899cda47SBarry Smith if (lda>A->rmap.n) { 148899cda47SBarry Smith m = A->rmap.n; 149899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 150a5ce6ee0Svictorle ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 151a5ce6ee0Svictorle } 152a5ce6ee0Svictorle } else { 153899cda47SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 15402cad45dSBarry Smith } 155a5ce6ee0Svictorle } 15602cad45dSBarry Smith newi->assembled = PETSC_TRUE; 15702cad45dSBarry Smith *newmat = newi; 1583a40ed3dSBarry Smith PetscFunctionReturn(0); 15902cad45dSBarry Smith } 16002cad45dSBarry Smith 1614a2ae208SSatish Balay #undef __FUNCT__ 1624a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 163dfbe8321SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat A,IS row,IS col,MatFactorInfo *info,Mat *fact) 164289bc588SBarry Smith { 165dfbe8321SBarry Smith PetscErrorCode ierr; 1663a40ed3dSBarry Smith 1673a40ed3dSBarry Smith PetscFunctionBegin; 1685609ef8eSBarry Smith ierr = MatDuplicate_SeqDense(A,MAT_DO_NOT_COPY_VALUES,fact);CHKERRQ(ierr); 1693a40ed3dSBarry Smith PetscFunctionReturn(0); 170289bc588SBarry Smith } 1716ee01492SSatish Balay 1724a2ae208SSatish Balay #undef __FUNCT__ 1734a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 174af281ebdSHong Zhang PetscErrorCode MatLUFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact) 175289bc588SBarry Smith { 17602cad45dSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)(*fact)->data; 1776849ba73SBarry Smith PetscErrorCode ierr; 178899cda47SBarry Smith PetscInt lda1=mat->lda,lda2=l->lda, m=A->rmap.n,n=A->cmap.n, j; 1794482741eSBarry Smith MatFactorInfo info; 1803a40ed3dSBarry Smith 1813a40ed3dSBarry Smith PetscFunctionBegin; 18202cad45dSBarry Smith /* copy the numerical values */ 1831b807ce4Svictorle if (lda1>m || lda2>m ) { 1841b807ce4Svictorle for (j=0; j<n; j++) { 1851b807ce4Svictorle ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1861b807ce4Svictorle } 1871b807ce4Svictorle } else { 188899cda47SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 1891b807ce4Svictorle } 19002cad45dSBarry Smith (*fact)->factor = 0; 1914482741eSBarry Smith ierr = MatLUFactor(*fact,0,0,&info);CHKERRQ(ierr); 1923a40ed3dSBarry Smith PetscFunctionReturn(0); 193289bc588SBarry Smith } 1946ee01492SSatish Balay 1954a2ae208SSatish Balay #undef __FUNCT__ 1964a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 197dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat A,IS row,MatFactorInfo *info,Mat *fact) 198289bc588SBarry Smith { 199dfbe8321SBarry Smith PetscErrorCode ierr; 2003a40ed3dSBarry Smith 2013a40ed3dSBarry Smith PetscFunctionBegin; 202ceb03754SKris Buschelman ierr = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,fact);CHKERRQ(ierr); 2033a40ed3dSBarry Smith PetscFunctionReturn(0); 204289bc588SBarry Smith } 2056ee01492SSatish Balay 2064a2ae208SSatish Balay #undef __FUNCT__ 2074a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactor_SeqDense" 208dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo) 209289bc588SBarry Smith { 210a07ab388SBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 211a07ab388SBarry Smith PetscFunctionBegin; 212a07ab388SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 213a07ab388SBarry Smith #else 214c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2156849ba73SBarry Smith PetscErrorCode ierr; 216899cda47SBarry Smith PetscBLASInt n = (PetscBLASInt)A->cmap.n,info; 21755659b69SBarry Smith 2183a40ed3dSBarry Smith PetscFunctionBegin; 219606d414cSSatish Balay ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 220752f5567SLois Curfman McInnes mat->pivots = 0; 22105b42c5fSBarry Smith 222899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 22371044d3cSBarry Smith LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info); 22477431f27SBarry Smith if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 225c0bbcb79SLois Curfman McInnes A->factor = FACTOR_CHOLESKY; 226899cda47SBarry Smith ierr = PetscLogFlops((A->cmap.n*A->cmap.n*A->cmap.n)/3);CHKERRQ(ierr); 227a07ab388SBarry Smith #endif 2283a40ed3dSBarry Smith PetscFunctionReturn(0); 229289bc588SBarry Smith } 230289bc588SBarry Smith 2314a2ae208SSatish Balay #undef __FUNCT__ 2320b4b3355SBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 233af281ebdSHong Zhang PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat A,MatFactorInfo *info_dummy,Mat *fact) 2340b4b3355SBarry Smith { 235dfbe8321SBarry Smith PetscErrorCode ierr; 23615e8a5b3SHong Zhang MatFactorInfo info; 2370b4b3355SBarry Smith 2380b4b3355SBarry Smith PetscFunctionBegin; 23915e8a5b3SHong Zhang info.fill = 1.0; 24015e8a5b3SHong Zhang ierr = MatCholeskyFactor_SeqDense(*fact,0,&info);CHKERRQ(ierr); 2410b4b3355SBarry Smith PetscFunctionReturn(0); 2420b4b3355SBarry Smith } 2430b4b3355SBarry Smith 2440b4b3355SBarry Smith #undef __FUNCT__ 2454a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 246dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 247289bc588SBarry Smith { 248c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2496849ba73SBarry Smith PetscErrorCode ierr; 250899cda47SBarry Smith PetscBLASInt m = (PetscBLASInt)A->rmap.n, one = 1,info; 25187828ca2SBarry Smith PetscScalar *x,*y; 25267e560aaSBarry Smith 2533a40ed3dSBarry Smith PetscFunctionBegin; 2541ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2551ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 256899cda47SBarry Smith ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 257c0bbcb79SLois Curfman McInnes if (A->factor == FACTOR_LU) { 258ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 25980444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 260ae7cfcebSSatish Balay #else 26171044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 2624ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve"); 263ae7cfcebSSatish Balay #endif 2647a97a34bSBarry Smith } else if (A->factor == FACTOR_CHOLESKY){ 265ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 26680444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 267ae7cfcebSSatish Balay #else 26871044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 2694ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve"); 270ae7cfcebSSatish Balay #endif 271289bc588SBarry Smith } 27229bbc08cSBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 2731ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2741ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 275899cda47SBarry Smith ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr); 2763a40ed3dSBarry Smith PetscFunctionReturn(0); 277289bc588SBarry Smith } 2786ee01492SSatish Balay 2794a2ae208SSatish Balay #undef __FUNCT__ 2804a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 281dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 282da3a660dSBarry Smith { 283c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 284dfbe8321SBarry Smith PetscErrorCode ierr; 285899cda47SBarry Smith PetscBLASInt m = (PetscBLASInt) A->rmap.n,one = 1,info; 28687828ca2SBarry Smith PetscScalar *x,*y; 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); 292752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 293da3a660dSBarry Smith if (mat->pivots) { 294ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 29580444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 296ae7cfcebSSatish Balay #else 29771044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 2984ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve"); 299ae7cfcebSSatish Balay #endif 3007a97a34bSBarry Smith } else { 301ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 30280444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 303ae7cfcebSSatish Balay #else 30471044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 3054ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve"); 306ae7cfcebSSatish Balay #endif 307da3a660dSBarry Smith } 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); 312da3a660dSBarry Smith } 3136ee01492SSatish Balay 3144a2ae208SSatish Balay #undef __FUNCT__ 3154a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 316dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 317da3a660dSBarry Smith { 318c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 319dfbe8321SBarry Smith PetscErrorCode ierr; 320899cda47SBarry Smith PetscBLASInt m = (PetscBLASInt)A->rmap.n,one = 1,info; 32187828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 322da3a660dSBarry Smith Vec tmp = 0; 32367e560aaSBarry Smith 3243a40ed3dSBarry Smith PetscFunctionBegin; 3251ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3261ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 327899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 328da3a660dSBarry Smith if (yy == zz) { 32978b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 33052e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 33178b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 332da3a660dSBarry Smith } 333899cda47SBarry Smith ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 334752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 335da3a660dSBarry Smith if (mat->pivots) { 336ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 33780444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 338ae7cfcebSSatish Balay #else 33971044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 3402ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 341ae7cfcebSSatish Balay #endif 342a8c6a408SBarry Smith } else { 343ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 34480444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 345ae7cfcebSSatish Balay #else 34671044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 3472ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 348ae7cfcebSSatish Balay #endif 349da3a660dSBarry Smith } 3502dcb1b2aSMatthew Knepley if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 3512dcb1b2aSMatthew Knepley else {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);} 3521ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3531ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 354899cda47SBarry Smith ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr); 3553a40ed3dSBarry Smith PetscFunctionReturn(0); 356da3a660dSBarry Smith } 35767e560aaSBarry Smith 3584a2ae208SSatish Balay #undef __FUNCT__ 3594a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 360dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 361da3a660dSBarry Smith { 362c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3636849ba73SBarry Smith PetscErrorCode ierr; 364899cda47SBarry Smith PetscBLASInt m = (PetscBLASInt)A->rmap.n,one = 1,info; 36587828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 366da3a660dSBarry Smith Vec tmp; 36767e560aaSBarry Smith 3683a40ed3dSBarry Smith PetscFunctionBegin; 369899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 3701ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3711ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 372da3a660dSBarry Smith if (yy == zz) { 37378b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 37452e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 37578b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 376da3a660dSBarry Smith } 377899cda47SBarry Smith ierr = PetscMemcpy(y,x,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 378752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 379da3a660dSBarry Smith if (mat->pivots) { 380ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 38180444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 382ae7cfcebSSatish Balay #else 38371044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 3842ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 385ae7cfcebSSatish Balay #endif 3863a40ed3dSBarry Smith } else { 387ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 38880444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 389ae7cfcebSSatish Balay #else 39071044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 3912ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 392ae7cfcebSSatish Balay #endif 393da3a660dSBarry Smith } 39490f02eecSBarry Smith if (tmp) { 3952dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 39690f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 3973a40ed3dSBarry Smith } else { 3982dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 39990f02eecSBarry Smith } 4001ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4011ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 402899cda47SBarry Smith ierr = PetscLogFlops(2*A->cmap.n*A->cmap.n);CHKERRQ(ierr); 4033a40ed3dSBarry Smith PetscFunctionReturn(0); 404da3a660dSBarry Smith } 405289bc588SBarry Smith /* ------------------------------------------------------------------*/ 4064a2ae208SSatish Balay #undef __FUNCT__ 4074a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense" 40813f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 409289bc588SBarry Smith { 410c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 41187828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 412dfbe8321SBarry Smith PetscErrorCode ierr; 413899cda47SBarry Smith PetscInt m = A->rmap.n,i; 414aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 4154ce68768SBarry Smith PetscBLASInt bm = (PetscBLASInt)m, o = 1; 416bc1b551cSSatish Balay #endif 417289bc588SBarry Smith 4183a40ed3dSBarry Smith PetscFunctionBegin; 419289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 42071044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 4212dcb1b2aSMatthew Knepley ierr = VecSet(xx,zero);CHKERRQ(ierr); 422289bc588SBarry Smith } 4231ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4241ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 425b965ef7fSBarry Smith its = its*lits; 42677431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 427289bc588SBarry Smith while (its--) { 428fccaa45eSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 429289bc588SBarry Smith for (i=0; i<m; i++) { 430aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 431f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 432f1747703SBarry Smith not happy about returning a double complex */ 43313f74950SBarry Smith PetscInt _i; 43487828ca2SBarry Smith PetscScalar sum = b[i]; 435f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4363f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 437f1747703SBarry Smith } 438f1747703SBarry Smith xt = sum; 439f1747703SBarry Smith #else 44071044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 441f1747703SBarry Smith #endif 44255a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 443289bc588SBarry Smith } 444289bc588SBarry Smith } 445fccaa45eSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 446289bc588SBarry Smith for (i=m-1; i>=0; i--) { 447aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 448f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 449f1747703SBarry Smith not happy about returning a double complex */ 45013f74950SBarry Smith PetscInt _i; 45187828ca2SBarry Smith PetscScalar sum = b[i]; 452f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4533f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 454f1747703SBarry Smith } 455f1747703SBarry Smith xt = sum; 456f1747703SBarry Smith #else 45771044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 458f1747703SBarry Smith #endif 45955a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 460289bc588SBarry Smith } 461289bc588SBarry Smith } 462289bc588SBarry Smith } 4631ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 4641ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4653a40ed3dSBarry Smith PetscFunctionReturn(0); 466289bc588SBarry Smith } 467289bc588SBarry Smith 468289bc588SBarry Smith /* -----------------------------------------------------------------*/ 4694a2ae208SSatish Balay #undef __FUNCT__ 4704a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 471dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 472289bc588SBarry Smith { 473c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 47487828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 475dfbe8321SBarry Smith PetscErrorCode ierr; 476899cda47SBarry Smith PetscBLASInt m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n,_One=1; 477ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 4783a40ed3dSBarry Smith 4793a40ed3dSBarry Smith PetscFunctionBegin; 480899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 4811ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4821ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 48371044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One); 4841ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4851ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 486899cda47SBarry Smith ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->cmap.n);CHKERRQ(ierr); 4873a40ed3dSBarry Smith PetscFunctionReturn(0); 488289bc588SBarry Smith } 4896ee01492SSatish Balay 4904a2ae208SSatish Balay #undef __FUNCT__ 4914a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 492dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 493289bc588SBarry Smith { 494c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 49587828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 496dfbe8321SBarry Smith PetscErrorCode ierr; 497899cda47SBarry Smith PetscBLASInt m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1; 4983a40ed3dSBarry Smith 4993a40ed3dSBarry Smith PetscFunctionBegin; 500899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 5011ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5021ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 50371044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 5041ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5051ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 506899cda47SBarry Smith ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n - A->rmap.n);CHKERRQ(ierr); 5073a40ed3dSBarry Smith PetscFunctionReturn(0); 508289bc588SBarry Smith } 5096ee01492SSatish Balay 5104a2ae208SSatish Balay #undef __FUNCT__ 5114a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 512dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 513289bc588SBarry Smith { 514c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 51587828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 516dfbe8321SBarry Smith PetscErrorCode ierr; 517899cda47SBarry Smith PetscBLASInt m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1; 5183a40ed3dSBarry Smith 5193a40ed3dSBarry Smith PetscFunctionBegin; 520899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 521600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5221ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5231ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 52471044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5251ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5261ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 527899cda47SBarry Smith ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr); 5283a40ed3dSBarry Smith PetscFunctionReturn(0); 529289bc588SBarry Smith } 5306ee01492SSatish Balay 5314a2ae208SSatish Balay #undef __FUNCT__ 5324a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 533dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 534289bc588SBarry Smith { 535c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 53687828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 537dfbe8321SBarry Smith PetscErrorCode ierr; 538899cda47SBarry Smith PetscBLASInt m = (PetscBLASInt)A->rmap.n, n = (PetscBLASInt)A->cmap.n, _One=1; 53987828ca2SBarry Smith PetscScalar _DOne=1.0; 5403a40ed3dSBarry Smith 5413a40ed3dSBarry Smith PetscFunctionBegin; 542899cda47SBarry Smith if (!A->rmap.n || !A->cmap.n) PetscFunctionReturn(0); 543600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5441ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5451ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 54671044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5471ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5481ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 549899cda47SBarry Smith ierr = PetscLogFlops(2*A->rmap.n*A->cmap.n);CHKERRQ(ierr); 5503a40ed3dSBarry Smith PetscFunctionReturn(0); 551289bc588SBarry Smith } 552289bc588SBarry Smith 553289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5544a2ae208SSatish Balay #undef __FUNCT__ 5554a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 55613f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 557289bc588SBarry Smith { 558c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 55987828ca2SBarry Smith PetscScalar *v; 5606849ba73SBarry Smith PetscErrorCode ierr; 56113f74950SBarry Smith PetscInt i; 56267e560aaSBarry Smith 5633a40ed3dSBarry Smith PetscFunctionBegin; 564899cda47SBarry Smith *ncols = A->cmap.n; 565289bc588SBarry Smith if (cols) { 566899cda47SBarry Smith ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr); 567899cda47SBarry Smith for (i=0; i<A->cmap.n; i++) (*cols)[i] = i; 568289bc588SBarry Smith } 569289bc588SBarry Smith if (vals) { 570899cda47SBarry Smith ierr = PetscMalloc((A->cmap.n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 571289bc588SBarry Smith v = mat->v + row; 572899cda47SBarry Smith for (i=0; i<A->cmap.n; i++) {(*vals)[i] = *v; v += mat->lda;} 573289bc588SBarry Smith } 5743a40ed3dSBarry Smith PetscFunctionReturn(0); 575289bc588SBarry Smith } 5766ee01492SSatish Balay 5774a2ae208SSatish Balay #undef __FUNCT__ 5784a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 57913f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 580289bc588SBarry Smith { 581dfbe8321SBarry Smith PetscErrorCode ierr; 582606d414cSSatish Balay PetscFunctionBegin; 583606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 584606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 5853a40ed3dSBarry Smith PetscFunctionReturn(0); 586289bc588SBarry Smith } 587289bc588SBarry Smith /* ----------------------------------------------------------------*/ 5884a2ae208SSatish Balay #undef __FUNCT__ 5894a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 59013f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 591289bc588SBarry Smith { 592c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 59313f74950SBarry Smith PetscInt i,j,idx=0; 594d6dfbf8fSBarry Smith 5953a40ed3dSBarry Smith PetscFunctionBegin; 596289bc588SBarry Smith if (!mat->roworiented) { 597dbb450caSBarry Smith if (addv == INSERT_VALUES) { 598289bc588SBarry Smith for (j=0; j<n; j++) { 599cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6002515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 601899cda47SBarry 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); 60258804f6dSBarry Smith #endif 603289bc588SBarry Smith for (i=0; i<m; i++) { 604cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6052515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 606899cda47SBarry 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); 60758804f6dSBarry Smith #endif 608cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 609289bc588SBarry Smith } 610289bc588SBarry Smith } 6113a40ed3dSBarry Smith } else { 612289bc588SBarry Smith for (j=0; j<n; j++) { 613cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6142515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 615899cda47SBarry 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); 61658804f6dSBarry Smith #endif 617289bc588SBarry Smith for (i=0; i<m; i++) { 618cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6192515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 620899cda47SBarry 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); 62158804f6dSBarry Smith #endif 622cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 623289bc588SBarry Smith } 624289bc588SBarry Smith } 625289bc588SBarry Smith } 6263a40ed3dSBarry Smith } else { 627dbb450caSBarry Smith if (addv == INSERT_VALUES) { 628e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 629cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 6302515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 631899cda47SBarry 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); 63258804f6dSBarry Smith #endif 633e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 634cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 6352515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 636899cda47SBarry 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); 63758804f6dSBarry Smith #endif 638cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 639e8d4e0b9SBarry Smith } 640e8d4e0b9SBarry Smith } 6413a40ed3dSBarry Smith } else { 642289bc588SBarry Smith for (i=0; i<m; i++) { 643cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 6442515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 645899cda47SBarry 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); 64658804f6dSBarry Smith #endif 647289bc588SBarry Smith for (j=0; j<n; j++) { 648cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 6492515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 650899cda47SBarry 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); 65158804f6dSBarry Smith #endif 652cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 653289bc588SBarry Smith } 654289bc588SBarry Smith } 655289bc588SBarry Smith } 656e8d4e0b9SBarry Smith } 6573a40ed3dSBarry Smith PetscFunctionReturn(0); 658289bc588SBarry Smith } 659e8d4e0b9SBarry Smith 6604a2ae208SSatish Balay #undef __FUNCT__ 6614a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 66213f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 663ae80bb75SLois Curfman McInnes { 664ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 66513f74950SBarry Smith PetscInt i,j; 666ae80bb75SLois Curfman McInnes 6673a40ed3dSBarry Smith PetscFunctionBegin; 668ae80bb75SLois Curfman McInnes /* row-oriented output */ 669ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 67097e567efSBarry Smith if (indexm[i] < 0) {v += n;continue;} 671ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 67297e567efSBarry Smith if (indexn[i] < 0) {v++; continue;} 67397e567efSBarry Smith *v++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 674ae80bb75SLois Curfman McInnes } 675ae80bb75SLois Curfman McInnes } 6763a40ed3dSBarry Smith PetscFunctionReturn(0); 677ae80bb75SLois Curfman McInnes } 678ae80bb75SLois Curfman McInnes 679289bc588SBarry Smith /* -----------------------------------------------------------------*/ 680289bc588SBarry Smith 681e090d566SSatish Balay #include "petscsys.h" 68288e32edaSLois Curfman McInnes 6834a2ae208SSatish Balay #undef __FUNCT__ 6844a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense" 685f69a0ea3SMatthew Knepley PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, MatType type,Mat *A) 68688e32edaSLois Curfman McInnes { 68788e32edaSLois Curfman McInnes Mat_SeqDense *a; 68888e32edaSLois Curfman McInnes Mat B; 6896849ba73SBarry Smith PetscErrorCode ierr; 69013f74950SBarry Smith PetscInt *scols,i,j,nz,header[4]; 69113f74950SBarry Smith int fd; 69213f74950SBarry Smith PetscMPIInt size; 69313f74950SBarry Smith PetscInt *rowlengths = 0,M,N,*cols; 69487828ca2SBarry Smith PetscScalar *vals,*svals,*v,*w; 69519bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 69688e32edaSLois Curfman McInnes 6973a40ed3dSBarry Smith PetscFunctionBegin; 698d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 69929bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 700b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 7010752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 702552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 70388e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 70488e32edaSLois Curfman McInnes 70590ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 706f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 707f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr); 708be5d1d56SKris Buschelman ierr = MatSetType(*A,type);CHKERRQ(ierr); 7095c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr); 71090ace30eSBarry Smith B = *A; 71190ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 7124a41a4d3SSatish Balay v = a->v; 7134a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 7144a41a4d3SSatish Balay from row major to column major */ 715b72907f3SBarry Smith ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 71690ace30eSBarry Smith /* read in nonzero values */ 7174a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 7184a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 7194a41a4d3SSatish Balay for (j=0; j<N; j++) { 7204a41a4d3SSatish Balay for (i=0; i<M; i++) { 7214a41a4d3SSatish Balay *v++ =w[i*N+j]; 7224a41a4d3SSatish Balay } 7234a41a4d3SSatish Balay } 724606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 7256d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7266d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 72790ace30eSBarry Smith } else { 72888e32edaSLois Curfman McInnes /* read row lengths */ 72913f74950SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 7300752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 73188e32edaSLois Curfman McInnes 73288e32edaSLois Curfman McInnes /* create our matrix */ 733f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 734f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr); 735be5d1d56SKris Buschelman ierr = MatSetType(*A,type);CHKERRQ(ierr); 7365c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr); 73788e32edaSLois Curfman McInnes B = *A; 73888e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 73988e32edaSLois Curfman McInnes v = a->v; 74088e32edaSLois Curfman McInnes 74188e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 74213f74950SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr); 743b0a32e0cSBarry Smith cols = scols; 7440752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 74587828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 746b0a32e0cSBarry Smith vals = svals; 7470752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 74888e32edaSLois Curfman McInnes 74988e32edaSLois Curfman McInnes /* insert into matrix */ 75088e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 75188e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 75288e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 75388e32edaSLois Curfman McInnes } 754606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 755606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 756606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 75788e32edaSLois Curfman McInnes 7586d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7596d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 76090ace30eSBarry Smith } 7613a40ed3dSBarry Smith PetscFunctionReturn(0); 76288e32edaSLois Curfman McInnes } 76388e32edaSLois Curfman McInnes 764e090d566SSatish Balay #include "petscsys.h" 765932b0c3eSLois Curfman McInnes 7664a2ae208SSatish Balay #undef __FUNCT__ 7674a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 7686849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 769289bc588SBarry Smith { 770932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 771dfbe8321SBarry Smith PetscErrorCode ierr; 77213f74950SBarry Smith PetscInt i,j; 7732dcb1b2aSMatthew Knepley const char *name; 77487828ca2SBarry Smith PetscScalar *v; 775f3ef73ceSBarry Smith PetscViewerFormat format; 7765f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX) 7775f481a85SSatish Balay PetscTruth allreal = PETSC_TRUE; 7785f481a85SSatish Balay #endif 779932b0c3eSLois Curfman McInnes 7803a40ed3dSBarry Smith PetscFunctionBegin; 781435da068SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 782b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 783456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 7843a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 785fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 786b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 787899cda47SBarry Smith for (i=0; i<A->rmap.n; i++) { 78844cd7ae7SLois Curfman McInnes v = a->v + i; 78977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 790899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 791aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 792329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 793a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 794329f5518SBarry Smith } else if (PetscRealPart(*v)) { 795a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr); 7966831982aSBarry Smith } 79780cd9d93SLois Curfman McInnes #else 7986831982aSBarry Smith if (*v) { 799a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr); 8006831982aSBarry Smith } 80180cd9d93SLois Curfman McInnes #endif 8021b807ce4Svictorle v += a->lda; 80380cd9d93SLois Curfman McInnes } 804b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 80580cd9d93SLois Curfman McInnes } 806b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 8073a40ed3dSBarry Smith } else { 808b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 809aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 81047989497SBarry Smith /* determine if matrix has all real values */ 81147989497SBarry Smith v = a->v; 812899cda47SBarry Smith for (i=0; i<A->rmap.n*A->cmap.n; i++) { 813ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 81447989497SBarry Smith } 81547989497SBarry Smith #endif 816fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 8173a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 818899cda47SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap.n,A->cmap.n);CHKERRQ(ierr); 819899cda47SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap.n,A->cmap.n);CHKERRQ(ierr); 820fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 821ffac6cdbSBarry Smith } 822ffac6cdbSBarry Smith 823899cda47SBarry Smith for (i=0; i<A->rmap.n; i++) { 824932b0c3eSLois Curfman McInnes v = a->v + i; 825899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 826aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 82747989497SBarry Smith if (allreal) { 828f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr); 82947989497SBarry Smith } else { 830f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 83147989497SBarry Smith } 832289bc588SBarry Smith #else 833f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr); 834289bc588SBarry Smith #endif 8351b807ce4Svictorle v += a->lda; 836289bc588SBarry Smith } 837b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 838289bc588SBarry Smith } 839fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 840b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 841ffac6cdbSBarry Smith } 842b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 843da3a660dSBarry Smith } 844b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 8453a40ed3dSBarry Smith PetscFunctionReturn(0); 846289bc588SBarry Smith } 847289bc588SBarry Smith 8484a2ae208SSatish Balay #undef __FUNCT__ 8494a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 8506849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 851932b0c3eSLois Curfman McInnes { 852932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 8536849ba73SBarry Smith PetscErrorCode ierr; 85413f74950SBarry Smith int fd; 855899cda47SBarry Smith PetscInt ict,j,n = A->cmap.n,m = A->rmap.n,i,*col_lens,nz = m*n; 85687828ca2SBarry Smith PetscScalar *v,*anonz,*vals; 857f3ef73ceSBarry Smith PetscViewerFormat format; 858932b0c3eSLois Curfman McInnes 8593a40ed3dSBarry Smith PetscFunctionBegin; 860b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 86190ace30eSBarry Smith 862b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 863fb9695e5SSatish Balay if (format == PETSC_VIEWER_BINARY_NATIVE) { 86490ace30eSBarry Smith /* store the matrix as a dense matrix */ 86513f74950SBarry Smith ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 866552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 86790ace30eSBarry Smith col_lens[1] = m; 86890ace30eSBarry Smith col_lens[2] = n; 86990ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 8706f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 871606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 87290ace30eSBarry Smith 87390ace30eSBarry Smith /* write out matrix, by rows */ 87487828ca2SBarry Smith ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 87590ace30eSBarry Smith v = a->v; 87690ace30eSBarry Smith for (i=0; i<m; i++) { 87790ace30eSBarry Smith for (j=0; j<n; j++) { 87890ace30eSBarry Smith vals[i + j*m] = *v++; 87990ace30eSBarry Smith } 88090ace30eSBarry Smith } 8816f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 882606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 88390ace30eSBarry Smith } else { 88413f74950SBarry Smith ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 885552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 886932b0c3eSLois Curfman McInnes col_lens[1] = m; 887932b0c3eSLois Curfman McInnes col_lens[2] = n; 888932b0c3eSLois Curfman McInnes col_lens[3] = nz; 889932b0c3eSLois Curfman McInnes 890932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 891932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 8926f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 893932b0c3eSLois Curfman McInnes 894932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 895932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 896932b0c3eSLois Curfman McInnes ict = 0; 897932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 898932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 899932b0c3eSLois Curfman McInnes } 9006f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 901606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 902932b0c3eSLois Curfman McInnes 903932b0c3eSLois Curfman McInnes /* store nonzero values */ 90487828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 905932b0c3eSLois Curfman McInnes ict = 0; 906932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 907932b0c3eSLois Curfman McInnes v = a->v + i; 908932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 9091b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 910932b0c3eSLois Curfman McInnes } 911932b0c3eSLois Curfman McInnes } 9126f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 913606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 91490ace30eSBarry Smith } 9153a40ed3dSBarry Smith PetscFunctionReturn(0); 916932b0c3eSLois Curfman McInnes } 917932b0c3eSLois Curfman McInnes 9184a2ae208SSatish Balay #undef __FUNCT__ 9194a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 920dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 921f1af5d2fSBarry Smith { 922f1af5d2fSBarry Smith Mat A = (Mat) Aa; 923f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9246849ba73SBarry Smith PetscErrorCode ierr; 925899cda47SBarry Smith PetscInt m = A->rmap.n,n = A->cmap.n,color,i,j; 92687828ca2SBarry Smith PetscScalar *v = a->v; 927b0a32e0cSBarry Smith PetscViewer viewer; 928b0a32e0cSBarry Smith PetscDraw popup; 929329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 930f3ef73ceSBarry Smith PetscViewerFormat format; 931f1af5d2fSBarry Smith 932f1af5d2fSBarry Smith PetscFunctionBegin; 933f1af5d2fSBarry Smith 934f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 935b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 936b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 937f1af5d2fSBarry Smith 938f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 939fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 940f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 941b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 942f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 943f1af5d2fSBarry Smith x_l = j; 944f1af5d2fSBarry Smith x_r = x_l + 1.0; 945f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 946f1af5d2fSBarry Smith y_l = m - i - 1.0; 947f1af5d2fSBarry Smith y_r = y_l + 1.0; 948f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 949329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 950b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 951329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 952b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 953f1af5d2fSBarry Smith } else { 954f1af5d2fSBarry Smith continue; 955f1af5d2fSBarry Smith } 956f1af5d2fSBarry Smith #else 957f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 958b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 959f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 960b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 961f1af5d2fSBarry Smith } else { 962f1af5d2fSBarry Smith continue; 963f1af5d2fSBarry Smith } 964f1af5d2fSBarry Smith #endif 965b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 966f1af5d2fSBarry Smith } 967f1af5d2fSBarry Smith } 968f1af5d2fSBarry Smith } else { 969f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 970f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 971f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 972f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 973f1af5d2fSBarry Smith } 974b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 975b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 976b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 977f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 978f1af5d2fSBarry Smith x_l = j; 979f1af5d2fSBarry Smith x_r = x_l + 1.0; 980f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 981f1af5d2fSBarry Smith y_l = m - i - 1.0; 982f1af5d2fSBarry Smith y_r = y_l + 1.0; 983b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 984b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 985f1af5d2fSBarry Smith } 986f1af5d2fSBarry Smith } 987f1af5d2fSBarry Smith } 988f1af5d2fSBarry Smith PetscFunctionReturn(0); 989f1af5d2fSBarry Smith } 990f1af5d2fSBarry Smith 9914a2ae208SSatish Balay #undef __FUNCT__ 9924a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 993dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 994f1af5d2fSBarry Smith { 995b0a32e0cSBarry Smith PetscDraw draw; 996f1af5d2fSBarry Smith PetscTruth isnull; 997329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 998dfbe8321SBarry Smith PetscErrorCode ierr; 999f1af5d2fSBarry Smith 1000f1af5d2fSBarry Smith PetscFunctionBegin; 1001b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1002b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1003abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 1004f1af5d2fSBarry Smith 1005f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1006899cda47SBarry Smith xr = A->cmap.n; yr = A->rmap.n; h = yr/10.0; w = xr/10.0; 1007f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 1008b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1009b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 1010f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 1011f1af5d2fSBarry Smith PetscFunctionReturn(0); 1012f1af5d2fSBarry Smith } 1013f1af5d2fSBarry Smith 10144a2ae208SSatish Balay #undef __FUNCT__ 10154a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 1016dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1017932b0c3eSLois Curfman McInnes { 1018dfbe8321SBarry Smith PetscErrorCode ierr; 10196805f65bSBarry Smith PetscTruth iascii,isbinary,isdraw; 1020932b0c3eSLois Curfman McInnes 10213a40ed3dSBarry Smith PetscFunctionBegin; 102232077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1023fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 1024fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 10250f5bd95cSBarry Smith 1026c45a1595SBarry Smith if (iascii) { 1027c45a1595SBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 10280f5bd95cSBarry Smith } else if (isbinary) { 10293a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1030f1af5d2fSBarry Smith } else if (isdraw) { 1031f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 10325cd90555SBarry Smith } else { 1033958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1034932b0c3eSLois Curfman McInnes } 10353a40ed3dSBarry Smith PetscFunctionReturn(0); 1036932b0c3eSLois Curfman McInnes } 1037289bc588SBarry Smith 10384a2ae208SSatish Balay #undef __FUNCT__ 10394a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1040dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1041289bc588SBarry Smith { 1042ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1043dfbe8321SBarry Smith PetscErrorCode ierr; 104490f02eecSBarry Smith 10453a40ed3dSBarry Smith PetscFunctionBegin; 1046aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1047899cda47SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap.n,mat->cmap.n); 1048a5a9c739SBarry Smith #endif 104905b42c5fSBarry Smith ierr = PetscFree(l->pivots);CHKERRQ(ierr); 10506857c123SSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1051606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 1052dbd8c25aSHong Zhang 1053dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 1054901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 10554ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 10564ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 10574ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 10583a40ed3dSBarry Smith PetscFunctionReturn(0); 1059289bc588SBarry Smith } 1060289bc588SBarry Smith 10614a2ae208SSatish Balay #undef __FUNCT__ 10624a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1063dfbe8321SBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,Mat *matout) 1064289bc588SBarry Smith { 1065c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 10666849ba73SBarry Smith PetscErrorCode ierr; 106713f74950SBarry Smith PetscInt k,j,m,n,M; 106887828ca2SBarry Smith PetscScalar *v,tmp; 106948b35521SBarry Smith 10703a40ed3dSBarry Smith PetscFunctionBegin; 1071899cda47SBarry Smith v = mat->v; m = A->rmap.n; M = mat->lda; n = A->cmap.n; 10727c922b88SBarry Smith if (!matout) { /* in place transpose */ 1073a5ce6ee0Svictorle if (m != n) { 1074634064b4SBarry Smith SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1075d64ed03dSBarry Smith } else { 1076d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1077289bc588SBarry Smith for (k=0; k<j; k++) { 10781b807ce4Svictorle tmp = v[j + k*M]; 10791b807ce4Svictorle v[j + k*M] = v[k + j*M]; 10801b807ce4Svictorle v[k + j*M] = tmp; 1081289bc588SBarry Smith } 1082289bc588SBarry Smith } 1083d64ed03dSBarry Smith } 10843a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1085d3e5ee88SLois Curfman McInnes Mat tmat; 1086ec8511deSBarry Smith Mat_SeqDense *tmatd; 108787828ca2SBarry Smith PetscScalar *v2; 1088ea709b57SSatish Balay 1089f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&tmat);CHKERRQ(ierr); 1090899cda47SBarry Smith ierr = MatSetSizes(tmat,A->cmap.n,A->rmap.n,A->cmap.n,A->rmap.n);CHKERRQ(ierr); 10915c5985e7SKris Buschelman ierr = MatSetType(tmat,A->type_name);CHKERRQ(ierr); 10925c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1093ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 10940de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1095d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 10961b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1097d3e5ee88SLois Curfman McInnes } 10986d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10996d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1100d3e5ee88SLois Curfman McInnes *matout = tmat; 110148b35521SBarry Smith } 11023a40ed3dSBarry Smith PetscFunctionReturn(0); 1103289bc588SBarry Smith } 1104289bc588SBarry Smith 11054a2ae208SSatish Balay #undef __FUNCT__ 11064a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1107dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1108289bc588SBarry Smith { 1109c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1110c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 111113f74950SBarry Smith PetscInt i,j; 111287828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 11139ea5d5aeSSatish Balay 11143a40ed3dSBarry Smith PetscFunctionBegin; 1115899cda47SBarry Smith if (A1->rmap.n != A2->rmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1116899cda47SBarry Smith if (A1->cmap.n != A2->cmap.n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1117899cda47SBarry Smith for (i=0; i<A1->rmap.n; i++) { 11181b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1119899cda47SBarry Smith for (j=0; j<A1->cmap.n; j++) { 11203a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 11211b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 11221b807ce4Svictorle } 1123289bc588SBarry Smith } 112477c4ece6SBarry Smith *flg = PETSC_TRUE; 11253a40ed3dSBarry Smith PetscFunctionReturn(0); 1126289bc588SBarry Smith } 1127289bc588SBarry Smith 11284a2ae208SSatish Balay #undef __FUNCT__ 11294a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1130dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1131289bc588SBarry Smith { 1132c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1133dfbe8321SBarry Smith PetscErrorCode ierr; 113413f74950SBarry Smith PetscInt i,n,len; 113587828ca2SBarry Smith PetscScalar *x,zero = 0.0; 113644cd7ae7SLois Curfman McInnes 11373a40ed3dSBarry Smith PetscFunctionBegin; 11382dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 11397a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 11401ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1141899cda47SBarry Smith len = PetscMin(A->rmap.n,A->cmap.n); 1142899cda47SBarry Smith if (n != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 114344cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 11441b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1145289bc588SBarry Smith } 11461ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 11473a40ed3dSBarry Smith PetscFunctionReturn(0); 1148289bc588SBarry Smith } 1149289bc588SBarry Smith 11504a2ae208SSatish Balay #undef __FUNCT__ 11514a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1152dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1153289bc588SBarry Smith { 1154c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 115587828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1156dfbe8321SBarry Smith PetscErrorCode ierr; 1157899cda47SBarry Smith PetscInt i,j,m = A->rmap.n,n = A->cmap.n; 115855659b69SBarry Smith 11593a40ed3dSBarry Smith PetscFunctionBegin; 116028988994SBarry Smith if (ll) { 11617a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 11621ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1163899cda47SBarry Smith if (m != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1164da3a660dSBarry Smith for (i=0; i<m; i++) { 1165da3a660dSBarry Smith x = l[i]; 1166da3a660dSBarry Smith v = mat->v + i; 1167da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1168da3a660dSBarry Smith } 11691ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1170efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1171da3a660dSBarry Smith } 117228988994SBarry Smith if (rr) { 11737a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 11741ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1175899cda47SBarry Smith if (n != A->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1176da3a660dSBarry Smith for (i=0; i<n; i++) { 1177da3a660dSBarry Smith x = r[i]; 1178da3a660dSBarry Smith v = mat->v + i*m; 1179da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1180da3a660dSBarry Smith } 11811ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1182efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1183da3a660dSBarry Smith } 11843a40ed3dSBarry Smith PetscFunctionReturn(0); 1185289bc588SBarry Smith } 1186289bc588SBarry Smith 11874a2ae208SSatish Balay #undef __FUNCT__ 11884a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1189dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1190289bc588SBarry Smith { 1191c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 119287828ca2SBarry Smith PetscScalar *v = mat->v; 1193329f5518SBarry Smith PetscReal sum = 0.0; 1194899cda47SBarry Smith PetscInt lda=mat->lda,m=A->rmap.n,i,j; 1195efee365bSSatish Balay PetscErrorCode ierr; 119655659b69SBarry Smith 11973a40ed3dSBarry Smith PetscFunctionBegin; 1198289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1199a5ce6ee0Svictorle if (lda>m) { 1200899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 1201a5ce6ee0Svictorle v = mat->v+j*lda; 1202a5ce6ee0Svictorle for (i=0; i<m; i++) { 1203a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1204a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1205a5ce6ee0Svictorle #else 1206a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1207a5ce6ee0Svictorle #endif 1208a5ce6ee0Svictorle } 1209a5ce6ee0Svictorle } 1210a5ce6ee0Svictorle } else { 1211899cda47SBarry Smith for (i=0; i<A->cmap.n*A->rmap.n; i++) { 1212aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1213329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1214289bc588SBarry Smith #else 1215289bc588SBarry Smith sum += (*v)*(*v); v++; 1216289bc588SBarry Smith #endif 1217289bc588SBarry Smith } 1218a5ce6ee0Svictorle } 1219064f8208SBarry Smith *nrm = sqrt(sum); 1220899cda47SBarry Smith ierr = PetscLogFlops(2*A->cmap.n*A->rmap.n);CHKERRQ(ierr); 12213a40ed3dSBarry Smith } else if (type == NORM_1) { 1222064f8208SBarry Smith *nrm = 0.0; 1223899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 12241b807ce4Svictorle v = mat->v + j*mat->lda; 1225289bc588SBarry Smith sum = 0.0; 1226899cda47SBarry Smith for (i=0; i<A->rmap.n; i++) { 122733a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1228289bc588SBarry Smith } 1229064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1230289bc588SBarry Smith } 1231899cda47SBarry Smith ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr); 12323a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1233064f8208SBarry Smith *nrm = 0.0; 1234899cda47SBarry Smith for (j=0; j<A->rmap.n; j++) { 1235289bc588SBarry Smith v = mat->v + j; 1236289bc588SBarry Smith sum = 0.0; 1237899cda47SBarry Smith for (i=0; i<A->cmap.n; i++) { 12381b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1239289bc588SBarry Smith } 1240064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1241289bc588SBarry Smith } 1242899cda47SBarry Smith ierr = PetscLogFlops(A->cmap.n*A->rmap.n);CHKERRQ(ierr); 12433a40ed3dSBarry Smith } else { 124429bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1245289bc588SBarry Smith } 12463a40ed3dSBarry Smith PetscFunctionReturn(0); 1247289bc588SBarry Smith } 1248289bc588SBarry Smith 12494a2ae208SSatish Balay #undef __FUNCT__ 12504a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 12514e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscTruth flg) 1252289bc588SBarry Smith { 1253c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 125463ba0a88SBarry Smith PetscErrorCode ierr; 125567e560aaSBarry Smith 12563a40ed3dSBarry Smith PetscFunctionBegin; 1257b5a2b587SKris Buschelman switch (op) { 1258b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 12594e0d8c25SBarry Smith aij->roworiented = flg; 1260b5a2b587SKris Buschelman break; 1261512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 1262b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 1263*3971808eSMatthew Knepley case MAT_NEW_NONZERO_ALLOCATION_ERR: 12644e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 1265b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1266b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 126777e54ba9SKris Buschelman case MAT_SYMMETRIC: 126877e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 12699a4540c5SBarry Smith case MAT_HERMITIAN: 12709a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1271290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 127277e54ba9SKris Buschelman break; 1273b5a2b587SKris Buschelman default: 1274ad86a440SBarry Smith SETERRQ1(PETSC_ERR_SUP,"unknown option %d",op); 12753a40ed3dSBarry Smith } 12763a40ed3dSBarry Smith PetscFunctionReturn(0); 1277289bc588SBarry Smith } 1278289bc588SBarry Smith 12794a2ae208SSatish Balay #undef __FUNCT__ 12804a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1281dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 12826f0a148fSBarry Smith { 1283ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 12846849ba73SBarry Smith PetscErrorCode ierr; 1285899cda47SBarry Smith PetscInt lda=l->lda,m=A->rmap.n,j; 12863a40ed3dSBarry Smith 12873a40ed3dSBarry Smith PetscFunctionBegin; 1288a5ce6ee0Svictorle if (lda>m) { 1289899cda47SBarry Smith for (j=0; j<A->cmap.n; j++) { 1290a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1291a5ce6ee0Svictorle } 1292a5ce6ee0Svictorle } else { 1293899cda47SBarry Smith ierr = PetscMemzero(l->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 1294a5ce6ee0Svictorle } 12953a40ed3dSBarry Smith PetscFunctionReturn(0); 12966f0a148fSBarry Smith } 12976f0a148fSBarry Smith 12984a2ae208SSatish Balay #undef __FUNCT__ 12994a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 1300f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag) 13016f0a148fSBarry Smith { 1302ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1303899cda47SBarry Smith PetscInt n = A->cmap.n,i,j; 130487828ca2SBarry Smith PetscScalar *slot; 130555659b69SBarry Smith 13063a40ed3dSBarry Smith PetscFunctionBegin; 13076f0a148fSBarry Smith for (i=0; i<N; i++) { 13086f0a148fSBarry Smith slot = l->v + rows[i]; 13096f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 13106f0a148fSBarry Smith } 1311f4df32b1SMatthew Knepley if (diag != 0.0) { 13126f0a148fSBarry Smith for (i=0; i<N; i++) { 13136f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1314f4df32b1SMatthew Knepley *slot = diag; 13156f0a148fSBarry Smith } 13166f0a148fSBarry Smith } 13173a40ed3dSBarry Smith PetscFunctionReturn(0); 13186f0a148fSBarry Smith } 1319557bce09SLois Curfman McInnes 13204a2ae208SSatish Balay #undef __FUNCT__ 13214a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 1322dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 132364e87e97SBarry Smith { 1324c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 13253a40ed3dSBarry Smith 13263a40ed3dSBarry Smith PetscFunctionBegin; 1327899cda47SBarry Smith if (mat->lda != A->rmap.n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows"); 132864e87e97SBarry Smith *array = mat->v; 13293a40ed3dSBarry Smith PetscFunctionReturn(0); 133064e87e97SBarry Smith } 13310754003eSLois Curfman McInnes 13324a2ae208SSatish Balay #undef __FUNCT__ 13334a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 1334dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1335ff14e315SSatish Balay { 13363a40ed3dSBarry Smith PetscFunctionBegin; 133709b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 13383a40ed3dSBarry Smith PetscFunctionReturn(0); 1339ff14e315SSatish Balay } 13400754003eSLois Curfman McInnes 13414a2ae208SSatish Balay #undef __FUNCT__ 13424a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 134313f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 13440754003eSLois Curfman McInnes { 1345c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 13466849ba73SBarry Smith PetscErrorCode ierr; 134721a2c019SBarry Smith PetscInt i,j,*irow,*icol,nrows,ncols; 134887828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 13490754003eSLois Curfman McInnes Mat newmat; 13500754003eSLois Curfman McInnes 13513a40ed3dSBarry Smith PetscFunctionBegin; 135278b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 135378b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1354e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1355e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 13560754003eSLois Curfman McInnes 1357182d2002SSatish Balay /* Check submatrixcall */ 1358182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 135913f74950SBarry Smith PetscInt n_cols,n_rows; 1360182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 136121a2c019SBarry Smith if (n_rows != nrows || n_cols != ncols) { 136221a2c019SBarry Smith /* resize the result result matrix to match number of requested rows/columns */ 136321a2c019SBarry Smith ierr = MatSetSizes(*B,nrows,nrows,nrows,nrows);CHKERRQ(ierr); 136421a2c019SBarry Smith } 1365182d2002SSatish Balay newmat = *B; 1366182d2002SSatish Balay } else { 13670754003eSLois Curfman McInnes /* Create and fill new matrix */ 1368f69a0ea3SMatthew Knepley ierr = MatCreate(A->comm,&newmat);CHKERRQ(ierr); 1369f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 13705c5985e7SKris Buschelman ierr = MatSetType(newmat,A->type_name);CHKERRQ(ierr); 13715c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1372182d2002SSatish Balay } 1373182d2002SSatish Balay 1374182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1375182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1376182d2002SSatish Balay 1377182d2002SSatish Balay for (i=0; i<ncols; i++) { 13786de62eeeSBarry Smith av = v + mat->lda*icol[i]; 1379182d2002SSatish Balay for (j=0; j<nrows; j++) { 1380182d2002SSatish Balay *bv++ = av[irow[j]]; 13810754003eSLois Curfman McInnes } 13820754003eSLois Curfman McInnes } 1383182d2002SSatish Balay 1384182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 13856d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13866d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13870754003eSLois Curfman McInnes 13880754003eSLois Curfman McInnes /* Free work space */ 138978b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 139078b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1391182d2002SSatish Balay *B = newmat; 13923a40ed3dSBarry Smith PetscFunctionReturn(0); 13930754003eSLois Curfman McInnes } 13940754003eSLois Curfman McInnes 13954a2ae208SSatish Balay #undef __FUNCT__ 13964a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 139713f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1398905e6a2fSBarry Smith { 13996849ba73SBarry Smith PetscErrorCode ierr; 140013f74950SBarry Smith PetscInt i; 1401905e6a2fSBarry Smith 14023a40ed3dSBarry Smith PetscFunctionBegin; 1403905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1404b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1405905e6a2fSBarry Smith } 1406905e6a2fSBarry Smith 1407905e6a2fSBarry Smith for (i=0; i<n; i++) { 14086a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1409905e6a2fSBarry Smith } 14103a40ed3dSBarry Smith PetscFunctionReturn(0); 1411905e6a2fSBarry Smith } 1412905e6a2fSBarry Smith 14134a2ae208SSatish Balay #undef __FUNCT__ 1414c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1415c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1416c0aa2d19SHong Zhang { 1417c0aa2d19SHong Zhang PetscFunctionBegin; 1418c0aa2d19SHong Zhang PetscFunctionReturn(0); 1419c0aa2d19SHong Zhang } 1420c0aa2d19SHong Zhang 1421c0aa2d19SHong Zhang #undef __FUNCT__ 1422c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1423c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1424c0aa2d19SHong Zhang { 1425c0aa2d19SHong Zhang PetscFunctionBegin; 1426c0aa2d19SHong Zhang PetscFunctionReturn(0); 1427c0aa2d19SHong Zhang } 1428c0aa2d19SHong Zhang 1429c0aa2d19SHong Zhang #undef __FUNCT__ 14304a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1431dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 14324b0e389bSBarry Smith { 14334b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 14346849ba73SBarry Smith PetscErrorCode ierr; 1435899cda47SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap.n,n=A->cmap.n, j; 14363a40ed3dSBarry Smith 14373a40ed3dSBarry Smith PetscFunctionBegin; 143833f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 143933f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1440cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 14413a40ed3dSBarry Smith PetscFunctionReturn(0); 14423a40ed3dSBarry Smith } 1443899cda47SBarry Smith if (m != B->rmap.n || n != B->cmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1444a5ce6ee0Svictorle if (lda1>m || lda2>m) { 14450dbb7854Svictorle for (j=0; j<n; j++) { 1446a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1447a5ce6ee0Svictorle } 1448a5ce6ee0Svictorle } else { 1449899cda47SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->rmap.n*A->cmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 1450a5ce6ee0Svictorle } 1451273d9f13SBarry Smith PetscFunctionReturn(0); 1452273d9f13SBarry Smith } 1453273d9f13SBarry Smith 14544a2ae208SSatish Balay #undef __FUNCT__ 14554a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1456dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A) 1457273d9f13SBarry Smith { 1458dfbe8321SBarry Smith PetscErrorCode ierr; 1459273d9f13SBarry Smith 1460273d9f13SBarry Smith PetscFunctionBegin; 1461273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 14623a40ed3dSBarry Smith PetscFunctionReturn(0); 14634b0e389bSBarry Smith } 14644b0e389bSBarry Smith 1465284134d9SBarry Smith #undef __FUNCT__ 1466284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense" 1467284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N) 1468284134d9SBarry Smith { 1469284134d9SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 147021a2c019SBarry Smith PetscErrorCode ierr; 1471284134d9SBarry Smith PetscFunctionBegin; 147221a2c019SBarry Smith /* this will not be called before lda, Mmax, and Nmax have been set */ 1473284134d9SBarry Smith m = PetscMax(m,M); 1474284134d9SBarry Smith n = PetscMax(n,N); 147521a2c019SBarry 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); 1476284134d9SBarry 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); 1477899cda47SBarry Smith A->rmap.n = A->rmap.n = m; 1478899cda47SBarry Smith A->cmap.n = A->cmap.N = n; 147921a2c019SBarry Smith if (a->changelda) a->lda = m; 148021a2c019SBarry Smith ierr = PetscMemzero(a->v,m*n*sizeof(PetscScalar));CHKERRQ(ierr); 1481284134d9SBarry Smith PetscFunctionReturn(0); 1482284134d9SBarry Smith } 1483284134d9SBarry Smith 1484a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/ 1485a9fe9ddaSSatish Balay #undef __FUNCT__ 1486a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1487a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1488a9fe9ddaSSatish Balay { 1489a9fe9ddaSSatish Balay PetscErrorCode ierr; 1490a9fe9ddaSSatish Balay 1491a9fe9ddaSSatish Balay PetscFunctionBegin; 1492a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1493a9fe9ddaSSatish Balay ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1494a9fe9ddaSSatish Balay } 1495a9fe9ddaSSatish Balay ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1496a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1497a9fe9ddaSSatish Balay } 1498a9fe9ddaSSatish Balay 1499a9fe9ddaSSatish Balay #undef __FUNCT__ 1500a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1501a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1502a9fe9ddaSSatish Balay { 1503ee16a9a1SHong Zhang PetscErrorCode ierr; 1504899cda47SBarry Smith PetscInt m=A->rmap.n,n=B->cmap.n; 1505ee16a9a1SHong Zhang Mat Cmat; 1506a9fe9ddaSSatish Balay 1507ee16a9a1SHong Zhang PetscFunctionBegin; 1508899cda47SBarry 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); 1509ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1510ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1511ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1512ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1513ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1514ee16a9a1SHong Zhang *C = Cmat; 1515ee16a9a1SHong Zhang PetscFunctionReturn(0); 1516ee16a9a1SHong Zhang } 1517a9fe9ddaSSatish Balay 1518a9fe9ddaSSatish Balay #undef __FUNCT__ 1519a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1520a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1521a9fe9ddaSSatish Balay { 1522a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1523a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1524a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 1525899cda47SBarry Smith PetscBLASInt m=(PetscBLASInt)A->rmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->cmap.n; 1526a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1527a9fe9ddaSSatish Balay 1528a9fe9ddaSSatish Balay PetscFunctionBegin; 1529a9fe9ddaSSatish Balay BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1530a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1531a9fe9ddaSSatish Balay } 1532a9fe9ddaSSatish Balay 1533a9fe9ddaSSatish Balay #undef __FUNCT__ 1534a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense" 1535a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1536a9fe9ddaSSatish Balay { 1537a9fe9ddaSSatish Balay PetscErrorCode ierr; 1538a9fe9ddaSSatish Balay 1539a9fe9ddaSSatish Balay PetscFunctionBegin; 1540a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1541a9fe9ddaSSatish Balay ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1542a9fe9ddaSSatish Balay } 1543a9fe9ddaSSatish Balay ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1544a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1545a9fe9ddaSSatish Balay } 1546a9fe9ddaSSatish Balay 1547a9fe9ddaSSatish Balay #undef __FUNCT__ 1548a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense" 1549a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1550a9fe9ddaSSatish Balay { 1551ee16a9a1SHong Zhang PetscErrorCode ierr; 1552899cda47SBarry Smith PetscInt m=A->cmap.n,n=B->cmap.n; 1553ee16a9a1SHong Zhang Mat Cmat; 1554a9fe9ddaSSatish Balay 1555ee16a9a1SHong Zhang PetscFunctionBegin; 1556899cda47SBarry 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); 1557ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1558ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1559ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1560ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1561ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1562ee16a9a1SHong Zhang *C = Cmat; 1563ee16a9a1SHong Zhang PetscFunctionReturn(0); 1564ee16a9a1SHong Zhang } 1565a9fe9ddaSSatish Balay 1566a9fe9ddaSSatish Balay #undef __FUNCT__ 1567a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense" 1568a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1569a9fe9ddaSSatish Balay { 1570a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1571a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1572a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 1573899cda47SBarry Smith PetscBLASInt m=(PetscBLASInt)A->cmap.n,n=(PetscBLASInt)B->cmap.n,k=(PetscBLASInt)A->rmap.n; 1574a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1575a9fe9ddaSSatish Balay 1576a9fe9ddaSSatish Balay PetscFunctionBegin; 1577a9fe9ddaSSatish Balay BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1578a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1579a9fe9ddaSSatish Balay } 1580985db425SBarry Smith 1581985db425SBarry Smith #undef __FUNCT__ 1582985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense" 1583985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 1584985db425SBarry Smith { 1585985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1586985db425SBarry Smith PetscErrorCode ierr; 1587985db425SBarry Smith PetscInt i,j,m = A->rmap.n,n = A->cmap.n,p; 1588985db425SBarry Smith PetscScalar *x; 1589985db425SBarry Smith MatScalar *aa = a->v; 1590985db425SBarry Smith 1591985db425SBarry Smith PetscFunctionBegin; 1592985db425SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1593985db425SBarry Smith 1594985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1595985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1596985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1597985db425SBarry Smith if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1598985db425SBarry Smith for (i=0; i<m; i++) { 1599985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1600985db425SBarry Smith for (j=1; j<n; j++){ 1601985db425SBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1602985db425SBarry Smith } 1603985db425SBarry Smith } 1604985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1605985db425SBarry Smith PetscFunctionReturn(0); 1606985db425SBarry Smith } 1607985db425SBarry Smith 1608985db425SBarry Smith #undef __FUNCT__ 1609985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense" 1610985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 1611985db425SBarry Smith { 1612985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1613985db425SBarry Smith PetscErrorCode ierr; 1614985db425SBarry Smith PetscInt i,j,m = A->rmap.n,n = A->cmap.n,p; 1615985db425SBarry Smith PetscScalar *x; 1616985db425SBarry Smith PetscReal atmp; 1617985db425SBarry Smith MatScalar *aa = a->v; 1618985db425SBarry Smith 1619985db425SBarry Smith PetscFunctionBegin; 1620985db425SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1621985db425SBarry Smith 1622985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1623985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1624985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1625985db425SBarry Smith if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1626985db425SBarry Smith for (i=0; i<m; i++) { 1627985db425SBarry Smith x[i] = PetscAbsScalar(aa[i]); if (idx) idx[i] = 0; 1628985db425SBarry Smith for (j=1; j<n; j++){ 1629985db425SBarry Smith atmp = PetscAbsScalar(aa[i+m*j]); 1630985db425SBarry Smith if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 1631985db425SBarry Smith } 1632985db425SBarry Smith } 1633985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1634985db425SBarry Smith PetscFunctionReturn(0); 1635985db425SBarry Smith } 1636985db425SBarry Smith 1637985db425SBarry Smith #undef __FUNCT__ 1638985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense" 1639985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 1640985db425SBarry Smith { 1641985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1642985db425SBarry Smith PetscErrorCode ierr; 1643985db425SBarry Smith PetscInt i,j,m = A->rmap.n,n = A->cmap.n,p; 1644985db425SBarry Smith PetscScalar *x; 1645985db425SBarry Smith MatScalar *aa = a->v; 1646985db425SBarry Smith 1647985db425SBarry Smith PetscFunctionBegin; 1648985db425SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1649985db425SBarry Smith 1650985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1651985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1652985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1653985db425SBarry Smith if (p != A->rmap.n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1654985db425SBarry Smith for (i=0; i<m; i++) { 1655985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1656985db425SBarry Smith for (j=1; j<n; j++){ 1657985db425SBarry Smith if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1658985db425SBarry Smith } 1659985db425SBarry Smith } 1660985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1661985db425SBarry Smith PetscFunctionReturn(0); 1662985db425SBarry Smith } 1663985db425SBarry Smith 16648d0534beSBarry Smith #undef __FUNCT__ 16658d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense" 16668d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 16678d0534beSBarry Smith { 16688d0534beSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 16698d0534beSBarry Smith PetscErrorCode ierr; 16708d0534beSBarry Smith PetscScalar *x; 16718d0534beSBarry Smith 16728d0534beSBarry Smith PetscFunctionBegin; 16738d0534beSBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 16748d0534beSBarry Smith 16758d0534beSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 16768d0534beSBarry Smith ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap.n*sizeof(PetscScalar));CHKERRQ(ierr); 16778d0534beSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 16788d0534beSBarry Smith PetscFunctionReturn(0); 16798d0534beSBarry Smith } 16808d0534beSBarry Smith 1681289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1682a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1683905e6a2fSBarry Smith MatGetRow_SeqDense, 1684905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1685905e6a2fSBarry Smith MatMult_SeqDense, 168697304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 16877c922b88SBarry Smith MatMultTranspose_SeqDense, 16887c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1689905e6a2fSBarry Smith MatSolve_SeqDense, 1690905e6a2fSBarry Smith MatSolveAdd_SeqDense, 16917c922b88SBarry Smith MatSolveTranspose_SeqDense, 169297304618SKris Buschelman /*10*/ MatSolveTransposeAdd_SeqDense, 1693905e6a2fSBarry Smith MatLUFactor_SeqDense, 1694905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1695ec8511deSBarry Smith MatRelax_SeqDense, 1696ec8511deSBarry Smith MatTranspose_SeqDense, 169797304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1698905e6a2fSBarry Smith MatEqual_SeqDense, 1699905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1700905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1701905e6a2fSBarry Smith MatNorm_SeqDense, 1702c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense, 1703c0aa2d19SHong Zhang MatAssemblyEnd_SeqDense, 1704905e6a2fSBarry Smith 0, 1705905e6a2fSBarry Smith MatSetOption_SeqDense, 1706905e6a2fSBarry Smith MatZeroEntries_SeqDense, 170797304618SKris Buschelman /*25*/ MatZeroRows_SeqDense, 1708905e6a2fSBarry Smith MatLUFactorSymbolic_SeqDense, 1709905e6a2fSBarry Smith MatLUFactorNumeric_SeqDense, 1710905e6a2fSBarry Smith MatCholeskyFactorSymbolic_SeqDense, 1711905e6a2fSBarry Smith MatCholeskyFactorNumeric_SeqDense, 171297304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense, 1713273d9f13SBarry Smith 0, 1714905e6a2fSBarry Smith 0, 1715905e6a2fSBarry Smith MatGetArray_SeqDense, 1716905e6a2fSBarry Smith MatRestoreArray_SeqDense, 171797304618SKris Buschelman /*35*/ MatDuplicate_SeqDense, 1718a5ae1ecdSBarry Smith 0, 1719a5ae1ecdSBarry Smith 0, 1720a5ae1ecdSBarry Smith 0, 1721a5ae1ecdSBarry Smith 0, 172297304618SKris Buschelman /*40*/ MatAXPY_SeqDense, 1723a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1724a5ae1ecdSBarry Smith 0, 17254b0e389bSBarry Smith MatGetValues_SeqDense, 1726a5ae1ecdSBarry Smith MatCopy_SeqDense, 1727985db425SBarry Smith /*45*/ MatGetRowMax_SeqDense, 1728a5ae1ecdSBarry Smith MatScale_SeqDense, 1729a5ae1ecdSBarry Smith 0, 1730a5ae1ecdSBarry Smith 0, 1731a5ae1ecdSBarry Smith 0, 1732521d7252SBarry Smith /*50*/ 0, 1733a5ae1ecdSBarry Smith 0, 1734a5ae1ecdSBarry Smith 0, 1735a5ae1ecdSBarry Smith 0, 1736a5ae1ecdSBarry Smith 0, 173797304618SKris Buschelman /*55*/ 0, 1738a5ae1ecdSBarry Smith 0, 1739a5ae1ecdSBarry Smith 0, 1740a5ae1ecdSBarry Smith 0, 1741a5ae1ecdSBarry Smith 0, 174297304618SKris Buschelman /*60*/ 0, 1743e03a110bSBarry Smith MatDestroy_SeqDense, 1744e03a110bSBarry Smith MatView_SeqDense, 1745357abbc8SBarry Smith 0, 174697304618SKris Buschelman 0, 174797304618SKris Buschelman /*65*/ 0, 174897304618SKris Buschelman 0, 174997304618SKris Buschelman 0, 175097304618SKris Buschelman 0, 175197304618SKris Buschelman 0, 1752985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqDense, 175397304618SKris Buschelman 0, 175497304618SKris Buschelman 0, 175597304618SKris Buschelman 0, 175697304618SKris Buschelman 0, 175797304618SKris Buschelman /*75*/ 0, 175897304618SKris Buschelman 0, 175997304618SKris Buschelman 0, 176097304618SKris Buschelman 0, 176197304618SKris Buschelman 0, 176297304618SKris Buschelman /*80*/ 0, 176397304618SKris Buschelman 0, 176497304618SKris Buschelman 0, 176597304618SKris Buschelman 0, 1766865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense, 1767865e5f61SKris Buschelman 0, 17681cbb95d3SBarry Smith MatIsHermitian_SeqDense, 1769865e5f61SKris Buschelman 0, 1770865e5f61SKris Buschelman 0, 1771865e5f61SKris Buschelman 0, 1772a9fe9ddaSSatish Balay /*90*/ MatMatMult_SeqDense_SeqDense, 1773a9fe9ddaSSatish Balay MatMatMultSymbolic_SeqDense_SeqDense, 1774a9fe9ddaSSatish Balay MatMatMultNumeric_SeqDense_SeqDense, 1775865e5f61SKris Buschelman 0, 1776865e5f61SKris Buschelman 0, 1777865e5f61SKris Buschelman /*95*/ 0, 1778a9fe9ddaSSatish Balay MatMatMultTranspose_SeqDense_SeqDense, 1779a9fe9ddaSSatish Balay MatMatMultTransposeSymbolic_SeqDense_SeqDense, 1780a9fe9ddaSSatish Balay MatMatMultTransposeNumeric_SeqDense_SeqDense, 1781284134d9SBarry Smith 0, 1782284134d9SBarry Smith /*100*/0, 1783284134d9SBarry Smith 0, 1784284134d9SBarry Smith 0, 1785284134d9SBarry Smith 0, 1786985db425SBarry Smith MatSetSizes_SeqDense, 1787985db425SBarry Smith 0, 1788985db425SBarry Smith 0, 1789985db425SBarry Smith 0, 1790985db425SBarry Smith 0, 1791985db425SBarry Smith 0, 1792985db425SBarry Smith /*110*/0, 1793985db425SBarry Smith 0, 17948d0534beSBarry Smith MatGetRowMin_SeqDense, 17958d0534beSBarry Smith MatGetColumnVector_SeqDense 1796985db425SBarry Smith }; 179790ace30eSBarry Smith 17984a2ae208SSatish Balay #undef __FUNCT__ 17994a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 18004b828684SBarry Smith /*@C 1801fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1802d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1803d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1804289bc588SBarry Smith 1805db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1806db81eaa0SLois Curfman McInnes 180720563c6bSBarry Smith Input Parameters: 1808db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 18090c775827SLois Curfman McInnes . m - number of rows 181018f449edSLois Curfman McInnes . n - number of columns 1811db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1812dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 181320563c6bSBarry Smith 181420563c6bSBarry Smith Output Parameter: 181544cd7ae7SLois Curfman McInnes . A - the matrix 181620563c6bSBarry Smith 1817b259b22eSLois Curfman McInnes Notes: 181818f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 181918f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1820b4fd4287SBarry Smith set data=PETSC_NULL. 182118f449edSLois Curfman McInnes 1822027ccd11SLois Curfman McInnes Level: intermediate 1823027ccd11SLois Curfman McInnes 1824dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1825d65003e9SLois Curfman McInnes 1826db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 182720563c6bSBarry Smith @*/ 1828be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 1829289bc588SBarry Smith { 1830dfbe8321SBarry Smith PetscErrorCode ierr; 18313b2fbd54SBarry Smith 18323a40ed3dSBarry Smith PetscFunctionBegin; 1833f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1834f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1835273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1836273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1837273d9f13SBarry Smith PetscFunctionReturn(0); 1838273d9f13SBarry Smith } 1839273d9f13SBarry Smith 18404a2ae208SSatish Balay #undef __FUNCT__ 1841afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation" 1842273d9f13SBarry Smith /*@C 1843273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1844273d9f13SBarry Smith 1845273d9f13SBarry Smith Collective on MPI_Comm 1846273d9f13SBarry Smith 1847273d9f13SBarry Smith Input Parameters: 1848273d9f13SBarry Smith + A - the matrix 1849273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1850273d9f13SBarry Smith 1851273d9f13SBarry Smith Notes: 1852273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1853273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1854284134d9SBarry Smith need not call this routine. 1855273d9f13SBarry Smith 1856273d9f13SBarry Smith Level: intermediate 1857273d9f13SBarry Smith 1858273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1859273d9f13SBarry Smith 1860273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1861273d9f13SBarry Smith @*/ 1862be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 1863273d9f13SBarry Smith { 1864dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar[]); 1865a23d5eceSKris Buschelman 1866a23d5eceSKris Buschelman PetscFunctionBegin; 1867a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1868a23d5eceSKris Buschelman if (f) { 1869a23d5eceSKris Buschelman ierr = (*f)(B,data);CHKERRQ(ierr); 1870a23d5eceSKris Buschelman } 1871a23d5eceSKris Buschelman PetscFunctionReturn(0); 1872a23d5eceSKris Buschelman } 1873a23d5eceSKris Buschelman 1874a23d5eceSKris Buschelman EXTERN_C_BEGIN 1875a23d5eceSKris Buschelman #undef __FUNCT__ 1876afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense" 1877be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 1878a23d5eceSKris Buschelman { 1879273d9f13SBarry Smith Mat_SeqDense *b; 1880dfbe8321SBarry Smith PetscErrorCode ierr; 1881273d9f13SBarry Smith 1882273d9f13SBarry Smith PetscFunctionBegin; 1883273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1884273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 1885273d9f13SBarry Smith if (!data) { 1886284134d9SBarry Smith ierr = PetscMalloc((b->lda*b->Nmax+1)*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 1887284134d9SBarry Smith ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 1888273d9f13SBarry Smith b->user_alloc = PETSC_FALSE; 1889284134d9SBarry Smith ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 1890273d9f13SBarry Smith } else { /* user-allocated storage */ 1891273d9f13SBarry Smith b->v = data; 1892273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 1893273d9f13SBarry Smith } 1894273d9f13SBarry Smith PetscFunctionReturn(0); 1895273d9f13SBarry Smith } 1896a23d5eceSKris Buschelman EXTERN_C_END 1897273d9f13SBarry Smith 18981b807ce4Svictorle #undef __FUNCT__ 18991b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 19001b807ce4Svictorle /*@C 19011b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 19021b807ce4Svictorle 19031b807ce4Svictorle Input parameter: 19041b807ce4Svictorle + A - the matrix 19051b807ce4Svictorle - lda - the leading dimension 19061b807ce4Svictorle 19071b807ce4Svictorle Notes: 19081b807ce4Svictorle This routine is to be used in conjunction with MatSeqDenseSetPreallocation; 19091b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 19101b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 19111b807ce4Svictorle 19121b807ce4Svictorle Level: intermediate 19131b807ce4Svictorle 19141b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 19151b807ce4Svictorle 1916284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 19171b807ce4Svictorle @*/ 1918be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda) 19191b807ce4Svictorle { 19201b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 192121a2c019SBarry Smith 19221b807ce4Svictorle PetscFunctionBegin; 1923899cda47SBarry Smith if (lda < B->rmap.n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap.n); 19241b807ce4Svictorle b->lda = lda; 192521a2c019SBarry Smith b->changelda = PETSC_FALSE; 192621a2c019SBarry Smith b->Mmax = PetscMax(b->Mmax,lda); 19271b807ce4Svictorle PetscFunctionReturn(0); 19281b807ce4Svictorle } 19291b807ce4Svictorle 19300bad9183SKris Buschelman /*MC 1931fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 19320bad9183SKris Buschelman 19330bad9183SKris Buschelman Options Database Keys: 19340bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 19350bad9183SKris Buschelman 19360bad9183SKris Buschelman Level: beginner 19370bad9183SKris Buschelman 193889665df3SBarry Smith .seealso: MatCreateSeqDense() 193989665df3SBarry Smith 19400bad9183SKris Buschelman M*/ 19410bad9183SKris Buschelman 1942273d9f13SBarry Smith EXTERN_C_BEGIN 19434a2ae208SSatish Balay #undef __FUNCT__ 19444a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 1945be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B) 1946273d9f13SBarry Smith { 1947273d9f13SBarry Smith Mat_SeqDense *b; 1948dfbe8321SBarry Smith PetscErrorCode ierr; 19497c334f02SBarry Smith PetscMPIInt size; 1950273d9f13SBarry Smith 1951273d9f13SBarry Smith PetscFunctionBegin; 1952273d9f13SBarry Smith ierr = MPI_Comm_size(B->comm,&size);CHKERRQ(ierr); 195329bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 195455659b69SBarry Smith 1955899cda47SBarry Smith B->rmap.bs = B->cmap.bs = 1; 19566148ca0dSBarry Smith ierr = PetscMapSetUp(&B->rmap);CHKERRQ(ierr); 19576148ca0dSBarry Smith ierr = PetscMapSetUp(&B->cmap);CHKERRQ(ierr); 1958273d9f13SBarry Smith 195938f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr); 1960549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 196144cd7ae7SLois Curfman McInnes B->factor = 0; 196290f02eecSBarry Smith B->mapping = 0; 196344cd7ae7SLois Curfman McInnes B->data = (void*)b; 196418f449edSLois Curfman McInnes 1965a5ae1ecdSBarry Smith 196644cd7ae7SLois Curfman McInnes b->pivots = 0; 1967273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 1968273d9f13SBarry Smith b->v = 0; 1969899cda47SBarry Smith b->lda = B->rmap.n; 197021a2c019SBarry Smith b->changelda = PETSC_FALSE; 1971899cda47SBarry Smith b->Mmax = B->rmap.n; 1972899cda47SBarry Smith b->Nmax = B->cmap.n; 19734e220ebcSLois Curfman McInnes 1974a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 1975a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 1976a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 19774ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C", 19784ae313f4SHong Zhang "MatMatMult_SeqAIJ_SeqDense", 19794ae313f4SHong Zhang MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 19804ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C", 19814ae313f4SHong Zhang "MatMatMultSymbolic_SeqAIJ_SeqDense", 19824ae313f4SHong Zhang MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 19834ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C", 19844ae313f4SHong Zhang "MatMatMultNumeric_SeqAIJ_SeqDense", 19854ae313f4SHong Zhang MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 198617667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 19873a40ed3dSBarry Smith PetscFunctionReturn(0); 1988289bc588SBarry Smith } 1989c0aa2d19SHong Zhang 1990c0aa2d19SHong Zhang 1991273d9f13SBarry Smith EXTERN_C_END 1992