1be1d678aSKris Buschelman 267e560aaSBarry Smith /* 367e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 467e560aaSBarry Smith */ 5289bc588SBarry Smith 6c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> 7c6db04a5SJed Brown #include <petscblaslapack.h> 8289bc588SBarry Smith 94a2ae208SSatish Balay #undef __FUNCT__ 104a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense" 11f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 121987afe7SBarry Smith { 131987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 14f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 1513f74950SBarry Smith PetscInt j; 160805154bSBarry Smith PetscBLASInt N,m,ldax,lday,one = 1; 17efee365bSSatish Balay PetscErrorCode ierr; 183a40ed3dSBarry Smith 193a40ed3dSBarry Smith PetscFunctionBegin; 20d0f46423SBarry Smith N = PetscBLASIntCast(X->rmap->n*X->cmap->n); 21d0f46423SBarry Smith m = PetscBLASIntCast(X->rmap->n); 220805154bSBarry Smith ldax = PetscBLASIntCast(x->lda); 230805154bSBarry Smith lday = PetscBLASIntCast(y->lda); 24a5ce6ee0Svictorle if (ldax>m || lday>m) { 25d0f46423SBarry Smith for (j=0; j<X->cmap->n; j++) { 26f4df32b1SMatthew Knepley BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one); 27a5ce6ee0Svictorle } 28a5ce6ee0Svictorle } else { 29f4df32b1SMatthew Knepley BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one); 30a5ce6ee0Svictorle } 310450473dSBarry Smith ierr = PetscLogFlops(PetscMax(2*N-1,0));CHKERRQ(ierr); 323a40ed3dSBarry Smith PetscFunctionReturn(0); 331987afe7SBarry Smith } 341987afe7SBarry Smith 354a2ae208SSatish Balay #undef __FUNCT__ 364a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 37dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 38289bc588SBarry Smith { 39d0f46423SBarry Smith PetscInt N = A->rmap->n*A->cmap->n; 403a40ed3dSBarry Smith 413a40ed3dSBarry Smith PetscFunctionBegin; 424e220ebcSLois Curfman McInnes info->block_size = 1.0; 434e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 446de62eeeSBarry Smith info->nz_used = (double)N; 456de62eeeSBarry Smith info->nz_unneeded = (double)0; 464e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 474e220ebcSLois Curfman McInnes info->mallocs = 0; 487adad957SLisandro Dalcin info->memory = ((PetscObject)A)->mem; 494e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 504e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 514e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 523a40ed3dSBarry Smith PetscFunctionReturn(0); 53289bc588SBarry Smith } 54289bc588SBarry Smith 554a2ae208SSatish Balay #undef __FUNCT__ 564a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense" 57f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha) 5880cd9d93SLois Curfman McInnes { 59273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 60f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 61efee365bSSatish Balay PetscErrorCode ierr; 620805154bSBarry Smith PetscBLASInt one = 1,j,nz,lda = PetscBLASIntCast(a->lda); 6380cd9d93SLois Curfman McInnes 643a40ed3dSBarry Smith PetscFunctionBegin; 65d0f46423SBarry Smith if (lda>A->rmap->n) { 66d0f46423SBarry Smith nz = PetscBLASIntCast(A->rmap->n); 67d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 68f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v+j*lda,&one); 69a5ce6ee0Svictorle } 70a5ce6ee0Svictorle } else { 71d0f46423SBarry Smith nz = PetscBLASIntCast(A->rmap->n*A->cmap->n); 72f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v,&one); 73a5ce6ee0Svictorle } 74efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 753a40ed3dSBarry Smith PetscFunctionReturn(0); 7680cd9d93SLois Curfman McInnes } 7780cd9d93SLois Curfman McInnes 781cbb95d3SBarry Smith #undef __FUNCT__ 791cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense" 80ace3abfcSBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool *fl) 811cbb95d3SBarry Smith { 821cbb95d3SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 83d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,N; 841cbb95d3SBarry Smith PetscScalar *v = a->v; 851cbb95d3SBarry Smith 861cbb95d3SBarry Smith PetscFunctionBegin; 871cbb95d3SBarry Smith *fl = PETSC_FALSE; 88d0f46423SBarry Smith if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0); 891cbb95d3SBarry Smith N = a->lda; 901cbb95d3SBarry Smith 911cbb95d3SBarry Smith for (i=0; i<m; i++) { 921cbb95d3SBarry Smith for (j=i+1; j<m; j++) { 931cbb95d3SBarry Smith if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0); 941cbb95d3SBarry Smith } 951cbb95d3SBarry Smith } 961cbb95d3SBarry Smith *fl = PETSC_TRUE; 971cbb95d3SBarry Smith PetscFunctionReturn(0); 981cbb95d3SBarry Smith } 991cbb95d3SBarry Smith 100b24902e0SBarry Smith #undef __FUNCT__ 101b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense" 102719d5645SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues) 103b24902e0SBarry Smith { 104b24902e0SBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 105b24902e0SBarry Smith PetscErrorCode ierr; 106b24902e0SBarry Smith PetscInt lda = (PetscInt)mat->lda,j,m; 107b24902e0SBarry Smith 108b24902e0SBarry Smith PetscFunctionBegin; 109*aa5ea44dSBarry Smith ierr = PetscLayoutReference(A->rmap,&newi->rmap);CHKERRQ(ierr); 110*aa5ea44dSBarry Smith ierr = PetscLayoutReference(A->cmap,&newi->cmap);CHKERRQ(ierr); 111719d5645SBarry Smith ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr); 112b24902e0SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 113b24902e0SBarry Smith l = (Mat_SeqDense*)newi->data; 114d0f46423SBarry Smith if (lda>A->rmap->n) { 115d0f46423SBarry Smith m = A->rmap->n; 116d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 117b24902e0SBarry Smith ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 118b24902e0SBarry Smith } 119b24902e0SBarry Smith } else { 120d0f46423SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 121b24902e0SBarry Smith } 122b24902e0SBarry Smith } 123b24902e0SBarry Smith newi->assembled = PETSC_TRUE; 124b24902e0SBarry Smith PetscFunctionReturn(0); 125b24902e0SBarry Smith } 126b24902e0SBarry Smith 1274a2ae208SSatish Balay #undef __FUNCT__ 1284a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 129dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 13002cad45dSBarry Smith { 1316849ba73SBarry Smith PetscErrorCode ierr; 13202cad45dSBarry Smith 1333a40ed3dSBarry Smith PetscFunctionBegin; 1345c9eb25fSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr); 135d0f46423SBarry Smith ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 1365c9eb25fSBarry Smith ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 137719d5645SBarry Smith ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr); 138b24902e0SBarry Smith PetscFunctionReturn(0); 139b24902e0SBarry Smith } 140b24902e0SBarry Smith 1416ee01492SSatish Balay 1420481f469SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*); 143719d5645SBarry Smith 1444a2ae208SSatish Balay #undef __FUNCT__ 1454a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 1460481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 147289bc588SBarry Smith { 1484482741eSBarry Smith MatFactorInfo info; 149a093e273SMatthew Knepley PetscErrorCode ierr; 1503a40ed3dSBarry Smith 1513a40ed3dSBarry Smith PetscFunctionBegin; 152c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 153719d5645SBarry Smith ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr); 1543a40ed3dSBarry Smith PetscFunctionReturn(0); 155289bc588SBarry Smith } 1566ee01492SSatish Balay 1570b4b3355SBarry Smith #undef __FUNCT__ 1584a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 159dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 160289bc588SBarry Smith { 161c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1626849ba73SBarry Smith PetscErrorCode ierr; 16387828ca2SBarry Smith PetscScalar *x,*y; 164d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 16567e560aaSBarry Smith 1663a40ed3dSBarry Smith PetscFunctionBegin; 1671ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1681ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 169d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 170d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_LU) { 171ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 172e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 173ae7cfcebSSatish Balay #else 17471044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 175e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 176ae7cfcebSSatish Balay #endif 177d5f3da31SBarry Smith } else if (A->factortype == MAT_FACTOR_CHOLESKY){ 178ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 179e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 180ae7cfcebSSatish Balay #else 18171044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 182e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 183ae7cfcebSSatish Balay #endif 184289bc588SBarry Smith } 185e32f2f54SBarry Smith else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 1861ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1871ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 188dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 1893a40ed3dSBarry Smith PetscFunctionReturn(0); 190289bc588SBarry Smith } 1916ee01492SSatish Balay 1924a2ae208SSatish Balay #undef __FUNCT__ 19385e2c93fSHong Zhang #define __FUNCT__ "MatMatSolve_SeqDense" 19485e2c93fSHong Zhang PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X) 19585e2c93fSHong Zhang { 19685e2c93fSHong Zhang Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 19785e2c93fSHong Zhang PetscErrorCode ierr; 19885e2c93fSHong Zhang PetscScalar *b,*x; 19985e2c93fSHong Zhang PetscBLASInt nrhs,info,m=PetscBLASIntCast(A->rmap->n); 20085e2c93fSHong Zhang 20185e2c93fSHong Zhang PetscFunctionBegin; 20285e2c93fSHong Zhang ierr = MatGetSize(B,PETSC_NULL,&nrhs);CHKERRQ(ierr); 20385e2c93fSHong Zhang ierr = MatGetArray(B,&b);CHKERRQ(ierr); 20485e2c93fSHong Zhang ierr = MatGetArray(X,&x);CHKERRQ(ierr); 20585e2c93fSHong Zhang 20685e2c93fSHong Zhang ierr = PetscMemcpy(x,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 20785e2c93fSHong Zhang 20885e2c93fSHong Zhang if (A->factortype == MAT_FACTOR_LU) { 20985e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_GETRS) 21085e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 21185e2c93fSHong Zhang #else 21285e2c93fSHong Zhang LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info); 21385e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 21485e2c93fSHong Zhang #endif 21585e2c93fSHong Zhang } else if (A->factortype == MAT_FACTOR_CHOLESKY){ 21685e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_POTRS) 21785e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 21885e2c93fSHong Zhang #else 21985e2c93fSHong Zhang LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info); 22085e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 22185e2c93fSHong Zhang #endif 22285e2c93fSHong Zhang } 22385e2c93fSHong Zhang else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 22485e2c93fSHong Zhang 22585e2c93fSHong Zhang ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 22685e2c93fSHong Zhang ierr = MatRestoreArray(X,&x);CHKERRQ(ierr); 22785e2c93fSHong Zhang ierr = PetscLogFlops(nrhs*(2.0*m*m- m));CHKERRQ(ierr); 22885e2c93fSHong Zhang PetscFunctionReturn(0); 22985e2c93fSHong Zhang } 23085e2c93fSHong Zhang 23185e2c93fSHong Zhang #undef __FUNCT__ 2324a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 233dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 234da3a660dSBarry Smith { 235c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 236dfbe8321SBarry Smith PetscErrorCode ierr; 23787828ca2SBarry Smith PetscScalar *x,*y; 238d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 23967e560aaSBarry Smith 2403a40ed3dSBarry Smith PetscFunctionBegin; 2411ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2421ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 243d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 244752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 245da3a660dSBarry Smith if (mat->pivots) { 246ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 247e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 248ae7cfcebSSatish Balay #else 24971044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 250e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 251ae7cfcebSSatish Balay #endif 2527a97a34bSBarry Smith } else { 253ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 254e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 255ae7cfcebSSatish Balay #else 25671044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 257e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 258ae7cfcebSSatish Balay #endif 259da3a660dSBarry Smith } 2601ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2611ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 262dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 2633a40ed3dSBarry Smith PetscFunctionReturn(0); 264da3a660dSBarry Smith } 2656ee01492SSatish Balay 2664a2ae208SSatish Balay #undef __FUNCT__ 2674a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 268dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 269da3a660dSBarry Smith { 270c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 271dfbe8321SBarry Smith PetscErrorCode ierr; 27287828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 273da3a660dSBarry Smith Vec tmp = 0; 274d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 27567e560aaSBarry Smith 2763a40ed3dSBarry Smith PetscFunctionBegin; 2771ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2781ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 279d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 280da3a660dSBarry Smith if (yy == zz) { 28178b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 28252e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 28378b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 284da3a660dSBarry Smith } 285d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 286752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 287da3a660dSBarry Smith if (mat->pivots) { 288ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 289e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 290ae7cfcebSSatish Balay #else 29171044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 292e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 293ae7cfcebSSatish Balay #endif 294a8c6a408SBarry Smith } else { 295ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 296e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 297ae7cfcebSSatish Balay #else 29871044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 299e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 300ae7cfcebSSatish Balay #endif 301da3a660dSBarry Smith } 3026bf464f9SBarry Smith if (tmp) { 3036bf464f9SBarry Smith ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 3046bf464f9SBarry Smith ierr = VecDestroy(&tmp);CHKERRQ(ierr); 3056bf464f9SBarry Smith } else { 3066bf464f9SBarry Smith ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 3076bf464f9SBarry Smith } 3081ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3091ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 310dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 3113a40ed3dSBarry Smith PetscFunctionReturn(0); 312da3a660dSBarry Smith } 31367e560aaSBarry Smith 3144a2ae208SSatish Balay #undef __FUNCT__ 3154a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 316dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 317da3a660dSBarry Smith { 318c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3196849ba73SBarry Smith PetscErrorCode ierr; 32087828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 321da3a660dSBarry Smith Vec tmp; 322d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 32367e560aaSBarry Smith 3243a40ed3dSBarry Smith PetscFunctionBegin; 325d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 3261ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3271ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 328da3a660dSBarry Smith if (yy == zz) { 32978b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 33052e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 33178b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 332da3a660dSBarry Smith } 333d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 334752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 335da3a660dSBarry Smith if (mat->pivots) { 336ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 337e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 338ae7cfcebSSatish Balay #else 33971044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 340e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 341ae7cfcebSSatish Balay #endif 3423a40ed3dSBarry Smith } else { 343ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 344e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 345ae7cfcebSSatish Balay #else 34671044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 347e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 348ae7cfcebSSatish Balay #endif 349da3a660dSBarry Smith } 35090f02eecSBarry Smith if (tmp) { 3512dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 3526bf464f9SBarry Smith ierr = VecDestroy(&tmp);CHKERRQ(ierr); 3533a40ed3dSBarry Smith } else { 3542dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 35590f02eecSBarry Smith } 3561ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3571ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 358dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 3593a40ed3dSBarry Smith PetscFunctionReturn(0); 360da3a660dSBarry Smith } 361db4efbfdSBarry Smith 362db4efbfdSBarry Smith /* ---------------------------------------------------------------*/ 363db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 364db4efbfdSBarry Smith rather than put it in the Mat->row slot.*/ 365db4efbfdSBarry Smith #undef __FUNCT__ 366db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense" 3670481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo) 368db4efbfdSBarry Smith { 369db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 370db4efbfdSBarry Smith PetscFunctionBegin; 371e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 372db4efbfdSBarry Smith #else 373db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 374db4efbfdSBarry Smith PetscErrorCode ierr; 375db4efbfdSBarry Smith PetscBLASInt n,m,info; 376db4efbfdSBarry Smith 377db4efbfdSBarry Smith PetscFunctionBegin; 378db4efbfdSBarry Smith n = PetscBLASIntCast(A->cmap->n); 379db4efbfdSBarry Smith m = PetscBLASIntCast(A->rmap->n); 380db4efbfdSBarry Smith if (!mat->pivots) { 381db4efbfdSBarry Smith ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr); 382db4efbfdSBarry Smith ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr); 383db4efbfdSBarry Smith } 384db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 385db4efbfdSBarry Smith LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info); 386e32f2f54SBarry Smith if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization"); 387e32f2f54SBarry Smith if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 388db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 389db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 390db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 391db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 392d5f3da31SBarry Smith A->factortype = MAT_FACTOR_LU; 393db4efbfdSBarry Smith 394dc0b31edSSatish Balay ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr); 395db4efbfdSBarry Smith #endif 396db4efbfdSBarry Smith PetscFunctionReturn(0); 397db4efbfdSBarry Smith } 398db4efbfdSBarry Smith 399db4efbfdSBarry Smith #undef __FUNCT__ 400db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense" 4010481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo) 402db4efbfdSBarry Smith { 403db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 404db4efbfdSBarry Smith PetscFunctionBegin; 405e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 406db4efbfdSBarry Smith #else 407db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 408db4efbfdSBarry Smith PetscErrorCode ierr; 409db4efbfdSBarry Smith PetscBLASInt info,n = PetscBLASIntCast(A->cmap->n); 410db4efbfdSBarry Smith 411db4efbfdSBarry Smith PetscFunctionBegin; 412db4efbfdSBarry Smith ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 413db4efbfdSBarry Smith 414db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 415db4efbfdSBarry Smith LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info); 416e32f2f54SBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 417db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 418db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 419db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 420db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 421d5f3da31SBarry Smith A->factortype = MAT_FACTOR_CHOLESKY; 422dc0b31edSSatish Balay ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr); 423db4efbfdSBarry Smith #endif 424db4efbfdSBarry Smith PetscFunctionReturn(0); 425db4efbfdSBarry Smith } 426db4efbfdSBarry Smith 427db4efbfdSBarry Smith 428db4efbfdSBarry Smith #undef __FUNCT__ 429db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 4300481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 431db4efbfdSBarry Smith { 432db4efbfdSBarry Smith PetscErrorCode ierr; 433db4efbfdSBarry Smith MatFactorInfo info; 434db4efbfdSBarry Smith 435db4efbfdSBarry Smith PetscFunctionBegin; 436db4efbfdSBarry Smith info.fill = 1.0; 437c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 438719d5645SBarry Smith ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr); 439db4efbfdSBarry Smith PetscFunctionReturn(0); 440db4efbfdSBarry Smith } 441db4efbfdSBarry Smith 442db4efbfdSBarry Smith #undef __FUNCT__ 443db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 4440481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info) 445db4efbfdSBarry Smith { 446db4efbfdSBarry Smith PetscFunctionBegin; 447c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 448719d5645SBarry Smith fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense; 449db4efbfdSBarry Smith PetscFunctionReturn(0); 450db4efbfdSBarry Smith } 451db4efbfdSBarry Smith 452db4efbfdSBarry Smith #undef __FUNCT__ 453db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 4540481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 455db4efbfdSBarry Smith { 456db4efbfdSBarry Smith PetscFunctionBegin; 457c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 458719d5645SBarry Smith fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense; 459db4efbfdSBarry Smith PetscFunctionReturn(0); 460db4efbfdSBarry Smith } 461db4efbfdSBarry Smith 462bb5747d9SMatthew Knepley EXTERN_C_BEGIN 463db4efbfdSBarry Smith #undef __FUNCT__ 464db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc" 465db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact) 466db4efbfdSBarry Smith { 467db4efbfdSBarry Smith PetscErrorCode ierr; 468db4efbfdSBarry Smith 469db4efbfdSBarry Smith PetscFunctionBegin; 470db4efbfdSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr); 471db4efbfdSBarry Smith ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 472db4efbfdSBarry Smith ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr); 473db4efbfdSBarry Smith if (ftype == MAT_FACTOR_LU){ 474db4efbfdSBarry Smith (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense; 475db4efbfdSBarry Smith } else { 476db4efbfdSBarry Smith (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense; 477db4efbfdSBarry Smith } 478d5f3da31SBarry Smith (*fact)->factortype = ftype; 479db4efbfdSBarry Smith PetscFunctionReturn(0); 480db4efbfdSBarry Smith } 481bb5747d9SMatthew Knepley EXTERN_C_END 482db4efbfdSBarry Smith 483289bc588SBarry Smith /* ------------------------------------------------------------------*/ 4844a2ae208SSatish Balay #undef __FUNCT__ 48541f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense" 48641f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 487289bc588SBarry Smith { 488c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 48987828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 490dfbe8321SBarry Smith PetscErrorCode ierr; 491d0f46423SBarry Smith PetscInt m = A->rmap->n,i; 492aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 4930805154bSBarry Smith PetscBLASInt o = 1,bm = PetscBLASIntCast(m); 494bc1b551cSSatish Balay #endif 495289bc588SBarry Smith 4963a40ed3dSBarry Smith PetscFunctionBegin; 497289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 49871044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 4992dcb1b2aSMatthew Knepley ierr = VecSet(xx,zero);CHKERRQ(ierr); 500289bc588SBarry Smith } 5011ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5021ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 503b965ef7fSBarry Smith its = its*lits; 504e32f2f54SBarry 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); 505289bc588SBarry Smith while (its--) { 506fccaa45eSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 507289bc588SBarry Smith for (i=0; i<m; i++) { 508aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 509f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 510f1747703SBarry Smith not happy about returning a double complex */ 51113f74950SBarry Smith PetscInt _i; 51287828ca2SBarry Smith PetscScalar sum = b[i]; 513f1747703SBarry Smith for (_i=0; _i<m; _i++) { 5143f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 515f1747703SBarry Smith } 516f1747703SBarry Smith xt = sum; 517f1747703SBarry Smith #else 51871044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 519f1747703SBarry Smith #endif 52055a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 521289bc588SBarry Smith } 522289bc588SBarry Smith } 523fccaa45eSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 524289bc588SBarry Smith for (i=m-1; i>=0; i--) { 525aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 526f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 527f1747703SBarry Smith not happy about returning a double complex */ 52813f74950SBarry Smith PetscInt _i; 52987828ca2SBarry Smith PetscScalar sum = b[i]; 530f1747703SBarry Smith for (_i=0; _i<m; _i++) { 5313f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 532f1747703SBarry Smith } 533f1747703SBarry Smith xt = sum; 534f1747703SBarry Smith #else 53571044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 536f1747703SBarry Smith #endif 53755a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 538289bc588SBarry Smith } 539289bc588SBarry Smith } 540289bc588SBarry Smith } 5411ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 5421ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5433a40ed3dSBarry Smith PetscFunctionReturn(0); 544289bc588SBarry Smith } 545289bc588SBarry Smith 546289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5474a2ae208SSatish Balay #undef __FUNCT__ 5484a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 549dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 550289bc588SBarry Smith { 551c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 55287828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 553dfbe8321SBarry Smith PetscErrorCode ierr; 5540805154bSBarry Smith PetscBLASInt m, n,_One=1; 555ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 5563a40ed3dSBarry Smith 5573a40ed3dSBarry Smith PetscFunctionBegin; 558d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 559d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 560d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5611ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5621ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 56371044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One); 5641ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5651ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 566dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 5673a40ed3dSBarry Smith PetscFunctionReturn(0); 568289bc588SBarry Smith } 569800995b7SMatthew Knepley 5704a2ae208SSatish Balay #undef __FUNCT__ 5714a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 572dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 573289bc588SBarry Smith { 574c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 57587828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 576dfbe8321SBarry Smith PetscErrorCode ierr; 5770805154bSBarry Smith PetscBLASInt m, n, _One=1; 5783a40ed3dSBarry Smith 5793a40ed3dSBarry Smith PetscFunctionBegin; 580d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 581d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 582d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5831ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5841ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 58571044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 5861ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5871ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 588dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr); 5893a40ed3dSBarry Smith PetscFunctionReturn(0); 590289bc588SBarry Smith } 5916ee01492SSatish Balay 5924a2ae208SSatish Balay #undef __FUNCT__ 5934a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 594dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 595289bc588SBarry Smith { 596c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 59787828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 598dfbe8321SBarry Smith PetscErrorCode ierr; 5990805154bSBarry Smith PetscBLASInt m, n, _One=1; 6003a40ed3dSBarry Smith 6013a40ed3dSBarry Smith PetscFunctionBegin; 602d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 603d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 604d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 605600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6061ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6071ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 60871044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 6091ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6101ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 611dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6123a40ed3dSBarry Smith PetscFunctionReturn(0); 613289bc588SBarry Smith } 6146ee01492SSatish Balay 6154a2ae208SSatish Balay #undef __FUNCT__ 6164a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 617dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 618289bc588SBarry Smith { 619c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 62087828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 621dfbe8321SBarry Smith PetscErrorCode ierr; 6220805154bSBarry Smith PetscBLASInt m, n, _One=1; 62387828ca2SBarry Smith PetscScalar _DOne=1.0; 6243a40ed3dSBarry Smith 6253a40ed3dSBarry Smith PetscFunctionBegin; 626d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 627d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 628d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 629600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6301ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6311ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 63271044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 6331ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6341ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 635dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6363a40ed3dSBarry Smith PetscFunctionReturn(0); 637289bc588SBarry Smith } 638289bc588SBarry Smith 639289bc588SBarry Smith /* -----------------------------------------------------------------*/ 6404a2ae208SSatish Balay #undef __FUNCT__ 6414a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 64213f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 643289bc588SBarry Smith { 644c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 64587828ca2SBarry Smith PetscScalar *v; 6466849ba73SBarry Smith PetscErrorCode ierr; 64713f74950SBarry Smith PetscInt i; 64867e560aaSBarry Smith 6493a40ed3dSBarry Smith PetscFunctionBegin; 650d0f46423SBarry Smith *ncols = A->cmap->n; 651289bc588SBarry Smith if (cols) { 652d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr); 653d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) (*cols)[i] = i; 654289bc588SBarry Smith } 655289bc588SBarry Smith if (vals) { 656d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 657289bc588SBarry Smith v = mat->v + row; 658d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;} 659289bc588SBarry Smith } 6603a40ed3dSBarry Smith PetscFunctionReturn(0); 661289bc588SBarry Smith } 6626ee01492SSatish Balay 6634a2ae208SSatish Balay #undef __FUNCT__ 6644a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 66513f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 666289bc588SBarry Smith { 667dfbe8321SBarry Smith PetscErrorCode ierr; 668606d414cSSatish Balay PetscFunctionBegin; 669606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 670606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 6713a40ed3dSBarry Smith PetscFunctionReturn(0); 672289bc588SBarry Smith } 673289bc588SBarry Smith /* ----------------------------------------------------------------*/ 6744a2ae208SSatish Balay #undef __FUNCT__ 6754a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 67613f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 677289bc588SBarry Smith { 678c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 67913f74950SBarry Smith PetscInt i,j,idx=0; 680d6dfbf8fSBarry Smith 6813a40ed3dSBarry Smith PetscFunctionBegin; 68271fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 683289bc588SBarry Smith if (!mat->roworiented) { 684dbb450caSBarry Smith if (addv == INSERT_VALUES) { 685289bc588SBarry Smith for (j=0; j<n; j++) { 686cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6872515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 688e32f2f54SBarry 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); 68958804f6dSBarry Smith #endif 690289bc588SBarry Smith for (i=0; i<m; i++) { 691cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6922515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 693e32f2f54SBarry 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); 69458804f6dSBarry Smith #endif 695cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 696289bc588SBarry Smith } 697289bc588SBarry Smith } 6983a40ed3dSBarry Smith } else { 699289bc588SBarry Smith for (j=0; j<n; j++) { 700cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 7012515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 702e32f2f54SBarry 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); 70358804f6dSBarry Smith #endif 704289bc588SBarry Smith for (i=0; i<m; i++) { 705cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 7062515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 707e32f2f54SBarry 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); 70858804f6dSBarry Smith #endif 709cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 710289bc588SBarry Smith } 711289bc588SBarry Smith } 712289bc588SBarry Smith } 7133a40ed3dSBarry Smith } else { 714dbb450caSBarry Smith if (addv == INSERT_VALUES) { 715e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 716cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7172515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 718e32f2f54SBarry 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); 71958804f6dSBarry Smith #endif 720e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 721cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7222515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 723e32f2f54SBarry 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); 72458804f6dSBarry Smith #endif 725cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 726e8d4e0b9SBarry Smith } 727e8d4e0b9SBarry Smith } 7283a40ed3dSBarry Smith } else { 729289bc588SBarry Smith for (i=0; i<m; i++) { 730cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7312515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 732e32f2f54SBarry 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); 73358804f6dSBarry Smith #endif 734289bc588SBarry Smith for (j=0; j<n; j++) { 735cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7362515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 737e32f2f54SBarry 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); 73858804f6dSBarry Smith #endif 739cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 740289bc588SBarry Smith } 741289bc588SBarry Smith } 742289bc588SBarry Smith } 743e8d4e0b9SBarry Smith } 7443a40ed3dSBarry Smith PetscFunctionReturn(0); 745289bc588SBarry Smith } 746e8d4e0b9SBarry Smith 7474a2ae208SSatish Balay #undef __FUNCT__ 7484a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 74913f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 750ae80bb75SLois Curfman McInnes { 751ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 75213f74950SBarry Smith PetscInt i,j; 753ae80bb75SLois Curfman McInnes 7543a40ed3dSBarry Smith PetscFunctionBegin; 755ae80bb75SLois Curfman McInnes /* row-oriented output */ 756ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 75797e567efSBarry Smith if (indexm[i] < 0) {v += n;continue;} 758e32f2f54SBarry 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); 759ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 7606f31f424SBarry Smith if (indexn[j] < 0) {v++; continue;} 761e32f2f54SBarry 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); 76297e567efSBarry Smith *v++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 763ae80bb75SLois Curfman McInnes } 764ae80bb75SLois Curfman McInnes } 7653a40ed3dSBarry Smith PetscFunctionReturn(0); 766ae80bb75SLois Curfman McInnes } 767ae80bb75SLois Curfman McInnes 768289bc588SBarry Smith /* -----------------------------------------------------------------*/ 769289bc588SBarry Smith 7704a2ae208SSatish Balay #undef __FUNCT__ 7715bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_SeqDense" 772112444f4SShri Abhyankar PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer) 773aabbc4fbSShri Abhyankar { 774aabbc4fbSShri Abhyankar Mat_SeqDense *a; 775aabbc4fbSShri Abhyankar PetscErrorCode ierr; 776aabbc4fbSShri Abhyankar PetscInt *scols,i,j,nz,header[4]; 777aabbc4fbSShri Abhyankar int fd; 778aabbc4fbSShri Abhyankar PetscMPIInt size; 779aabbc4fbSShri Abhyankar PetscInt *rowlengths = 0,M,N,*cols,grows,gcols; 780aabbc4fbSShri Abhyankar PetscScalar *vals,*svals,*v,*w; 781aabbc4fbSShri Abhyankar MPI_Comm comm = ((PetscObject)viewer)->comm; 782aabbc4fbSShri Abhyankar 783aabbc4fbSShri Abhyankar PetscFunctionBegin; 784aabbc4fbSShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 785aabbc4fbSShri Abhyankar if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor"); 786aabbc4fbSShri Abhyankar ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 787aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 788aabbc4fbSShri Abhyankar if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 789aabbc4fbSShri Abhyankar M = header[1]; N = header[2]; nz = header[3]; 790aabbc4fbSShri Abhyankar 791aabbc4fbSShri Abhyankar /* set global size if not set already*/ 792aabbc4fbSShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 793aabbc4fbSShri Abhyankar ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr); 794aabbc4fbSShri Abhyankar } else { 795aabbc4fbSShri Abhyankar /* if sizes and type are already set, check if the vector global sizes are correct */ 796aabbc4fbSShri Abhyankar ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr); 797aabbc4fbSShri 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); 798aabbc4fbSShri Abhyankar } 799aabbc4fbSShri Abhyankar ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 800aabbc4fbSShri Abhyankar 801aabbc4fbSShri Abhyankar if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 802aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 803aabbc4fbSShri Abhyankar v = a->v; 804aabbc4fbSShri Abhyankar /* Allocate some temp space to read in the values and then flip them 805aabbc4fbSShri Abhyankar from row major to column major */ 806aabbc4fbSShri Abhyankar ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 807aabbc4fbSShri Abhyankar /* read in nonzero values */ 808aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 809aabbc4fbSShri Abhyankar /* now flip the values and store them in the matrix*/ 810aabbc4fbSShri Abhyankar for (j=0; j<N; j++) { 811aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 812aabbc4fbSShri Abhyankar *v++ =w[i*N+j]; 813aabbc4fbSShri Abhyankar } 814aabbc4fbSShri Abhyankar } 815aabbc4fbSShri Abhyankar ierr = PetscFree(w);CHKERRQ(ierr); 816aabbc4fbSShri Abhyankar } else { 817aabbc4fbSShri Abhyankar /* read row lengths */ 818aabbc4fbSShri Abhyankar ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 819aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 820aabbc4fbSShri Abhyankar 821aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 822aabbc4fbSShri Abhyankar v = a->v; 823aabbc4fbSShri Abhyankar 824aabbc4fbSShri Abhyankar /* read column indices and nonzeros */ 825aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr); 826aabbc4fbSShri Abhyankar cols = scols; 827aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 828aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 829aabbc4fbSShri Abhyankar vals = svals; 830aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 831aabbc4fbSShri Abhyankar 832aabbc4fbSShri Abhyankar /* insert into matrix */ 833aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 834aabbc4fbSShri Abhyankar for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 835aabbc4fbSShri Abhyankar svals += rowlengths[i]; scols += rowlengths[i]; 836aabbc4fbSShri Abhyankar } 837aabbc4fbSShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 838aabbc4fbSShri Abhyankar ierr = PetscFree(cols);CHKERRQ(ierr); 839aabbc4fbSShri Abhyankar ierr = PetscFree(rowlengths);CHKERRQ(ierr); 840aabbc4fbSShri Abhyankar } 841aabbc4fbSShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 842aabbc4fbSShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 843aabbc4fbSShri Abhyankar 844aabbc4fbSShri Abhyankar PetscFunctionReturn(0); 845aabbc4fbSShri Abhyankar } 846aabbc4fbSShri Abhyankar 847aabbc4fbSShri Abhyankar #undef __FUNCT__ 8484a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 8496849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 850289bc588SBarry Smith { 851932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 852dfbe8321SBarry Smith PetscErrorCode ierr; 85313f74950SBarry Smith PetscInt i,j; 8542dcb1b2aSMatthew Knepley const char *name; 85587828ca2SBarry Smith PetscScalar *v; 856f3ef73ceSBarry Smith PetscViewerFormat format; 8575f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX) 858ace3abfcSBarry Smith PetscBool allreal = PETSC_TRUE; 8595f481a85SSatish Balay #endif 860932b0c3eSLois Curfman McInnes 8613a40ed3dSBarry Smith PetscFunctionBegin; 862b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 863456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 8643a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 865fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 866d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 8677566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 868d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 86944cd7ae7SLois Curfman McInnes v = a->v + i; 87077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 871d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 872aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 873329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 874a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 875329f5518SBarry Smith } else if (PetscRealPart(*v)) { 876a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr); 8776831982aSBarry Smith } 87880cd9d93SLois Curfman McInnes #else 8796831982aSBarry Smith if (*v) { 880a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr); 8816831982aSBarry Smith } 88280cd9d93SLois Curfman McInnes #endif 8831b807ce4Svictorle v += a->lda; 88480cd9d93SLois Curfman McInnes } 885b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 88680cd9d93SLois Curfman McInnes } 887d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 8883a40ed3dSBarry Smith } else { 889d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 890aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 89147989497SBarry Smith /* determine if matrix has all real values */ 89247989497SBarry Smith v = a->v; 893d0f46423SBarry Smith for (i=0; i<A->rmap->n*A->cmap->n; i++) { 894ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 89547989497SBarry Smith } 89647989497SBarry Smith #endif 897fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 8983a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 899d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr); 900d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 901fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 9027566de4bSShri Abhyankar } else { 9037566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 904ffac6cdbSBarry Smith } 905ffac6cdbSBarry Smith 906d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 907932b0c3eSLois Curfman McInnes v = a->v + i; 908d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 909aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 91047989497SBarry Smith if (allreal) { 911f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr); 91247989497SBarry Smith } else { 913f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 91447989497SBarry Smith } 915289bc588SBarry Smith #else 916f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr); 917289bc588SBarry Smith #endif 9181b807ce4Svictorle v += a->lda; 919289bc588SBarry Smith } 920b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 921289bc588SBarry Smith } 922fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 923b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 924ffac6cdbSBarry Smith } 925d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 926da3a660dSBarry Smith } 927b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 9283a40ed3dSBarry Smith PetscFunctionReturn(0); 929289bc588SBarry Smith } 930289bc588SBarry Smith 9314a2ae208SSatish Balay #undef __FUNCT__ 9324a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 9336849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 934932b0c3eSLois Curfman McInnes { 935932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9366849ba73SBarry Smith PetscErrorCode ierr; 93713f74950SBarry Smith int fd; 938d0f46423SBarry Smith PetscInt ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n; 939f4403165SShri Abhyankar PetscScalar *v,*anonz,*vals; 940f4403165SShri Abhyankar PetscViewerFormat format; 941932b0c3eSLois Curfman McInnes 9423a40ed3dSBarry Smith PetscFunctionBegin; 943b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 94490ace30eSBarry Smith 945f4403165SShri Abhyankar ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 946f4403165SShri Abhyankar if (format == PETSC_VIEWER_NATIVE) { 947f4403165SShri Abhyankar /* store the matrix as a dense matrix */ 948f4403165SShri Abhyankar ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 949f4403165SShri Abhyankar col_lens[0] = MAT_FILE_CLASSID; 950f4403165SShri Abhyankar col_lens[1] = m; 951f4403165SShri Abhyankar col_lens[2] = n; 952f4403165SShri Abhyankar col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 953f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 954f4403165SShri Abhyankar ierr = PetscFree(col_lens);CHKERRQ(ierr); 955f4403165SShri Abhyankar 956f4403165SShri Abhyankar /* write out matrix, by rows */ 957f4403165SShri Abhyankar ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 958f4403165SShri Abhyankar v = a->v; 959f4403165SShri Abhyankar for (j=0; j<n; j++) { 960f4403165SShri Abhyankar for (i=0; i<m; i++) { 961f4403165SShri Abhyankar vals[j + i*n] = *v++; 962f4403165SShri Abhyankar } 963f4403165SShri Abhyankar } 964f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 965f4403165SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 966f4403165SShri Abhyankar } else { 96713f74950SBarry Smith ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 9680700a824SBarry Smith col_lens[0] = MAT_FILE_CLASSID; 969932b0c3eSLois Curfman McInnes col_lens[1] = m; 970932b0c3eSLois Curfman McInnes col_lens[2] = n; 971932b0c3eSLois Curfman McInnes col_lens[3] = nz; 972932b0c3eSLois Curfman McInnes 973932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 974932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 9756f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 976932b0c3eSLois Curfman McInnes 977932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 978932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 979932b0c3eSLois Curfman McInnes ict = 0; 980932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 981932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 982932b0c3eSLois Curfman McInnes } 9836f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 984606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 985932b0c3eSLois Curfman McInnes 986932b0c3eSLois Curfman McInnes /* store nonzero values */ 98787828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 988932b0c3eSLois Curfman McInnes ict = 0; 989932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 990932b0c3eSLois Curfman McInnes v = a->v + i; 991932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 9921b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 993932b0c3eSLois Curfman McInnes } 994932b0c3eSLois Curfman McInnes } 9956f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 996606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 997f4403165SShri Abhyankar } 9983a40ed3dSBarry Smith PetscFunctionReturn(0); 999932b0c3eSLois Curfman McInnes } 1000932b0c3eSLois Curfman McInnes 10014a2ae208SSatish Balay #undef __FUNCT__ 10024a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 1003dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 1004f1af5d2fSBarry Smith { 1005f1af5d2fSBarry Smith Mat A = (Mat) Aa; 1006f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 10076849ba73SBarry Smith PetscErrorCode ierr; 1008d0f46423SBarry Smith PetscInt m = A->rmap->n,n = A->cmap->n,color,i,j; 100987828ca2SBarry Smith PetscScalar *v = a->v; 1010b0a32e0cSBarry Smith PetscViewer viewer; 1011b0a32e0cSBarry Smith PetscDraw popup; 1012329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 1013f3ef73ceSBarry Smith PetscViewerFormat format; 1014f1af5d2fSBarry Smith 1015f1af5d2fSBarry Smith PetscFunctionBegin; 1016f1af5d2fSBarry Smith 1017f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 1018b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1019b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1020f1af5d2fSBarry Smith 1021f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 1022fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 1023f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 1024b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1025f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1026f1af5d2fSBarry Smith x_l = j; 1027f1af5d2fSBarry Smith x_r = x_l + 1.0; 1028f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1029f1af5d2fSBarry Smith y_l = m - i - 1.0; 1030f1af5d2fSBarry Smith y_r = y_l + 1.0; 1031f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 1032329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 1033b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1034329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 1035b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1036f1af5d2fSBarry Smith } else { 1037f1af5d2fSBarry Smith continue; 1038f1af5d2fSBarry Smith } 1039f1af5d2fSBarry Smith #else 1040f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 1041b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1042f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 1043b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1044f1af5d2fSBarry Smith } else { 1045f1af5d2fSBarry Smith continue; 1046f1af5d2fSBarry Smith } 1047f1af5d2fSBarry Smith #endif 1048b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1049f1af5d2fSBarry Smith } 1050f1af5d2fSBarry Smith } 1051f1af5d2fSBarry Smith } else { 1052f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 1053f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 1054f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 1055f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 1056f1af5d2fSBarry Smith } 1057b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 1058b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 1059b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 1060f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1061f1af5d2fSBarry Smith x_l = j; 1062f1af5d2fSBarry Smith x_r = x_l + 1.0; 1063f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1064f1af5d2fSBarry Smith y_l = m - i - 1.0; 1065f1af5d2fSBarry Smith y_r = y_l + 1.0; 1066b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 1067b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1068f1af5d2fSBarry Smith } 1069f1af5d2fSBarry Smith } 1070f1af5d2fSBarry Smith } 1071f1af5d2fSBarry Smith PetscFunctionReturn(0); 1072f1af5d2fSBarry Smith } 1073f1af5d2fSBarry Smith 10744a2ae208SSatish Balay #undef __FUNCT__ 10754a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 1076dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 1077f1af5d2fSBarry Smith { 1078b0a32e0cSBarry Smith PetscDraw draw; 1079ace3abfcSBarry Smith PetscBool isnull; 1080329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 1081dfbe8321SBarry Smith PetscErrorCode ierr; 1082f1af5d2fSBarry Smith 1083f1af5d2fSBarry Smith PetscFunctionBegin; 1084b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1085b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1086abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 1087f1af5d2fSBarry Smith 1088f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1089d0f46423SBarry Smith xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0; 1090f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 1091b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1092b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 1093f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 1094f1af5d2fSBarry Smith PetscFunctionReturn(0); 1095f1af5d2fSBarry Smith } 1096f1af5d2fSBarry Smith 10974a2ae208SSatish Balay #undef __FUNCT__ 10984a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 1099dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1100932b0c3eSLois Curfman McInnes { 1101dfbe8321SBarry Smith PetscErrorCode ierr; 1102ace3abfcSBarry Smith PetscBool iascii,isbinary,isdraw; 1103932b0c3eSLois Curfman McInnes 11043a40ed3dSBarry Smith PetscFunctionBegin; 11052692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 11062692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 11072692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 11080f5bd95cSBarry Smith 1109c45a1595SBarry Smith if (iascii) { 1110c45a1595SBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 11110f5bd95cSBarry Smith } else if (isbinary) { 11123a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1113f1af5d2fSBarry Smith } else if (isdraw) { 1114f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 11155cd90555SBarry Smith } else { 1116e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1117932b0c3eSLois Curfman McInnes } 11183a40ed3dSBarry Smith PetscFunctionReturn(0); 1119932b0c3eSLois Curfman McInnes } 1120289bc588SBarry Smith 11214a2ae208SSatish Balay #undef __FUNCT__ 11224a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1123dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1124289bc588SBarry Smith { 1125ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1126dfbe8321SBarry Smith PetscErrorCode ierr; 112790f02eecSBarry Smith 11283a40ed3dSBarry Smith PetscFunctionBegin; 1129aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1130d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n); 1131a5a9c739SBarry Smith #endif 113205b42c5fSBarry Smith ierr = PetscFree(l->pivots);CHKERRQ(ierr); 11336857c123SSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1134bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 1135dbd8c25aSHong Zhang 1136dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 1137901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 11384ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11394ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11404ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11413a40ed3dSBarry Smith PetscFunctionReturn(0); 1142289bc588SBarry Smith } 1143289bc588SBarry Smith 11444a2ae208SSatish Balay #undef __FUNCT__ 11454a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1146fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout) 1147289bc588SBarry Smith { 1148c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11496849ba73SBarry Smith PetscErrorCode ierr; 115013f74950SBarry Smith PetscInt k,j,m,n,M; 115187828ca2SBarry Smith PetscScalar *v,tmp; 115248b35521SBarry Smith 11533a40ed3dSBarry Smith PetscFunctionBegin; 1154d0f46423SBarry Smith v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n; 1155e9695a30SBarry Smith if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */ 1156e7e72b3dSBarry Smith if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1157e7e72b3dSBarry Smith else { 1158d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1159289bc588SBarry Smith for (k=0; k<j; k++) { 11601b807ce4Svictorle tmp = v[j + k*M]; 11611b807ce4Svictorle v[j + k*M] = v[k + j*M]; 11621b807ce4Svictorle v[k + j*M] = tmp; 1163289bc588SBarry Smith } 1164289bc588SBarry Smith } 1165d64ed03dSBarry Smith } 11663a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1167d3e5ee88SLois Curfman McInnes Mat tmat; 1168ec8511deSBarry Smith Mat_SeqDense *tmatd; 116987828ca2SBarry Smith PetscScalar *v2; 1170ea709b57SSatish Balay 1171fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 11727adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr); 1173d0f46423SBarry Smith ierr = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr); 11747adad957SLisandro Dalcin ierr = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 11755c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1176fc4dec0aSBarry Smith } else { 1177fc4dec0aSBarry Smith tmat = *matout; 1178fc4dec0aSBarry Smith } 1179ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 11800de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1181d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 11821b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1183d3e5ee88SLois Curfman McInnes } 11846d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11856d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1186d3e5ee88SLois Curfman McInnes *matout = tmat; 118748b35521SBarry Smith } 11883a40ed3dSBarry Smith PetscFunctionReturn(0); 1189289bc588SBarry Smith } 1190289bc588SBarry Smith 11914a2ae208SSatish Balay #undef __FUNCT__ 11924a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1193ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool *flg) 1194289bc588SBarry Smith { 1195c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1196c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 119713f74950SBarry Smith PetscInt i,j; 119887828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 11999ea5d5aeSSatish Balay 12003a40ed3dSBarry Smith PetscFunctionBegin; 1201d0f46423SBarry Smith if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1202d0f46423SBarry Smith if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1203d0f46423SBarry Smith for (i=0; i<A1->rmap->n; i++) { 12041b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1205d0f46423SBarry Smith for (j=0; j<A1->cmap->n; j++) { 12063a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 12071b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 12081b807ce4Svictorle } 1209289bc588SBarry Smith } 121077c4ece6SBarry Smith *flg = PETSC_TRUE; 12113a40ed3dSBarry Smith PetscFunctionReturn(0); 1212289bc588SBarry Smith } 1213289bc588SBarry Smith 12144a2ae208SSatish Balay #undef __FUNCT__ 12154a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1216dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1217289bc588SBarry Smith { 1218c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1219dfbe8321SBarry Smith PetscErrorCode ierr; 122013f74950SBarry Smith PetscInt i,n,len; 122187828ca2SBarry Smith PetscScalar *x,zero = 0.0; 122244cd7ae7SLois Curfman McInnes 12233a40ed3dSBarry Smith PetscFunctionBegin; 12242dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 12257a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 12261ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1227d0f46423SBarry Smith len = PetscMin(A->rmap->n,A->cmap->n); 1228e32f2f54SBarry Smith if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 122944cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 12301b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1231289bc588SBarry Smith } 12321ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 12333a40ed3dSBarry Smith PetscFunctionReturn(0); 1234289bc588SBarry Smith } 1235289bc588SBarry Smith 12364a2ae208SSatish Balay #undef __FUNCT__ 12374a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1238dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1239289bc588SBarry Smith { 1240c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 124187828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1242dfbe8321SBarry Smith PetscErrorCode ierr; 1243d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n; 124455659b69SBarry Smith 12453a40ed3dSBarry Smith PetscFunctionBegin; 124628988994SBarry Smith if (ll) { 12477a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 12481ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1249e32f2f54SBarry Smith if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1250da3a660dSBarry Smith for (i=0; i<m; i++) { 1251da3a660dSBarry Smith x = l[i]; 1252da3a660dSBarry Smith v = mat->v + i; 1253da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1254da3a660dSBarry Smith } 12551ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1256efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1257da3a660dSBarry Smith } 125828988994SBarry Smith if (rr) { 12597a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 12601ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1261e32f2f54SBarry Smith if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1262da3a660dSBarry Smith for (i=0; i<n; i++) { 1263da3a660dSBarry Smith x = r[i]; 1264da3a660dSBarry Smith v = mat->v + i*m; 1265da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1266da3a660dSBarry Smith } 12671ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1268efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1269da3a660dSBarry Smith } 12703a40ed3dSBarry Smith PetscFunctionReturn(0); 1271289bc588SBarry Smith } 1272289bc588SBarry Smith 12734a2ae208SSatish Balay #undef __FUNCT__ 12744a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1275dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1276289bc588SBarry Smith { 1277c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 127887828ca2SBarry Smith PetscScalar *v = mat->v; 1279329f5518SBarry Smith PetscReal sum = 0.0; 1280d0f46423SBarry Smith PetscInt lda=mat->lda,m=A->rmap->n,i,j; 1281efee365bSSatish Balay PetscErrorCode ierr; 128255659b69SBarry Smith 12833a40ed3dSBarry Smith PetscFunctionBegin; 1284289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1285a5ce6ee0Svictorle if (lda>m) { 1286d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1287a5ce6ee0Svictorle v = mat->v+j*lda; 1288a5ce6ee0Svictorle for (i=0; i<m; i++) { 1289a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1290a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1291a5ce6ee0Svictorle #else 1292a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1293a5ce6ee0Svictorle #endif 1294a5ce6ee0Svictorle } 1295a5ce6ee0Svictorle } 1296a5ce6ee0Svictorle } else { 1297d0f46423SBarry Smith for (i=0; i<A->cmap->n*A->rmap->n; i++) { 1298aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1299329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1300289bc588SBarry Smith #else 1301289bc588SBarry Smith sum += (*v)*(*v); v++; 1302289bc588SBarry Smith #endif 1303289bc588SBarry Smith } 1304a5ce6ee0Svictorle } 1305064f8208SBarry Smith *nrm = sqrt(sum); 1306dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13073a40ed3dSBarry Smith } else if (type == NORM_1) { 1308064f8208SBarry Smith *nrm = 0.0; 1309d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 13101b807ce4Svictorle v = mat->v + j*mat->lda; 1311289bc588SBarry Smith sum = 0.0; 1312d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 131333a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1314289bc588SBarry Smith } 1315064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1316289bc588SBarry Smith } 1317d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13183a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1319064f8208SBarry Smith *nrm = 0.0; 1320d0f46423SBarry Smith for (j=0; j<A->rmap->n; j++) { 1321289bc588SBarry Smith v = mat->v + j; 1322289bc588SBarry Smith sum = 0.0; 1323d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) { 13241b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1325289bc588SBarry Smith } 1326064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1327289bc588SBarry Smith } 1328d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 1329e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 13303a40ed3dSBarry Smith PetscFunctionReturn(0); 1331289bc588SBarry Smith } 1332289bc588SBarry Smith 13334a2ae208SSatish Balay #undef __FUNCT__ 13344a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 1335ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg) 1336289bc588SBarry Smith { 1337c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 133863ba0a88SBarry Smith PetscErrorCode ierr; 133967e560aaSBarry Smith 13403a40ed3dSBarry Smith PetscFunctionBegin; 1341b5a2b587SKris Buschelman switch (op) { 1342b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 13434e0d8c25SBarry Smith aij->roworiented = flg; 1344b5a2b587SKris Buschelman break; 1345512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 1346b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 13473971808eSMatthew Knepley case MAT_NEW_NONZERO_ALLOCATION_ERR: 13484e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 1349b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1350b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 135177e54ba9SKris Buschelman case MAT_SYMMETRIC: 135277e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 13539a4540c5SBarry Smith case MAT_HERMITIAN: 13549a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1355600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 1356290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 135777e54ba9SKris Buschelman break; 1358b5a2b587SKris Buschelman default: 1359e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 13603a40ed3dSBarry Smith } 13613a40ed3dSBarry Smith PetscFunctionReturn(0); 1362289bc588SBarry Smith } 1363289bc588SBarry Smith 13644a2ae208SSatish Balay #undef __FUNCT__ 13654a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1366dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 13676f0a148fSBarry Smith { 1368ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 13696849ba73SBarry Smith PetscErrorCode ierr; 1370d0f46423SBarry Smith PetscInt lda=l->lda,m=A->rmap->n,j; 13713a40ed3dSBarry Smith 13723a40ed3dSBarry Smith PetscFunctionBegin; 1373a5ce6ee0Svictorle if (lda>m) { 1374d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1375a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1376a5ce6ee0Svictorle } 1377a5ce6ee0Svictorle } else { 1378d0f46423SBarry Smith ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1379a5ce6ee0Svictorle } 13803a40ed3dSBarry Smith PetscFunctionReturn(0); 13816f0a148fSBarry Smith } 13826f0a148fSBarry Smith 13834a2ae208SSatish Balay #undef __FUNCT__ 13844a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 13852b40b63fSBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 13866f0a148fSBarry Smith { 138797b48c8fSBarry Smith PetscErrorCode ierr; 1388ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1389d0f46423SBarry Smith PetscInt n = A->cmap->n,i,j; 139097b48c8fSBarry Smith PetscScalar *slot,*bb; 139197b48c8fSBarry Smith const PetscScalar *xx; 139255659b69SBarry Smith 13933a40ed3dSBarry Smith PetscFunctionBegin; 139497b48c8fSBarry Smith /* fix right hand side if needed */ 139597b48c8fSBarry Smith if (x && b) { 139697b48c8fSBarry Smith ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 139797b48c8fSBarry Smith ierr = VecGetArray(b,&bb);CHKERRQ(ierr); 139897b48c8fSBarry Smith for (i=0; i<N; i++) { 139997b48c8fSBarry Smith bb[rows[i]] = diag*xx[rows[i]]; 140097b48c8fSBarry Smith } 140197b48c8fSBarry Smith ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 140297b48c8fSBarry Smith ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr); 140397b48c8fSBarry Smith } 140497b48c8fSBarry Smith 14056f0a148fSBarry Smith for (i=0; i<N; i++) { 14066f0a148fSBarry Smith slot = l->v + rows[i]; 14076f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 14086f0a148fSBarry Smith } 1409f4df32b1SMatthew Knepley if (diag != 0.0) { 14106f0a148fSBarry Smith for (i=0; i<N; i++) { 14116f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1412f4df32b1SMatthew Knepley *slot = diag; 14136f0a148fSBarry Smith } 14146f0a148fSBarry Smith } 14153a40ed3dSBarry Smith PetscFunctionReturn(0); 14166f0a148fSBarry Smith } 1417557bce09SLois Curfman McInnes 14184a2ae208SSatish Balay #undef __FUNCT__ 14194a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 1420dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 142164e87e97SBarry Smith { 1422c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14233a40ed3dSBarry Smith 14243a40ed3dSBarry Smith PetscFunctionBegin; 1425e32f2f54SBarry 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"); 142664e87e97SBarry Smith *array = mat->v; 14273a40ed3dSBarry Smith PetscFunctionReturn(0); 142864e87e97SBarry Smith } 14290754003eSLois Curfman McInnes 14304a2ae208SSatish Balay #undef __FUNCT__ 14314a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 1432dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1433ff14e315SSatish Balay { 14343a40ed3dSBarry Smith PetscFunctionBegin; 143509b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 14363a40ed3dSBarry Smith PetscFunctionReturn(0); 1437ff14e315SSatish Balay } 14380754003eSLois Curfman McInnes 14394a2ae208SSatish Balay #undef __FUNCT__ 14404a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 144113f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 14420754003eSLois Curfman McInnes { 1443c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14446849ba73SBarry Smith PetscErrorCode ierr; 14455d0c19d7SBarry Smith PetscInt i,j,nrows,ncols; 14465d0c19d7SBarry Smith const PetscInt *irow,*icol; 144787828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 14480754003eSLois Curfman McInnes Mat newmat; 14490754003eSLois Curfman McInnes 14503a40ed3dSBarry Smith PetscFunctionBegin; 145178b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 145278b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1453e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1454e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 14550754003eSLois Curfman McInnes 1456182d2002SSatish Balay /* Check submatrixcall */ 1457182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 145813f74950SBarry Smith PetscInt n_cols,n_rows; 1459182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 146021a2c019SBarry Smith if (n_rows != nrows || n_cols != ncols) { 1461f746d493SDmitry Karpeev /* resize the result matrix to match number of requested rows/columns */ 1462c61587bbSBarry Smith ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 146321a2c019SBarry Smith } 1464182d2002SSatish Balay newmat = *B; 1465182d2002SSatish Balay } else { 14660754003eSLois Curfman McInnes /* Create and fill new matrix */ 14677adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr); 1468f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 14697adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 14705c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1471182d2002SSatish Balay } 1472182d2002SSatish Balay 1473182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1474182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1475182d2002SSatish Balay 1476182d2002SSatish Balay for (i=0; i<ncols; i++) { 14776de62eeeSBarry Smith av = v + mat->lda*icol[i]; 1478182d2002SSatish Balay for (j=0; j<nrows; j++) { 1479182d2002SSatish Balay *bv++ = av[irow[j]]; 14800754003eSLois Curfman McInnes } 14810754003eSLois Curfman McInnes } 1482182d2002SSatish Balay 1483182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 14846d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14856d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14860754003eSLois Curfman McInnes 14870754003eSLois Curfman McInnes /* Free work space */ 148878b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 148978b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1490182d2002SSatish Balay *B = newmat; 14913a40ed3dSBarry Smith PetscFunctionReturn(0); 14920754003eSLois Curfman McInnes } 14930754003eSLois Curfman McInnes 14944a2ae208SSatish Balay #undef __FUNCT__ 14954a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 149613f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1497905e6a2fSBarry Smith { 14986849ba73SBarry Smith PetscErrorCode ierr; 149913f74950SBarry Smith PetscInt i; 1500905e6a2fSBarry Smith 15013a40ed3dSBarry Smith PetscFunctionBegin; 1502905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1503b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1504905e6a2fSBarry Smith } 1505905e6a2fSBarry Smith 1506905e6a2fSBarry Smith for (i=0; i<n; i++) { 15076a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1508905e6a2fSBarry Smith } 15093a40ed3dSBarry Smith PetscFunctionReturn(0); 1510905e6a2fSBarry Smith } 1511905e6a2fSBarry Smith 15124a2ae208SSatish Balay #undef __FUNCT__ 1513c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1514c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1515c0aa2d19SHong Zhang { 1516c0aa2d19SHong Zhang PetscFunctionBegin; 1517c0aa2d19SHong Zhang PetscFunctionReturn(0); 1518c0aa2d19SHong Zhang } 1519c0aa2d19SHong Zhang 1520c0aa2d19SHong Zhang #undef __FUNCT__ 1521c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1522c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1523c0aa2d19SHong Zhang { 1524c0aa2d19SHong Zhang PetscFunctionBegin; 1525c0aa2d19SHong Zhang PetscFunctionReturn(0); 1526c0aa2d19SHong Zhang } 1527c0aa2d19SHong Zhang 1528c0aa2d19SHong Zhang #undef __FUNCT__ 15294a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1530dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 15314b0e389bSBarry Smith { 15324b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 15336849ba73SBarry Smith PetscErrorCode ierr; 1534d0f46423SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j; 15353a40ed3dSBarry Smith 15363a40ed3dSBarry Smith PetscFunctionBegin; 153733f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 153833f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1539cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 15403a40ed3dSBarry Smith PetscFunctionReturn(0); 15413a40ed3dSBarry Smith } 1542e32f2f54SBarry Smith if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1543a5ce6ee0Svictorle if (lda1>m || lda2>m) { 15440dbb7854Svictorle for (j=0; j<n; j++) { 1545a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1546a5ce6ee0Svictorle } 1547a5ce6ee0Svictorle } else { 1548d0f46423SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1549a5ce6ee0Svictorle } 1550273d9f13SBarry Smith PetscFunctionReturn(0); 1551273d9f13SBarry Smith } 1552273d9f13SBarry Smith 15534a2ae208SSatish Balay #undef __FUNCT__ 15544a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1555dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A) 1556273d9f13SBarry Smith { 1557dfbe8321SBarry Smith PetscErrorCode ierr; 1558273d9f13SBarry Smith 1559273d9f13SBarry Smith PetscFunctionBegin; 1560273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 15613a40ed3dSBarry Smith PetscFunctionReturn(0); 15624b0e389bSBarry Smith } 15634b0e389bSBarry Smith 1564284134d9SBarry Smith #undef __FUNCT__ 1565284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense" 1566284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N) 1567284134d9SBarry Smith { 1568284134d9SBarry Smith PetscFunctionBegin; 156921a2c019SBarry Smith /* this will not be called before lda, Mmax, and Nmax have been set */ 1570284134d9SBarry Smith m = PetscMax(m,M); 1571284134d9SBarry Smith n = PetscMax(n,N); 1572a868139aSShri Abhyankar 157386d161a7SShri 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); 157486d161a7SShri 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); 157586d161a7SShri Abhyankar */ 1576dc5cefdeSJed Brown A->rmap->n = A->rmap->N = m; 1577d0f46423SBarry Smith A->cmap->n = A->cmap->N = n; 1578284134d9SBarry Smith PetscFunctionReturn(0); 1579284134d9SBarry Smith } 1580170fe5c8SBarry Smith 1581ba337c44SJed Brown #undef __FUNCT__ 1582ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense" 1583ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A) 1584ba337c44SJed Brown { 1585ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1586ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1587ba337c44SJed Brown PetscScalar *aa = a->v; 1588ba337c44SJed Brown 1589ba337c44SJed Brown PetscFunctionBegin; 1590ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]); 1591ba337c44SJed Brown PetscFunctionReturn(0); 1592ba337c44SJed Brown } 1593ba337c44SJed Brown 1594ba337c44SJed Brown #undef __FUNCT__ 1595ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense" 1596ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A) 1597ba337c44SJed Brown { 1598ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1599ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1600ba337c44SJed Brown PetscScalar *aa = a->v; 1601ba337c44SJed Brown 1602ba337c44SJed Brown PetscFunctionBegin; 1603ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 1604ba337c44SJed Brown PetscFunctionReturn(0); 1605ba337c44SJed Brown } 1606ba337c44SJed Brown 1607ba337c44SJed Brown #undef __FUNCT__ 1608ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense" 1609ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A) 1610ba337c44SJed Brown { 1611ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1612ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1613ba337c44SJed Brown PetscScalar *aa = a->v; 1614ba337c44SJed Brown 1615ba337c44SJed Brown PetscFunctionBegin; 1616ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 1617ba337c44SJed Brown PetscFunctionReturn(0); 1618ba337c44SJed Brown } 1619284134d9SBarry Smith 1620a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/ 1621a9fe9ddaSSatish Balay #undef __FUNCT__ 1622a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1623a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1624a9fe9ddaSSatish Balay { 1625a9fe9ddaSSatish Balay PetscErrorCode ierr; 1626a9fe9ddaSSatish Balay 1627a9fe9ddaSSatish Balay PetscFunctionBegin; 1628a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1629a9fe9ddaSSatish Balay ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1630a9fe9ddaSSatish Balay } 1631a9fe9ddaSSatish Balay ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1632a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1633a9fe9ddaSSatish Balay } 1634a9fe9ddaSSatish Balay 1635a9fe9ddaSSatish Balay #undef __FUNCT__ 1636a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1637a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1638a9fe9ddaSSatish Balay { 1639ee16a9a1SHong Zhang PetscErrorCode ierr; 1640d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 1641ee16a9a1SHong Zhang Mat Cmat; 1642a9fe9ddaSSatish Balay 1643ee16a9a1SHong Zhang PetscFunctionBegin; 1644e32f2f54SBarry 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); 1645ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1646ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1647ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1648ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1649ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1650ee16a9a1SHong Zhang *C = Cmat; 1651ee16a9a1SHong Zhang PetscFunctionReturn(0); 1652ee16a9a1SHong Zhang } 1653a9fe9ddaSSatish Balay 165498a3b096SSatish Balay #undef __FUNCT__ 1655a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1656a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1657a9fe9ddaSSatish Balay { 1658a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1659a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1660a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 16610805154bSBarry Smith PetscBLASInt m,n,k; 1662a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1663a9fe9ddaSSatish Balay 1664a9fe9ddaSSatish Balay PetscFunctionBegin; 1665d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 1666d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1667d0f46423SBarry Smith k = PetscBLASIntCast(A->cmap->n); 1668a9fe9ddaSSatish Balay BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1669a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1670a9fe9ddaSSatish Balay } 1671a9fe9ddaSSatish Balay 1672a9fe9ddaSSatish Balay #undef __FUNCT__ 1673a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense" 1674a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1675a9fe9ddaSSatish Balay { 1676a9fe9ddaSSatish Balay PetscErrorCode ierr; 1677a9fe9ddaSSatish Balay 1678a9fe9ddaSSatish Balay PetscFunctionBegin; 1679a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1680a9fe9ddaSSatish Balay ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1681a9fe9ddaSSatish Balay } 1682a9fe9ddaSSatish Balay ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1683a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1684a9fe9ddaSSatish Balay } 1685a9fe9ddaSSatish Balay 1686a9fe9ddaSSatish Balay #undef __FUNCT__ 1687a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense" 1688a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1689a9fe9ddaSSatish Balay { 1690ee16a9a1SHong Zhang PetscErrorCode ierr; 1691d0f46423SBarry Smith PetscInt m=A->cmap->n,n=B->cmap->n; 1692ee16a9a1SHong Zhang Mat Cmat; 1693a9fe9ddaSSatish Balay 1694ee16a9a1SHong Zhang PetscFunctionBegin; 1695e32f2f54SBarry 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); 1696ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1697ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1698ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1699ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1700ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1701ee16a9a1SHong Zhang *C = Cmat; 1702ee16a9a1SHong Zhang PetscFunctionReturn(0); 1703ee16a9a1SHong Zhang } 1704a9fe9ddaSSatish Balay 1705a9fe9ddaSSatish Balay #undef __FUNCT__ 1706a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense" 1707a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1708a9fe9ddaSSatish Balay { 1709a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1710a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1711a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 17120805154bSBarry Smith PetscBLASInt m,n,k; 1713a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1714a9fe9ddaSSatish Balay 1715a9fe9ddaSSatish Balay PetscFunctionBegin; 1716d0f46423SBarry Smith m = PetscBLASIntCast(A->cmap->n); 1717d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1718d0f46423SBarry Smith k = PetscBLASIntCast(A->rmap->n); 17192fbe02b9SBarry Smith /* 17202fbe02b9SBarry Smith Note the m and n arguments below are the number rows and columns of A', not A! 17212fbe02b9SBarry Smith */ 1722a9fe9ddaSSatish Balay BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1723a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1724a9fe9ddaSSatish Balay } 1725985db425SBarry Smith 1726985db425SBarry Smith #undef __FUNCT__ 1727985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense" 1728985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 1729985db425SBarry Smith { 1730985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1731985db425SBarry Smith PetscErrorCode ierr; 1732d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1733985db425SBarry Smith PetscScalar *x; 1734985db425SBarry Smith MatScalar *aa = a->v; 1735985db425SBarry Smith 1736985db425SBarry Smith PetscFunctionBegin; 1737e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1738985db425SBarry Smith 1739985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1740985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1741985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1742e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1743985db425SBarry Smith for (i=0; i<m; i++) { 1744985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1745985db425SBarry Smith for (j=1; j<n; j++){ 1746985db425SBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1747985db425SBarry Smith } 1748985db425SBarry Smith } 1749985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1750985db425SBarry Smith PetscFunctionReturn(0); 1751985db425SBarry Smith } 1752985db425SBarry Smith 1753985db425SBarry Smith #undef __FUNCT__ 1754985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense" 1755985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 1756985db425SBarry Smith { 1757985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1758985db425SBarry Smith PetscErrorCode ierr; 1759d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1760985db425SBarry Smith PetscScalar *x; 1761985db425SBarry Smith PetscReal atmp; 1762985db425SBarry Smith MatScalar *aa = a->v; 1763985db425SBarry Smith 1764985db425SBarry Smith PetscFunctionBegin; 1765e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1766985db425SBarry Smith 1767985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1768985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1769985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1770e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1771985db425SBarry Smith for (i=0; i<m; i++) { 17729189402eSHong Zhang x[i] = PetscAbsScalar(aa[i]); 1773985db425SBarry Smith for (j=1; j<n; j++){ 1774985db425SBarry Smith atmp = PetscAbsScalar(aa[i+m*j]); 1775985db425SBarry Smith if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 1776985db425SBarry Smith } 1777985db425SBarry Smith } 1778985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1779985db425SBarry Smith PetscFunctionReturn(0); 1780985db425SBarry Smith } 1781985db425SBarry Smith 1782985db425SBarry Smith #undef __FUNCT__ 1783985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense" 1784985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 1785985db425SBarry Smith { 1786985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1787985db425SBarry Smith PetscErrorCode ierr; 1788d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1789985db425SBarry Smith PetscScalar *x; 1790985db425SBarry Smith MatScalar *aa = a->v; 1791985db425SBarry Smith 1792985db425SBarry Smith PetscFunctionBegin; 1793e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1794985db425SBarry Smith 1795985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1796985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1797985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1798e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1799985db425SBarry Smith for (i=0; i<m; i++) { 1800985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1801985db425SBarry Smith for (j=1; j<n; j++){ 1802985db425SBarry Smith if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1803985db425SBarry Smith } 1804985db425SBarry Smith } 1805985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1806985db425SBarry Smith PetscFunctionReturn(0); 1807985db425SBarry Smith } 1808985db425SBarry Smith 18098d0534beSBarry Smith #undef __FUNCT__ 18108d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense" 18118d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 18128d0534beSBarry Smith { 18138d0534beSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 18148d0534beSBarry Smith PetscErrorCode ierr; 18158d0534beSBarry Smith PetscScalar *x; 18168d0534beSBarry Smith 18178d0534beSBarry Smith PetscFunctionBegin; 1818e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 18198d0534beSBarry Smith 18208d0534beSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1821d0f46423SBarry Smith ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 18228d0534beSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 18238d0534beSBarry Smith PetscFunctionReturn(0); 18248d0534beSBarry Smith } 18258d0534beSBarry Smith 18260716a85fSBarry Smith 18270716a85fSBarry Smith #undef __FUNCT__ 18280716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqDense" 18290716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms) 18300716a85fSBarry Smith { 18310716a85fSBarry Smith PetscErrorCode ierr; 18320716a85fSBarry Smith PetscInt i,j,m,n; 18330716a85fSBarry Smith PetscScalar *a; 18340716a85fSBarry Smith 18350716a85fSBarry Smith PetscFunctionBegin; 18360716a85fSBarry Smith ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 18370716a85fSBarry Smith ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr); 18380716a85fSBarry Smith ierr = MatGetArray(A,&a);CHKERRQ(ierr); 18390716a85fSBarry Smith if (type == NORM_2) { 18400716a85fSBarry Smith for (i=0; i<n; i++ ){ 18410716a85fSBarry Smith for (j=0; j<m; j++) { 18420716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]*a[j]); 18430716a85fSBarry Smith } 18440716a85fSBarry Smith a += m; 18450716a85fSBarry Smith } 18460716a85fSBarry Smith } else if (type == NORM_1) { 18470716a85fSBarry Smith for (i=0; i<n; i++ ){ 18480716a85fSBarry Smith for (j=0; j<m; j++) { 18490716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]); 18500716a85fSBarry Smith } 18510716a85fSBarry Smith a += m; 18520716a85fSBarry Smith } 18530716a85fSBarry Smith } else if (type == NORM_INFINITY) { 18540716a85fSBarry Smith for (i=0; i<n; i++ ){ 18550716a85fSBarry Smith for (j=0; j<m; j++) { 18560716a85fSBarry Smith norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]); 18570716a85fSBarry Smith } 18580716a85fSBarry Smith a += m; 18590716a85fSBarry Smith } 18600716a85fSBarry Smith } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Unknown NormType"); 18610716a85fSBarry Smith if (type == NORM_2) { 18620716a85fSBarry Smith for (i=0; i<n; i++) norms[i] = sqrt(norms[i]); 18630716a85fSBarry Smith } 18640716a85fSBarry Smith PetscFunctionReturn(0); 18650716a85fSBarry Smith } 18660716a85fSBarry Smith 1867289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1868a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1869905e6a2fSBarry Smith MatGetRow_SeqDense, 1870905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1871905e6a2fSBarry Smith MatMult_SeqDense, 187297304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 18737c922b88SBarry Smith MatMultTranspose_SeqDense, 18747c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1875db4efbfdSBarry Smith 0, 1876db4efbfdSBarry Smith 0, 1877db4efbfdSBarry Smith 0, 1878db4efbfdSBarry Smith /*10*/ 0, 1879905e6a2fSBarry Smith MatLUFactor_SeqDense, 1880905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 188141f059aeSBarry Smith MatSOR_SeqDense, 1882ec8511deSBarry Smith MatTranspose_SeqDense, 188397304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1884905e6a2fSBarry Smith MatEqual_SeqDense, 1885905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1886905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1887905e6a2fSBarry Smith MatNorm_SeqDense, 1888c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense, 1889c0aa2d19SHong Zhang MatAssemblyEnd_SeqDense, 1890905e6a2fSBarry Smith MatSetOption_SeqDense, 1891905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1892d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense, 1893db4efbfdSBarry Smith 0, 1894db4efbfdSBarry Smith 0, 1895db4efbfdSBarry Smith 0, 1896db4efbfdSBarry Smith 0, 1897d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqDense, 1898273d9f13SBarry Smith 0, 1899905e6a2fSBarry Smith 0, 1900905e6a2fSBarry Smith MatGetArray_SeqDense, 1901905e6a2fSBarry Smith MatRestoreArray_SeqDense, 1902d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense, 1903a5ae1ecdSBarry Smith 0, 1904a5ae1ecdSBarry Smith 0, 1905a5ae1ecdSBarry Smith 0, 1906a5ae1ecdSBarry Smith 0, 1907d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense, 1908a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1909a5ae1ecdSBarry Smith 0, 19104b0e389bSBarry Smith MatGetValues_SeqDense, 1911a5ae1ecdSBarry Smith MatCopy_SeqDense, 1912d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense, 1913a5ae1ecdSBarry Smith MatScale_SeqDense, 1914a5ae1ecdSBarry Smith 0, 1915a5ae1ecdSBarry Smith 0, 1916a5ae1ecdSBarry Smith 0, 1917d519adbfSMatthew Knepley /*49*/ 0, 1918a5ae1ecdSBarry Smith 0, 1919a5ae1ecdSBarry Smith 0, 1920a5ae1ecdSBarry Smith 0, 1921a5ae1ecdSBarry Smith 0, 1922d519adbfSMatthew Knepley /*54*/ 0, 1923a5ae1ecdSBarry Smith 0, 1924a5ae1ecdSBarry Smith 0, 1925a5ae1ecdSBarry Smith 0, 1926a5ae1ecdSBarry Smith 0, 1927d519adbfSMatthew Knepley /*59*/ 0, 1928e03a110bSBarry Smith MatDestroy_SeqDense, 1929e03a110bSBarry Smith MatView_SeqDense, 1930357abbc8SBarry Smith 0, 193197304618SKris Buschelman 0, 1932d519adbfSMatthew Knepley /*64*/ 0, 193397304618SKris Buschelman 0, 193497304618SKris Buschelman 0, 193597304618SKris Buschelman 0, 193697304618SKris Buschelman 0, 1937d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense, 193897304618SKris Buschelman 0, 193997304618SKris Buschelman 0, 194097304618SKris Buschelman 0, 194197304618SKris Buschelman 0, 1942d519adbfSMatthew Knepley /*74*/ 0, 194397304618SKris Buschelman 0, 194497304618SKris Buschelman 0, 194597304618SKris Buschelman 0, 194697304618SKris Buschelman 0, 1947d519adbfSMatthew Knepley /*79*/ 0, 194897304618SKris Buschelman 0, 194997304618SKris Buschelman 0, 195097304618SKris Buschelman 0, 19515bba2384SShri Abhyankar /*83*/ MatLoad_SeqDense, 1952865e5f61SKris Buschelman 0, 19531cbb95d3SBarry Smith MatIsHermitian_SeqDense, 1954865e5f61SKris Buschelman 0, 1955865e5f61SKris Buschelman 0, 1956865e5f61SKris Buschelman 0, 1957d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense, 1958a9fe9ddaSSatish Balay MatMatMultSymbolic_SeqDense_SeqDense, 1959a9fe9ddaSSatish Balay MatMatMultNumeric_SeqDense_SeqDense, 1960865e5f61SKris Buschelman 0, 1961865e5f61SKris Buschelman 0, 1962d519adbfSMatthew Knepley /*94*/ 0, 1963a9fe9ddaSSatish Balay MatMatMultTranspose_SeqDense_SeqDense, 1964a9fe9ddaSSatish Balay MatMatMultTransposeSymbolic_SeqDense_SeqDense, 1965a9fe9ddaSSatish Balay MatMatMultTransposeNumeric_SeqDense_SeqDense, 1966284134d9SBarry Smith 0, 1967d519adbfSMatthew Knepley /*99*/ 0, 1968284134d9SBarry Smith 0, 1969284134d9SBarry Smith 0, 1970ba337c44SJed Brown MatConjugate_SeqDense, 1971985db425SBarry Smith MatSetSizes_SeqDense, 1972ba337c44SJed Brown /*104*/0, 1973ba337c44SJed Brown MatRealPart_SeqDense, 1974ba337c44SJed Brown MatImaginaryPart_SeqDense, 1975985db425SBarry Smith 0, 1976985db425SBarry Smith 0, 197785e2c93fSHong Zhang /*109*/MatMatSolve_SeqDense, 1978985db425SBarry Smith 0, 19798d0534beSBarry Smith MatGetRowMin_SeqDense, 1980aabbc4fbSShri Abhyankar MatGetColumnVector_SeqDense, 1981aabbc4fbSShri Abhyankar 0, 1982aabbc4fbSShri Abhyankar /*114*/0, 1983aabbc4fbSShri Abhyankar 0, 1984aabbc4fbSShri Abhyankar 0, 1985aabbc4fbSShri Abhyankar 0, 1986aabbc4fbSShri Abhyankar 0, 1987aabbc4fbSShri Abhyankar /*119*/0, 1988aabbc4fbSShri Abhyankar 0, 1989aabbc4fbSShri Abhyankar 0, 19900716a85fSBarry Smith 0, 19910716a85fSBarry Smith 0, 19920716a85fSBarry Smith /*124*/0, 19930716a85fSBarry Smith MatGetColumnNorms_SeqDense 1994985db425SBarry Smith }; 199590ace30eSBarry Smith 19964a2ae208SSatish Balay #undef __FUNCT__ 19974a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 19984b828684SBarry Smith /*@C 1999fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 2000d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 2001d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 2002289bc588SBarry Smith 2003db81eaa0SLois Curfman McInnes Collective on MPI_Comm 2004db81eaa0SLois Curfman McInnes 200520563c6bSBarry Smith Input Parameters: 2006db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 20070c775827SLois Curfman McInnes . m - number of rows 200818f449edSLois Curfman McInnes . n - number of columns 2009c0235b3cSMatthew Knepley - data - optional location of matrix data in column major order. Set data=PETSC_NULL for PETSc 2010dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 201120563c6bSBarry Smith 201220563c6bSBarry Smith Output Parameter: 201344cd7ae7SLois Curfman McInnes . A - the matrix 201420563c6bSBarry Smith 2015b259b22eSLois Curfman McInnes Notes: 201618f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 201718f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 2018b4fd4287SBarry Smith set data=PETSC_NULL. 201918f449edSLois Curfman McInnes 2020027ccd11SLois Curfman McInnes Level: intermediate 2021027ccd11SLois Curfman McInnes 2022dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 2023d65003e9SLois Curfman McInnes 2024db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 202520563c6bSBarry Smith @*/ 20267087cfbeSBarry Smith PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 2027289bc588SBarry Smith { 2028dfbe8321SBarry Smith PetscErrorCode ierr; 20293b2fbd54SBarry Smith 20303a40ed3dSBarry Smith PetscFunctionBegin; 2031f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 2032f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 2033273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 2034273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 2035273d9f13SBarry Smith PetscFunctionReturn(0); 2036273d9f13SBarry Smith } 2037273d9f13SBarry Smith 20384a2ae208SSatish Balay #undef __FUNCT__ 2039afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation" 2040273d9f13SBarry Smith /*@C 2041273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 2042273d9f13SBarry Smith 2043273d9f13SBarry Smith Collective on MPI_Comm 2044273d9f13SBarry Smith 2045273d9f13SBarry Smith Input Parameters: 2046273d9f13SBarry Smith + A - the matrix 2047273d9f13SBarry Smith - data - the array (or PETSC_NULL) 2048273d9f13SBarry Smith 2049273d9f13SBarry Smith Notes: 2050273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 2051273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 2052284134d9SBarry Smith need not call this routine. 2053273d9f13SBarry Smith 2054273d9f13SBarry Smith Level: intermediate 2055273d9f13SBarry Smith 2056273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 2057273d9f13SBarry Smith 2058867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA() 2059867c911aSBarry Smith 2060273d9f13SBarry Smith @*/ 20617087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 2062273d9f13SBarry Smith { 20634ac538c5SBarry Smith PetscErrorCode ierr; 2064a23d5eceSKris Buschelman 2065a23d5eceSKris Buschelman PetscFunctionBegin; 20664ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));CHKERRQ(ierr); 2067a23d5eceSKris Buschelman PetscFunctionReturn(0); 2068a23d5eceSKris Buschelman } 2069a23d5eceSKris Buschelman 2070a23d5eceSKris Buschelman EXTERN_C_BEGIN 2071a23d5eceSKris Buschelman #undef __FUNCT__ 2072afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense" 20737087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 2074a23d5eceSKris Buschelman { 2075273d9f13SBarry Smith Mat_SeqDense *b; 2076dfbe8321SBarry Smith PetscErrorCode ierr; 2077273d9f13SBarry Smith 2078273d9f13SBarry Smith PetscFunctionBegin; 2079273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 2080a868139aSShri Abhyankar 208134ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr); 208234ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr); 208334ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 208434ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 208534ef9618SShri Abhyankar 2086273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 208786d161a7SShri Abhyankar b->Mmax = B->rmap->n; 208886d161a7SShri Abhyankar b->Nmax = B->cmap->n; 208986d161a7SShri Abhyankar if(b->lda <= 0 || b->changelda) b->lda = B->rmap->n; 209086d161a7SShri Abhyankar 20919e8f95c4SLisandro Dalcin if (!data) { /* petsc-allocated storage */ 20929e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 20935afd5e0cSBarry Smith ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 2094284134d9SBarry Smith ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 2095284134d9SBarry Smith ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 20969e8f95c4SLisandro Dalcin b->user_alloc = PETSC_FALSE; 2097273d9f13SBarry Smith } else { /* user-allocated storage */ 20989e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 2099273d9f13SBarry Smith b->v = data; 2100273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 2101273d9f13SBarry Smith } 21020450473dSBarry Smith B->assembled = PETSC_TRUE; 2103273d9f13SBarry Smith PetscFunctionReturn(0); 2104273d9f13SBarry Smith } 2105a23d5eceSKris Buschelman EXTERN_C_END 2106273d9f13SBarry Smith 21071b807ce4Svictorle #undef __FUNCT__ 21081b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 21091b807ce4Svictorle /*@C 21101b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 21111b807ce4Svictorle 21121b807ce4Svictorle Input parameter: 21131b807ce4Svictorle + A - the matrix 21141b807ce4Svictorle - lda - the leading dimension 21151b807ce4Svictorle 21161b807ce4Svictorle Notes: 2117867c911aSBarry Smith This routine is to be used in conjunction with MatSeqDenseSetPreallocation(); 21181b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 21191b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 21201b807ce4Svictorle 21211b807ce4Svictorle Level: intermediate 21221b807ce4Svictorle 21231b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 21241b807ce4Svictorle 2125284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 2126867c911aSBarry Smith 21271b807ce4Svictorle @*/ 21287087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetLDA(Mat B,PetscInt lda) 21291b807ce4Svictorle { 21301b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 213121a2c019SBarry Smith 21321b807ce4Svictorle PetscFunctionBegin; 2133e32f2f54SBarry 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); 21341b807ce4Svictorle b->lda = lda; 213521a2c019SBarry Smith b->changelda = PETSC_FALSE; 213621a2c019SBarry Smith b->Mmax = PetscMax(b->Mmax,lda); 21371b807ce4Svictorle PetscFunctionReturn(0); 21381b807ce4Svictorle } 21391b807ce4Svictorle 21400bad9183SKris Buschelman /*MC 2141fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 21420bad9183SKris Buschelman 21430bad9183SKris Buschelman Options Database Keys: 21440bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 21450bad9183SKris Buschelman 21460bad9183SKris Buschelman Level: beginner 21470bad9183SKris Buschelman 214889665df3SBarry Smith .seealso: MatCreateSeqDense() 214989665df3SBarry Smith 21500bad9183SKris Buschelman M*/ 21510bad9183SKris Buschelman 2152273d9f13SBarry Smith EXTERN_C_BEGIN 21534a2ae208SSatish Balay #undef __FUNCT__ 21544a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 21557087cfbeSBarry Smith PetscErrorCode MatCreate_SeqDense(Mat B) 2156273d9f13SBarry Smith { 2157273d9f13SBarry Smith Mat_SeqDense *b; 2158dfbe8321SBarry Smith PetscErrorCode ierr; 21597c334f02SBarry Smith PetscMPIInt size; 2160273d9f13SBarry Smith 2161273d9f13SBarry Smith PetscFunctionBegin; 21627adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr); 2163e32f2f54SBarry Smith if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 216455659b69SBarry Smith 216538f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr); 2166549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 216744cd7ae7SLois Curfman McInnes B->data = (void*)b; 216818f449edSLois Curfman McInnes 216944cd7ae7SLois Curfman McInnes b->pivots = 0; 2170273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 2171273d9f13SBarry Smith b->v = 0; 217221a2c019SBarry Smith b->changelda = PETSC_FALSE; 21734e220ebcSLois Curfman McInnes 2174b24902e0SBarry Smith 2175ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C", 2176b24902e0SBarry Smith "MatGetFactor_seqdense_petsc", 2177b24902e0SBarry Smith MatGetFactor_seqdense_petsc);CHKERRQ(ierr); 2178a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 2179a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 2180a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 21814ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C", 21824ae313f4SHong Zhang "MatMatMult_SeqAIJ_SeqDense", 21834ae313f4SHong Zhang MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 21844ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C", 21854ae313f4SHong Zhang "MatMatMultSymbolic_SeqAIJ_SeqDense", 21864ae313f4SHong Zhang MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 21874ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C", 21884ae313f4SHong Zhang "MatMatMultNumeric_SeqAIJ_SeqDense", 21894ae313f4SHong Zhang MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 219017667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 21913a40ed3dSBarry Smith PetscFunctionReturn(0); 2192289bc588SBarry Smith } 2193273d9f13SBarry Smith EXTERN_C_END 2194