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" 81ace3abfcSBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool *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__ 192*85e2c93fSHong Zhang #define __FUNCT__ "MatMatSolve_SeqDense" 193*85e2c93fSHong Zhang PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X) 194*85e2c93fSHong Zhang { 195*85e2c93fSHong Zhang Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 196*85e2c93fSHong Zhang PetscErrorCode ierr; 197*85e2c93fSHong Zhang PetscScalar *b,*x; 198*85e2c93fSHong Zhang PetscBLASInt nrhs,info,m=PetscBLASIntCast(A->rmap->n); 199*85e2c93fSHong Zhang 200*85e2c93fSHong Zhang PetscFunctionBegin; 201*85e2c93fSHong Zhang ierr = MatGetSize(B,PETSC_NULL,&nrhs);CHKERRQ(ierr); 202*85e2c93fSHong Zhang ierr = MatGetArray(B,&b);CHKERRQ(ierr); 203*85e2c93fSHong Zhang ierr = MatGetArray(X,&x);CHKERRQ(ierr); 204*85e2c93fSHong Zhang 205*85e2c93fSHong Zhang ierr = PetscMemcpy(x,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 206*85e2c93fSHong Zhang 207*85e2c93fSHong Zhang if (A->factortype == MAT_FACTOR_LU) { 208*85e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_GETRS) 209*85e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 210*85e2c93fSHong Zhang #else 211*85e2c93fSHong Zhang LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info); 212*85e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 213*85e2c93fSHong Zhang #endif 214*85e2c93fSHong Zhang } else if (A->factortype == MAT_FACTOR_CHOLESKY){ 215*85e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_POTRS) 216*85e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 217*85e2c93fSHong Zhang #else 218*85e2c93fSHong Zhang LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info); 219*85e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 220*85e2c93fSHong Zhang #endif 221*85e2c93fSHong Zhang } 222*85e2c93fSHong Zhang else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 223*85e2c93fSHong Zhang 224*85e2c93fSHong Zhang ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 225*85e2c93fSHong Zhang ierr = MatRestoreArray(X,&x);CHKERRQ(ierr); 226*85e2c93fSHong Zhang ierr = PetscLogFlops(nrhs*(2.0*m*m- m));CHKERRQ(ierr); 227*85e2c93fSHong Zhang PetscFunctionReturn(0); 228*85e2c93fSHong Zhang } 229*85e2c93fSHong Zhang 230*85e2c93fSHong Zhang #undef __FUNCT__ 2314a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 232dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 233da3a660dSBarry Smith { 234c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 235dfbe8321SBarry Smith PetscErrorCode ierr; 23687828ca2SBarry Smith PetscScalar *x,*y; 237d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 23867e560aaSBarry Smith 2393a40ed3dSBarry Smith PetscFunctionBegin; 2401ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2411ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 242d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 243752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 244da3a660dSBarry Smith if (mat->pivots) { 245ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 246e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 247ae7cfcebSSatish Balay #else 24871044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 249e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 250ae7cfcebSSatish Balay #endif 2517a97a34bSBarry Smith } else { 252ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 253e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 254ae7cfcebSSatish Balay #else 25571044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 256e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 257ae7cfcebSSatish Balay #endif 258da3a660dSBarry Smith } 2591ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2601ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 261dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 2623a40ed3dSBarry Smith PetscFunctionReturn(0); 263da3a660dSBarry Smith } 2646ee01492SSatish Balay 2654a2ae208SSatish Balay #undef __FUNCT__ 2664a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 267dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 268da3a660dSBarry Smith { 269c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 270dfbe8321SBarry Smith PetscErrorCode ierr; 27187828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 272da3a660dSBarry Smith Vec tmp = 0; 273d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 27467e560aaSBarry Smith 2753a40ed3dSBarry Smith PetscFunctionBegin; 2761ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2771ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 278d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 279da3a660dSBarry Smith if (yy == zz) { 28078b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 28152e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 28278b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 283da3a660dSBarry Smith } 284d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 285752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 286da3a660dSBarry Smith if (mat->pivots) { 287ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 288e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 289ae7cfcebSSatish Balay #else 29071044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 291e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 292ae7cfcebSSatish Balay #endif 293a8c6a408SBarry Smith } else { 294ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 295e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 296ae7cfcebSSatish Balay #else 29771044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 298e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 299ae7cfcebSSatish Balay #endif 300da3a660dSBarry Smith } 3012dcb1b2aSMatthew Knepley if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 3022dcb1b2aSMatthew Knepley else {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);} 3031ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3041ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 305dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 3063a40ed3dSBarry Smith PetscFunctionReturn(0); 307da3a660dSBarry Smith } 30867e560aaSBarry Smith 3094a2ae208SSatish Balay #undef __FUNCT__ 3104a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 311dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 312da3a660dSBarry Smith { 313c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3146849ba73SBarry Smith PetscErrorCode ierr; 31587828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 316da3a660dSBarry Smith Vec tmp; 317d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 31867e560aaSBarry Smith 3193a40ed3dSBarry Smith PetscFunctionBegin; 320d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 3211ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3221ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 323da3a660dSBarry Smith if (yy == zz) { 32478b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 32552e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 32678b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 327da3a660dSBarry Smith } 328d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 329752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 330da3a660dSBarry Smith if (mat->pivots) { 331ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 332e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 333ae7cfcebSSatish Balay #else 33471044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 335e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 336ae7cfcebSSatish Balay #endif 3373a40ed3dSBarry Smith } else { 338ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 339e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 340ae7cfcebSSatish Balay #else 34171044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 342e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 343ae7cfcebSSatish Balay #endif 344da3a660dSBarry Smith } 34590f02eecSBarry Smith if (tmp) { 3462dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 34790f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 3483a40ed3dSBarry Smith } else { 3492dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 35090f02eecSBarry Smith } 3511ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3521ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 353dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 3543a40ed3dSBarry Smith PetscFunctionReturn(0); 355da3a660dSBarry Smith } 356db4efbfdSBarry Smith 357db4efbfdSBarry Smith /* ---------------------------------------------------------------*/ 358db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 359db4efbfdSBarry Smith rather than put it in the Mat->row slot.*/ 360db4efbfdSBarry Smith #undef __FUNCT__ 361db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense" 3620481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo) 363db4efbfdSBarry Smith { 364db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 365db4efbfdSBarry Smith PetscFunctionBegin; 366e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 367db4efbfdSBarry Smith #else 368db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 369db4efbfdSBarry Smith PetscErrorCode ierr; 370db4efbfdSBarry Smith PetscBLASInt n,m,info; 371db4efbfdSBarry Smith 372db4efbfdSBarry Smith PetscFunctionBegin; 373db4efbfdSBarry Smith n = PetscBLASIntCast(A->cmap->n); 374db4efbfdSBarry Smith m = PetscBLASIntCast(A->rmap->n); 375db4efbfdSBarry Smith if (!mat->pivots) { 376db4efbfdSBarry Smith ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr); 377db4efbfdSBarry Smith ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr); 378db4efbfdSBarry Smith } 379db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 380db4efbfdSBarry Smith LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info); 381e32f2f54SBarry Smith if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization"); 382e32f2f54SBarry Smith if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 383db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 384db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 385db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 386db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 387d5f3da31SBarry Smith A->factortype = MAT_FACTOR_LU; 388db4efbfdSBarry Smith 389dc0b31edSSatish Balay ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr); 390db4efbfdSBarry Smith #endif 391db4efbfdSBarry Smith PetscFunctionReturn(0); 392db4efbfdSBarry Smith } 393db4efbfdSBarry Smith 394db4efbfdSBarry Smith #undef __FUNCT__ 395db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense" 3960481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo) 397db4efbfdSBarry Smith { 398db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 399db4efbfdSBarry Smith PetscFunctionBegin; 400e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 401db4efbfdSBarry Smith #else 402db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 403db4efbfdSBarry Smith PetscErrorCode ierr; 404db4efbfdSBarry Smith PetscBLASInt info,n = PetscBLASIntCast(A->cmap->n); 405db4efbfdSBarry Smith 406db4efbfdSBarry Smith PetscFunctionBegin; 407db4efbfdSBarry Smith ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 408db4efbfdSBarry Smith 409db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 410db4efbfdSBarry Smith LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info); 411e32f2f54SBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 412db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 413db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 414db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 415db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 416d5f3da31SBarry Smith A->factortype = MAT_FACTOR_CHOLESKY; 417dc0b31edSSatish Balay ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr); 418db4efbfdSBarry Smith #endif 419db4efbfdSBarry Smith PetscFunctionReturn(0); 420db4efbfdSBarry Smith } 421db4efbfdSBarry Smith 422db4efbfdSBarry Smith 423db4efbfdSBarry Smith #undef __FUNCT__ 424db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 4250481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 426db4efbfdSBarry Smith { 427db4efbfdSBarry Smith PetscErrorCode ierr; 428db4efbfdSBarry Smith MatFactorInfo info; 429db4efbfdSBarry Smith 430db4efbfdSBarry Smith PetscFunctionBegin; 431db4efbfdSBarry Smith info.fill = 1.0; 432c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 433719d5645SBarry Smith ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr); 434db4efbfdSBarry Smith PetscFunctionReturn(0); 435db4efbfdSBarry Smith } 436db4efbfdSBarry Smith 437db4efbfdSBarry Smith #undef __FUNCT__ 438db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 4390481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info) 440db4efbfdSBarry Smith { 441db4efbfdSBarry Smith PetscFunctionBegin; 442c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 443719d5645SBarry Smith fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense; 444db4efbfdSBarry Smith PetscFunctionReturn(0); 445db4efbfdSBarry Smith } 446db4efbfdSBarry Smith 447db4efbfdSBarry Smith #undef __FUNCT__ 448db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 4490481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 450db4efbfdSBarry Smith { 451db4efbfdSBarry Smith PetscFunctionBegin; 452c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 453719d5645SBarry Smith fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense; 454db4efbfdSBarry Smith PetscFunctionReturn(0); 455db4efbfdSBarry Smith } 456db4efbfdSBarry Smith 457bb5747d9SMatthew Knepley EXTERN_C_BEGIN 458db4efbfdSBarry Smith #undef __FUNCT__ 459db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc" 460db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact) 461db4efbfdSBarry Smith { 462db4efbfdSBarry Smith PetscErrorCode ierr; 463db4efbfdSBarry Smith 464db4efbfdSBarry Smith PetscFunctionBegin; 465db4efbfdSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr); 466db4efbfdSBarry Smith ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 467db4efbfdSBarry Smith ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr); 468db4efbfdSBarry Smith if (ftype == MAT_FACTOR_LU){ 469db4efbfdSBarry Smith (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense; 470db4efbfdSBarry Smith } else { 471db4efbfdSBarry Smith (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense; 472db4efbfdSBarry Smith } 473d5f3da31SBarry Smith (*fact)->factortype = ftype; 474db4efbfdSBarry Smith PetscFunctionReturn(0); 475db4efbfdSBarry Smith } 476bb5747d9SMatthew Knepley EXTERN_C_END 477db4efbfdSBarry Smith 478289bc588SBarry Smith /* ------------------------------------------------------------------*/ 4794a2ae208SSatish Balay #undef __FUNCT__ 48041f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense" 48141f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 482289bc588SBarry Smith { 483c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 48487828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 485dfbe8321SBarry Smith PetscErrorCode ierr; 486d0f46423SBarry Smith PetscInt m = A->rmap->n,i; 487aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 4880805154bSBarry Smith PetscBLASInt o = 1,bm = PetscBLASIntCast(m); 489bc1b551cSSatish Balay #endif 490289bc588SBarry Smith 4913a40ed3dSBarry Smith PetscFunctionBegin; 492289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 49371044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 4942dcb1b2aSMatthew Knepley ierr = VecSet(xx,zero);CHKERRQ(ierr); 495289bc588SBarry Smith } 4961ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4971ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 498b965ef7fSBarry Smith its = its*lits; 499e32f2f54SBarry 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); 500289bc588SBarry Smith while (its--) { 501fccaa45eSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 502289bc588SBarry Smith for (i=0; i<m; i++) { 503aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 504f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 505f1747703SBarry Smith not happy about returning a double complex */ 50613f74950SBarry Smith PetscInt _i; 50787828ca2SBarry Smith PetscScalar sum = b[i]; 508f1747703SBarry Smith for (_i=0; _i<m; _i++) { 5093f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 510f1747703SBarry Smith } 511f1747703SBarry Smith xt = sum; 512f1747703SBarry Smith #else 51371044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 514f1747703SBarry Smith #endif 51555a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 516289bc588SBarry Smith } 517289bc588SBarry Smith } 518fccaa45eSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 519289bc588SBarry Smith for (i=m-1; i>=0; i--) { 520aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 521f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 522f1747703SBarry Smith not happy about returning a double complex */ 52313f74950SBarry Smith PetscInt _i; 52487828ca2SBarry Smith PetscScalar sum = b[i]; 525f1747703SBarry Smith for (_i=0; _i<m; _i++) { 5263f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 527f1747703SBarry Smith } 528f1747703SBarry Smith xt = sum; 529f1747703SBarry Smith #else 53071044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 531f1747703SBarry Smith #endif 53255a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 533289bc588SBarry Smith } 534289bc588SBarry Smith } 535289bc588SBarry Smith } 5361ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 5371ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5383a40ed3dSBarry Smith PetscFunctionReturn(0); 539289bc588SBarry Smith } 540289bc588SBarry Smith 541289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5424a2ae208SSatish Balay #undef __FUNCT__ 5434a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 544dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 545289bc588SBarry Smith { 546c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 54787828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 548dfbe8321SBarry Smith PetscErrorCode ierr; 5490805154bSBarry Smith PetscBLASInt m, n,_One=1; 550ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 5513a40ed3dSBarry Smith 5523a40ed3dSBarry Smith PetscFunctionBegin; 553d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 554d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 555d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5561ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5571ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 55871044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One); 5591ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5601ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 561dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 5623a40ed3dSBarry Smith PetscFunctionReturn(0); 563289bc588SBarry Smith } 564800995b7SMatthew Knepley 5654a2ae208SSatish Balay #undef __FUNCT__ 5664a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 567dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 568289bc588SBarry Smith { 569c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 57087828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 571dfbe8321SBarry Smith PetscErrorCode ierr; 5720805154bSBarry Smith PetscBLASInt m, n, _One=1; 5733a40ed3dSBarry Smith 5743a40ed3dSBarry Smith PetscFunctionBegin; 575d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 576d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 577d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5781ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5791ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 58071044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 5811ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5821ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 583dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr); 5843a40ed3dSBarry Smith PetscFunctionReturn(0); 585289bc588SBarry Smith } 5866ee01492SSatish Balay 5874a2ae208SSatish Balay #undef __FUNCT__ 5884a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 589dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 590289bc588SBarry Smith { 591c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 59287828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 593dfbe8321SBarry Smith PetscErrorCode ierr; 5940805154bSBarry Smith PetscBLASInt m, n, _One=1; 5953a40ed3dSBarry Smith 5963a40ed3dSBarry Smith PetscFunctionBegin; 597d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 598d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 599d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 600600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6011ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6021ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 60371044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 6041ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6051ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 606dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6073a40ed3dSBarry Smith PetscFunctionReturn(0); 608289bc588SBarry Smith } 6096ee01492SSatish Balay 6104a2ae208SSatish Balay #undef __FUNCT__ 6114a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 612dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 613289bc588SBarry Smith { 614c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 61587828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 616dfbe8321SBarry Smith PetscErrorCode ierr; 6170805154bSBarry Smith PetscBLASInt m, n, _One=1; 61887828ca2SBarry Smith PetscScalar _DOne=1.0; 6193a40ed3dSBarry Smith 6203a40ed3dSBarry Smith PetscFunctionBegin; 621d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 622d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 623d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 624600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6251ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6261ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 62771044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 6281ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6291ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 630dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6313a40ed3dSBarry Smith PetscFunctionReturn(0); 632289bc588SBarry Smith } 633289bc588SBarry Smith 634289bc588SBarry Smith /* -----------------------------------------------------------------*/ 6354a2ae208SSatish Balay #undef __FUNCT__ 6364a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 63713f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 638289bc588SBarry Smith { 639c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 64087828ca2SBarry Smith PetscScalar *v; 6416849ba73SBarry Smith PetscErrorCode ierr; 64213f74950SBarry Smith PetscInt i; 64367e560aaSBarry Smith 6443a40ed3dSBarry Smith PetscFunctionBegin; 645d0f46423SBarry Smith *ncols = A->cmap->n; 646289bc588SBarry Smith if (cols) { 647d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr); 648d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) (*cols)[i] = i; 649289bc588SBarry Smith } 650289bc588SBarry Smith if (vals) { 651d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 652289bc588SBarry Smith v = mat->v + row; 653d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;} 654289bc588SBarry Smith } 6553a40ed3dSBarry Smith PetscFunctionReturn(0); 656289bc588SBarry Smith } 6576ee01492SSatish Balay 6584a2ae208SSatish Balay #undef __FUNCT__ 6594a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 66013f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 661289bc588SBarry Smith { 662dfbe8321SBarry Smith PetscErrorCode ierr; 663606d414cSSatish Balay PetscFunctionBegin; 664606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 665606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 6663a40ed3dSBarry Smith PetscFunctionReturn(0); 667289bc588SBarry Smith } 668289bc588SBarry Smith /* ----------------------------------------------------------------*/ 6694a2ae208SSatish Balay #undef __FUNCT__ 6704a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 67113f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 672289bc588SBarry Smith { 673c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 67413f74950SBarry Smith PetscInt i,j,idx=0; 675d6dfbf8fSBarry Smith 6763a40ed3dSBarry Smith PetscFunctionBegin; 67771fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 678289bc588SBarry Smith if (!mat->roworiented) { 679dbb450caSBarry Smith if (addv == INSERT_VALUES) { 680289bc588SBarry Smith for (j=0; j<n; j++) { 681cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6822515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 683e32f2f54SBarry 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); 68458804f6dSBarry Smith #endif 685289bc588SBarry Smith for (i=0; i<m; i++) { 686cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6872515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 688e32f2f54SBarry 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); 68958804f6dSBarry Smith #endif 690cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 691289bc588SBarry Smith } 692289bc588SBarry Smith } 6933a40ed3dSBarry Smith } else { 694289bc588SBarry Smith for (j=0; j<n; j++) { 695cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6962515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 697e32f2f54SBarry 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); 69858804f6dSBarry Smith #endif 699289bc588SBarry Smith for (i=0; i<m; i++) { 700cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 7012515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 702e32f2f54SBarry 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); 70358804f6dSBarry Smith #endif 704cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 705289bc588SBarry Smith } 706289bc588SBarry Smith } 707289bc588SBarry Smith } 7083a40ed3dSBarry Smith } else { 709dbb450caSBarry Smith if (addv == INSERT_VALUES) { 710e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 711cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7122515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 713e32f2f54SBarry 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); 71458804f6dSBarry Smith #endif 715e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 716cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7172515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 718e32f2f54SBarry 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); 71958804f6dSBarry Smith #endif 720cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 721e8d4e0b9SBarry Smith } 722e8d4e0b9SBarry Smith } 7233a40ed3dSBarry Smith } else { 724289bc588SBarry Smith for (i=0; i<m; i++) { 725cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7262515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 727e32f2f54SBarry 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); 72858804f6dSBarry Smith #endif 729289bc588SBarry Smith for (j=0; j<n; j++) { 730cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7312515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 732e32f2f54SBarry 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); 73358804f6dSBarry Smith #endif 734cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 735289bc588SBarry Smith } 736289bc588SBarry Smith } 737289bc588SBarry Smith } 738e8d4e0b9SBarry Smith } 7393a40ed3dSBarry Smith PetscFunctionReturn(0); 740289bc588SBarry Smith } 741e8d4e0b9SBarry Smith 7424a2ae208SSatish Balay #undef __FUNCT__ 7434a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 74413f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 745ae80bb75SLois Curfman McInnes { 746ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 74713f74950SBarry Smith PetscInt i,j; 748ae80bb75SLois Curfman McInnes 7493a40ed3dSBarry Smith PetscFunctionBegin; 750ae80bb75SLois Curfman McInnes /* row-oriented output */ 751ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 75297e567efSBarry Smith if (indexm[i] < 0) {v += n;continue;} 753e32f2f54SBarry 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); 754ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 7556f31f424SBarry Smith if (indexn[j] < 0) {v++; continue;} 756e32f2f54SBarry 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); 75797e567efSBarry Smith *v++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 758ae80bb75SLois Curfman McInnes } 759ae80bb75SLois Curfman McInnes } 7603a40ed3dSBarry Smith PetscFunctionReturn(0); 761ae80bb75SLois Curfman McInnes } 762ae80bb75SLois Curfman McInnes 763289bc588SBarry Smith /* -----------------------------------------------------------------*/ 764289bc588SBarry Smith 7654a2ae208SSatish Balay #undef __FUNCT__ 7665bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_SeqDense" 767112444f4SShri Abhyankar PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer) 768aabbc4fbSShri Abhyankar { 769aabbc4fbSShri Abhyankar Mat_SeqDense *a; 770aabbc4fbSShri Abhyankar PetscErrorCode ierr; 771aabbc4fbSShri Abhyankar PetscInt *scols,i,j,nz,header[4]; 772aabbc4fbSShri Abhyankar int fd; 773aabbc4fbSShri Abhyankar PetscMPIInt size; 774aabbc4fbSShri Abhyankar PetscInt *rowlengths = 0,M,N,*cols,grows,gcols; 775aabbc4fbSShri Abhyankar PetscScalar *vals,*svals,*v,*w; 776aabbc4fbSShri Abhyankar MPI_Comm comm = ((PetscObject)viewer)->comm; 777aabbc4fbSShri Abhyankar 778aabbc4fbSShri Abhyankar PetscFunctionBegin; 779aabbc4fbSShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 780aabbc4fbSShri Abhyankar if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor"); 781aabbc4fbSShri Abhyankar ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 782aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 783aabbc4fbSShri Abhyankar if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 784aabbc4fbSShri Abhyankar M = header[1]; N = header[2]; nz = header[3]; 785aabbc4fbSShri Abhyankar 786aabbc4fbSShri Abhyankar /* set global size if not set already*/ 787aabbc4fbSShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 788aabbc4fbSShri Abhyankar ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr); 789aabbc4fbSShri Abhyankar } else { 790aabbc4fbSShri Abhyankar /* if sizes and type are already set, check if the vector global sizes are correct */ 791aabbc4fbSShri Abhyankar ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr); 792aabbc4fbSShri 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); 793aabbc4fbSShri Abhyankar } 794aabbc4fbSShri Abhyankar ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 795aabbc4fbSShri Abhyankar 796aabbc4fbSShri Abhyankar if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 797aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 798aabbc4fbSShri Abhyankar v = a->v; 799aabbc4fbSShri Abhyankar /* Allocate some temp space to read in the values and then flip them 800aabbc4fbSShri Abhyankar from row major to column major */ 801aabbc4fbSShri Abhyankar ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 802aabbc4fbSShri Abhyankar /* read in nonzero values */ 803aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 804aabbc4fbSShri Abhyankar /* now flip the values and store them in the matrix*/ 805aabbc4fbSShri Abhyankar for (j=0; j<N; j++) { 806aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 807aabbc4fbSShri Abhyankar *v++ =w[i*N+j]; 808aabbc4fbSShri Abhyankar } 809aabbc4fbSShri Abhyankar } 810aabbc4fbSShri Abhyankar ierr = PetscFree(w);CHKERRQ(ierr); 811aabbc4fbSShri Abhyankar } else { 812aabbc4fbSShri Abhyankar /* read row lengths */ 813aabbc4fbSShri Abhyankar ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 814aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 815aabbc4fbSShri Abhyankar 816aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 817aabbc4fbSShri Abhyankar v = a->v; 818aabbc4fbSShri Abhyankar 819aabbc4fbSShri Abhyankar /* read column indices and nonzeros */ 820aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr); 821aabbc4fbSShri Abhyankar cols = scols; 822aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 823aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 824aabbc4fbSShri Abhyankar vals = svals; 825aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 826aabbc4fbSShri Abhyankar 827aabbc4fbSShri Abhyankar /* insert into matrix */ 828aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 829aabbc4fbSShri Abhyankar for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 830aabbc4fbSShri Abhyankar svals += rowlengths[i]; scols += rowlengths[i]; 831aabbc4fbSShri Abhyankar } 832aabbc4fbSShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 833aabbc4fbSShri Abhyankar ierr = PetscFree(cols);CHKERRQ(ierr); 834aabbc4fbSShri Abhyankar ierr = PetscFree(rowlengths);CHKERRQ(ierr); 835aabbc4fbSShri Abhyankar } 836aabbc4fbSShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 837aabbc4fbSShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 838aabbc4fbSShri Abhyankar 839aabbc4fbSShri Abhyankar PetscFunctionReturn(0); 840aabbc4fbSShri Abhyankar } 841aabbc4fbSShri Abhyankar 842aabbc4fbSShri Abhyankar #undef __FUNCT__ 8434a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 8446849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 845289bc588SBarry Smith { 846932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 847dfbe8321SBarry Smith PetscErrorCode ierr; 84813f74950SBarry Smith PetscInt i,j; 8492dcb1b2aSMatthew Knepley const char *name; 85087828ca2SBarry Smith PetscScalar *v; 851f3ef73ceSBarry Smith PetscViewerFormat format; 8525f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX) 853ace3abfcSBarry Smith PetscBool allreal = PETSC_TRUE; 8545f481a85SSatish Balay #endif 855932b0c3eSLois Curfman McInnes 8563a40ed3dSBarry Smith PetscFunctionBegin; 857b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 858456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 8593a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 860fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 861d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 8627566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 863d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 86444cd7ae7SLois Curfman McInnes v = a->v + i; 86577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 866d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 867aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 868329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 869a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 870329f5518SBarry Smith } else if (PetscRealPart(*v)) { 871a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr); 8726831982aSBarry Smith } 87380cd9d93SLois Curfman McInnes #else 8746831982aSBarry Smith if (*v) { 875a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr); 8766831982aSBarry Smith } 87780cd9d93SLois Curfman McInnes #endif 8781b807ce4Svictorle v += a->lda; 87980cd9d93SLois Curfman McInnes } 880b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 88180cd9d93SLois Curfman McInnes } 882d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 8833a40ed3dSBarry Smith } else { 884d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 885aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 88647989497SBarry Smith /* determine if matrix has all real values */ 88747989497SBarry Smith v = a->v; 888d0f46423SBarry Smith for (i=0; i<A->rmap->n*A->cmap->n; i++) { 889ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 89047989497SBarry Smith } 89147989497SBarry Smith #endif 892fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 8933a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 894d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr); 895d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 896fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 8977566de4bSShri Abhyankar } else { 8987566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 899ffac6cdbSBarry Smith } 900ffac6cdbSBarry Smith 901d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 902932b0c3eSLois Curfman McInnes v = a->v + i; 903d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 904aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 90547989497SBarry Smith if (allreal) { 906f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr); 90747989497SBarry Smith } else { 908f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 90947989497SBarry Smith } 910289bc588SBarry Smith #else 911f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr); 912289bc588SBarry Smith #endif 9131b807ce4Svictorle v += a->lda; 914289bc588SBarry Smith } 915b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 916289bc588SBarry Smith } 917fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 918b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 919ffac6cdbSBarry Smith } 920d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 921da3a660dSBarry Smith } 922b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 9233a40ed3dSBarry Smith PetscFunctionReturn(0); 924289bc588SBarry Smith } 925289bc588SBarry Smith 9264a2ae208SSatish Balay #undef __FUNCT__ 9274a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 9286849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 929932b0c3eSLois Curfman McInnes { 930932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9316849ba73SBarry Smith PetscErrorCode ierr; 93213f74950SBarry Smith int fd; 933d0f46423SBarry Smith PetscInt ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n; 934f4403165SShri Abhyankar PetscScalar *v,*anonz,*vals; 935f4403165SShri Abhyankar PetscViewerFormat format; 936932b0c3eSLois Curfman McInnes 9373a40ed3dSBarry Smith PetscFunctionBegin; 938b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 93990ace30eSBarry Smith 940f4403165SShri Abhyankar ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 941f4403165SShri Abhyankar if (format == PETSC_VIEWER_NATIVE) { 942f4403165SShri Abhyankar /* store the matrix as a dense matrix */ 943f4403165SShri Abhyankar ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 944f4403165SShri Abhyankar col_lens[0] = MAT_FILE_CLASSID; 945f4403165SShri Abhyankar col_lens[1] = m; 946f4403165SShri Abhyankar col_lens[2] = n; 947f4403165SShri Abhyankar col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 948f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 949f4403165SShri Abhyankar ierr = PetscFree(col_lens);CHKERRQ(ierr); 950f4403165SShri Abhyankar 951f4403165SShri Abhyankar /* write out matrix, by rows */ 952f4403165SShri Abhyankar ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 953f4403165SShri Abhyankar v = a->v; 954f4403165SShri Abhyankar for (j=0; j<n; j++) { 955f4403165SShri Abhyankar for (i=0; i<m; i++) { 956f4403165SShri Abhyankar vals[j + i*n] = *v++; 957f4403165SShri Abhyankar } 958f4403165SShri Abhyankar } 959f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 960f4403165SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 961f4403165SShri Abhyankar } else { 96213f74950SBarry Smith ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 9630700a824SBarry Smith col_lens[0] = MAT_FILE_CLASSID; 964932b0c3eSLois Curfman McInnes col_lens[1] = m; 965932b0c3eSLois Curfman McInnes col_lens[2] = n; 966932b0c3eSLois Curfman McInnes col_lens[3] = nz; 967932b0c3eSLois Curfman McInnes 968932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 969932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 9706f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 971932b0c3eSLois Curfman McInnes 972932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 973932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 974932b0c3eSLois Curfman McInnes ict = 0; 975932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 976932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 977932b0c3eSLois Curfman McInnes } 9786f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 979606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 980932b0c3eSLois Curfman McInnes 981932b0c3eSLois Curfman McInnes /* store nonzero values */ 98287828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 983932b0c3eSLois Curfman McInnes ict = 0; 984932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 985932b0c3eSLois Curfman McInnes v = a->v + i; 986932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 9871b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 988932b0c3eSLois Curfman McInnes } 989932b0c3eSLois Curfman McInnes } 9906f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 991606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 992f4403165SShri Abhyankar } 9933a40ed3dSBarry Smith PetscFunctionReturn(0); 994932b0c3eSLois Curfman McInnes } 995932b0c3eSLois Curfman McInnes 9964a2ae208SSatish Balay #undef __FUNCT__ 9974a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 998dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 999f1af5d2fSBarry Smith { 1000f1af5d2fSBarry Smith Mat A = (Mat) Aa; 1001f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 10026849ba73SBarry Smith PetscErrorCode ierr; 1003d0f46423SBarry Smith PetscInt m = A->rmap->n,n = A->cmap->n,color,i,j; 100487828ca2SBarry Smith PetscScalar *v = a->v; 1005b0a32e0cSBarry Smith PetscViewer viewer; 1006b0a32e0cSBarry Smith PetscDraw popup; 1007329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 1008f3ef73ceSBarry Smith PetscViewerFormat format; 1009f1af5d2fSBarry Smith 1010f1af5d2fSBarry Smith PetscFunctionBegin; 1011f1af5d2fSBarry Smith 1012f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 1013b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1014b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1015f1af5d2fSBarry Smith 1016f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 1017fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 1018f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 1019b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1020f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1021f1af5d2fSBarry Smith x_l = j; 1022f1af5d2fSBarry Smith x_r = x_l + 1.0; 1023f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1024f1af5d2fSBarry Smith y_l = m - i - 1.0; 1025f1af5d2fSBarry Smith y_r = y_l + 1.0; 1026f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 1027329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 1028b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1029329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 1030b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1031f1af5d2fSBarry Smith } else { 1032f1af5d2fSBarry Smith continue; 1033f1af5d2fSBarry Smith } 1034f1af5d2fSBarry Smith #else 1035f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 1036b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1037f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 1038b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1039f1af5d2fSBarry Smith } else { 1040f1af5d2fSBarry Smith continue; 1041f1af5d2fSBarry Smith } 1042f1af5d2fSBarry Smith #endif 1043b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1044f1af5d2fSBarry Smith } 1045f1af5d2fSBarry Smith } 1046f1af5d2fSBarry Smith } else { 1047f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 1048f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 1049f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 1050f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 1051f1af5d2fSBarry Smith } 1052b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 1053b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 1054b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 1055f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1056f1af5d2fSBarry Smith x_l = j; 1057f1af5d2fSBarry Smith x_r = x_l + 1.0; 1058f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1059f1af5d2fSBarry Smith y_l = m - i - 1.0; 1060f1af5d2fSBarry Smith y_r = y_l + 1.0; 1061b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 1062b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1063f1af5d2fSBarry Smith } 1064f1af5d2fSBarry Smith } 1065f1af5d2fSBarry Smith } 1066f1af5d2fSBarry Smith PetscFunctionReturn(0); 1067f1af5d2fSBarry Smith } 1068f1af5d2fSBarry Smith 10694a2ae208SSatish Balay #undef __FUNCT__ 10704a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 1071dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 1072f1af5d2fSBarry Smith { 1073b0a32e0cSBarry Smith PetscDraw draw; 1074ace3abfcSBarry Smith PetscBool isnull; 1075329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 1076dfbe8321SBarry Smith PetscErrorCode ierr; 1077f1af5d2fSBarry Smith 1078f1af5d2fSBarry Smith PetscFunctionBegin; 1079b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1080b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1081abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 1082f1af5d2fSBarry Smith 1083f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1084d0f46423SBarry Smith xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0; 1085f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 1086b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1087b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 1088f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 1089f1af5d2fSBarry Smith PetscFunctionReturn(0); 1090f1af5d2fSBarry Smith } 1091f1af5d2fSBarry Smith 10924a2ae208SSatish Balay #undef __FUNCT__ 10934a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 1094dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1095932b0c3eSLois Curfman McInnes { 1096dfbe8321SBarry Smith PetscErrorCode ierr; 1097ace3abfcSBarry Smith PetscBool iascii,isbinary,isdraw; 1098932b0c3eSLois Curfman McInnes 10993a40ed3dSBarry Smith PetscFunctionBegin; 11002692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 11012692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 11022692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 11030f5bd95cSBarry Smith 1104c45a1595SBarry Smith if (iascii) { 1105c45a1595SBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 11060f5bd95cSBarry Smith } else if (isbinary) { 11073a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1108f1af5d2fSBarry Smith } else if (isdraw) { 1109f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 11105cd90555SBarry Smith } else { 1111e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1112932b0c3eSLois Curfman McInnes } 11133a40ed3dSBarry Smith PetscFunctionReturn(0); 1114932b0c3eSLois Curfman McInnes } 1115289bc588SBarry Smith 11164a2ae208SSatish Balay #undef __FUNCT__ 11174a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1118dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1119289bc588SBarry Smith { 1120ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1121dfbe8321SBarry Smith PetscErrorCode ierr; 112290f02eecSBarry Smith 11233a40ed3dSBarry Smith PetscFunctionBegin; 1124aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1125d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n); 1126a5a9c739SBarry Smith #endif 112705b42c5fSBarry Smith ierr = PetscFree(l->pivots);CHKERRQ(ierr); 11286857c123SSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1129606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 1130dbd8c25aSHong Zhang 1131dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 1132901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 11334ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11344ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11354ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11363a40ed3dSBarry Smith PetscFunctionReturn(0); 1137289bc588SBarry Smith } 1138289bc588SBarry Smith 11394a2ae208SSatish Balay #undef __FUNCT__ 11404a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1141fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout) 1142289bc588SBarry Smith { 1143c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11446849ba73SBarry Smith PetscErrorCode ierr; 114513f74950SBarry Smith PetscInt k,j,m,n,M; 114687828ca2SBarry Smith PetscScalar *v,tmp; 114748b35521SBarry Smith 11483a40ed3dSBarry Smith PetscFunctionBegin; 1149d0f46423SBarry Smith v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n; 1150e9695a30SBarry Smith if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */ 1151e7e72b3dSBarry Smith if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1152e7e72b3dSBarry Smith else { 1153d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1154289bc588SBarry Smith for (k=0; k<j; k++) { 11551b807ce4Svictorle tmp = v[j + k*M]; 11561b807ce4Svictorle v[j + k*M] = v[k + j*M]; 11571b807ce4Svictorle v[k + j*M] = tmp; 1158289bc588SBarry Smith } 1159289bc588SBarry Smith } 1160d64ed03dSBarry Smith } 11613a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1162d3e5ee88SLois Curfman McInnes Mat tmat; 1163ec8511deSBarry Smith Mat_SeqDense *tmatd; 116487828ca2SBarry Smith PetscScalar *v2; 1165ea709b57SSatish Balay 1166fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 11677adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr); 1168d0f46423SBarry Smith ierr = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr); 11697adad957SLisandro Dalcin ierr = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 11705c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1171fc4dec0aSBarry Smith } else { 1172fc4dec0aSBarry Smith tmat = *matout; 1173fc4dec0aSBarry Smith } 1174ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 11750de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1176d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 11771b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1178d3e5ee88SLois Curfman McInnes } 11796d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11806d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1181d3e5ee88SLois Curfman McInnes *matout = tmat; 118248b35521SBarry Smith } 11833a40ed3dSBarry Smith PetscFunctionReturn(0); 1184289bc588SBarry Smith } 1185289bc588SBarry Smith 11864a2ae208SSatish Balay #undef __FUNCT__ 11874a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1188ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool *flg) 1189289bc588SBarry Smith { 1190c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1191c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 119213f74950SBarry Smith PetscInt i,j; 119387828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 11949ea5d5aeSSatish Balay 11953a40ed3dSBarry Smith PetscFunctionBegin; 1196d0f46423SBarry Smith if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1197d0f46423SBarry Smith if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1198d0f46423SBarry Smith for (i=0; i<A1->rmap->n; i++) { 11991b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1200d0f46423SBarry Smith for (j=0; j<A1->cmap->n; j++) { 12013a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 12021b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 12031b807ce4Svictorle } 1204289bc588SBarry Smith } 120577c4ece6SBarry Smith *flg = PETSC_TRUE; 12063a40ed3dSBarry Smith PetscFunctionReturn(0); 1207289bc588SBarry Smith } 1208289bc588SBarry Smith 12094a2ae208SSatish Balay #undef __FUNCT__ 12104a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1211dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1212289bc588SBarry Smith { 1213c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1214dfbe8321SBarry Smith PetscErrorCode ierr; 121513f74950SBarry Smith PetscInt i,n,len; 121687828ca2SBarry Smith PetscScalar *x,zero = 0.0; 121744cd7ae7SLois Curfman McInnes 12183a40ed3dSBarry Smith PetscFunctionBegin; 12192dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 12207a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 12211ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1222d0f46423SBarry Smith len = PetscMin(A->rmap->n,A->cmap->n); 1223e32f2f54SBarry Smith if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 122444cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 12251b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1226289bc588SBarry Smith } 12271ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 12283a40ed3dSBarry Smith PetscFunctionReturn(0); 1229289bc588SBarry Smith } 1230289bc588SBarry Smith 12314a2ae208SSatish Balay #undef __FUNCT__ 12324a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1233dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1234289bc588SBarry Smith { 1235c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 123687828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1237dfbe8321SBarry Smith PetscErrorCode ierr; 1238d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n; 123955659b69SBarry Smith 12403a40ed3dSBarry Smith PetscFunctionBegin; 124128988994SBarry Smith if (ll) { 12427a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 12431ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1244e32f2f54SBarry Smith if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1245da3a660dSBarry Smith for (i=0; i<m; i++) { 1246da3a660dSBarry Smith x = l[i]; 1247da3a660dSBarry Smith v = mat->v + i; 1248da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1249da3a660dSBarry Smith } 12501ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1251efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1252da3a660dSBarry Smith } 125328988994SBarry Smith if (rr) { 12547a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 12551ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1256e32f2f54SBarry Smith if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1257da3a660dSBarry Smith for (i=0; i<n; i++) { 1258da3a660dSBarry Smith x = r[i]; 1259da3a660dSBarry Smith v = mat->v + i*m; 1260da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1261da3a660dSBarry Smith } 12621ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1263efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1264da3a660dSBarry Smith } 12653a40ed3dSBarry Smith PetscFunctionReturn(0); 1266289bc588SBarry Smith } 1267289bc588SBarry Smith 12684a2ae208SSatish Balay #undef __FUNCT__ 12694a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1270dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1271289bc588SBarry Smith { 1272c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 127387828ca2SBarry Smith PetscScalar *v = mat->v; 1274329f5518SBarry Smith PetscReal sum = 0.0; 1275d0f46423SBarry Smith PetscInt lda=mat->lda,m=A->rmap->n,i,j; 1276efee365bSSatish Balay PetscErrorCode ierr; 127755659b69SBarry Smith 12783a40ed3dSBarry Smith PetscFunctionBegin; 1279289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1280a5ce6ee0Svictorle if (lda>m) { 1281d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1282a5ce6ee0Svictorle v = mat->v+j*lda; 1283a5ce6ee0Svictorle for (i=0; i<m; i++) { 1284a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1285a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1286a5ce6ee0Svictorle #else 1287a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1288a5ce6ee0Svictorle #endif 1289a5ce6ee0Svictorle } 1290a5ce6ee0Svictorle } 1291a5ce6ee0Svictorle } else { 1292d0f46423SBarry Smith for (i=0; i<A->cmap->n*A->rmap->n; i++) { 1293aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1294329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1295289bc588SBarry Smith #else 1296289bc588SBarry Smith sum += (*v)*(*v); v++; 1297289bc588SBarry Smith #endif 1298289bc588SBarry Smith } 1299a5ce6ee0Svictorle } 1300064f8208SBarry Smith *nrm = sqrt(sum); 1301dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13023a40ed3dSBarry Smith } else if (type == NORM_1) { 1303064f8208SBarry Smith *nrm = 0.0; 1304d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 13051b807ce4Svictorle v = mat->v + j*mat->lda; 1306289bc588SBarry Smith sum = 0.0; 1307d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 130833a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1309289bc588SBarry Smith } 1310064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1311289bc588SBarry Smith } 1312d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13133a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1314064f8208SBarry Smith *nrm = 0.0; 1315d0f46423SBarry Smith for (j=0; j<A->rmap->n; j++) { 1316289bc588SBarry Smith v = mat->v + j; 1317289bc588SBarry Smith sum = 0.0; 1318d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) { 13191b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1320289bc588SBarry Smith } 1321064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1322289bc588SBarry Smith } 1323d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 1324e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 13253a40ed3dSBarry Smith PetscFunctionReturn(0); 1326289bc588SBarry Smith } 1327289bc588SBarry Smith 13284a2ae208SSatish Balay #undef __FUNCT__ 13294a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 1330ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg) 1331289bc588SBarry Smith { 1332c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 133363ba0a88SBarry Smith PetscErrorCode ierr; 133467e560aaSBarry Smith 13353a40ed3dSBarry Smith PetscFunctionBegin; 1336b5a2b587SKris Buschelman switch (op) { 1337b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 13384e0d8c25SBarry Smith aij->roworiented = flg; 1339b5a2b587SKris Buschelman break; 1340512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 1341b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 13423971808eSMatthew Knepley case MAT_NEW_NONZERO_ALLOCATION_ERR: 13434e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 1344b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1345b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 134677e54ba9SKris Buschelman case MAT_SYMMETRIC: 134777e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 13489a4540c5SBarry Smith case MAT_HERMITIAN: 13499a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1350600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 1351290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 135277e54ba9SKris Buschelman break; 1353b5a2b587SKris Buschelman default: 1354e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 13553a40ed3dSBarry Smith } 13563a40ed3dSBarry Smith PetscFunctionReturn(0); 1357289bc588SBarry Smith } 1358289bc588SBarry Smith 13594a2ae208SSatish Balay #undef __FUNCT__ 13604a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1361dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 13626f0a148fSBarry Smith { 1363ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 13646849ba73SBarry Smith PetscErrorCode ierr; 1365d0f46423SBarry Smith PetscInt lda=l->lda,m=A->rmap->n,j; 13663a40ed3dSBarry Smith 13673a40ed3dSBarry Smith PetscFunctionBegin; 1368a5ce6ee0Svictorle if (lda>m) { 1369d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1370a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1371a5ce6ee0Svictorle } 1372a5ce6ee0Svictorle } else { 1373d0f46423SBarry Smith ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1374a5ce6ee0Svictorle } 13753a40ed3dSBarry Smith PetscFunctionReturn(0); 13766f0a148fSBarry Smith } 13776f0a148fSBarry Smith 13784a2ae208SSatish Balay #undef __FUNCT__ 13794a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 13802b40b63fSBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 13816f0a148fSBarry Smith { 138297b48c8fSBarry Smith PetscErrorCode ierr; 1383ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1384d0f46423SBarry Smith PetscInt n = A->cmap->n,i,j; 138597b48c8fSBarry Smith PetscScalar *slot,*bb; 138697b48c8fSBarry Smith const PetscScalar *xx; 138755659b69SBarry Smith 13883a40ed3dSBarry Smith PetscFunctionBegin; 138997b48c8fSBarry Smith /* fix right hand side if needed */ 139097b48c8fSBarry Smith if (x && b) { 139197b48c8fSBarry Smith ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 139297b48c8fSBarry Smith ierr = VecGetArray(b,&bb);CHKERRQ(ierr); 139397b48c8fSBarry Smith for (i=0; i<N; i++) { 139497b48c8fSBarry Smith bb[rows[i]] = diag*xx[rows[i]]; 139597b48c8fSBarry Smith } 139697b48c8fSBarry Smith ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 139797b48c8fSBarry Smith ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr); 139897b48c8fSBarry Smith } 139997b48c8fSBarry Smith 14006f0a148fSBarry Smith for (i=0; i<N; i++) { 14016f0a148fSBarry Smith slot = l->v + rows[i]; 14026f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 14036f0a148fSBarry Smith } 1404f4df32b1SMatthew Knepley if (diag != 0.0) { 14056f0a148fSBarry Smith for (i=0; i<N; i++) { 14066f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1407f4df32b1SMatthew Knepley *slot = diag; 14086f0a148fSBarry Smith } 14096f0a148fSBarry Smith } 14103a40ed3dSBarry Smith PetscFunctionReturn(0); 14116f0a148fSBarry Smith } 1412557bce09SLois Curfman McInnes 14134a2ae208SSatish Balay #undef __FUNCT__ 14144a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 1415dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 141664e87e97SBarry Smith { 1417c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14183a40ed3dSBarry Smith 14193a40ed3dSBarry Smith PetscFunctionBegin; 1420e32f2f54SBarry 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"); 142164e87e97SBarry Smith *array = mat->v; 14223a40ed3dSBarry Smith PetscFunctionReturn(0); 142364e87e97SBarry Smith } 14240754003eSLois Curfman McInnes 14254a2ae208SSatish Balay #undef __FUNCT__ 14264a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 1427dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1428ff14e315SSatish Balay { 14293a40ed3dSBarry Smith PetscFunctionBegin; 143009b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 14313a40ed3dSBarry Smith PetscFunctionReturn(0); 1432ff14e315SSatish Balay } 14330754003eSLois Curfman McInnes 14344a2ae208SSatish Balay #undef __FUNCT__ 14354a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 143613f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 14370754003eSLois Curfman McInnes { 1438c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14396849ba73SBarry Smith PetscErrorCode ierr; 14405d0c19d7SBarry Smith PetscInt i,j,nrows,ncols; 14415d0c19d7SBarry Smith const PetscInt *irow,*icol; 144287828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 14430754003eSLois Curfman McInnes Mat newmat; 14440754003eSLois Curfman McInnes 14453a40ed3dSBarry Smith PetscFunctionBegin; 144678b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 144778b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1448e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1449e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 14500754003eSLois Curfman McInnes 1451182d2002SSatish Balay /* Check submatrixcall */ 1452182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 145313f74950SBarry Smith PetscInt n_cols,n_rows; 1454182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 145521a2c019SBarry Smith if (n_rows != nrows || n_cols != ncols) { 1456f746d493SDmitry Karpeev /* resize the result matrix to match number of requested rows/columns */ 1457c61587bbSBarry Smith ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 145821a2c019SBarry Smith } 1459182d2002SSatish Balay newmat = *B; 1460182d2002SSatish Balay } else { 14610754003eSLois Curfman McInnes /* Create and fill new matrix */ 14627adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr); 1463f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 14647adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 14655c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1466182d2002SSatish Balay } 1467182d2002SSatish Balay 1468182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1469182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1470182d2002SSatish Balay 1471182d2002SSatish Balay for (i=0; i<ncols; i++) { 14726de62eeeSBarry Smith av = v + mat->lda*icol[i]; 1473182d2002SSatish Balay for (j=0; j<nrows; j++) { 1474182d2002SSatish Balay *bv++ = av[irow[j]]; 14750754003eSLois Curfman McInnes } 14760754003eSLois Curfman McInnes } 1477182d2002SSatish Balay 1478182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 14796d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14806d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14810754003eSLois Curfman McInnes 14820754003eSLois Curfman McInnes /* Free work space */ 148378b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 148478b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1485182d2002SSatish Balay *B = newmat; 14863a40ed3dSBarry Smith PetscFunctionReturn(0); 14870754003eSLois Curfman McInnes } 14880754003eSLois Curfman McInnes 14894a2ae208SSatish Balay #undef __FUNCT__ 14904a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 149113f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1492905e6a2fSBarry Smith { 14936849ba73SBarry Smith PetscErrorCode ierr; 149413f74950SBarry Smith PetscInt i; 1495905e6a2fSBarry Smith 14963a40ed3dSBarry Smith PetscFunctionBegin; 1497905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1498b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1499905e6a2fSBarry Smith } 1500905e6a2fSBarry Smith 1501905e6a2fSBarry Smith for (i=0; i<n; i++) { 15026a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1503905e6a2fSBarry Smith } 15043a40ed3dSBarry Smith PetscFunctionReturn(0); 1505905e6a2fSBarry Smith } 1506905e6a2fSBarry Smith 15074a2ae208SSatish Balay #undef __FUNCT__ 1508c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1509c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1510c0aa2d19SHong Zhang { 1511c0aa2d19SHong Zhang PetscFunctionBegin; 1512c0aa2d19SHong Zhang PetscFunctionReturn(0); 1513c0aa2d19SHong Zhang } 1514c0aa2d19SHong Zhang 1515c0aa2d19SHong Zhang #undef __FUNCT__ 1516c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1517c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1518c0aa2d19SHong Zhang { 1519c0aa2d19SHong Zhang PetscFunctionBegin; 1520c0aa2d19SHong Zhang PetscFunctionReturn(0); 1521c0aa2d19SHong Zhang } 1522c0aa2d19SHong Zhang 1523c0aa2d19SHong Zhang #undef __FUNCT__ 15244a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1525dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 15264b0e389bSBarry Smith { 15274b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 15286849ba73SBarry Smith PetscErrorCode ierr; 1529d0f46423SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j; 15303a40ed3dSBarry Smith 15313a40ed3dSBarry Smith PetscFunctionBegin; 153233f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 153333f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1534cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 15353a40ed3dSBarry Smith PetscFunctionReturn(0); 15363a40ed3dSBarry Smith } 1537e32f2f54SBarry Smith if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1538a5ce6ee0Svictorle if (lda1>m || lda2>m) { 15390dbb7854Svictorle for (j=0; j<n; j++) { 1540a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1541a5ce6ee0Svictorle } 1542a5ce6ee0Svictorle } else { 1543d0f46423SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1544a5ce6ee0Svictorle } 1545273d9f13SBarry Smith PetscFunctionReturn(0); 1546273d9f13SBarry Smith } 1547273d9f13SBarry Smith 15484a2ae208SSatish Balay #undef __FUNCT__ 15494a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1550dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A) 1551273d9f13SBarry Smith { 1552dfbe8321SBarry Smith PetscErrorCode ierr; 1553273d9f13SBarry Smith 1554273d9f13SBarry Smith PetscFunctionBegin; 1555273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 15563a40ed3dSBarry Smith PetscFunctionReturn(0); 15574b0e389bSBarry Smith } 15584b0e389bSBarry Smith 1559284134d9SBarry Smith #undef __FUNCT__ 1560284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense" 1561284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N) 1562284134d9SBarry Smith { 1563284134d9SBarry Smith PetscFunctionBegin; 156421a2c019SBarry Smith /* this will not be called before lda, Mmax, and Nmax have been set */ 1565284134d9SBarry Smith m = PetscMax(m,M); 1566284134d9SBarry Smith n = PetscMax(n,N); 1567a868139aSShri Abhyankar 156886d161a7SShri 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); 156986d161a7SShri 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); 157086d161a7SShri Abhyankar */ 1571dc5cefdeSJed Brown A->rmap->n = A->rmap->N = m; 1572d0f46423SBarry Smith A->cmap->n = A->cmap->N = n; 1573284134d9SBarry Smith PetscFunctionReturn(0); 1574284134d9SBarry Smith } 1575170fe5c8SBarry Smith 1576ba337c44SJed Brown #undef __FUNCT__ 1577ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense" 1578ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A) 1579ba337c44SJed Brown { 1580ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1581ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1582ba337c44SJed Brown PetscScalar *aa = a->v; 1583ba337c44SJed Brown 1584ba337c44SJed Brown PetscFunctionBegin; 1585ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]); 1586ba337c44SJed Brown PetscFunctionReturn(0); 1587ba337c44SJed Brown } 1588ba337c44SJed Brown 1589ba337c44SJed Brown #undef __FUNCT__ 1590ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense" 1591ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A) 1592ba337c44SJed Brown { 1593ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1594ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1595ba337c44SJed Brown PetscScalar *aa = a->v; 1596ba337c44SJed Brown 1597ba337c44SJed Brown PetscFunctionBegin; 1598ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 1599ba337c44SJed Brown PetscFunctionReturn(0); 1600ba337c44SJed Brown } 1601ba337c44SJed Brown 1602ba337c44SJed Brown #undef __FUNCT__ 1603ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense" 1604ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A) 1605ba337c44SJed Brown { 1606ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1607ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1608ba337c44SJed Brown PetscScalar *aa = a->v; 1609ba337c44SJed Brown 1610ba337c44SJed Brown PetscFunctionBegin; 1611ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 1612ba337c44SJed Brown PetscFunctionReturn(0); 1613ba337c44SJed Brown } 1614284134d9SBarry Smith 1615a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/ 1616a9fe9ddaSSatish Balay #undef __FUNCT__ 1617a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1618a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1619a9fe9ddaSSatish Balay { 1620a9fe9ddaSSatish Balay PetscErrorCode ierr; 1621a9fe9ddaSSatish Balay 1622a9fe9ddaSSatish Balay PetscFunctionBegin; 1623a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1624a9fe9ddaSSatish Balay ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1625a9fe9ddaSSatish Balay } 1626a9fe9ddaSSatish Balay ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1627a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1628a9fe9ddaSSatish Balay } 1629a9fe9ddaSSatish Balay 1630a9fe9ddaSSatish Balay #undef __FUNCT__ 1631a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1632a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1633a9fe9ddaSSatish Balay { 1634ee16a9a1SHong Zhang PetscErrorCode ierr; 1635d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 1636ee16a9a1SHong Zhang Mat Cmat; 1637a9fe9ddaSSatish Balay 1638ee16a9a1SHong Zhang PetscFunctionBegin; 1639e32f2f54SBarry 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); 1640ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1641ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1642ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1643ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1644ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1645ee16a9a1SHong Zhang *C = Cmat; 1646ee16a9a1SHong Zhang PetscFunctionReturn(0); 1647ee16a9a1SHong Zhang } 1648a9fe9ddaSSatish Balay 164998a3b096SSatish Balay #undef __FUNCT__ 1650a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1651a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1652a9fe9ddaSSatish Balay { 1653a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1654a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1655a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 16560805154bSBarry Smith PetscBLASInt m,n,k; 1657a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1658a9fe9ddaSSatish Balay 1659a9fe9ddaSSatish Balay PetscFunctionBegin; 1660d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 1661d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1662d0f46423SBarry Smith k = PetscBLASIntCast(A->cmap->n); 1663a9fe9ddaSSatish Balay BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1664a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1665a9fe9ddaSSatish Balay } 1666a9fe9ddaSSatish Balay 1667a9fe9ddaSSatish Balay #undef __FUNCT__ 1668a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense" 1669a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1670a9fe9ddaSSatish Balay { 1671a9fe9ddaSSatish Balay PetscErrorCode ierr; 1672a9fe9ddaSSatish Balay 1673a9fe9ddaSSatish Balay PetscFunctionBegin; 1674a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1675a9fe9ddaSSatish Balay ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1676a9fe9ddaSSatish Balay } 1677a9fe9ddaSSatish Balay ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1678a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1679a9fe9ddaSSatish Balay } 1680a9fe9ddaSSatish Balay 1681a9fe9ddaSSatish Balay #undef __FUNCT__ 1682a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense" 1683a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1684a9fe9ddaSSatish Balay { 1685ee16a9a1SHong Zhang PetscErrorCode ierr; 1686d0f46423SBarry Smith PetscInt m=A->cmap->n,n=B->cmap->n; 1687ee16a9a1SHong Zhang Mat Cmat; 1688a9fe9ddaSSatish Balay 1689ee16a9a1SHong Zhang PetscFunctionBegin; 1690e32f2f54SBarry 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); 1691ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1692ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1693ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1694ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1695ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1696ee16a9a1SHong Zhang *C = Cmat; 1697ee16a9a1SHong Zhang PetscFunctionReturn(0); 1698ee16a9a1SHong Zhang } 1699a9fe9ddaSSatish Balay 1700a9fe9ddaSSatish Balay #undef __FUNCT__ 1701a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense" 1702a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1703a9fe9ddaSSatish Balay { 1704a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1705a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1706a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 17070805154bSBarry Smith PetscBLASInt m,n,k; 1708a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1709a9fe9ddaSSatish Balay 1710a9fe9ddaSSatish Balay PetscFunctionBegin; 1711d0f46423SBarry Smith m = PetscBLASIntCast(A->cmap->n); 1712d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1713d0f46423SBarry Smith k = PetscBLASIntCast(A->rmap->n); 17142fbe02b9SBarry Smith /* 17152fbe02b9SBarry Smith Note the m and n arguments below are the number rows and columns of A', not A! 17162fbe02b9SBarry Smith */ 1717a9fe9ddaSSatish Balay BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1718a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1719a9fe9ddaSSatish Balay } 1720985db425SBarry Smith 1721985db425SBarry Smith #undef __FUNCT__ 1722985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense" 1723985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 1724985db425SBarry Smith { 1725985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1726985db425SBarry Smith PetscErrorCode ierr; 1727d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1728985db425SBarry Smith PetscScalar *x; 1729985db425SBarry Smith MatScalar *aa = a->v; 1730985db425SBarry Smith 1731985db425SBarry Smith PetscFunctionBegin; 1732e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1733985db425SBarry Smith 1734985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1735985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1736985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1737e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1738985db425SBarry Smith for (i=0; i<m; i++) { 1739985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1740985db425SBarry Smith for (j=1; j<n; j++){ 1741985db425SBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1742985db425SBarry Smith } 1743985db425SBarry Smith } 1744985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1745985db425SBarry Smith PetscFunctionReturn(0); 1746985db425SBarry Smith } 1747985db425SBarry Smith 1748985db425SBarry Smith #undef __FUNCT__ 1749985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense" 1750985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 1751985db425SBarry Smith { 1752985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1753985db425SBarry Smith PetscErrorCode ierr; 1754d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1755985db425SBarry Smith PetscScalar *x; 1756985db425SBarry Smith PetscReal atmp; 1757985db425SBarry Smith MatScalar *aa = a->v; 1758985db425SBarry Smith 1759985db425SBarry Smith PetscFunctionBegin; 1760e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1761985db425SBarry Smith 1762985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1763985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1764985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1765e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1766985db425SBarry Smith for (i=0; i<m; i++) { 17679189402eSHong Zhang x[i] = PetscAbsScalar(aa[i]); 1768985db425SBarry Smith for (j=1; j<n; j++){ 1769985db425SBarry Smith atmp = PetscAbsScalar(aa[i+m*j]); 1770985db425SBarry Smith if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 1771985db425SBarry Smith } 1772985db425SBarry Smith } 1773985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1774985db425SBarry Smith PetscFunctionReturn(0); 1775985db425SBarry Smith } 1776985db425SBarry Smith 1777985db425SBarry Smith #undef __FUNCT__ 1778985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense" 1779985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 1780985db425SBarry Smith { 1781985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1782985db425SBarry Smith PetscErrorCode ierr; 1783d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1784985db425SBarry Smith PetscScalar *x; 1785985db425SBarry Smith MatScalar *aa = a->v; 1786985db425SBarry Smith 1787985db425SBarry Smith PetscFunctionBegin; 1788e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1789985db425SBarry Smith 1790985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1791985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1792985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1793e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1794985db425SBarry Smith for (i=0; i<m; i++) { 1795985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1796985db425SBarry Smith for (j=1; j<n; j++){ 1797985db425SBarry Smith if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1798985db425SBarry Smith } 1799985db425SBarry Smith } 1800985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1801985db425SBarry Smith PetscFunctionReturn(0); 1802985db425SBarry Smith } 1803985db425SBarry Smith 18048d0534beSBarry Smith #undef __FUNCT__ 18058d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense" 18068d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 18078d0534beSBarry Smith { 18088d0534beSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 18098d0534beSBarry Smith PetscErrorCode ierr; 18108d0534beSBarry Smith PetscScalar *x; 18118d0534beSBarry Smith 18128d0534beSBarry Smith PetscFunctionBegin; 1813e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 18148d0534beSBarry Smith 18158d0534beSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1816d0f46423SBarry Smith ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 18178d0534beSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 18188d0534beSBarry Smith PetscFunctionReturn(0); 18198d0534beSBarry Smith } 18208d0534beSBarry Smith 1821289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1822a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1823905e6a2fSBarry Smith MatGetRow_SeqDense, 1824905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1825905e6a2fSBarry Smith MatMult_SeqDense, 182697304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 18277c922b88SBarry Smith MatMultTranspose_SeqDense, 18287c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1829db4efbfdSBarry Smith 0, 1830db4efbfdSBarry Smith 0, 1831db4efbfdSBarry Smith 0, 1832db4efbfdSBarry Smith /*10*/ 0, 1833905e6a2fSBarry Smith MatLUFactor_SeqDense, 1834905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 183541f059aeSBarry Smith MatSOR_SeqDense, 1836ec8511deSBarry Smith MatTranspose_SeqDense, 183797304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1838905e6a2fSBarry Smith MatEqual_SeqDense, 1839905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1840905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1841905e6a2fSBarry Smith MatNorm_SeqDense, 1842c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense, 1843c0aa2d19SHong Zhang MatAssemblyEnd_SeqDense, 1844905e6a2fSBarry Smith MatSetOption_SeqDense, 1845905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1846d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense, 1847db4efbfdSBarry Smith 0, 1848db4efbfdSBarry Smith 0, 1849db4efbfdSBarry Smith 0, 1850db4efbfdSBarry Smith 0, 1851d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqDense, 1852273d9f13SBarry Smith 0, 1853905e6a2fSBarry Smith 0, 1854905e6a2fSBarry Smith MatGetArray_SeqDense, 1855905e6a2fSBarry Smith MatRestoreArray_SeqDense, 1856d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense, 1857a5ae1ecdSBarry Smith 0, 1858a5ae1ecdSBarry Smith 0, 1859a5ae1ecdSBarry Smith 0, 1860a5ae1ecdSBarry Smith 0, 1861d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense, 1862a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1863a5ae1ecdSBarry Smith 0, 18644b0e389bSBarry Smith MatGetValues_SeqDense, 1865a5ae1ecdSBarry Smith MatCopy_SeqDense, 1866d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense, 1867a5ae1ecdSBarry Smith MatScale_SeqDense, 1868a5ae1ecdSBarry Smith 0, 1869a5ae1ecdSBarry Smith 0, 1870a5ae1ecdSBarry Smith 0, 1871d519adbfSMatthew Knepley /*49*/ 0, 1872a5ae1ecdSBarry Smith 0, 1873a5ae1ecdSBarry Smith 0, 1874a5ae1ecdSBarry Smith 0, 1875a5ae1ecdSBarry Smith 0, 1876d519adbfSMatthew Knepley /*54*/ 0, 1877a5ae1ecdSBarry Smith 0, 1878a5ae1ecdSBarry Smith 0, 1879a5ae1ecdSBarry Smith 0, 1880a5ae1ecdSBarry Smith 0, 1881d519adbfSMatthew Knepley /*59*/ 0, 1882e03a110bSBarry Smith MatDestroy_SeqDense, 1883e03a110bSBarry Smith MatView_SeqDense, 1884357abbc8SBarry Smith 0, 188597304618SKris Buschelman 0, 1886d519adbfSMatthew Knepley /*64*/ 0, 188797304618SKris Buschelman 0, 188897304618SKris Buschelman 0, 188997304618SKris Buschelman 0, 189097304618SKris Buschelman 0, 1891d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense, 189297304618SKris Buschelman 0, 189397304618SKris Buschelman 0, 189497304618SKris Buschelman 0, 189597304618SKris Buschelman 0, 1896d519adbfSMatthew Knepley /*74*/ 0, 189797304618SKris Buschelman 0, 189897304618SKris Buschelman 0, 189997304618SKris Buschelman 0, 190097304618SKris Buschelman 0, 1901d519adbfSMatthew Knepley /*79*/ 0, 190297304618SKris Buschelman 0, 190397304618SKris Buschelman 0, 190497304618SKris Buschelman 0, 19055bba2384SShri Abhyankar /*83*/ MatLoad_SeqDense, 1906865e5f61SKris Buschelman 0, 19071cbb95d3SBarry Smith MatIsHermitian_SeqDense, 1908865e5f61SKris Buschelman 0, 1909865e5f61SKris Buschelman 0, 1910865e5f61SKris Buschelman 0, 1911d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense, 1912a9fe9ddaSSatish Balay MatMatMultSymbolic_SeqDense_SeqDense, 1913a9fe9ddaSSatish Balay MatMatMultNumeric_SeqDense_SeqDense, 1914865e5f61SKris Buschelman 0, 1915865e5f61SKris Buschelman 0, 1916d519adbfSMatthew Knepley /*94*/ 0, 1917a9fe9ddaSSatish Balay MatMatMultTranspose_SeqDense_SeqDense, 1918a9fe9ddaSSatish Balay MatMatMultTransposeSymbolic_SeqDense_SeqDense, 1919a9fe9ddaSSatish Balay MatMatMultTransposeNumeric_SeqDense_SeqDense, 1920284134d9SBarry Smith 0, 1921d519adbfSMatthew Knepley /*99*/ 0, 1922284134d9SBarry Smith 0, 1923284134d9SBarry Smith 0, 1924ba337c44SJed Brown MatConjugate_SeqDense, 1925985db425SBarry Smith MatSetSizes_SeqDense, 1926ba337c44SJed Brown /*104*/0, 1927ba337c44SJed Brown MatRealPart_SeqDense, 1928ba337c44SJed Brown MatImaginaryPart_SeqDense, 1929985db425SBarry Smith 0, 1930985db425SBarry Smith 0, 1931*85e2c93fSHong Zhang /*109*/MatMatSolve_SeqDense, 1932985db425SBarry Smith 0, 19338d0534beSBarry Smith MatGetRowMin_SeqDense, 1934aabbc4fbSShri Abhyankar MatGetColumnVector_SeqDense, 1935aabbc4fbSShri Abhyankar 0, 1936aabbc4fbSShri Abhyankar /*114*/0, 1937aabbc4fbSShri Abhyankar 0, 1938aabbc4fbSShri Abhyankar 0, 1939aabbc4fbSShri Abhyankar 0, 1940aabbc4fbSShri Abhyankar 0, 1941aabbc4fbSShri Abhyankar /*119*/0, 1942aabbc4fbSShri Abhyankar 0, 1943aabbc4fbSShri Abhyankar 0, 19445bba2384SShri Abhyankar 0 1945985db425SBarry Smith }; 194690ace30eSBarry Smith 19474a2ae208SSatish Balay #undef __FUNCT__ 19484a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 19494b828684SBarry Smith /*@C 1950fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1951d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1952d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1953289bc588SBarry Smith 1954db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1955db81eaa0SLois Curfman McInnes 195620563c6bSBarry Smith Input Parameters: 1957db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 19580c775827SLois Curfman McInnes . m - number of rows 195918f449edSLois Curfman McInnes . n - number of columns 1960c0235b3cSMatthew Knepley - data - optional location of matrix data in column major order. Set data=PETSC_NULL for PETSc 1961dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 196220563c6bSBarry Smith 196320563c6bSBarry Smith Output Parameter: 196444cd7ae7SLois Curfman McInnes . A - the matrix 196520563c6bSBarry Smith 1966b259b22eSLois Curfman McInnes Notes: 196718f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 196818f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1969b4fd4287SBarry Smith set data=PETSC_NULL. 197018f449edSLois Curfman McInnes 1971027ccd11SLois Curfman McInnes Level: intermediate 1972027ccd11SLois Curfman McInnes 1973dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1974d65003e9SLois Curfman McInnes 1975db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 197620563c6bSBarry Smith @*/ 19777087cfbeSBarry Smith PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 1978289bc588SBarry Smith { 1979dfbe8321SBarry Smith PetscErrorCode ierr; 19803b2fbd54SBarry Smith 19813a40ed3dSBarry Smith PetscFunctionBegin; 1982f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1983f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1984273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1985273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1986273d9f13SBarry Smith PetscFunctionReturn(0); 1987273d9f13SBarry Smith } 1988273d9f13SBarry Smith 19894a2ae208SSatish Balay #undef __FUNCT__ 1990afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation" 1991273d9f13SBarry Smith /*@C 1992273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1993273d9f13SBarry Smith 1994273d9f13SBarry Smith Collective on MPI_Comm 1995273d9f13SBarry Smith 1996273d9f13SBarry Smith Input Parameters: 1997273d9f13SBarry Smith + A - the matrix 1998273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1999273d9f13SBarry Smith 2000273d9f13SBarry Smith Notes: 2001273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 2002273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 2003284134d9SBarry Smith need not call this routine. 2004273d9f13SBarry Smith 2005273d9f13SBarry Smith Level: intermediate 2006273d9f13SBarry Smith 2007273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 2008273d9f13SBarry Smith 2009867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA() 2010867c911aSBarry Smith 2011273d9f13SBarry Smith @*/ 20127087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 2013273d9f13SBarry Smith { 20144ac538c5SBarry Smith PetscErrorCode ierr; 2015a23d5eceSKris Buschelman 2016a23d5eceSKris Buschelman PetscFunctionBegin; 20174ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));CHKERRQ(ierr); 2018a23d5eceSKris Buschelman PetscFunctionReturn(0); 2019a23d5eceSKris Buschelman } 2020a23d5eceSKris Buschelman 2021a23d5eceSKris Buschelman EXTERN_C_BEGIN 2022a23d5eceSKris Buschelman #undef __FUNCT__ 2023afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense" 20247087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 2025a23d5eceSKris Buschelman { 2026273d9f13SBarry Smith Mat_SeqDense *b; 2027dfbe8321SBarry Smith PetscErrorCode ierr; 2028273d9f13SBarry Smith 2029273d9f13SBarry Smith PetscFunctionBegin; 2030273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 2031a868139aSShri Abhyankar 203234ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr); 203334ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr); 203434ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 203534ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 203634ef9618SShri Abhyankar 2037273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 203886d161a7SShri Abhyankar b->Mmax = B->rmap->n; 203986d161a7SShri Abhyankar b->Nmax = B->cmap->n; 204086d161a7SShri Abhyankar if(b->lda <= 0 || b->changelda) b->lda = B->rmap->n; 204186d161a7SShri Abhyankar 20429e8f95c4SLisandro Dalcin if (!data) { /* petsc-allocated storage */ 20439e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 20445afd5e0cSBarry Smith ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 2045284134d9SBarry Smith ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 2046284134d9SBarry Smith ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 20479e8f95c4SLisandro Dalcin b->user_alloc = PETSC_FALSE; 2048273d9f13SBarry Smith } else { /* user-allocated storage */ 20499e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 2050273d9f13SBarry Smith b->v = data; 2051273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 2052273d9f13SBarry Smith } 20530450473dSBarry Smith B->assembled = PETSC_TRUE; 2054273d9f13SBarry Smith PetscFunctionReturn(0); 2055273d9f13SBarry Smith } 2056a23d5eceSKris Buschelman EXTERN_C_END 2057273d9f13SBarry Smith 20581b807ce4Svictorle #undef __FUNCT__ 20591b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 20601b807ce4Svictorle /*@C 20611b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 20621b807ce4Svictorle 20631b807ce4Svictorle Input parameter: 20641b807ce4Svictorle + A - the matrix 20651b807ce4Svictorle - lda - the leading dimension 20661b807ce4Svictorle 20671b807ce4Svictorle Notes: 2068867c911aSBarry Smith This routine is to be used in conjunction with MatSeqDenseSetPreallocation(); 20691b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 20701b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 20711b807ce4Svictorle 20721b807ce4Svictorle Level: intermediate 20731b807ce4Svictorle 20741b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 20751b807ce4Svictorle 2076284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 2077867c911aSBarry Smith 20781b807ce4Svictorle @*/ 20797087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetLDA(Mat B,PetscInt lda) 20801b807ce4Svictorle { 20811b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 208221a2c019SBarry Smith 20831b807ce4Svictorle PetscFunctionBegin; 2084e32f2f54SBarry 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); 20851b807ce4Svictorle b->lda = lda; 208621a2c019SBarry Smith b->changelda = PETSC_FALSE; 208721a2c019SBarry Smith b->Mmax = PetscMax(b->Mmax,lda); 20881b807ce4Svictorle PetscFunctionReturn(0); 20891b807ce4Svictorle } 20901b807ce4Svictorle 20910bad9183SKris Buschelman /*MC 2092fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 20930bad9183SKris Buschelman 20940bad9183SKris Buschelman Options Database Keys: 20950bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 20960bad9183SKris Buschelman 20970bad9183SKris Buschelman Level: beginner 20980bad9183SKris Buschelman 209989665df3SBarry Smith .seealso: MatCreateSeqDense() 210089665df3SBarry Smith 21010bad9183SKris Buschelman M*/ 21020bad9183SKris Buschelman 2103273d9f13SBarry Smith EXTERN_C_BEGIN 21044a2ae208SSatish Balay #undef __FUNCT__ 21054a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 21067087cfbeSBarry Smith PetscErrorCode MatCreate_SeqDense(Mat B) 2107273d9f13SBarry Smith { 2108273d9f13SBarry Smith Mat_SeqDense *b; 2109dfbe8321SBarry Smith PetscErrorCode ierr; 21107c334f02SBarry Smith PetscMPIInt size; 2111273d9f13SBarry Smith 2112273d9f13SBarry Smith PetscFunctionBegin; 21137adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr); 2114e32f2f54SBarry Smith if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 211555659b69SBarry Smith 211638f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr); 2117549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 211844cd7ae7SLois Curfman McInnes B->data = (void*)b; 211918f449edSLois Curfman McInnes 212044cd7ae7SLois Curfman McInnes b->pivots = 0; 2121273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 2122273d9f13SBarry Smith b->v = 0; 212321a2c019SBarry Smith b->changelda = PETSC_FALSE; 21244e220ebcSLois Curfman McInnes 2125b24902e0SBarry Smith 2126ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C", 2127b24902e0SBarry Smith "MatGetFactor_seqdense_petsc", 2128b24902e0SBarry Smith MatGetFactor_seqdense_petsc);CHKERRQ(ierr); 2129a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 2130a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 2131a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 21324ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C", 21334ae313f4SHong Zhang "MatMatMult_SeqAIJ_SeqDense", 21344ae313f4SHong Zhang MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 21354ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C", 21364ae313f4SHong Zhang "MatMatMultSymbolic_SeqAIJ_SeqDense", 21374ae313f4SHong Zhang MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 21384ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C", 21394ae313f4SHong Zhang "MatMatMultNumeric_SeqAIJ_SeqDense", 21404ae313f4SHong Zhang MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 214117667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 21423a40ed3dSBarry Smith PetscFunctionReturn(0); 2143289bc588SBarry Smith } 2144273d9f13SBarry Smith EXTERN_C_END 2145