1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 367e560aaSBarry Smith /* 467e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 567e560aaSBarry Smith */ 6289bc588SBarry Smith 77c4f633dSBarry Smith #include "../src/mat/impls/dense/seq/dense.h" 8f3da1532SBarry Smith #include "petscblaslapack.h" 9289bc588SBarry Smith 104a2ae208SSatish Balay #undef __FUNCT__ 114a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense" 12f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 131987afe7SBarry Smith { 141987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 15f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 1613f74950SBarry Smith PetscInt j; 170805154bSBarry Smith PetscBLASInt N,m,ldax,lday,one = 1; 18efee365bSSatish Balay PetscErrorCode ierr; 193a40ed3dSBarry Smith 203a40ed3dSBarry Smith PetscFunctionBegin; 21d0f46423SBarry Smith N = PetscBLASIntCast(X->rmap->n*X->cmap->n); 22d0f46423SBarry Smith m = PetscBLASIntCast(X->rmap->n); 230805154bSBarry Smith ldax = PetscBLASIntCast(x->lda); 240805154bSBarry Smith lday = PetscBLASIntCast(y->lda); 25a5ce6ee0Svictorle if (ldax>m || lday>m) { 26d0f46423SBarry Smith for (j=0; j<X->cmap->n; j++) { 27f4df32b1SMatthew Knepley BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one); 28a5ce6ee0Svictorle } 29a5ce6ee0Svictorle } else { 30f4df32b1SMatthew Knepley BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one); 31a5ce6ee0Svictorle } 320450473dSBarry Smith ierr = PetscLogFlops(PetscMax(2*N-1,0));CHKERRQ(ierr); 333a40ed3dSBarry Smith PetscFunctionReturn(0); 341987afe7SBarry Smith } 351987afe7SBarry Smith 364a2ae208SSatish Balay #undef __FUNCT__ 374a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 38dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 39289bc588SBarry Smith { 40d0f46423SBarry Smith PetscInt N = A->rmap->n*A->cmap->n; 413a40ed3dSBarry Smith 423a40ed3dSBarry Smith PetscFunctionBegin; 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; 497adad957SLisandro Dalcin info->memory = ((PetscObject)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; 62efee365bSSatish Balay PetscErrorCode ierr; 630805154bSBarry Smith PetscBLASInt one = 1,j,nz,lda = PetscBLASIntCast(a->lda); 6480cd9d93SLois Curfman McInnes 653a40ed3dSBarry Smith PetscFunctionBegin; 66d0f46423SBarry Smith if (lda>A->rmap->n) { 67d0f46423SBarry Smith nz = PetscBLASIntCast(A->rmap->n); 68d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 69f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v+j*lda,&one); 70a5ce6ee0Svictorle } 71a5ce6ee0Svictorle } else { 72d0f46423SBarry Smith nz = PetscBLASIntCast(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; 84d0f46423SBarry 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; 89d0f46423SBarry 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 101b24902e0SBarry Smith #undef __FUNCT__ 102b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense" 103719d5645SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues) 104b24902e0SBarry Smith { 105b24902e0SBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 106b24902e0SBarry Smith PetscErrorCode ierr; 107b24902e0SBarry Smith PetscInt lda = (PetscInt)mat->lda,j,m; 108b24902e0SBarry Smith 109b24902e0SBarry Smith PetscFunctionBegin; 110719d5645SBarry Smith ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr); 111b24902e0SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 112b24902e0SBarry Smith l = (Mat_SeqDense*)newi->data; 113d0f46423SBarry Smith if (lda>A->rmap->n) { 114d0f46423SBarry Smith m = A->rmap->n; 115d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 116b24902e0SBarry Smith ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 117b24902e0SBarry Smith } 118b24902e0SBarry Smith } else { 119d0f46423SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 120b24902e0SBarry Smith } 121b24902e0SBarry Smith } 122b24902e0SBarry Smith newi->assembled = PETSC_TRUE; 123b24902e0SBarry Smith PetscFunctionReturn(0); 124b24902e0SBarry Smith } 125b24902e0SBarry Smith 1264a2ae208SSatish Balay #undef __FUNCT__ 1274a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 128dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 12902cad45dSBarry Smith { 1306849ba73SBarry Smith PetscErrorCode ierr; 13102cad45dSBarry Smith 1323a40ed3dSBarry Smith PetscFunctionBegin; 1335c9eb25fSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr); 134d0f46423SBarry Smith ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 1355c9eb25fSBarry Smith ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 136719d5645SBarry Smith ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr); 137b24902e0SBarry Smith PetscFunctionReturn(0); 138b24902e0SBarry Smith } 139b24902e0SBarry Smith 1406ee01492SSatish Balay 1410481f469SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*); 142719d5645SBarry Smith 1434a2ae208SSatish Balay #undef __FUNCT__ 1444a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 1450481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 146289bc588SBarry Smith { 1474482741eSBarry Smith MatFactorInfo info; 148a093e273SMatthew Knepley PetscErrorCode ierr; 1493a40ed3dSBarry Smith 1503a40ed3dSBarry Smith PetscFunctionBegin; 151c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 152719d5645SBarry Smith ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr); 1533a40ed3dSBarry Smith PetscFunctionReturn(0); 154289bc588SBarry Smith } 1556ee01492SSatish Balay 1560b4b3355SBarry Smith #undef __FUNCT__ 1574a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 158dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 159289bc588SBarry Smith { 160c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1616849ba73SBarry Smith PetscErrorCode ierr; 16287828ca2SBarry Smith PetscScalar *x,*y; 163d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 16467e560aaSBarry Smith 1653a40ed3dSBarry Smith PetscFunctionBegin; 1661ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1671ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 168d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 169d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_LU) { 170ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 171e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 172ae7cfcebSSatish Balay #else 17371044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 174e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 175ae7cfcebSSatish Balay #endif 176d5f3da31SBarry Smith } else if (A->factortype == MAT_FACTOR_CHOLESKY){ 177ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 178e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 179ae7cfcebSSatish Balay #else 18071044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 181e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 182ae7cfcebSSatish Balay #endif 183289bc588SBarry Smith } 184e32f2f54SBarry Smith else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 1851ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1861ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 187dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 1883a40ed3dSBarry Smith PetscFunctionReturn(0); 189289bc588SBarry Smith } 1906ee01492SSatish Balay 1914a2ae208SSatish Balay #undef __FUNCT__ 1924a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 193dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 194da3a660dSBarry Smith { 195c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 196dfbe8321SBarry Smith PetscErrorCode ierr; 19787828ca2SBarry Smith PetscScalar *x,*y; 198d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 19967e560aaSBarry Smith 2003a40ed3dSBarry Smith PetscFunctionBegin; 2011ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2021ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 203d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 204752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 205da3a660dSBarry Smith if (mat->pivots) { 206ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 207e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 208ae7cfcebSSatish Balay #else 20971044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 210e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 211ae7cfcebSSatish Balay #endif 2127a97a34bSBarry Smith } else { 213ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 214e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 215ae7cfcebSSatish Balay #else 21671044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 217e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 218ae7cfcebSSatish Balay #endif 219da3a660dSBarry Smith } 2201ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2211ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 222dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 2233a40ed3dSBarry Smith PetscFunctionReturn(0); 224da3a660dSBarry Smith } 2256ee01492SSatish Balay 2264a2ae208SSatish Balay #undef __FUNCT__ 2274a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 228dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 229da3a660dSBarry Smith { 230c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 231dfbe8321SBarry Smith PetscErrorCode ierr; 23287828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 233da3a660dSBarry Smith Vec tmp = 0; 234d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 23567e560aaSBarry Smith 2363a40ed3dSBarry Smith PetscFunctionBegin; 2371ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2381ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 239d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 240da3a660dSBarry Smith if (yy == zz) { 24178b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 24252e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 24378b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 244da3a660dSBarry Smith } 245d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 246752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 247da3a660dSBarry Smith if (mat->pivots) { 248ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 249e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 250ae7cfcebSSatish Balay #else 25171044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 252e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 253ae7cfcebSSatish Balay #endif 254a8c6a408SBarry Smith } else { 255ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 256e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 257ae7cfcebSSatish Balay #else 25871044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 259e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 260ae7cfcebSSatish Balay #endif 261da3a660dSBarry Smith } 2622dcb1b2aSMatthew Knepley if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 2632dcb1b2aSMatthew Knepley else {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);} 2641ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2651ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 266dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 2673a40ed3dSBarry Smith PetscFunctionReturn(0); 268da3a660dSBarry Smith } 26967e560aaSBarry Smith 2704a2ae208SSatish Balay #undef __FUNCT__ 2714a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 272dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 273da3a660dSBarry Smith { 274c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2756849ba73SBarry Smith PetscErrorCode ierr; 27687828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 277da3a660dSBarry Smith Vec tmp; 278d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 27967e560aaSBarry Smith 2803a40ed3dSBarry Smith PetscFunctionBegin; 281d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 2821ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2831ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 284da3a660dSBarry Smith if (yy == zz) { 28578b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 28652e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 28778b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 288da3a660dSBarry Smith } 289d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 290752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 291da3a660dSBarry Smith if (mat->pivots) { 292ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 293e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 294ae7cfcebSSatish Balay #else 29571044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 296e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 297ae7cfcebSSatish Balay #endif 2983a40ed3dSBarry Smith } else { 299ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 300e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 301ae7cfcebSSatish Balay #else 30271044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 303e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 304ae7cfcebSSatish Balay #endif 305da3a660dSBarry Smith } 30690f02eecSBarry Smith if (tmp) { 3072dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 30890f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 3093a40ed3dSBarry Smith } else { 3102dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 31190f02eecSBarry Smith } 3121ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3131ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 314dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 3153a40ed3dSBarry Smith PetscFunctionReturn(0); 316da3a660dSBarry Smith } 317db4efbfdSBarry Smith 318db4efbfdSBarry Smith /* ---------------------------------------------------------------*/ 319db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 320db4efbfdSBarry Smith rather than put it in the Mat->row slot.*/ 321db4efbfdSBarry Smith #undef __FUNCT__ 322db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense" 3230481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo) 324db4efbfdSBarry Smith { 325db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 326db4efbfdSBarry Smith PetscFunctionBegin; 327e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 328db4efbfdSBarry Smith #else 329db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 330db4efbfdSBarry Smith PetscErrorCode ierr; 331db4efbfdSBarry Smith PetscBLASInt n,m,info; 332db4efbfdSBarry Smith 333db4efbfdSBarry Smith PetscFunctionBegin; 334db4efbfdSBarry Smith n = PetscBLASIntCast(A->cmap->n); 335db4efbfdSBarry Smith m = PetscBLASIntCast(A->rmap->n); 336db4efbfdSBarry Smith if (!mat->pivots) { 337db4efbfdSBarry Smith ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr); 338db4efbfdSBarry Smith ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr); 339db4efbfdSBarry Smith } 340db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 341db4efbfdSBarry Smith LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info); 342e32f2f54SBarry Smith if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization"); 343e32f2f54SBarry Smith if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 344db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 345db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 346db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 347db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 348d5f3da31SBarry Smith A->factortype = MAT_FACTOR_LU; 349db4efbfdSBarry Smith 350dc0b31edSSatish Balay ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr); 351db4efbfdSBarry Smith #endif 352db4efbfdSBarry Smith PetscFunctionReturn(0); 353db4efbfdSBarry Smith } 354db4efbfdSBarry Smith 355db4efbfdSBarry Smith #undef __FUNCT__ 356db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense" 3570481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo) 358db4efbfdSBarry Smith { 359db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 360db4efbfdSBarry Smith PetscFunctionBegin; 361e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 362db4efbfdSBarry Smith #else 363db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 364db4efbfdSBarry Smith PetscErrorCode ierr; 365db4efbfdSBarry Smith PetscBLASInt info,n = PetscBLASIntCast(A->cmap->n); 366db4efbfdSBarry Smith 367db4efbfdSBarry Smith PetscFunctionBegin; 368db4efbfdSBarry Smith ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 369db4efbfdSBarry Smith mat->pivots = 0; 370db4efbfdSBarry Smith 371db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 372db4efbfdSBarry Smith LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info); 373e32f2f54SBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 374db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 375db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 376db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 377db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 378d5f3da31SBarry Smith A->factortype = MAT_FACTOR_CHOLESKY; 379dc0b31edSSatish Balay ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr); 380db4efbfdSBarry Smith #endif 381db4efbfdSBarry Smith PetscFunctionReturn(0); 382db4efbfdSBarry Smith } 383db4efbfdSBarry Smith 384db4efbfdSBarry Smith 385db4efbfdSBarry Smith #undef __FUNCT__ 386db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 3870481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 388db4efbfdSBarry Smith { 389db4efbfdSBarry Smith PetscErrorCode ierr; 390db4efbfdSBarry Smith MatFactorInfo info; 391db4efbfdSBarry Smith 392db4efbfdSBarry Smith PetscFunctionBegin; 393db4efbfdSBarry Smith info.fill = 1.0; 394c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 395719d5645SBarry Smith ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr); 396db4efbfdSBarry Smith PetscFunctionReturn(0); 397db4efbfdSBarry Smith } 398db4efbfdSBarry Smith 399db4efbfdSBarry Smith #undef __FUNCT__ 400db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 4010481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info) 402db4efbfdSBarry Smith { 403db4efbfdSBarry Smith PetscFunctionBegin; 404c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 405719d5645SBarry Smith fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense; 406db4efbfdSBarry Smith PetscFunctionReturn(0); 407db4efbfdSBarry Smith } 408db4efbfdSBarry Smith 409db4efbfdSBarry Smith #undef __FUNCT__ 410db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 4110481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 412db4efbfdSBarry Smith { 413db4efbfdSBarry Smith PetscFunctionBegin; 414c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 415719d5645SBarry Smith fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense; 416db4efbfdSBarry Smith PetscFunctionReturn(0); 417db4efbfdSBarry Smith } 418db4efbfdSBarry Smith 419bb5747d9SMatthew Knepley EXTERN_C_BEGIN 420db4efbfdSBarry Smith #undef __FUNCT__ 421db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc" 422db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact) 423db4efbfdSBarry Smith { 424db4efbfdSBarry Smith PetscErrorCode ierr; 425db4efbfdSBarry Smith 426db4efbfdSBarry Smith PetscFunctionBegin; 427db4efbfdSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr); 428db4efbfdSBarry Smith ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 429db4efbfdSBarry Smith ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr); 430db4efbfdSBarry Smith if (ftype == MAT_FACTOR_LU){ 431db4efbfdSBarry Smith (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense; 432db4efbfdSBarry Smith } else { 433db4efbfdSBarry Smith (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense; 434db4efbfdSBarry Smith } 435d5f3da31SBarry Smith (*fact)->factortype = ftype; 436db4efbfdSBarry Smith PetscFunctionReturn(0); 437db4efbfdSBarry Smith } 438bb5747d9SMatthew Knepley EXTERN_C_END 439db4efbfdSBarry Smith 440289bc588SBarry Smith /* ------------------------------------------------------------------*/ 4414a2ae208SSatish Balay #undef __FUNCT__ 44241f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense" 44341f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 444289bc588SBarry Smith { 445c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 44687828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 447dfbe8321SBarry Smith PetscErrorCode ierr; 448d0f46423SBarry Smith PetscInt m = A->rmap->n,i; 449aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 4500805154bSBarry Smith PetscBLASInt o = 1,bm = PetscBLASIntCast(m); 451bc1b551cSSatish Balay #endif 452289bc588SBarry Smith 4533a40ed3dSBarry Smith PetscFunctionBegin; 454289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 45571044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 4562dcb1b2aSMatthew Knepley ierr = VecSet(xx,zero);CHKERRQ(ierr); 457289bc588SBarry Smith } 4581ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4591ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 460b965ef7fSBarry Smith its = its*lits; 461e32f2f54SBarry Smith if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 462289bc588SBarry Smith while (its--) { 463fccaa45eSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 464289bc588SBarry Smith for (i=0; i<m; i++) { 465aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 466f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 467f1747703SBarry Smith not happy about returning a double complex */ 46813f74950SBarry Smith PetscInt _i; 46987828ca2SBarry Smith PetscScalar sum = b[i]; 470f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4713f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 472f1747703SBarry Smith } 473f1747703SBarry Smith xt = sum; 474f1747703SBarry Smith #else 47571044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 476f1747703SBarry Smith #endif 47755a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 478289bc588SBarry Smith } 479289bc588SBarry Smith } 480fccaa45eSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 481289bc588SBarry Smith for (i=m-1; i>=0; i--) { 482aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 483f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 484f1747703SBarry Smith not happy about returning a double complex */ 48513f74950SBarry Smith PetscInt _i; 48687828ca2SBarry Smith PetscScalar sum = b[i]; 487f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4883f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 489f1747703SBarry Smith } 490f1747703SBarry Smith xt = sum; 491f1747703SBarry Smith #else 49271044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 493f1747703SBarry Smith #endif 49455a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 495289bc588SBarry Smith } 496289bc588SBarry Smith } 497289bc588SBarry Smith } 4981ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 4991ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5003a40ed3dSBarry Smith PetscFunctionReturn(0); 501289bc588SBarry Smith } 502289bc588SBarry Smith 503289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5044a2ae208SSatish Balay #undef __FUNCT__ 5054a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 506dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 507289bc588SBarry Smith { 508c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 50987828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 510dfbe8321SBarry Smith PetscErrorCode ierr; 5110805154bSBarry Smith PetscBLASInt m, n,_One=1; 512ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 5133a40ed3dSBarry Smith 5143a40ed3dSBarry Smith PetscFunctionBegin; 515d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 516d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 517d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5181ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5191ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 52071044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One); 5211ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5221ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 523dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 5243a40ed3dSBarry Smith PetscFunctionReturn(0); 525289bc588SBarry Smith } 526800995b7SMatthew Knepley 5274a2ae208SSatish Balay #undef __FUNCT__ 5284a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 529dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 530289bc588SBarry Smith { 531c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 53287828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 533dfbe8321SBarry Smith PetscErrorCode ierr; 5340805154bSBarry Smith PetscBLASInt m, n, _One=1; 5353a40ed3dSBarry Smith 5363a40ed3dSBarry Smith PetscFunctionBegin; 537d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 538d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 539d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5401ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5411ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 54271044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 5431ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5441ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 545dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr); 5463a40ed3dSBarry Smith PetscFunctionReturn(0); 547289bc588SBarry Smith } 5486ee01492SSatish Balay 5494a2ae208SSatish Balay #undef __FUNCT__ 5504a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 551dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 552289bc588SBarry Smith { 553c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 55487828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 555dfbe8321SBarry Smith PetscErrorCode ierr; 5560805154bSBarry Smith PetscBLASInt m, n, _One=1; 5573a40ed3dSBarry Smith 5583a40ed3dSBarry Smith PetscFunctionBegin; 559d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 560d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 561d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 562600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5631ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5641ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 56571044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5661ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5671ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 568dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 5693a40ed3dSBarry Smith PetscFunctionReturn(0); 570289bc588SBarry Smith } 5716ee01492SSatish Balay 5724a2ae208SSatish Balay #undef __FUNCT__ 5734a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 574dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 575289bc588SBarry Smith { 576c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 57787828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 578dfbe8321SBarry Smith PetscErrorCode ierr; 5790805154bSBarry Smith PetscBLASInt m, n, _One=1; 58087828ca2SBarry Smith PetscScalar _DOne=1.0; 5813a40ed3dSBarry Smith 5823a40ed3dSBarry Smith PetscFunctionBegin; 583d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 584d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 585d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 586600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5871ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5881ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 58971044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5901ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5911ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 592dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 5933a40ed3dSBarry Smith PetscFunctionReturn(0); 594289bc588SBarry Smith } 595289bc588SBarry Smith 596289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5974a2ae208SSatish Balay #undef __FUNCT__ 5984a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 59913f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 600289bc588SBarry Smith { 601c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 60287828ca2SBarry Smith PetscScalar *v; 6036849ba73SBarry Smith PetscErrorCode ierr; 60413f74950SBarry Smith PetscInt i; 60567e560aaSBarry Smith 6063a40ed3dSBarry Smith PetscFunctionBegin; 607d0f46423SBarry Smith *ncols = A->cmap->n; 608289bc588SBarry Smith if (cols) { 609d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr); 610d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) (*cols)[i] = i; 611289bc588SBarry Smith } 612289bc588SBarry Smith if (vals) { 613d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 614289bc588SBarry Smith v = mat->v + row; 615d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;} 616289bc588SBarry Smith } 6173a40ed3dSBarry Smith PetscFunctionReturn(0); 618289bc588SBarry Smith } 6196ee01492SSatish Balay 6204a2ae208SSatish Balay #undef __FUNCT__ 6214a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 62213f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 623289bc588SBarry Smith { 624dfbe8321SBarry Smith PetscErrorCode ierr; 625606d414cSSatish Balay PetscFunctionBegin; 626606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 627606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 6283a40ed3dSBarry Smith PetscFunctionReturn(0); 629289bc588SBarry Smith } 630289bc588SBarry Smith /* ----------------------------------------------------------------*/ 6314a2ae208SSatish Balay #undef __FUNCT__ 6324a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 63313f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 634289bc588SBarry Smith { 635c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 63613f74950SBarry Smith PetscInt i,j,idx=0; 637d6dfbf8fSBarry Smith 6383a40ed3dSBarry Smith PetscFunctionBegin; 63971fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 640289bc588SBarry Smith if (!mat->roworiented) { 641dbb450caSBarry Smith if (addv == INSERT_VALUES) { 642289bc588SBarry Smith for (j=0; j<n; j++) { 643cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6442515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 645e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 64658804f6dSBarry Smith #endif 647289bc588SBarry Smith for (i=0; i<m; i++) { 648cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6492515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 650e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 65158804f6dSBarry Smith #endif 652cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 653289bc588SBarry Smith } 654289bc588SBarry Smith } 6553a40ed3dSBarry Smith } else { 656289bc588SBarry Smith for (j=0; j<n; j++) { 657cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6582515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 659e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 66058804f6dSBarry Smith #endif 661289bc588SBarry Smith for (i=0; i<m; i++) { 662cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6632515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 664e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 66558804f6dSBarry Smith #endif 666cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 667289bc588SBarry Smith } 668289bc588SBarry Smith } 669289bc588SBarry Smith } 6703a40ed3dSBarry Smith } else { 671dbb450caSBarry Smith if (addv == INSERT_VALUES) { 672e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 673cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 6742515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 675e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 67658804f6dSBarry Smith #endif 677e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 678cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 6792515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 680e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 68158804f6dSBarry Smith #endif 682cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 683e8d4e0b9SBarry Smith } 684e8d4e0b9SBarry Smith } 6853a40ed3dSBarry Smith } else { 686289bc588SBarry Smith for (i=0; i<m; i++) { 687cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 6882515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 689e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 69058804f6dSBarry Smith #endif 691289bc588SBarry Smith for (j=0; j<n; j++) { 692cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 6932515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 694e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 69558804f6dSBarry Smith #endif 696cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 697289bc588SBarry Smith } 698289bc588SBarry Smith } 699289bc588SBarry Smith } 700e8d4e0b9SBarry Smith } 7013a40ed3dSBarry Smith PetscFunctionReturn(0); 702289bc588SBarry Smith } 703e8d4e0b9SBarry Smith 7044a2ae208SSatish Balay #undef __FUNCT__ 7054a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 70613f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 707ae80bb75SLois Curfman McInnes { 708ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 70913f74950SBarry Smith PetscInt i,j; 710ae80bb75SLois Curfman McInnes 7113a40ed3dSBarry Smith PetscFunctionBegin; 712ae80bb75SLois Curfman McInnes /* row-oriented output */ 713ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 71497e567efSBarry Smith if (indexm[i] < 0) {v += n;continue;} 715e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap->n); 716ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 7176f31f424SBarry Smith if (indexn[j] < 0) {v++; continue;} 718e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap->n); 71997e567efSBarry Smith *v++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 720ae80bb75SLois Curfman McInnes } 721ae80bb75SLois Curfman McInnes } 7223a40ed3dSBarry Smith PetscFunctionReturn(0); 723ae80bb75SLois Curfman McInnes } 724ae80bb75SLois Curfman McInnes 725289bc588SBarry Smith /* -----------------------------------------------------------------*/ 726289bc588SBarry Smith 7274a2ae208SSatish Balay #undef __FUNCT__ 728*5bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_SeqDense" 729*5bba2384SShri Abhyankar PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, Mat newmat) 730aabbc4fbSShri Abhyankar { 731aabbc4fbSShri Abhyankar Mat_SeqDense *a; 732aabbc4fbSShri Abhyankar PetscErrorCode ierr; 733aabbc4fbSShri Abhyankar PetscInt *scols,i,j,nz,header[4]; 734aabbc4fbSShri Abhyankar int fd; 735aabbc4fbSShri Abhyankar PetscMPIInt size; 736aabbc4fbSShri Abhyankar PetscInt *rowlengths = 0,M,N,*cols,grows,gcols; 737aabbc4fbSShri Abhyankar PetscScalar *vals,*svals,*v,*w; 738aabbc4fbSShri Abhyankar MPI_Comm comm = ((PetscObject)viewer)->comm; 739aabbc4fbSShri Abhyankar 740aabbc4fbSShri Abhyankar PetscFunctionBegin; 741aabbc4fbSShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 742aabbc4fbSShri Abhyankar if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor"); 743aabbc4fbSShri Abhyankar ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 744aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 745aabbc4fbSShri Abhyankar if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 746aabbc4fbSShri Abhyankar M = header[1]; N = header[2]; nz = header[3]; 747aabbc4fbSShri Abhyankar 748aabbc4fbSShri Abhyankar /* set global size if not set already*/ 749aabbc4fbSShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 750aabbc4fbSShri Abhyankar ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr); 751aabbc4fbSShri Abhyankar } else { 752aabbc4fbSShri Abhyankar /* if sizes and type are already set, check if the vector global sizes are correct */ 753aabbc4fbSShri Abhyankar ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr); 754aabbc4fbSShri Abhyankar if (M != grows || N != gcols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,grows,gcols); 755aabbc4fbSShri Abhyankar } 756aabbc4fbSShri Abhyankar ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 757aabbc4fbSShri Abhyankar 758aabbc4fbSShri Abhyankar if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 759aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 760aabbc4fbSShri Abhyankar v = a->v; 761aabbc4fbSShri Abhyankar /* Allocate some temp space to read in the values and then flip them 762aabbc4fbSShri Abhyankar from row major to column major */ 763aabbc4fbSShri Abhyankar ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 764aabbc4fbSShri Abhyankar /* read in nonzero values */ 765aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 766aabbc4fbSShri Abhyankar /* now flip the values and store them in the matrix*/ 767aabbc4fbSShri Abhyankar for (j=0; j<N; j++) { 768aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 769aabbc4fbSShri Abhyankar *v++ =w[i*N+j]; 770aabbc4fbSShri Abhyankar } 771aabbc4fbSShri Abhyankar } 772aabbc4fbSShri Abhyankar ierr = PetscFree(w);CHKERRQ(ierr); 773aabbc4fbSShri Abhyankar } else { 774aabbc4fbSShri Abhyankar /* read row lengths */ 775aabbc4fbSShri Abhyankar ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 776aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 777aabbc4fbSShri Abhyankar 778aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 779aabbc4fbSShri Abhyankar v = a->v; 780aabbc4fbSShri Abhyankar 781aabbc4fbSShri Abhyankar /* read column indices and nonzeros */ 782aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr); 783aabbc4fbSShri Abhyankar cols = scols; 784aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 785aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 786aabbc4fbSShri Abhyankar vals = svals; 787aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 788aabbc4fbSShri Abhyankar 789aabbc4fbSShri Abhyankar /* insert into matrix */ 790aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 791aabbc4fbSShri Abhyankar for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 792aabbc4fbSShri Abhyankar svals += rowlengths[i]; scols += rowlengths[i]; 793aabbc4fbSShri Abhyankar } 794aabbc4fbSShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 795aabbc4fbSShri Abhyankar ierr = PetscFree(cols);CHKERRQ(ierr); 796aabbc4fbSShri Abhyankar ierr = PetscFree(rowlengths);CHKERRQ(ierr); 797aabbc4fbSShri Abhyankar } 798aabbc4fbSShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 799aabbc4fbSShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 800aabbc4fbSShri Abhyankar 801aabbc4fbSShri Abhyankar PetscFunctionReturn(0); 802aabbc4fbSShri Abhyankar } 803aabbc4fbSShri Abhyankar 804aabbc4fbSShri Abhyankar #undef __FUNCT__ 8054a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 8066849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 807289bc588SBarry Smith { 808932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 809dfbe8321SBarry Smith PetscErrorCode ierr; 81013f74950SBarry Smith PetscInt i,j; 8112dcb1b2aSMatthew Knepley const char *name; 81287828ca2SBarry Smith PetscScalar *v; 813f3ef73ceSBarry Smith PetscViewerFormat format; 8145f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX) 8155f481a85SSatish Balay PetscTruth allreal = PETSC_TRUE; 8165f481a85SSatish Balay #endif 817932b0c3eSLois Curfman McInnes 8183a40ed3dSBarry Smith PetscFunctionBegin; 819435da068SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 820b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 821456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 8223a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 823fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 824b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 825d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 82644cd7ae7SLois Curfman McInnes v = a->v + i; 82777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 828d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 829aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 830329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 831a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 832329f5518SBarry Smith } else if (PetscRealPart(*v)) { 833a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr); 8346831982aSBarry Smith } 83580cd9d93SLois Curfman McInnes #else 8366831982aSBarry Smith if (*v) { 837a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr); 8386831982aSBarry Smith } 83980cd9d93SLois Curfman McInnes #endif 8401b807ce4Svictorle v += a->lda; 84180cd9d93SLois Curfman McInnes } 842b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 84380cd9d93SLois Curfman McInnes } 844b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 8453a40ed3dSBarry Smith } else { 846b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 847aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 84847989497SBarry Smith /* determine if matrix has all real values */ 84947989497SBarry Smith v = a->v; 850d0f46423SBarry Smith for (i=0; i<A->rmap->n*A->cmap->n; i++) { 851ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 85247989497SBarry Smith } 85347989497SBarry Smith #endif 854fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 8553a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 856d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr); 857d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 858fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 859ffac6cdbSBarry Smith } 860ffac6cdbSBarry Smith 861d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 862932b0c3eSLois Curfman McInnes v = a->v + i; 863d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 864aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 86547989497SBarry Smith if (allreal) { 866f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr); 86747989497SBarry Smith } else { 868f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 86947989497SBarry Smith } 870289bc588SBarry Smith #else 871f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr); 872289bc588SBarry Smith #endif 8731b807ce4Svictorle v += a->lda; 874289bc588SBarry Smith } 875b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 876289bc588SBarry Smith } 877fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 878b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 879ffac6cdbSBarry Smith } 880b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 881da3a660dSBarry Smith } 882b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 8833a40ed3dSBarry Smith PetscFunctionReturn(0); 884289bc588SBarry Smith } 885289bc588SBarry Smith 8864a2ae208SSatish Balay #undef __FUNCT__ 8874a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 8886849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 889932b0c3eSLois Curfman McInnes { 890932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 8916849ba73SBarry Smith PetscErrorCode ierr; 89213f74950SBarry Smith int fd; 893d0f46423SBarry Smith PetscInt ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n; 89487828ca2SBarry Smith PetscScalar *v,*anonz,*vals; 895f3ef73ceSBarry Smith PetscViewerFormat format; 896932b0c3eSLois Curfman McInnes 8973a40ed3dSBarry Smith PetscFunctionBegin; 898b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 89990ace30eSBarry Smith 900b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 9017973ad04SBarry Smith if (format == PETSC_VIEWER_NATIVE) { 90290ace30eSBarry Smith /* store the matrix as a dense matrix */ 90313f74950SBarry Smith ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 9040700a824SBarry Smith col_lens[0] = MAT_FILE_CLASSID; 90590ace30eSBarry Smith col_lens[1] = m; 90690ace30eSBarry Smith col_lens[2] = n; 90790ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 9086f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 909606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 91090ace30eSBarry Smith 91190ace30eSBarry Smith /* write out matrix, by rows */ 91287828ca2SBarry Smith ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 91390ace30eSBarry Smith v = a->v; 91490ace30eSBarry Smith for (j=0; j<n; j++) { 915578230a0SSatish Balay for (i=0; i<m; i++) { 916578230a0SSatish Balay vals[j + i*n] = *v++; 91790ace30eSBarry Smith } 91890ace30eSBarry Smith } 9196f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 920606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 92190ace30eSBarry Smith } else { 92213f74950SBarry Smith ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 9230700a824SBarry Smith col_lens[0] = MAT_FILE_CLASSID; 924932b0c3eSLois Curfman McInnes col_lens[1] = m; 925932b0c3eSLois Curfman McInnes col_lens[2] = n; 926932b0c3eSLois Curfman McInnes col_lens[3] = nz; 927932b0c3eSLois Curfman McInnes 928932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 929932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 9306f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 931932b0c3eSLois Curfman McInnes 932932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 933932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 934932b0c3eSLois Curfman McInnes ict = 0; 935932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 936932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 937932b0c3eSLois Curfman McInnes } 9386f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 939606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 940932b0c3eSLois Curfman McInnes 941932b0c3eSLois Curfman McInnes /* store nonzero values */ 94287828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 943932b0c3eSLois Curfman McInnes ict = 0; 944932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 945932b0c3eSLois Curfman McInnes v = a->v + i; 946932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 9471b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 948932b0c3eSLois Curfman McInnes } 949932b0c3eSLois Curfman McInnes } 9506f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 951606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 95290ace30eSBarry Smith } 9533a40ed3dSBarry Smith PetscFunctionReturn(0); 954932b0c3eSLois Curfman McInnes } 955932b0c3eSLois Curfman McInnes 9564a2ae208SSatish Balay #undef __FUNCT__ 9574a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 958dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 959f1af5d2fSBarry Smith { 960f1af5d2fSBarry Smith Mat A = (Mat) Aa; 961f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9626849ba73SBarry Smith PetscErrorCode ierr; 963d0f46423SBarry Smith PetscInt m = A->rmap->n,n = A->cmap->n,color,i,j; 96487828ca2SBarry Smith PetscScalar *v = a->v; 965b0a32e0cSBarry Smith PetscViewer viewer; 966b0a32e0cSBarry Smith PetscDraw popup; 967329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 968f3ef73ceSBarry Smith PetscViewerFormat format; 969f1af5d2fSBarry Smith 970f1af5d2fSBarry Smith PetscFunctionBegin; 971f1af5d2fSBarry Smith 972f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 973b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 974b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 975f1af5d2fSBarry Smith 976f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 977fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 978f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 979b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 980f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 981f1af5d2fSBarry Smith x_l = j; 982f1af5d2fSBarry Smith x_r = x_l + 1.0; 983f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 984f1af5d2fSBarry Smith y_l = m - i - 1.0; 985f1af5d2fSBarry Smith y_r = y_l + 1.0; 986f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 987329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 988b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 989329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 990b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 991f1af5d2fSBarry Smith } else { 992f1af5d2fSBarry Smith continue; 993f1af5d2fSBarry Smith } 994f1af5d2fSBarry Smith #else 995f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 996b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 997f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 998b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 999f1af5d2fSBarry Smith } else { 1000f1af5d2fSBarry Smith continue; 1001f1af5d2fSBarry Smith } 1002f1af5d2fSBarry Smith #endif 1003b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1004f1af5d2fSBarry Smith } 1005f1af5d2fSBarry Smith } 1006f1af5d2fSBarry Smith } else { 1007f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 1008f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 1009f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 1010f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 1011f1af5d2fSBarry Smith } 1012b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 1013b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 1014b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 1015f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1016f1af5d2fSBarry Smith x_l = j; 1017f1af5d2fSBarry Smith x_r = x_l + 1.0; 1018f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1019f1af5d2fSBarry Smith y_l = m - i - 1.0; 1020f1af5d2fSBarry Smith y_r = y_l + 1.0; 1021b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 1022b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1023f1af5d2fSBarry Smith } 1024f1af5d2fSBarry Smith } 1025f1af5d2fSBarry Smith } 1026f1af5d2fSBarry Smith PetscFunctionReturn(0); 1027f1af5d2fSBarry Smith } 1028f1af5d2fSBarry Smith 10294a2ae208SSatish Balay #undef __FUNCT__ 10304a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 1031dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 1032f1af5d2fSBarry Smith { 1033b0a32e0cSBarry Smith PetscDraw draw; 1034f1af5d2fSBarry Smith PetscTruth isnull; 1035329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 1036dfbe8321SBarry Smith PetscErrorCode ierr; 1037f1af5d2fSBarry Smith 1038f1af5d2fSBarry Smith PetscFunctionBegin; 1039b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1040b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1041abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 1042f1af5d2fSBarry Smith 1043f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1044d0f46423SBarry Smith xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0; 1045f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 1046b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1047b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 1048f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 1049f1af5d2fSBarry Smith PetscFunctionReturn(0); 1050f1af5d2fSBarry Smith } 1051f1af5d2fSBarry Smith 10524a2ae208SSatish Balay #undef __FUNCT__ 10534a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 1054dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1055932b0c3eSLois Curfman McInnes { 1056dfbe8321SBarry Smith PetscErrorCode ierr; 10576805f65bSBarry Smith PetscTruth iascii,isbinary,isdraw; 1058932b0c3eSLois Curfman McInnes 10593a40ed3dSBarry Smith PetscFunctionBegin; 10602692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 10612692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 10622692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 10630f5bd95cSBarry Smith 1064c45a1595SBarry Smith if (iascii) { 1065c45a1595SBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 10660f5bd95cSBarry Smith } else if (isbinary) { 10673a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1068f1af5d2fSBarry Smith } else if (isdraw) { 1069f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 10705cd90555SBarry Smith } else { 1071e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1072932b0c3eSLois Curfman McInnes } 10733a40ed3dSBarry Smith PetscFunctionReturn(0); 1074932b0c3eSLois Curfman McInnes } 1075289bc588SBarry Smith 10764a2ae208SSatish Balay #undef __FUNCT__ 10774a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1078dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1079289bc588SBarry Smith { 1080ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1081dfbe8321SBarry Smith PetscErrorCode ierr; 108290f02eecSBarry Smith 10833a40ed3dSBarry Smith PetscFunctionBegin; 1084aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1085d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n); 1086a5a9c739SBarry Smith #endif 108705b42c5fSBarry Smith ierr = PetscFree(l->pivots);CHKERRQ(ierr); 10886857c123SSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1089606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 1090dbd8c25aSHong Zhang 1091dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 1092901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 10934ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 10944ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 10954ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 10963a40ed3dSBarry Smith PetscFunctionReturn(0); 1097289bc588SBarry Smith } 1098289bc588SBarry Smith 10994a2ae208SSatish Balay #undef __FUNCT__ 11004a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1101fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout) 1102289bc588SBarry Smith { 1103c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11046849ba73SBarry Smith PetscErrorCode ierr; 110513f74950SBarry Smith PetscInt k,j,m,n,M; 110687828ca2SBarry Smith PetscScalar *v,tmp; 110748b35521SBarry Smith 11083a40ed3dSBarry Smith PetscFunctionBegin; 1109d0f46423SBarry Smith v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n; 1110e9695a30SBarry Smith if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */ 1111e7e72b3dSBarry Smith if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1112e7e72b3dSBarry Smith else { 1113d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1114289bc588SBarry Smith for (k=0; k<j; k++) { 11151b807ce4Svictorle tmp = v[j + k*M]; 11161b807ce4Svictorle v[j + k*M] = v[k + j*M]; 11171b807ce4Svictorle v[k + j*M] = tmp; 1118289bc588SBarry Smith } 1119289bc588SBarry Smith } 1120d64ed03dSBarry Smith } 11213a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1122d3e5ee88SLois Curfman McInnes Mat tmat; 1123ec8511deSBarry Smith Mat_SeqDense *tmatd; 112487828ca2SBarry Smith PetscScalar *v2; 1125ea709b57SSatish Balay 1126fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 11277adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr); 1128d0f46423SBarry Smith ierr = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr); 11297adad957SLisandro Dalcin ierr = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 11305c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1131fc4dec0aSBarry Smith } else { 1132fc4dec0aSBarry Smith tmat = *matout; 1133fc4dec0aSBarry Smith } 1134ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 11350de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1136d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 11371b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1138d3e5ee88SLois Curfman McInnes } 11396d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11406d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1141d3e5ee88SLois Curfman McInnes *matout = tmat; 114248b35521SBarry Smith } 11433a40ed3dSBarry Smith PetscFunctionReturn(0); 1144289bc588SBarry Smith } 1145289bc588SBarry Smith 11464a2ae208SSatish Balay #undef __FUNCT__ 11474a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1148dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1149289bc588SBarry Smith { 1150c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1151c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 115213f74950SBarry Smith PetscInt i,j; 115387828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 11549ea5d5aeSSatish Balay 11553a40ed3dSBarry Smith PetscFunctionBegin; 1156d0f46423SBarry Smith if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1157d0f46423SBarry Smith if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1158d0f46423SBarry Smith for (i=0; i<A1->rmap->n; i++) { 11591b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1160d0f46423SBarry Smith for (j=0; j<A1->cmap->n; j++) { 11613a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 11621b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 11631b807ce4Svictorle } 1164289bc588SBarry Smith } 116577c4ece6SBarry Smith *flg = PETSC_TRUE; 11663a40ed3dSBarry Smith PetscFunctionReturn(0); 1167289bc588SBarry Smith } 1168289bc588SBarry Smith 11694a2ae208SSatish Balay #undef __FUNCT__ 11704a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1171dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1172289bc588SBarry Smith { 1173c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1174dfbe8321SBarry Smith PetscErrorCode ierr; 117513f74950SBarry Smith PetscInt i,n,len; 117687828ca2SBarry Smith PetscScalar *x,zero = 0.0; 117744cd7ae7SLois Curfman McInnes 11783a40ed3dSBarry Smith PetscFunctionBegin; 11792dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 11807a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 11811ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1182d0f46423SBarry Smith len = PetscMin(A->rmap->n,A->cmap->n); 1183e32f2f54SBarry Smith if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 118444cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 11851b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1186289bc588SBarry Smith } 11871ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 11883a40ed3dSBarry Smith PetscFunctionReturn(0); 1189289bc588SBarry Smith } 1190289bc588SBarry Smith 11914a2ae208SSatish Balay #undef __FUNCT__ 11924a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1193dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1194289bc588SBarry Smith { 1195c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 119687828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1197dfbe8321SBarry Smith PetscErrorCode ierr; 1198d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n; 119955659b69SBarry Smith 12003a40ed3dSBarry Smith PetscFunctionBegin; 120128988994SBarry Smith if (ll) { 12027a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 12031ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1204e32f2f54SBarry Smith if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1205da3a660dSBarry Smith for (i=0; i<m; i++) { 1206da3a660dSBarry Smith x = l[i]; 1207da3a660dSBarry Smith v = mat->v + i; 1208da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1209da3a660dSBarry Smith } 12101ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1211efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1212da3a660dSBarry Smith } 121328988994SBarry Smith if (rr) { 12147a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 12151ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1216e32f2f54SBarry Smith if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1217da3a660dSBarry Smith for (i=0; i<n; i++) { 1218da3a660dSBarry Smith x = r[i]; 1219da3a660dSBarry Smith v = mat->v + i*m; 1220da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1221da3a660dSBarry Smith } 12221ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1223efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1224da3a660dSBarry Smith } 12253a40ed3dSBarry Smith PetscFunctionReturn(0); 1226289bc588SBarry Smith } 1227289bc588SBarry Smith 12284a2ae208SSatish Balay #undef __FUNCT__ 12294a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1230dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1231289bc588SBarry Smith { 1232c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 123387828ca2SBarry Smith PetscScalar *v = mat->v; 1234329f5518SBarry Smith PetscReal sum = 0.0; 1235d0f46423SBarry Smith PetscInt lda=mat->lda,m=A->rmap->n,i,j; 1236efee365bSSatish Balay PetscErrorCode ierr; 123755659b69SBarry Smith 12383a40ed3dSBarry Smith PetscFunctionBegin; 1239289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1240a5ce6ee0Svictorle if (lda>m) { 1241d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1242a5ce6ee0Svictorle v = mat->v+j*lda; 1243a5ce6ee0Svictorle for (i=0; i<m; i++) { 1244a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1245a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1246a5ce6ee0Svictorle #else 1247a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1248a5ce6ee0Svictorle #endif 1249a5ce6ee0Svictorle } 1250a5ce6ee0Svictorle } 1251a5ce6ee0Svictorle } else { 1252d0f46423SBarry Smith for (i=0; i<A->cmap->n*A->rmap->n; i++) { 1253aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1254329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1255289bc588SBarry Smith #else 1256289bc588SBarry Smith sum += (*v)*(*v); v++; 1257289bc588SBarry Smith #endif 1258289bc588SBarry Smith } 1259a5ce6ee0Svictorle } 1260064f8208SBarry Smith *nrm = sqrt(sum); 1261dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr); 12623a40ed3dSBarry Smith } else if (type == NORM_1) { 1263064f8208SBarry Smith *nrm = 0.0; 1264d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 12651b807ce4Svictorle v = mat->v + j*mat->lda; 1266289bc588SBarry Smith sum = 0.0; 1267d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 126833a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1269289bc588SBarry Smith } 1270064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1271289bc588SBarry Smith } 1272d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 12733a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1274064f8208SBarry Smith *nrm = 0.0; 1275d0f46423SBarry Smith for (j=0; j<A->rmap->n; j++) { 1276289bc588SBarry Smith v = mat->v + j; 1277289bc588SBarry Smith sum = 0.0; 1278d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) { 12791b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1280289bc588SBarry Smith } 1281064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1282289bc588SBarry Smith } 1283d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 1284e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 12853a40ed3dSBarry Smith PetscFunctionReturn(0); 1286289bc588SBarry Smith } 1287289bc588SBarry Smith 12884a2ae208SSatish Balay #undef __FUNCT__ 12894a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 12904e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscTruth flg) 1291289bc588SBarry Smith { 1292c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 129363ba0a88SBarry Smith PetscErrorCode ierr; 129467e560aaSBarry Smith 12953a40ed3dSBarry Smith PetscFunctionBegin; 1296b5a2b587SKris Buschelman switch (op) { 1297b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 12984e0d8c25SBarry Smith aij->roworiented = flg; 1299b5a2b587SKris Buschelman break; 1300512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 1301b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 13023971808eSMatthew Knepley case MAT_NEW_NONZERO_ALLOCATION_ERR: 13034e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 1304b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1305b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 130677e54ba9SKris Buschelman case MAT_SYMMETRIC: 130777e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 13089a4540c5SBarry Smith case MAT_HERMITIAN: 13099a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1310600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 1311290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 131277e54ba9SKris Buschelman break; 1313b5a2b587SKris Buschelman default: 1314e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 13153a40ed3dSBarry Smith } 13163a40ed3dSBarry Smith PetscFunctionReturn(0); 1317289bc588SBarry Smith } 1318289bc588SBarry Smith 13194a2ae208SSatish Balay #undef __FUNCT__ 13204a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1321dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 13226f0a148fSBarry Smith { 1323ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 13246849ba73SBarry Smith PetscErrorCode ierr; 1325d0f46423SBarry Smith PetscInt lda=l->lda,m=A->rmap->n,j; 13263a40ed3dSBarry Smith 13273a40ed3dSBarry Smith PetscFunctionBegin; 1328a5ce6ee0Svictorle if (lda>m) { 1329d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1330a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1331a5ce6ee0Svictorle } 1332a5ce6ee0Svictorle } else { 1333d0f46423SBarry Smith ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1334a5ce6ee0Svictorle } 13353a40ed3dSBarry Smith PetscFunctionReturn(0); 13366f0a148fSBarry Smith } 13376f0a148fSBarry Smith 13384a2ae208SSatish Balay #undef __FUNCT__ 13394a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 1340f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag) 13416f0a148fSBarry Smith { 1342ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1343d0f46423SBarry Smith PetscInt n = A->cmap->n,i,j; 134487828ca2SBarry Smith PetscScalar *slot; 134555659b69SBarry Smith 13463a40ed3dSBarry Smith PetscFunctionBegin; 13476f0a148fSBarry Smith for (i=0; i<N; i++) { 13486f0a148fSBarry Smith slot = l->v + rows[i]; 13496f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 13506f0a148fSBarry Smith } 1351f4df32b1SMatthew Knepley if (diag != 0.0) { 13526f0a148fSBarry Smith for (i=0; i<N; i++) { 13536f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1354f4df32b1SMatthew Knepley *slot = diag; 13556f0a148fSBarry Smith } 13566f0a148fSBarry Smith } 13573a40ed3dSBarry Smith PetscFunctionReturn(0); 13586f0a148fSBarry Smith } 1359557bce09SLois Curfman McInnes 13604a2ae208SSatish Balay #undef __FUNCT__ 13614a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 1362dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 136364e87e97SBarry Smith { 1364c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 13653a40ed3dSBarry Smith 13663a40ed3dSBarry Smith PetscFunctionBegin; 1367e32f2f54SBarry Smith if (mat->lda != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows"); 136864e87e97SBarry Smith *array = mat->v; 13693a40ed3dSBarry Smith PetscFunctionReturn(0); 137064e87e97SBarry Smith } 13710754003eSLois Curfman McInnes 13724a2ae208SSatish Balay #undef __FUNCT__ 13734a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 1374dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1375ff14e315SSatish Balay { 13763a40ed3dSBarry Smith PetscFunctionBegin; 137709b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 13783a40ed3dSBarry Smith PetscFunctionReturn(0); 1379ff14e315SSatish Balay } 13800754003eSLois Curfman McInnes 13814a2ae208SSatish Balay #undef __FUNCT__ 13824a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 138313f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 13840754003eSLois Curfman McInnes { 1385c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 13866849ba73SBarry Smith PetscErrorCode ierr; 13875d0c19d7SBarry Smith PetscInt i,j,nrows,ncols; 13885d0c19d7SBarry Smith const PetscInt *irow,*icol; 138987828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 13900754003eSLois Curfman McInnes Mat newmat; 13910754003eSLois Curfman McInnes 13923a40ed3dSBarry Smith PetscFunctionBegin; 139378b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 139478b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1395e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1396e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 13970754003eSLois Curfman McInnes 1398182d2002SSatish Balay /* Check submatrixcall */ 1399182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 140013f74950SBarry Smith PetscInt n_cols,n_rows; 1401182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 140221a2c019SBarry Smith if (n_rows != nrows || n_cols != ncols) { 140321a2c019SBarry Smith /* resize the result result matrix to match number of requested rows/columns */ 1404c61587bbSBarry Smith ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 140521a2c019SBarry Smith } 1406182d2002SSatish Balay newmat = *B; 1407182d2002SSatish Balay } else { 14080754003eSLois Curfman McInnes /* Create and fill new matrix */ 14097adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr); 1410f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 14117adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 14125c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1413182d2002SSatish Balay } 1414182d2002SSatish Balay 1415182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1416182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1417182d2002SSatish Balay 1418182d2002SSatish Balay for (i=0; i<ncols; i++) { 14196de62eeeSBarry Smith av = v + mat->lda*icol[i]; 1420182d2002SSatish Balay for (j=0; j<nrows; j++) { 1421182d2002SSatish Balay *bv++ = av[irow[j]]; 14220754003eSLois Curfman McInnes } 14230754003eSLois Curfman McInnes } 1424182d2002SSatish Balay 1425182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 14266d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14276d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14280754003eSLois Curfman McInnes 14290754003eSLois Curfman McInnes /* Free work space */ 143078b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 143178b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1432182d2002SSatish Balay *B = newmat; 14333a40ed3dSBarry Smith PetscFunctionReturn(0); 14340754003eSLois Curfman McInnes } 14350754003eSLois Curfman McInnes 14364a2ae208SSatish Balay #undef __FUNCT__ 14374a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 143813f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1439905e6a2fSBarry Smith { 14406849ba73SBarry Smith PetscErrorCode ierr; 144113f74950SBarry Smith PetscInt i; 1442905e6a2fSBarry Smith 14433a40ed3dSBarry Smith PetscFunctionBegin; 1444905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1445b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1446905e6a2fSBarry Smith } 1447905e6a2fSBarry Smith 1448905e6a2fSBarry Smith for (i=0; i<n; i++) { 14496a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1450905e6a2fSBarry Smith } 14513a40ed3dSBarry Smith PetscFunctionReturn(0); 1452905e6a2fSBarry Smith } 1453905e6a2fSBarry Smith 14544a2ae208SSatish Balay #undef __FUNCT__ 1455c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1456c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1457c0aa2d19SHong Zhang { 1458c0aa2d19SHong Zhang PetscFunctionBegin; 1459c0aa2d19SHong Zhang PetscFunctionReturn(0); 1460c0aa2d19SHong Zhang } 1461c0aa2d19SHong Zhang 1462c0aa2d19SHong Zhang #undef __FUNCT__ 1463c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1464c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1465c0aa2d19SHong Zhang { 1466c0aa2d19SHong Zhang PetscFunctionBegin; 1467c0aa2d19SHong Zhang PetscFunctionReturn(0); 1468c0aa2d19SHong Zhang } 1469c0aa2d19SHong Zhang 1470c0aa2d19SHong Zhang #undef __FUNCT__ 14714a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1472dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 14734b0e389bSBarry Smith { 14744b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 14756849ba73SBarry Smith PetscErrorCode ierr; 1476d0f46423SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j; 14773a40ed3dSBarry Smith 14783a40ed3dSBarry Smith PetscFunctionBegin; 147933f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 148033f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1481cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 14823a40ed3dSBarry Smith PetscFunctionReturn(0); 14833a40ed3dSBarry Smith } 1484e32f2f54SBarry Smith if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1485a5ce6ee0Svictorle if (lda1>m || lda2>m) { 14860dbb7854Svictorle for (j=0; j<n; j++) { 1487a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1488a5ce6ee0Svictorle } 1489a5ce6ee0Svictorle } else { 1490d0f46423SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1491a5ce6ee0Svictorle } 1492273d9f13SBarry Smith PetscFunctionReturn(0); 1493273d9f13SBarry Smith } 1494273d9f13SBarry Smith 14954a2ae208SSatish Balay #undef __FUNCT__ 14964a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1497dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A) 1498273d9f13SBarry Smith { 1499dfbe8321SBarry Smith PetscErrorCode ierr; 1500273d9f13SBarry Smith 1501273d9f13SBarry Smith PetscFunctionBegin; 1502273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 15033a40ed3dSBarry Smith PetscFunctionReturn(0); 15044b0e389bSBarry Smith } 15054b0e389bSBarry Smith 1506284134d9SBarry Smith #undef __FUNCT__ 1507284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense" 1508284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N) 1509284134d9SBarry Smith { 1510284134d9SBarry Smith PetscFunctionBegin; 151121a2c019SBarry Smith /* this will not be called before lda, Mmax, and Nmax have been set */ 1512284134d9SBarry Smith m = PetscMax(m,M); 1513284134d9SBarry Smith n = PetscMax(n,N); 1514a868139aSShri Abhyankar 151586d161a7SShri Abhyankar /* if (m > a->Mmax) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot yet resize number rows of dense matrix larger then its initial size %d, requested %d",a->lda,(int)m); 151686d161a7SShri Abhyankar if (n > a->Nmax) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot yet resize number columns of dense matrix larger then its initial size %d, requested %d",a->Nmax,(int)n); 151786d161a7SShri Abhyankar */ 1518dc5cefdeSJed Brown A->rmap->n = A->rmap->N = m; 1519d0f46423SBarry Smith A->cmap->n = A->cmap->N = n; 1520284134d9SBarry Smith PetscFunctionReturn(0); 1521284134d9SBarry Smith } 1522170fe5c8SBarry Smith 1523ba337c44SJed Brown #undef __FUNCT__ 1524ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense" 1525ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A) 1526ba337c44SJed Brown { 1527ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1528ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1529ba337c44SJed Brown PetscScalar *aa = a->v; 1530ba337c44SJed Brown 1531ba337c44SJed Brown PetscFunctionBegin; 1532ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]); 1533ba337c44SJed Brown PetscFunctionReturn(0); 1534ba337c44SJed Brown } 1535ba337c44SJed Brown 1536ba337c44SJed Brown #undef __FUNCT__ 1537ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense" 1538ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A) 1539ba337c44SJed Brown { 1540ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1541ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1542ba337c44SJed Brown PetscScalar *aa = a->v; 1543ba337c44SJed Brown 1544ba337c44SJed Brown PetscFunctionBegin; 1545ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 1546ba337c44SJed Brown PetscFunctionReturn(0); 1547ba337c44SJed Brown } 1548ba337c44SJed Brown 1549ba337c44SJed Brown #undef __FUNCT__ 1550ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense" 1551ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A) 1552ba337c44SJed Brown { 1553ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1554ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1555ba337c44SJed Brown PetscScalar *aa = a->v; 1556ba337c44SJed Brown 1557ba337c44SJed Brown PetscFunctionBegin; 1558ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 1559ba337c44SJed Brown PetscFunctionReturn(0); 1560ba337c44SJed Brown } 1561284134d9SBarry Smith 1562a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/ 1563a9fe9ddaSSatish Balay #undef __FUNCT__ 1564a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1565a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1566a9fe9ddaSSatish Balay { 1567a9fe9ddaSSatish Balay PetscErrorCode ierr; 1568a9fe9ddaSSatish Balay 1569a9fe9ddaSSatish Balay PetscFunctionBegin; 1570a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1571a9fe9ddaSSatish Balay ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1572a9fe9ddaSSatish Balay } 1573a9fe9ddaSSatish Balay ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1574a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1575a9fe9ddaSSatish Balay } 1576a9fe9ddaSSatish Balay 1577a9fe9ddaSSatish Balay #undef __FUNCT__ 1578a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1579a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1580a9fe9ddaSSatish Balay { 1581ee16a9a1SHong Zhang PetscErrorCode ierr; 1582d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 1583ee16a9a1SHong Zhang Mat Cmat; 1584a9fe9ddaSSatish Balay 1585ee16a9a1SHong Zhang PetscFunctionBegin; 1586e32f2f54SBarry Smith if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->rmap->n %d\n",A->cmap->n,B->rmap->n); 1587ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1588ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1589ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1590ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1591ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1592ee16a9a1SHong Zhang *C = Cmat; 1593ee16a9a1SHong Zhang PetscFunctionReturn(0); 1594ee16a9a1SHong Zhang } 1595a9fe9ddaSSatish Balay 159698a3b096SSatish Balay #undef __FUNCT__ 1597a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1598a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1599a9fe9ddaSSatish Balay { 1600a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1601a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1602a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 16030805154bSBarry Smith PetscBLASInt m,n,k; 1604a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1605a9fe9ddaSSatish Balay 1606a9fe9ddaSSatish Balay PetscFunctionBegin; 1607d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 1608d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1609d0f46423SBarry Smith k = PetscBLASIntCast(A->cmap->n); 1610a9fe9ddaSSatish Balay BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1611a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1612a9fe9ddaSSatish Balay } 1613a9fe9ddaSSatish Balay 1614a9fe9ddaSSatish Balay #undef __FUNCT__ 1615a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense" 1616a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1617a9fe9ddaSSatish Balay { 1618a9fe9ddaSSatish Balay PetscErrorCode ierr; 1619a9fe9ddaSSatish Balay 1620a9fe9ddaSSatish Balay PetscFunctionBegin; 1621a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1622a9fe9ddaSSatish Balay ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1623a9fe9ddaSSatish Balay } 1624a9fe9ddaSSatish Balay ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1625a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1626a9fe9ddaSSatish Balay } 1627a9fe9ddaSSatish Balay 1628a9fe9ddaSSatish Balay #undef __FUNCT__ 1629a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense" 1630a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1631a9fe9ddaSSatish Balay { 1632ee16a9a1SHong Zhang PetscErrorCode ierr; 1633d0f46423SBarry Smith PetscInt m=A->cmap->n,n=B->cmap->n; 1634ee16a9a1SHong Zhang Mat Cmat; 1635a9fe9ddaSSatish Balay 1636ee16a9a1SHong Zhang PetscFunctionBegin; 1637e32f2f54SBarry Smith if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->rmap->n %d != B->rmap->n %d\n",A->rmap->n,B->rmap->n); 1638ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1639ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1640ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1641ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1642ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1643ee16a9a1SHong Zhang *C = Cmat; 1644ee16a9a1SHong Zhang PetscFunctionReturn(0); 1645ee16a9a1SHong Zhang } 1646a9fe9ddaSSatish Balay 1647a9fe9ddaSSatish Balay #undef __FUNCT__ 1648a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense" 1649a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1650a9fe9ddaSSatish Balay { 1651a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1652a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1653a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 16540805154bSBarry Smith PetscBLASInt m,n,k; 1655a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1656a9fe9ddaSSatish Balay 1657a9fe9ddaSSatish Balay PetscFunctionBegin; 1658d0f46423SBarry Smith m = PetscBLASIntCast(A->cmap->n); 1659d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1660d0f46423SBarry Smith k = PetscBLASIntCast(A->rmap->n); 16612fbe02b9SBarry Smith /* 16622fbe02b9SBarry Smith Note the m and n arguments below are the number rows and columns of A', not A! 16632fbe02b9SBarry Smith */ 1664a9fe9ddaSSatish Balay BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1665a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1666a9fe9ddaSSatish Balay } 1667985db425SBarry Smith 1668985db425SBarry Smith #undef __FUNCT__ 1669985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense" 1670985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 1671985db425SBarry Smith { 1672985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1673985db425SBarry Smith PetscErrorCode ierr; 1674d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1675985db425SBarry Smith PetscScalar *x; 1676985db425SBarry Smith MatScalar *aa = a->v; 1677985db425SBarry Smith 1678985db425SBarry Smith PetscFunctionBegin; 1679e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1680985db425SBarry Smith 1681985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1682985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1683985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1684e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1685985db425SBarry Smith for (i=0; i<m; i++) { 1686985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1687985db425SBarry Smith for (j=1; j<n; j++){ 1688985db425SBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1689985db425SBarry Smith } 1690985db425SBarry Smith } 1691985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1692985db425SBarry Smith PetscFunctionReturn(0); 1693985db425SBarry Smith } 1694985db425SBarry Smith 1695985db425SBarry Smith #undef __FUNCT__ 1696985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense" 1697985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 1698985db425SBarry Smith { 1699985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1700985db425SBarry Smith PetscErrorCode ierr; 1701d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1702985db425SBarry Smith PetscScalar *x; 1703985db425SBarry Smith PetscReal atmp; 1704985db425SBarry Smith MatScalar *aa = a->v; 1705985db425SBarry Smith 1706985db425SBarry Smith PetscFunctionBegin; 1707e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1708985db425SBarry Smith 1709985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1710985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1711985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1712e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1713985db425SBarry Smith for (i=0; i<m; i++) { 17149189402eSHong Zhang x[i] = PetscAbsScalar(aa[i]); 1715985db425SBarry Smith for (j=1; j<n; j++){ 1716985db425SBarry Smith atmp = PetscAbsScalar(aa[i+m*j]); 1717985db425SBarry Smith if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 1718985db425SBarry Smith } 1719985db425SBarry Smith } 1720985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1721985db425SBarry Smith PetscFunctionReturn(0); 1722985db425SBarry Smith } 1723985db425SBarry Smith 1724985db425SBarry Smith #undef __FUNCT__ 1725985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense" 1726985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 1727985db425SBarry Smith { 1728985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1729985db425SBarry Smith PetscErrorCode ierr; 1730d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1731985db425SBarry Smith PetscScalar *x; 1732985db425SBarry Smith MatScalar *aa = a->v; 1733985db425SBarry Smith 1734985db425SBarry Smith PetscFunctionBegin; 1735e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1736985db425SBarry Smith 1737985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1738985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1739985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1740e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1741985db425SBarry Smith for (i=0; i<m; i++) { 1742985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1743985db425SBarry Smith for (j=1; j<n; j++){ 1744985db425SBarry Smith if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1745985db425SBarry Smith } 1746985db425SBarry Smith } 1747985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1748985db425SBarry Smith PetscFunctionReturn(0); 1749985db425SBarry Smith } 1750985db425SBarry Smith 17518d0534beSBarry Smith #undef __FUNCT__ 17528d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense" 17538d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 17548d0534beSBarry Smith { 17558d0534beSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 17568d0534beSBarry Smith PetscErrorCode ierr; 17578d0534beSBarry Smith PetscScalar *x; 17588d0534beSBarry Smith 17598d0534beSBarry Smith PetscFunctionBegin; 1760e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 17618d0534beSBarry Smith 17628d0534beSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1763d0f46423SBarry Smith ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 17648d0534beSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 17658d0534beSBarry Smith PetscFunctionReturn(0); 17668d0534beSBarry Smith } 17678d0534beSBarry Smith 1768289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1769a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1770905e6a2fSBarry Smith MatGetRow_SeqDense, 1771905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1772905e6a2fSBarry Smith MatMult_SeqDense, 177397304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 17747c922b88SBarry Smith MatMultTranspose_SeqDense, 17757c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1776db4efbfdSBarry Smith 0, 1777db4efbfdSBarry Smith 0, 1778db4efbfdSBarry Smith 0, 1779db4efbfdSBarry Smith /*10*/ 0, 1780905e6a2fSBarry Smith MatLUFactor_SeqDense, 1781905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 178241f059aeSBarry Smith MatSOR_SeqDense, 1783ec8511deSBarry Smith MatTranspose_SeqDense, 178497304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1785905e6a2fSBarry Smith MatEqual_SeqDense, 1786905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1787905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1788905e6a2fSBarry Smith MatNorm_SeqDense, 1789c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense, 1790c0aa2d19SHong Zhang MatAssemblyEnd_SeqDense, 1791905e6a2fSBarry Smith MatSetOption_SeqDense, 1792905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1793d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense, 1794db4efbfdSBarry Smith 0, 1795db4efbfdSBarry Smith 0, 1796db4efbfdSBarry Smith 0, 1797db4efbfdSBarry Smith 0, 1798d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqDense, 1799273d9f13SBarry Smith 0, 1800905e6a2fSBarry Smith 0, 1801905e6a2fSBarry Smith MatGetArray_SeqDense, 1802905e6a2fSBarry Smith MatRestoreArray_SeqDense, 1803d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense, 1804a5ae1ecdSBarry Smith 0, 1805a5ae1ecdSBarry Smith 0, 1806a5ae1ecdSBarry Smith 0, 1807a5ae1ecdSBarry Smith 0, 1808d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense, 1809a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1810a5ae1ecdSBarry Smith 0, 18114b0e389bSBarry Smith MatGetValues_SeqDense, 1812a5ae1ecdSBarry Smith MatCopy_SeqDense, 1813d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense, 1814a5ae1ecdSBarry Smith MatScale_SeqDense, 1815a5ae1ecdSBarry Smith 0, 1816a5ae1ecdSBarry Smith 0, 1817a5ae1ecdSBarry Smith 0, 1818d519adbfSMatthew Knepley /*49*/ 0, 1819a5ae1ecdSBarry Smith 0, 1820a5ae1ecdSBarry Smith 0, 1821a5ae1ecdSBarry Smith 0, 1822a5ae1ecdSBarry Smith 0, 1823d519adbfSMatthew Knepley /*54*/ 0, 1824a5ae1ecdSBarry Smith 0, 1825a5ae1ecdSBarry Smith 0, 1826a5ae1ecdSBarry Smith 0, 1827a5ae1ecdSBarry Smith 0, 1828d519adbfSMatthew Knepley /*59*/ 0, 1829e03a110bSBarry Smith MatDestroy_SeqDense, 1830e03a110bSBarry Smith MatView_SeqDense, 1831357abbc8SBarry Smith 0, 183297304618SKris Buschelman 0, 1833d519adbfSMatthew Knepley /*64*/ 0, 183497304618SKris Buschelman 0, 183597304618SKris Buschelman 0, 183697304618SKris Buschelman 0, 183797304618SKris Buschelman 0, 1838d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense, 183997304618SKris Buschelman 0, 184097304618SKris Buschelman 0, 184197304618SKris Buschelman 0, 184297304618SKris Buschelman 0, 1843d519adbfSMatthew Knepley /*74*/ 0, 184497304618SKris Buschelman 0, 184597304618SKris Buschelman 0, 184697304618SKris Buschelman 0, 184797304618SKris Buschelman 0, 1848d519adbfSMatthew Knepley /*79*/ 0, 184997304618SKris Buschelman 0, 185097304618SKris Buschelman 0, 185197304618SKris Buschelman 0, 1852*5bba2384SShri Abhyankar /*83*/ MatLoad_SeqDense, 1853865e5f61SKris Buschelman 0, 18541cbb95d3SBarry Smith MatIsHermitian_SeqDense, 1855865e5f61SKris Buschelman 0, 1856865e5f61SKris Buschelman 0, 1857865e5f61SKris Buschelman 0, 1858d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense, 1859a9fe9ddaSSatish Balay MatMatMultSymbolic_SeqDense_SeqDense, 1860a9fe9ddaSSatish Balay MatMatMultNumeric_SeqDense_SeqDense, 1861865e5f61SKris Buschelman 0, 1862865e5f61SKris Buschelman 0, 1863d519adbfSMatthew Knepley /*94*/ 0, 1864a9fe9ddaSSatish Balay MatMatMultTranspose_SeqDense_SeqDense, 1865a9fe9ddaSSatish Balay MatMatMultTransposeSymbolic_SeqDense_SeqDense, 1866a9fe9ddaSSatish Balay MatMatMultTransposeNumeric_SeqDense_SeqDense, 1867284134d9SBarry Smith 0, 1868d519adbfSMatthew Knepley /*99*/ 0, 1869284134d9SBarry Smith 0, 1870284134d9SBarry Smith 0, 1871ba337c44SJed Brown MatConjugate_SeqDense, 1872985db425SBarry Smith MatSetSizes_SeqDense, 1873ba337c44SJed Brown /*104*/0, 1874ba337c44SJed Brown MatRealPart_SeqDense, 1875ba337c44SJed Brown MatImaginaryPart_SeqDense, 1876985db425SBarry Smith 0, 1877985db425SBarry Smith 0, 1878d519adbfSMatthew Knepley /*109*/0, 1879985db425SBarry Smith 0, 18808d0534beSBarry Smith MatGetRowMin_SeqDense, 1881aabbc4fbSShri Abhyankar MatGetColumnVector_SeqDense, 1882aabbc4fbSShri Abhyankar 0, 1883aabbc4fbSShri Abhyankar /*114*/0, 1884aabbc4fbSShri Abhyankar 0, 1885aabbc4fbSShri Abhyankar 0, 1886aabbc4fbSShri Abhyankar 0, 1887aabbc4fbSShri Abhyankar 0, 1888aabbc4fbSShri Abhyankar /*119*/0, 1889aabbc4fbSShri Abhyankar 0, 1890aabbc4fbSShri Abhyankar 0, 1891*5bba2384SShri Abhyankar 0 1892985db425SBarry Smith }; 189390ace30eSBarry Smith 18944a2ae208SSatish Balay #undef __FUNCT__ 18954a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 18964b828684SBarry Smith /*@C 1897fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1898d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1899d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1900289bc588SBarry Smith 1901db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1902db81eaa0SLois Curfman McInnes 190320563c6bSBarry Smith Input Parameters: 1904db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 19050c775827SLois Curfman McInnes . m - number of rows 190618f449edSLois Curfman McInnes . n - number of columns 1907c0235b3cSMatthew Knepley - data - optional location of matrix data in column major order. Set data=PETSC_NULL for PETSc 1908dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 190920563c6bSBarry Smith 191020563c6bSBarry Smith Output Parameter: 191144cd7ae7SLois Curfman McInnes . A - the matrix 191220563c6bSBarry Smith 1913b259b22eSLois Curfman McInnes Notes: 191418f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 191518f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1916b4fd4287SBarry Smith set data=PETSC_NULL. 191718f449edSLois Curfman McInnes 1918027ccd11SLois Curfman McInnes Level: intermediate 1919027ccd11SLois Curfman McInnes 1920dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1921d65003e9SLois Curfman McInnes 1922db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 192320563c6bSBarry Smith @*/ 1924be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 1925289bc588SBarry Smith { 1926dfbe8321SBarry Smith PetscErrorCode ierr; 19273b2fbd54SBarry Smith 19283a40ed3dSBarry Smith PetscFunctionBegin; 1929f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1930f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1931273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1932273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1933273d9f13SBarry Smith PetscFunctionReturn(0); 1934273d9f13SBarry Smith } 1935273d9f13SBarry Smith 19364a2ae208SSatish Balay #undef __FUNCT__ 1937afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation" 1938273d9f13SBarry Smith /*@C 1939273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1940273d9f13SBarry Smith 1941273d9f13SBarry Smith Collective on MPI_Comm 1942273d9f13SBarry Smith 1943273d9f13SBarry Smith Input Parameters: 1944273d9f13SBarry Smith + A - the matrix 1945273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1946273d9f13SBarry Smith 1947273d9f13SBarry Smith Notes: 1948273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1949273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1950284134d9SBarry Smith need not call this routine. 1951273d9f13SBarry Smith 1952273d9f13SBarry Smith Level: intermediate 1953273d9f13SBarry Smith 1954273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1955273d9f13SBarry Smith 1956867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA() 1957867c911aSBarry Smith 1958273d9f13SBarry Smith @*/ 1959be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 1960273d9f13SBarry Smith { 1961dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar[]); 1962a23d5eceSKris Buschelman 1963a23d5eceSKris Buschelman PetscFunctionBegin; 1964a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1965a23d5eceSKris Buschelman if (f) { 1966a23d5eceSKris Buschelman ierr = (*f)(B,data);CHKERRQ(ierr); 1967a23d5eceSKris Buschelman } 1968a23d5eceSKris Buschelman PetscFunctionReturn(0); 1969a23d5eceSKris Buschelman } 1970a23d5eceSKris Buschelman 1971a23d5eceSKris Buschelman EXTERN_C_BEGIN 1972a23d5eceSKris Buschelman #undef __FUNCT__ 1973afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense" 1974be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 1975a23d5eceSKris Buschelman { 1976273d9f13SBarry Smith Mat_SeqDense *b; 1977dfbe8321SBarry Smith PetscErrorCode ierr; 1978273d9f13SBarry Smith 1979273d9f13SBarry Smith PetscFunctionBegin; 1980273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1981a868139aSShri Abhyankar 198234ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr); 198334ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr); 198434ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 198534ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 198634ef9618SShri Abhyankar 1987273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 198886d161a7SShri Abhyankar b->Mmax = B->rmap->n; 198986d161a7SShri Abhyankar b->Nmax = B->cmap->n; 199086d161a7SShri Abhyankar if(b->lda <= 0 || b->changelda) b->lda = B->rmap->n; 199186d161a7SShri Abhyankar 19929e8f95c4SLisandro Dalcin if (!data) { /* petsc-allocated storage */ 19939e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 19945afd5e0cSBarry Smith ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 1995284134d9SBarry Smith ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 1996284134d9SBarry Smith ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 19979e8f95c4SLisandro Dalcin b->user_alloc = PETSC_FALSE; 1998273d9f13SBarry Smith } else { /* user-allocated storage */ 19999e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 2000273d9f13SBarry Smith b->v = data; 2001273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 2002273d9f13SBarry Smith } 20030450473dSBarry Smith B->assembled = PETSC_TRUE; 2004273d9f13SBarry Smith PetscFunctionReturn(0); 2005273d9f13SBarry Smith } 2006a23d5eceSKris Buschelman EXTERN_C_END 2007273d9f13SBarry Smith 20081b807ce4Svictorle #undef __FUNCT__ 20091b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 20101b807ce4Svictorle /*@C 20111b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 20121b807ce4Svictorle 20131b807ce4Svictorle Input parameter: 20141b807ce4Svictorle + A - the matrix 20151b807ce4Svictorle - lda - the leading dimension 20161b807ce4Svictorle 20171b807ce4Svictorle Notes: 2018867c911aSBarry Smith This routine is to be used in conjunction with MatSeqDenseSetPreallocation(); 20191b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 20201b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 20211b807ce4Svictorle 20221b807ce4Svictorle Level: intermediate 20231b807ce4Svictorle 20241b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 20251b807ce4Svictorle 2026284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 2027867c911aSBarry Smith 20281b807ce4Svictorle @*/ 2029be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda) 20301b807ce4Svictorle { 20311b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 203221a2c019SBarry Smith 20331b807ce4Svictorle PetscFunctionBegin; 2034e32f2f54SBarry Smith if (lda < B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n); 20351b807ce4Svictorle b->lda = lda; 203621a2c019SBarry Smith b->changelda = PETSC_FALSE; 203721a2c019SBarry Smith b->Mmax = PetscMax(b->Mmax,lda); 20381b807ce4Svictorle PetscFunctionReturn(0); 20391b807ce4Svictorle } 20401b807ce4Svictorle 20410bad9183SKris Buschelman /*MC 2042fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 20430bad9183SKris Buschelman 20440bad9183SKris Buschelman Options Database Keys: 20450bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 20460bad9183SKris Buschelman 20470bad9183SKris Buschelman Level: beginner 20480bad9183SKris Buschelman 204989665df3SBarry Smith .seealso: MatCreateSeqDense() 205089665df3SBarry Smith 20510bad9183SKris Buschelman M*/ 20520bad9183SKris Buschelman 2053273d9f13SBarry Smith EXTERN_C_BEGIN 20544a2ae208SSatish Balay #undef __FUNCT__ 20554a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 2056be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B) 2057273d9f13SBarry Smith { 2058273d9f13SBarry Smith Mat_SeqDense *b; 2059dfbe8321SBarry Smith PetscErrorCode ierr; 20607c334f02SBarry Smith PetscMPIInt size; 2061273d9f13SBarry Smith 2062273d9f13SBarry Smith PetscFunctionBegin; 20637adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr); 2064e32f2f54SBarry Smith if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 206555659b69SBarry Smith 206638f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr); 2067549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 206890f02eecSBarry Smith B->mapping = 0; 206944cd7ae7SLois Curfman McInnes B->data = (void*)b; 207018f449edSLois Curfman McInnes 207144cd7ae7SLois Curfman McInnes b->pivots = 0; 2072273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 2073273d9f13SBarry Smith b->v = 0; 207421a2c019SBarry Smith b->changelda = PETSC_FALSE; 20754e220ebcSLois Curfman McInnes 2076b24902e0SBarry Smith 2077ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C", 2078b24902e0SBarry Smith "MatGetFactor_seqdense_petsc", 2079b24902e0SBarry Smith MatGetFactor_seqdense_petsc);CHKERRQ(ierr); 2080a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 2081a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 2082a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 20834ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C", 20844ae313f4SHong Zhang "MatMatMult_SeqAIJ_SeqDense", 20854ae313f4SHong Zhang MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 20864ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C", 20874ae313f4SHong Zhang "MatMatMultSymbolic_SeqAIJ_SeqDense", 20884ae313f4SHong Zhang MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 20894ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C", 20904ae313f4SHong Zhang "MatMatMultNumeric_SeqAIJ_SeqDense", 20914ae313f4SHong Zhang MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 209217667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 20933a40ed3dSBarry Smith PetscFunctionReturn(0); 2094289bc588SBarry Smith } 2095273d9f13SBarry Smith EXTERN_C_END 2096