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; 109aa5ea44dSBarry Smith ierr = PetscLayoutReference(A->rmap,&newi->rmap);CHKERRQ(ierr); 110aa5ea44dSBarry 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; 199*efb80c78SLisandro Dalcin PetscInt n; 20085e2c93fSHong Zhang PetscBLASInt nrhs,info,m=PetscBLASIntCast(A->rmap->n); 20185e2c93fSHong Zhang 20285e2c93fSHong Zhang PetscFunctionBegin; 203*efb80c78SLisandro Dalcin ierr = MatGetSize(B,PETSC_NULL,&n);CHKERRQ(ierr); 204*efb80c78SLisandro Dalcin nrhs = PetscBLASIntCast(n); 20585e2c93fSHong Zhang ierr = MatGetArray(B,&b);CHKERRQ(ierr); 20685e2c93fSHong Zhang ierr = MatGetArray(X,&x);CHKERRQ(ierr); 20785e2c93fSHong Zhang 20885e2c93fSHong Zhang ierr = PetscMemcpy(x,b,m*sizeof(PetscScalar));CHKERRQ(ierr); 20985e2c93fSHong Zhang 21085e2c93fSHong Zhang if (A->factortype == MAT_FACTOR_LU) { 21185e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_GETRS) 21285e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 21385e2c93fSHong Zhang #else 21485e2c93fSHong Zhang LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info); 21585e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 21685e2c93fSHong Zhang #endif 21785e2c93fSHong Zhang } else if (A->factortype == MAT_FACTOR_CHOLESKY){ 21885e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_POTRS) 21985e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 22085e2c93fSHong Zhang #else 22185e2c93fSHong Zhang LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info); 22285e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 22385e2c93fSHong Zhang #endif 22485e2c93fSHong Zhang } 22585e2c93fSHong Zhang else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 22685e2c93fSHong Zhang 22785e2c93fSHong Zhang ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 22885e2c93fSHong Zhang ierr = MatRestoreArray(X,&x);CHKERRQ(ierr); 22985e2c93fSHong Zhang ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr); 23085e2c93fSHong Zhang PetscFunctionReturn(0); 23185e2c93fSHong Zhang } 23285e2c93fSHong Zhang 23385e2c93fSHong Zhang #undef __FUNCT__ 2344a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 235dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 236da3a660dSBarry Smith { 237c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 238dfbe8321SBarry Smith PetscErrorCode ierr; 23987828ca2SBarry Smith PetscScalar *x,*y; 240d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 24167e560aaSBarry Smith 2423a40ed3dSBarry Smith PetscFunctionBegin; 2431ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2441ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 245d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 246752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 247da3a660dSBarry Smith if (mat->pivots) { 248ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 249e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 250ae7cfcebSSatish Balay #else 25171044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 252e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 253ae7cfcebSSatish Balay #endif 2547a97a34bSBarry Smith } else { 255ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 256e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 257ae7cfcebSSatish Balay #else 25871044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 259e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 260ae7cfcebSSatish Balay #endif 261da3a660dSBarry Smith } 2621ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2631ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 264dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 2653a40ed3dSBarry Smith PetscFunctionReturn(0); 266da3a660dSBarry Smith } 2676ee01492SSatish Balay 2684a2ae208SSatish Balay #undef __FUNCT__ 2694a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 270dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 271da3a660dSBarry Smith { 272c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 273dfbe8321SBarry Smith PetscErrorCode ierr; 27487828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 275da3a660dSBarry Smith Vec tmp = 0; 276d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 27767e560aaSBarry Smith 2783a40ed3dSBarry Smith PetscFunctionBegin; 2791ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2801ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 281d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 282da3a660dSBarry Smith if (yy == zz) { 28378b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 28452e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 28578b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 286da3a660dSBarry Smith } 287d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 288752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 289da3a660dSBarry Smith if (mat->pivots) { 290ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 291e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 292ae7cfcebSSatish Balay #else 29371044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 294e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 295ae7cfcebSSatish Balay #endif 296a8c6a408SBarry Smith } else { 297ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 298e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 299ae7cfcebSSatish Balay #else 30071044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 301e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 302ae7cfcebSSatish Balay #endif 303da3a660dSBarry Smith } 3046bf464f9SBarry Smith if (tmp) { 3056bf464f9SBarry Smith ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 3066bf464f9SBarry Smith ierr = VecDestroy(&tmp);CHKERRQ(ierr); 3076bf464f9SBarry Smith } else { 3086bf464f9SBarry Smith ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 3096bf464f9SBarry Smith } 3101ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3111ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 312dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 3133a40ed3dSBarry Smith PetscFunctionReturn(0); 314da3a660dSBarry Smith } 31567e560aaSBarry Smith 3164a2ae208SSatish Balay #undef __FUNCT__ 3174a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 318dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 319da3a660dSBarry Smith { 320c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3216849ba73SBarry Smith PetscErrorCode ierr; 32287828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 323da3a660dSBarry Smith Vec tmp; 324d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 32567e560aaSBarry Smith 3263a40ed3dSBarry Smith PetscFunctionBegin; 327d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 3281ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3291ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 330da3a660dSBarry Smith if (yy == zz) { 33178b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 33252e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 33378b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 334da3a660dSBarry Smith } 335d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 336752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 337da3a660dSBarry Smith if (mat->pivots) { 338ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 339e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 340ae7cfcebSSatish Balay #else 34171044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 342e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 343ae7cfcebSSatish Balay #endif 3443a40ed3dSBarry Smith } else { 345ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 346e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 347ae7cfcebSSatish Balay #else 34871044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 349e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 350ae7cfcebSSatish Balay #endif 351da3a660dSBarry Smith } 35290f02eecSBarry Smith if (tmp) { 3532dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 3546bf464f9SBarry Smith ierr = VecDestroy(&tmp);CHKERRQ(ierr); 3553a40ed3dSBarry Smith } else { 3562dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 35790f02eecSBarry Smith } 3581ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3591ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 360dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 3613a40ed3dSBarry Smith PetscFunctionReturn(0); 362da3a660dSBarry Smith } 363db4efbfdSBarry Smith 364db4efbfdSBarry Smith /* ---------------------------------------------------------------*/ 365db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 366db4efbfdSBarry Smith rather than put it in the Mat->row slot.*/ 367db4efbfdSBarry Smith #undef __FUNCT__ 368db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense" 3690481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo) 370db4efbfdSBarry Smith { 371db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 372db4efbfdSBarry Smith PetscFunctionBegin; 373e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 374db4efbfdSBarry Smith #else 375db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 376db4efbfdSBarry Smith PetscErrorCode ierr; 377db4efbfdSBarry Smith PetscBLASInt n,m,info; 378db4efbfdSBarry Smith 379db4efbfdSBarry Smith PetscFunctionBegin; 380db4efbfdSBarry Smith n = PetscBLASIntCast(A->cmap->n); 381db4efbfdSBarry Smith m = PetscBLASIntCast(A->rmap->n); 382db4efbfdSBarry Smith if (!mat->pivots) { 383db4efbfdSBarry Smith ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr); 384db4efbfdSBarry Smith ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr); 385db4efbfdSBarry Smith } 386db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 387db4efbfdSBarry Smith LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info); 388e32f2f54SBarry Smith if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization"); 389e32f2f54SBarry Smith if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 390db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 391db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 392db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 393db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 394d5f3da31SBarry Smith A->factortype = MAT_FACTOR_LU; 395db4efbfdSBarry Smith 396dc0b31edSSatish Balay ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr); 397db4efbfdSBarry Smith #endif 398db4efbfdSBarry Smith PetscFunctionReturn(0); 399db4efbfdSBarry Smith } 400db4efbfdSBarry Smith 401db4efbfdSBarry Smith #undef __FUNCT__ 402db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense" 4030481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo) 404db4efbfdSBarry Smith { 405db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 406db4efbfdSBarry Smith PetscFunctionBegin; 407e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 408db4efbfdSBarry Smith #else 409db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 410db4efbfdSBarry Smith PetscErrorCode ierr; 411db4efbfdSBarry Smith PetscBLASInt info,n = PetscBLASIntCast(A->cmap->n); 412db4efbfdSBarry Smith 413db4efbfdSBarry Smith PetscFunctionBegin; 414db4efbfdSBarry Smith ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 415db4efbfdSBarry Smith 416db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 417db4efbfdSBarry Smith LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info); 418e32f2f54SBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 419db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 420db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 421db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 422db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 423d5f3da31SBarry Smith A->factortype = MAT_FACTOR_CHOLESKY; 424dc0b31edSSatish Balay ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr); 425db4efbfdSBarry Smith #endif 426db4efbfdSBarry Smith PetscFunctionReturn(0); 427db4efbfdSBarry Smith } 428db4efbfdSBarry Smith 429db4efbfdSBarry Smith 430db4efbfdSBarry Smith #undef __FUNCT__ 431db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 4320481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 433db4efbfdSBarry Smith { 434db4efbfdSBarry Smith PetscErrorCode ierr; 435db4efbfdSBarry Smith MatFactorInfo info; 436db4efbfdSBarry Smith 437db4efbfdSBarry Smith PetscFunctionBegin; 438db4efbfdSBarry Smith info.fill = 1.0; 439c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 440719d5645SBarry Smith ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr); 441db4efbfdSBarry Smith PetscFunctionReturn(0); 442db4efbfdSBarry Smith } 443db4efbfdSBarry Smith 444db4efbfdSBarry Smith #undef __FUNCT__ 445db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 4460481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info) 447db4efbfdSBarry Smith { 448db4efbfdSBarry Smith PetscFunctionBegin; 449c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 450719d5645SBarry Smith fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense; 451db4efbfdSBarry Smith PetscFunctionReturn(0); 452db4efbfdSBarry Smith } 453db4efbfdSBarry Smith 454db4efbfdSBarry Smith #undef __FUNCT__ 455db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 4560481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 457db4efbfdSBarry Smith { 458db4efbfdSBarry Smith PetscFunctionBegin; 459c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 460719d5645SBarry Smith fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense; 461db4efbfdSBarry Smith PetscFunctionReturn(0); 462db4efbfdSBarry Smith } 463db4efbfdSBarry Smith 464bb5747d9SMatthew Knepley EXTERN_C_BEGIN 465db4efbfdSBarry Smith #undef __FUNCT__ 466db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc" 467db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact) 468db4efbfdSBarry Smith { 469db4efbfdSBarry Smith PetscErrorCode ierr; 470db4efbfdSBarry Smith 471db4efbfdSBarry Smith PetscFunctionBegin; 472db4efbfdSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr); 473db4efbfdSBarry Smith ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 474db4efbfdSBarry Smith ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr); 475db4efbfdSBarry Smith if (ftype == MAT_FACTOR_LU){ 476db4efbfdSBarry Smith (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense; 477db4efbfdSBarry Smith } else { 478db4efbfdSBarry Smith (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense; 479db4efbfdSBarry Smith } 480d5f3da31SBarry Smith (*fact)->factortype = ftype; 481db4efbfdSBarry Smith PetscFunctionReturn(0); 482db4efbfdSBarry Smith } 483bb5747d9SMatthew Knepley EXTERN_C_END 484db4efbfdSBarry Smith 485289bc588SBarry Smith /* ------------------------------------------------------------------*/ 4864a2ae208SSatish Balay #undef __FUNCT__ 48741f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense" 48841f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 489289bc588SBarry Smith { 490c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 49187828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 492dfbe8321SBarry Smith PetscErrorCode ierr; 493d0f46423SBarry Smith PetscInt m = A->rmap->n,i; 494aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 4950805154bSBarry Smith PetscBLASInt o = 1,bm = PetscBLASIntCast(m); 496bc1b551cSSatish Balay #endif 497289bc588SBarry Smith 4983a40ed3dSBarry Smith PetscFunctionBegin; 499289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 50071044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 5012dcb1b2aSMatthew Knepley ierr = VecSet(xx,zero);CHKERRQ(ierr); 502289bc588SBarry Smith } 5031ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5041ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 505b965ef7fSBarry Smith its = its*lits; 506e32f2f54SBarry 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); 507289bc588SBarry Smith while (its--) { 508fccaa45eSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 509289bc588SBarry Smith for (i=0; i<m; i++) { 510aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 511f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 512f1747703SBarry Smith not happy about returning a double complex */ 51313f74950SBarry Smith PetscInt _i; 51487828ca2SBarry Smith PetscScalar sum = b[i]; 515f1747703SBarry Smith for (_i=0; _i<m; _i++) { 5163f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 517f1747703SBarry Smith } 518f1747703SBarry Smith xt = sum; 519f1747703SBarry Smith #else 52071044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 521f1747703SBarry Smith #endif 52255a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 523289bc588SBarry Smith } 524289bc588SBarry Smith } 525fccaa45eSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 526289bc588SBarry Smith for (i=m-1; i>=0; i--) { 527aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 528f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 529f1747703SBarry Smith not happy about returning a double complex */ 53013f74950SBarry Smith PetscInt _i; 53187828ca2SBarry Smith PetscScalar sum = b[i]; 532f1747703SBarry Smith for (_i=0; _i<m; _i++) { 5333f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 534f1747703SBarry Smith } 535f1747703SBarry Smith xt = sum; 536f1747703SBarry Smith #else 53771044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 538f1747703SBarry Smith #endif 53955a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 540289bc588SBarry Smith } 541289bc588SBarry Smith } 542289bc588SBarry Smith } 5431ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 5441ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5453a40ed3dSBarry Smith PetscFunctionReturn(0); 546289bc588SBarry Smith } 547289bc588SBarry Smith 548289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5494a2ae208SSatish Balay #undef __FUNCT__ 5504a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 551dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 552289bc588SBarry Smith { 553c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 55487828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 555dfbe8321SBarry Smith PetscErrorCode ierr; 5560805154bSBarry Smith PetscBLASInt m, n,_One=1; 557ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 5583a40ed3dSBarry Smith 5593a40ed3dSBarry Smith PetscFunctionBegin; 560d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 561d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 562d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5631ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5641ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 56571044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One); 5661ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5671ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 568dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 5693a40ed3dSBarry Smith PetscFunctionReturn(0); 570289bc588SBarry Smith } 571800995b7SMatthew Knepley 5724a2ae208SSatish Balay #undef __FUNCT__ 5734a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 574dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 575289bc588SBarry Smith { 576c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 57787828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 578dfbe8321SBarry Smith PetscErrorCode ierr; 5790805154bSBarry Smith PetscBLASInt m, n, _One=1; 5803a40ed3dSBarry Smith 5813a40ed3dSBarry Smith PetscFunctionBegin; 582d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 583d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 584d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5851ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5861ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 58771044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 5881ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5891ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 590dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr); 5913a40ed3dSBarry Smith PetscFunctionReturn(0); 592289bc588SBarry Smith } 5936ee01492SSatish Balay 5944a2ae208SSatish Balay #undef __FUNCT__ 5954a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 596dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 597289bc588SBarry Smith { 598c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 59987828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 600dfbe8321SBarry Smith PetscErrorCode ierr; 6010805154bSBarry Smith PetscBLASInt m, n, _One=1; 6023a40ed3dSBarry Smith 6033a40ed3dSBarry Smith PetscFunctionBegin; 604d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 605d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 606d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 607600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6081ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6091ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 61071044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 6111ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6121ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 613dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6143a40ed3dSBarry Smith PetscFunctionReturn(0); 615289bc588SBarry Smith } 6166ee01492SSatish Balay 6174a2ae208SSatish Balay #undef __FUNCT__ 6184a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 619dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 620289bc588SBarry Smith { 621c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 62287828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 623dfbe8321SBarry Smith PetscErrorCode ierr; 6240805154bSBarry Smith PetscBLASInt m, n, _One=1; 62587828ca2SBarry Smith PetscScalar _DOne=1.0; 6263a40ed3dSBarry Smith 6273a40ed3dSBarry Smith PetscFunctionBegin; 628d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 629d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 630d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 631600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6321ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6331ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 63471044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 6351ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6361ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 637dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6383a40ed3dSBarry Smith PetscFunctionReturn(0); 639289bc588SBarry Smith } 640289bc588SBarry Smith 641289bc588SBarry Smith /* -----------------------------------------------------------------*/ 6424a2ae208SSatish Balay #undef __FUNCT__ 6434a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 64413f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 645289bc588SBarry Smith { 646c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 64787828ca2SBarry Smith PetscScalar *v; 6486849ba73SBarry Smith PetscErrorCode ierr; 64913f74950SBarry Smith PetscInt i; 65067e560aaSBarry Smith 6513a40ed3dSBarry Smith PetscFunctionBegin; 652d0f46423SBarry Smith *ncols = A->cmap->n; 653289bc588SBarry Smith if (cols) { 654d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr); 655d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) (*cols)[i] = i; 656289bc588SBarry Smith } 657289bc588SBarry Smith if (vals) { 658d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 659289bc588SBarry Smith v = mat->v + row; 660d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;} 661289bc588SBarry Smith } 6623a40ed3dSBarry Smith PetscFunctionReturn(0); 663289bc588SBarry Smith } 6646ee01492SSatish Balay 6654a2ae208SSatish Balay #undef __FUNCT__ 6664a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 66713f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 668289bc588SBarry Smith { 669dfbe8321SBarry Smith PetscErrorCode ierr; 670606d414cSSatish Balay PetscFunctionBegin; 671606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 672606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 6733a40ed3dSBarry Smith PetscFunctionReturn(0); 674289bc588SBarry Smith } 675289bc588SBarry Smith /* ----------------------------------------------------------------*/ 6764a2ae208SSatish Balay #undef __FUNCT__ 6774a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 67813f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 679289bc588SBarry Smith { 680c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 68113f74950SBarry Smith PetscInt i,j,idx=0; 682d6dfbf8fSBarry Smith 6833a40ed3dSBarry Smith PetscFunctionBegin; 68471fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 685289bc588SBarry Smith if (!mat->roworiented) { 686dbb450caSBarry Smith if (addv == INSERT_VALUES) { 687289bc588SBarry Smith for (j=0; j<n; j++) { 688cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6892515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 690e32f2f54SBarry 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); 69158804f6dSBarry Smith #endif 692289bc588SBarry Smith for (i=0; i<m; i++) { 693cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6942515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 695e32f2f54SBarry 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); 69658804f6dSBarry Smith #endif 697cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 698289bc588SBarry Smith } 699289bc588SBarry Smith } 7003a40ed3dSBarry Smith } else { 701289bc588SBarry Smith for (j=0; j<n; j++) { 702cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 7032515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 704e32f2f54SBarry 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); 70558804f6dSBarry Smith #endif 706289bc588SBarry Smith for (i=0; i<m; i++) { 707cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 7082515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 709e32f2f54SBarry 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); 71058804f6dSBarry Smith #endif 711cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 712289bc588SBarry Smith } 713289bc588SBarry Smith } 714289bc588SBarry Smith } 7153a40ed3dSBarry Smith } else { 716dbb450caSBarry Smith if (addv == INSERT_VALUES) { 717e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 718cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7192515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 720e32f2f54SBarry 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); 72158804f6dSBarry Smith #endif 722e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 723cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7242515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 725e32f2f54SBarry 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); 72658804f6dSBarry Smith #endif 727cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 728e8d4e0b9SBarry Smith } 729e8d4e0b9SBarry Smith } 7303a40ed3dSBarry Smith } else { 731289bc588SBarry Smith for (i=0; i<m; i++) { 732cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7332515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 734e32f2f54SBarry 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); 73558804f6dSBarry Smith #endif 736289bc588SBarry Smith for (j=0; j<n; j++) { 737cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7382515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 739e32f2f54SBarry 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); 74058804f6dSBarry Smith #endif 741cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 742289bc588SBarry Smith } 743289bc588SBarry Smith } 744289bc588SBarry Smith } 745e8d4e0b9SBarry Smith } 7463a40ed3dSBarry Smith PetscFunctionReturn(0); 747289bc588SBarry Smith } 748e8d4e0b9SBarry Smith 7494a2ae208SSatish Balay #undef __FUNCT__ 7504a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 75113f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 752ae80bb75SLois Curfman McInnes { 753ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 75413f74950SBarry Smith PetscInt i,j; 755ae80bb75SLois Curfman McInnes 7563a40ed3dSBarry Smith PetscFunctionBegin; 757ae80bb75SLois Curfman McInnes /* row-oriented output */ 758ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 75997e567efSBarry Smith if (indexm[i] < 0) {v += n;continue;} 760e32f2f54SBarry 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); 761ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 7626f31f424SBarry Smith if (indexn[j] < 0) {v++; continue;} 763e32f2f54SBarry 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); 76497e567efSBarry Smith *v++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 765ae80bb75SLois Curfman McInnes } 766ae80bb75SLois Curfman McInnes } 7673a40ed3dSBarry Smith PetscFunctionReturn(0); 768ae80bb75SLois Curfman McInnes } 769ae80bb75SLois Curfman McInnes 770289bc588SBarry Smith /* -----------------------------------------------------------------*/ 771289bc588SBarry Smith 7724a2ae208SSatish Balay #undef __FUNCT__ 7735bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_SeqDense" 774112444f4SShri Abhyankar PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer) 775aabbc4fbSShri Abhyankar { 776aabbc4fbSShri Abhyankar Mat_SeqDense *a; 777aabbc4fbSShri Abhyankar PetscErrorCode ierr; 778aabbc4fbSShri Abhyankar PetscInt *scols,i,j,nz,header[4]; 779aabbc4fbSShri Abhyankar int fd; 780aabbc4fbSShri Abhyankar PetscMPIInt size; 781aabbc4fbSShri Abhyankar PetscInt *rowlengths = 0,M,N,*cols,grows,gcols; 782aabbc4fbSShri Abhyankar PetscScalar *vals,*svals,*v,*w; 783aabbc4fbSShri Abhyankar MPI_Comm comm = ((PetscObject)viewer)->comm; 784aabbc4fbSShri Abhyankar 785aabbc4fbSShri Abhyankar PetscFunctionBegin; 786aabbc4fbSShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 787aabbc4fbSShri Abhyankar if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor"); 788aabbc4fbSShri Abhyankar ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 789aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 790aabbc4fbSShri Abhyankar if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 791aabbc4fbSShri Abhyankar M = header[1]; N = header[2]; nz = header[3]; 792aabbc4fbSShri Abhyankar 793aabbc4fbSShri Abhyankar /* set global size if not set already*/ 794aabbc4fbSShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 795aabbc4fbSShri Abhyankar ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr); 796aabbc4fbSShri Abhyankar } else { 797aabbc4fbSShri Abhyankar /* if sizes and type are already set, check if the vector global sizes are correct */ 798aabbc4fbSShri Abhyankar ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr); 799aabbc4fbSShri 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); 800aabbc4fbSShri Abhyankar } 801aabbc4fbSShri Abhyankar ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 802aabbc4fbSShri Abhyankar 803aabbc4fbSShri Abhyankar if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 804aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 805aabbc4fbSShri Abhyankar v = a->v; 806aabbc4fbSShri Abhyankar /* Allocate some temp space to read in the values and then flip them 807aabbc4fbSShri Abhyankar from row major to column major */ 808aabbc4fbSShri Abhyankar ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 809aabbc4fbSShri Abhyankar /* read in nonzero values */ 810aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 811aabbc4fbSShri Abhyankar /* now flip the values and store them in the matrix*/ 812aabbc4fbSShri Abhyankar for (j=0; j<N; j++) { 813aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 814aabbc4fbSShri Abhyankar *v++ =w[i*N+j]; 815aabbc4fbSShri Abhyankar } 816aabbc4fbSShri Abhyankar } 817aabbc4fbSShri Abhyankar ierr = PetscFree(w);CHKERRQ(ierr); 818aabbc4fbSShri Abhyankar } else { 819aabbc4fbSShri Abhyankar /* read row lengths */ 820aabbc4fbSShri Abhyankar ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 821aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 822aabbc4fbSShri Abhyankar 823aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 824aabbc4fbSShri Abhyankar v = a->v; 825aabbc4fbSShri Abhyankar 826aabbc4fbSShri Abhyankar /* read column indices and nonzeros */ 827aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr); 828aabbc4fbSShri Abhyankar cols = scols; 829aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 830aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 831aabbc4fbSShri Abhyankar vals = svals; 832aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 833aabbc4fbSShri Abhyankar 834aabbc4fbSShri Abhyankar /* insert into matrix */ 835aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 836aabbc4fbSShri Abhyankar for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 837aabbc4fbSShri Abhyankar svals += rowlengths[i]; scols += rowlengths[i]; 838aabbc4fbSShri Abhyankar } 839aabbc4fbSShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 840aabbc4fbSShri Abhyankar ierr = PetscFree(cols);CHKERRQ(ierr); 841aabbc4fbSShri Abhyankar ierr = PetscFree(rowlengths);CHKERRQ(ierr); 842aabbc4fbSShri Abhyankar } 843aabbc4fbSShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 844aabbc4fbSShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 845aabbc4fbSShri Abhyankar 846aabbc4fbSShri Abhyankar PetscFunctionReturn(0); 847aabbc4fbSShri Abhyankar } 848aabbc4fbSShri Abhyankar 849aabbc4fbSShri Abhyankar #undef __FUNCT__ 8504a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 8516849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 852289bc588SBarry Smith { 853932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 854dfbe8321SBarry Smith PetscErrorCode ierr; 85513f74950SBarry Smith PetscInt i,j; 8562dcb1b2aSMatthew Knepley const char *name; 85787828ca2SBarry Smith PetscScalar *v; 858f3ef73ceSBarry Smith PetscViewerFormat format; 8595f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX) 860ace3abfcSBarry Smith PetscBool allreal = PETSC_TRUE; 8615f481a85SSatish Balay #endif 862932b0c3eSLois Curfman McInnes 8633a40ed3dSBarry Smith PetscFunctionBegin; 864b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 865456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 8663a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 867fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 868d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 8697566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 870d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 87144cd7ae7SLois Curfman McInnes v = a->v + i; 87277431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 873d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 874aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 875329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 876a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 877329f5518SBarry Smith } else if (PetscRealPart(*v)) { 878a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr); 8796831982aSBarry Smith } 88080cd9d93SLois Curfman McInnes #else 8816831982aSBarry Smith if (*v) { 882a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr); 8836831982aSBarry Smith } 88480cd9d93SLois Curfman McInnes #endif 8851b807ce4Svictorle v += a->lda; 88680cd9d93SLois Curfman McInnes } 887b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 88880cd9d93SLois Curfman McInnes } 889d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 8903a40ed3dSBarry Smith } else { 891d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 892aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 89347989497SBarry Smith /* determine if matrix has all real values */ 89447989497SBarry Smith v = a->v; 895d0f46423SBarry Smith for (i=0; i<A->rmap->n*A->cmap->n; i++) { 896ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 89747989497SBarry Smith } 89847989497SBarry Smith #endif 899fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 9003a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 901d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr); 902d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 903fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 9047566de4bSShri Abhyankar } else { 9057566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 906ffac6cdbSBarry Smith } 907ffac6cdbSBarry Smith 908d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 909932b0c3eSLois Curfman McInnes v = a->v + i; 910d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 911aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 91247989497SBarry Smith if (allreal) { 913f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr); 91447989497SBarry Smith } else { 915f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 91647989497SBarry Smith } 917289bc588SBarry Smith #else 918f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr); 919289bc588SBarry Smith #endif 9201b807ce4Svictorle v += a->lda; 921289bc588SBarry Smith } 922b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 923289bc588SBarry Smith } 924fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 925b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 926ffac6cdbSBarry Smith } 927d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 928da3a660dSBarry Smith } 929b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 9303a40ed3dSBarry Smith PetscFunctionReturn(0); 931289bc588SBarry Smith } 932289bc588SBarry Smith 9334a2ae208SSatish Balay #undef __FUNCT__ 9344a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 9356849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 936932b0c3eSLois Curfman McInnes { 937932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9386849ba73SBarry Smith PetscErrorCode ierr; 93913f74950SBarry Smith int fd; 940d0f46423SBarry Smith PetscInt ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n; 941f4403165SShri Abhyankar PetscScalar *v,*anonz,*vals; 942f4403165SShri Abhyankar PetscViewerFormat format; 943932b0c3eSLois Curfman McInnes 9443a40ed3dSBarry Smith PetscFunctionBegin; 945b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 94690ace30eSBarry Smith 947f4403165SShri Abhyankar ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 948f4403165SShri Abhyankar if (format == PETSC_VIEWER_NATIVE) { 949f4403165SShri Abhyankar /* store the matrix as a dense matrix */ 950f4403165SShri Abhyankar ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 951f4403165SShri Abhyankar col_lens[0] = MAT_FILE_CLASSID; 952f4403165SShri Abhyankar col_lens[1] = m; 953f4403165SShri Abhyankar col_lens[2] = n; 954f4403165SShri Abhyankar col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 955f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 956f4403165SShri Abhyankar ierr = PetscFree(col_lens);CHKERRQ(ierr); 957f4403165SShri Abhyankar 958f4403165SShri Abhyankar /* write out matrix, by rows */ 959f4403165SShri Abhyankar ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 960f4403165SShri Abhyankar v = a->v; 961f4403165SShri Abhyankar for (j=0; j<n; j++) { 962f4403165SShri Abhyankar for (i=0; i<m; i++) { 963f4403165SShri Abhyankar vals[j + i*n] = *v++; 964f4403165SShri Abhyankar } 965f4403165SShri Abhyankar } 966f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 967f4403165SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 968f4403165SShri Abhyankar } else { 96913f74950SBarry Smith ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 9700700a824SBarry Smith col_lens[0] = MAT_FILE_CLASSID; 971932b0c3eSLois Curfman McInnes col_lens[1] = m; 972932b0c3eSLois Curfman McInnes col_lens[2] = n; 973932b0c3eSLois Curfman McInnes col_lens[3] = nz; 974932b0c3eSLois Curfman McInnes 975932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 976932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 9776f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 978932b0c3eSLois Curfman McInnes 979932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 980932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 981932b0c3eSLois Curfman McInnes ict = 0; 982932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 983932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 984932b0c3eSLois Curfman McInnes } 9856f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 986606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 987932b0c3eSLois Curfman McInnes 988932b0c3eSLois Curfman McInnes /* store nonzero values */ 98987828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 990932b0c3eSLois Curfman McInnes ict = 0; 991932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 992932b0c3eSLois Curfman McInnes v = a->v + i; 993932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 9941b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 995932b0c3eSLois Curfman McInnes } 996932b0c3eSLois Curfman McInnes } 9976f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 998606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 999f4403165SShri Abhyankar } 10003a40ed3dSBarry Smith PetscFunctionReturn(0); 1001932b0c3eSLois Curfman McInnes } 1002932b0c3eSLois Curfman McInnes 10034a2ae208SSatish Balay #undef __FUNCT__ 10044a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 1005dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 1006f1af5d2fSBarry Smith { 1007f1af5d2fSBarry Smith Mat A = (Mat) Aa; 1008f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 10096849ba73SBarry Smith PetscErrorCode ierr; 1010d0f46423SBarry Smith PetscInt m = A->rmap->n,n = A->cmap->n,color,i,j; 101187828ca2SBarry Smith PetscScalar *v = a->v; 1012b0a32e0cSBarry Smith PetscViewer viewer; 1013b0a32e0cSBarry Smith PetscDraw popup; 1014329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 1015f3ef73ceSBarry Smith PetscViewerFormat format; 1016f1af5d2fSBarry Smith 1017f1af5d2fSBarry Smith PetscFunctionBegin; 1018f1af5d2fSBarry Smith 1019f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 1020b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1021b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1022f1af5d2fSBarry Smith 1023f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 1024fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 1025f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 1026b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1027f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1028f1af5d2fSBarry Smith x_l = j; 1029f1af5d2fSBarry Smith x_r = x_l + 1.0; 1030f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1031f1af5d2fSBarry Smith y_l = m - i - 1.0; 1032f1af5d2fSBarry Smith y_r = y_l + 1.0; 1033f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 1034329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 1035b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1036329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 1037b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1038f1af5d2fSBarry Smith } else { 1039f1af5d2fSBarry Smith continue; 1040f1af5d2fSBarry Smith } 1041f1af5d2fSBarry Smith #else 1042f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 1043b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1044f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 1045b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1046f1af5d2fSBarry Smith } else { 1047f1af5d2fSBarry Smith continue; 1048f1af5d2fSBarry Smith } 1049f1af5d2fSBarry Smith #endif 1050b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1051f1af5d2fSBarry Smith } 1052f1af5d2fSBarry Smith } 1053f1af5d2fSBarry Smith } else { 1054f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 1055f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 1056f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 1057f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 1058f1af5d2fSBarry Smith } 1059b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 1060b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 1061b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 1062f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1063f1af5d2fSBarry Smith x_l = j; 1064f1af5d2fSBarry Smith x_r = x_l + 1.0; 1065f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1066f1af5d2fSBarry Smith y_l = m - i - 1.0; 1067f1af5d2fSBarry Smith y_r = y_l + 1.0; 1068b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 1069b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1070f1af5d2fSBarry Smith } 1071f1af5d2fSBarry Smith } 1072f1af5d2fSBarry Smith } 1073f1af5d2fSBarry Smith PetscFunctionReturn(0); 1074f1af5d2fSBarry Smith } 1075f1af5d2fSBarry Smith 10764a2ae208SSatish Balay #undef __FUNCT__ 10774a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 1078dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 1079f1af5d2fSBarry Smith { 1080b0a32e0cSBarry Smith PetscDraw draw; 1081ace3abfcSBarry Smith PetscBool isnull; 1082329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 1083dfbe8321SBarry Smith PetscErrorCode ierr; 1084f1af5d2fSBarry Smith 1085f1af5d2fSBarry Smith PetscFunctionBegin; 1086b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1087b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1088abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 1089f1af5d2fSBarry Smith 1090f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1091d0f46423SBarry Smith xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0; 1092f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 1093b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1094b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 1095f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 1096f1af5d2fSBarry Smith PetscFunctionReturn(0); 1097f1af5d2fSBarry Smith } 1098f1af5d2fSBarry Smith 10994a2ae208SSatish Balay #undef __FUNCT__ 11004a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 1101dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1102932b0c3eSLois Curfman McInnes { 1103dfbe8321SBarry Smith PetscErrorCode ierr; 1104ace3abfcSBarry Smith PetscBool iascii,isbinary,isdraw; 1105932b0c3eSLois Curfman McInnes 11063a40ed3dSBarry Smith PetscFunctionBegin; 11072692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 11082692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 11092692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 11100f5bd95cSBarry Smith 1111c45a1595SBarry Smith if (iascii) { 1112c45a1595SBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 11130f5bd95cSBarry Smith } else if (isbinary) { 11143a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1115f1af5d2fSBarry Smith } else if (isdraw) { 1116f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 11175cd90555SBarry Smith } else { 1118e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1119932b0c3eSLois Curfman McInnes } 11203a40ed3dSBarry Smith PetscFunctionReturn(0); 1121932b0c3eSLois Curfman McInnes } 1122289bc588SBarry Smith 11234a2ae208SSatish Balay #undef __FUNCT__ 11244a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1125dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1126289bc588SBarry Smith { 1127ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1128dfbe8321SBarry Smith PetscErrorCode ierr; 112990f02eecSBarry Smith 11303a40ed3dSBarry Smith PetscFunctionBegin; 1131aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1132d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n); 1133a5a9c739SBarry Smith #endif 113405b42c5fSBarry Smith ierr = PetscFree(l->pivots);CHKERRQ(ierr); 11356857c123SSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1136bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 1137dbd8c25aSHong Zhang 1138dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 1139901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 11404ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11414ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11424ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11433a40ed3dSBarry Smith PetscFunctionReturn(0); 1144289bc588SBarry Smith } 1145289bc588SBarry Smith 11464a2ae208SSatish Balay #undef __FUNCT__ 11474a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1148fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout) 1149289bc588SBarry Smith { 1150c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11516849ba73SBarry Smith PetscErrorCode ierr; 115213f74950SBarry Smith PetscInt k,j,m,n,M; 115387828ca2SBarry Smith PetscScalar *v,tmp; 115448b35521SBarry Smith 11553a40ed3dSBarry Smith PetscFunctionBegin; 1156d0f46423SBarry Smith v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n; 1157e9695a30SBarry Smith if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */ 1158e7e72b3dSBarry Smith if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1159e7e72b3dSBarry Smith else { 1160d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1161289bc588SBarry Smith for (k=0; k<j; k++) { 11621b807ce4Svictorle tmp = v[j + k*M]; 11631b807ce4Svictorle v[j + k*M] = v[k + j*M]; 11641b807ce4Svictorle v[k + j*M] = tmp; 1165289bc588SBarry Smith } 1166289bc588SBarry Smith } 1167d64ed03dSBarry Smith } 11683a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1169d3e5ee88SLois Curfman McInnes Mat tmat; 1170ec8511deSBarry Smith Mat_SeqDense *tmatd; 117187828ca2SBarry Smith PetscScalar *v2; 1172ea709b57SSatish Balay 1173fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 11747adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr); 1175d0f46423SBarry Smith ierr = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr); 11767adad957SLisandro Dalcin ierr = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 11775c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1178fc4dec0aSBarry Smith } else { 1179fc4dec0aSBarry Smith tmat = *matout; 1180fc4dec0aSBarry Smith } 1181ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 11820de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1183d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 11841b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1185d3e5ee88SLois Curfman McInnes } 11866d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11876d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1188d3e5ee88SLois Curfman McInnes *matout = tmat; 118948b35521SBarry Smith } 11903a40ed3dSBarry Smith PetscFunctionReturn(0); 1191289bc588SBarry Smith } 1192289bc588SBarry Smith 11934a2ae208SSatish Balay #undef __FUNCT__ 11944a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1195ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool *flg) 1196289bc588SBarry Smith { 1197c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1198c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 119913f74950SBarry Smith PetscInt i,j; 120087828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 12019ea5d5aeSSatish Balay 12023a40ed3dSBarry Smith PetscFunctionBegin; 1203d0f46423SBarry Smith if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1204d0f46423SBarry Smith if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1205d0f46423SBarry Smith for (i=0; i<A1->rmap->n; i++) { 12061b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1207d0f46423SBarry Smith for (j=0; j<A1->cmap->n; j++) { 12083a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 12091b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 12101b807ce4Svictorle } 1211289bc588SBarry Smith } 121277c4ece6SBarry Smith *flg = PETSC_TRUE; 12133a40ed3dSBarry Smith PetscFunctionReturn(0); 1214289bc588SBarry Smith } 1215289bc588SBarry Smith 12164a2ae208SSatish Balay #undef __FUNCT__ 12174a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1218dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1219289bc588SBarry Smith { 1220c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1221dfbe8321SBarry Smith PetscErrorCode ierr; 122213f74950SBarry Smith PetscInt i,n,len; 122387828ca2SBarry Smith PetscScalar *x,zero = 0.0; 122444cd7ae7SLois Curfman McInnes 12253a40ed3dSBarry Smith PetscFunctionBegin; 12262dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 12277a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 12281ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1229d0f46423SBarry Smith len = PetscMin(A->rmap->n,A->cmap->n); 1230e32f2f54SBarry Smith if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 123144cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 12321b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1233289bc588SBarry Smith } 12341ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 12353a40ed3dSBarry Smith PetscFunctionReturn(0); 1236289bc588SBarry Smith } 1237289bc588SBarry Smith 12384a2ae208SSatish Balay #undef __FUNCT__ 12394a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1240dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1241289bc588SBarry Smith { 1242c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 124387828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1244dfbe8321SBarry Smith PetscErrorCode ierr; 1245d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n; 124655659b69SBarry Smith 12473a40ed3dSBarry Smith PetscFunctionBegin; 124828988994SBarry Smith if (ll) { 12497a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 12501ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1251e32f2f54SBarry Smith if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1252da3a660dSBarry Smith for (i=0; i<m; i++) { 1253da3a660dSBarry Smith x = l[i]; 1254da3a660dSBarry Smith v = mat->v + i; 1255da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1256da3a660dSBarry Smith } 12571ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1258efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1259da3a660dSBarry Smith } 126028988994SBarry Smith if (rr) { 12617a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 12621ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1263e32f2f54SBarry Smith if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1264da3a660dSBarry Smith for (i=0; i<n; i++) { 1265da3a660dSBarry Smith x = r[i]; 1266da3a660dSBarry Smith v = mat->v + i*m; 1267da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1268da3a660dSBarry Smith } 12691ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1270efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1271da3a660dSBarry Smith } 12723a40ed3dSBarry Smith PetscFunctionReturn(0); 1273289bc588SBarry Smith } 1274289bc588SBarry Smith 12754a2ae208SSatish Balay #undef __FUNCT__ 12764a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1277dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1278289bc588SBarry Smith { 1279c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 128087828ca2SBarry Smith PetscScalar *v = mat->v; 1281329f5518SBarry Smith PetscReal sum = 0.0; 1282d0f46423SBarry Smith PetscInt lda=mat->lda,m=A->rmap->n,i,j; 1283efee365bSSatish Balay PetscErrorCode ierr; 128455659b69SBarry Smith 12853a40ed3dSBarry Smith PetscFunctionBegin; 1286289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1287a5ce6ee0Svictorle if (lda>m) { 1288d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1289a5ce6ee0Svictorle v = mat->v+j*lda; 1290a5ce6ee0Svictorle for (i=0; i<m; i++) { 1291a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1292a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1293a5ce6ee0Svictorle #else 1294a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1295a5ce6ee0Svictorle #endif 1296a5ce6ee0Svictorle } 1297a5ce6ee0Svictorle } 1298a5ce6ee0Svictorle } else { 1299d0f46423SBarry Smith for (i=0; i<A->cmap->n*A->rmap->n; i++) { 1300aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1301329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1302289bc588SBarry Smith #else 1303289bc588SBarry Smith sum += (*v)*(*v); v++; 1304289bc588SBarry Smith #endif 1305289bc588SBarry Smith } 1306a5ce6ee0Svictorle } 1307064f8208SBarry Smith *nrm = sqrt(sum); 1308dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13093a40ed3dSBarry Smith } else if (type == NORM_1) { 1310064f8208SBarry Smith *nrm = 0.0; 1311d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 13121b807ce4Svictorle v = mat->v + j*mat->lda; 1313289bc588SBarry Smith sum = 0.0; 1314d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 131533a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1316289bc588SBarry Smith } 1317064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1318289bc588SBarry Smith } 1319d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13203a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1321064f8208SBarry Smith *nrm = 0.0; 1322d0f46423SBarry Smith for (j=0; j<A->rmap->n; j++) { 1323289bc588SBarry Smith v = mat->v + j; 1324289bc588SBarry Smith sum = 0.0; 1325d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) { 13261b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1327289bc588SBarry Smith } 1328064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1329289bc588SBarry Smith } 1330d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 1331e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 13323a40ed3dSBarry Smith PetscFunctionReturn(0); 1333289bc588SBarry Smith } 1334289bc588SBarry Smith 13354a2ae208SSatish Balay #undef __FUNCT__ 13364a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 1337ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg) 1338289bc588SBarry Smith { 1339c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 134063ba0a88SBarry Smith PetscErrorCode ierr; 134167e560aaSBarry Smith 13423a40ed3dSBarry Smith PetscFunctionBegin; 1343b5a2b587SKris Buschelman switch (op) { 1344b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 13454e0d8c25SBarry Smith aij->roworiented = flg; 1346b5a2b587SKris Buschelman break; 1347512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 1348b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 13493971808eSMatthew Knepley case MAT_NEW_NONZERO_ALLOCATION_ERR: 13504e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 1351b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1352b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 135377e54ba9SKris Buschelman case MAT_SYMMETRIC: 135477e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 13559a4540c5SBarry Smith case MAT_HERMITIAN: 13569a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1357600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 1358290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 135977e54ba9SKris Buschelman break; 1360b5a2b587SKris Buschelman default: 1361e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 13623a40ed3dSBarry Smith } 13633a40ed3dSBarry Smith PetscFunctionReturn(0); 1364289bc588SBarry Smith } 1365289bc588SBarry Smith 13664a2ae208SSatish Balay #undef __FUNCT__ 13674a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1368dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 13696f0a148fSBarry Smith { 1370ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 13716849ba73SBarry Smith PetscErrorCode ierr; 1372d0f46423SBarry Smith PetscInt lda=l->lda,m=A->rmap->n,j; 13733a40ed3dSBarry Smith 13743a40ed3dSBarry Smith PetscFunctionBegin; 1375a5ce6ee0Svictorle if (lda>m) { 1376d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1377a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1378a5ce6ee0Svictorle } 1379a5ce6ee0Svictorle } else { 1380d0f46423SBarry Smith ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1381a5ce6ee0Svictorle } 13823a40ed3dSBarry Smith PetscFunctionReturn(0); 13836f0a148fSBarry Smith } 13846f0a148fSBarry Smith 13854a2ae208SSatish Balay #undef __FUNCT__ 13864a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 13872b40b63fSBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 13886f0a148fSBarry Smith { 138997b48c8fSBarry Smith PetscErrorCode ierr; 1390ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1391d0f46423SBarry Smith PetscInt n = A->cmap->n,i,j; 139297b48c8fSBarry Smith PetscScalar *slot,*bb; 139397b48c8fSBarry Smith const PetscScalar *xx; 139455659b69SBarry Smith 13953a40ed3dSBarry Smith PetscFunctionBegin; 139697b48c8fSBarry Smith /* fix right hand side if needed */ 139797b48c8fSBarry Smith if (x && b) { 139897b48c8fSBarry Smith ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 139997b48c8fSBarry Smith ierr = VecGetArray(b,&bb);CHKERRQ(ierr); 140097b48c8fSBarry Smith for (i=0; i<N; i++) { 140197b48c8fSBarry Smith bb[rows[i]] = diag*xx[rows[i]]; 140297b48c8fSBarry Smith } 140397b48c8fSBarry Smith ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 140497b48c8fSBarry Smith ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr); 140597b48c8fSBarry Smith } 140697b48c8fSBarry Smith 14076f0a148fSBarry Smith for (i=0; i<N; i++) { 14086f0a148fSBarry Smith slot = l->v + rows[i]; 14096f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 14106f0a148fSBarry Smith } 1411f4df32b1SMatthew Knepley if (diag != 0.0) { 14126f0a148fSBarry Smith for (i=0; i<N; i++) { 14136f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1414f4df32b1SMatthew Knepley *slot = diag; 14156f0a148fSBarry Smith } 14166f0a148fSBarry Smith } 14173a40ed3dSBarry Smith PetscFunctionReturn(0); 14186f0a148fSBarry Smith } 1419557bce09SLois Curfman McInnes 14204a2ae208SSatish Balay #undef __FUNCT__ 14214a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 1422dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 142364e87e97SBarry Smith { 1424c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14253a40ed3dSBarry Smith 14263a40ed3dSBarry Smith PetscFunctionBegin; 1427e32f2f54SBarry 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"); 142864e87e97SBarry Smith *array = mat->v; 14293a40ed3dSBarry Smith PetscFunctionReturn(0); 143064e87e97SBarry Smith } 14310754003eSLois Curfman McInnes 14324a2ae208SSatish Balay #undef __FUNCT__ 14334a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 1434dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1435ff14e315SSatish Balay { 14363a40ed3dSBarry Smith PetscFunctionBegin; 143709b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 14383a40ed3dSBarry Smith PetscFunctionReturn(0); 1439ff14e315SSatish Balay } 14400754003eSLois Curfman McInnes 14414a2ae208SSatish Balay #undef __FUNCT__ 14424a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 144313f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 14440754003eSLois Curfman McInnes { 1445c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14466849ba73SBarry Smith PetscErrorCode ierr; 14475d0c19d7SBarry Smith PetscInt i,j,nrows,ncols; 14485d0c19d7SBarry Smith const PetscInt *irow,*icol; 144987828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 14500754003eSLois Curfman McInnes Mat newmat; 14510754003eSLois Curfman McInnes 14523a40ed3dSBarry Smith PetscFunctionBegin; 145378b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 145478b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1455e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1456e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 14570754003eSLois Curfman McInnes 1458182d2002SSatish Balay /* Check submatrixcall */ 1459182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 146013f74950SBarry Smith PetscInt n_cols,n_rows; 1461182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 146221a2c019SBarry Smith if (n_rows != nrows || n_cols != ncols) { 1463f746d493SDmitry Karpeev /* resize the result matrix to match number of requested rows/columns */ 1464c61587bbSBarry Smith ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 146521a2c019SBarry Smith } 1466182d2002SSatish Balay newmat = *B; 1467182d2002SSatish Balay } else { 14680754003eSLois Curfman McInnes /* Create and fill new matrix */ 14697adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr); 1470f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 14717adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 14725c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1473182d2002SSatish Balay } 1474182d2002SSatish Balay 1475182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1476182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1477182d2002SSatish Balay 1478182d2002SSatish Balay for (i=0; i<ncols; i++) { 14796de62eeeSBarry Smith av = v + mat->lda*icol[i]; 1480182d2002SSatish Balay for (j=0; j<nrows; j++) { 1481182d2002SSatish Balay *bv++ = av[irow[j]]; 14820754003eSLois Curfman McInnes } 14830754003eSLois Curfman McInnes } 1484182d2002SSatish Balay 1485182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 14866d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14876d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14880754003eSLois Curfman McInnes 14890754003eSLois Curfman McInnes /* Free work space */ 149078b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 149178b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1492182d2002SSatish Balay *B = newmat; 14933a40ed3dSBarry Smith PetscFunctionReturn(0); 14940754003eSLois Curfman McInnes } 14950754003eSLois Curfman McInnes 14964a2ae208SSatish Balay #undef __FUNCT__ 14974a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 149813f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1499905e6a2fSBarry Smith { 15006849ba73SBarry Smith PetscErrorCode ierr; 150113f74950SBarry Smith PetscInt i; 1502905e6a2fSBarry Smith 15033a40ed3dSBarry Smith PetscFunctionBegin; 1504905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1505b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1506905e6a2fSBarry Smith } 1507905e6a2fSBarry Smith 1508905e6a2fSBarry Smith for (i=0; i<n; i++) { 15096a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1510905e6a2fSBarry Smith } 15113a40ed3dSBarry Smith PetscFunctionReturn(0); 1512905e6a2fSBarry Smith } 1513905e6a2fSBarry Smith 15144a2ae208SSatish Balay #undef __FUNCT__ 1515c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1516c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1517c0aa2d19SHong Zhang { 1518c0aa2d19SHong Zhang PetscFunctionBegin; 1519c0aa2d19SHong Zhang PetscFunctionReturn(0); 1520c0aa2d19SHong Zhang } 1521c0aa2d19SHong Zhang 1522c0aa2d19SHong Zhang #undef __FUNCT__ 1523c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1524c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1525c0aa2d19SHong Zhang { 1526c0aa2d19SHong Zhang PetscFunctionBegin; 1527c0aa2d19SHong Zhang PetscFunctionReturn(0); 1528c0aa2d19SHong Zhang } 1529c0aa2d19SHong Zhang 1530c0aa2d19SHong Zhang #undef __FUNCT__ 15314a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1532dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 15334b0e389bSBarry Smith { 15344b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 15356849ba73SBarry Smith PetscErrorCode ierr; 1536d0f46423SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j; 15373a40ed3dSBarry Smith 15383a40ed3dSBarry Smith PetscFunctionBegin; 153933f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 154033f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1541cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 15423a40ed3dSBarry Smith PetscFunctionReturn(0); 15433a40ed3dSBarry Smith } 1544e32f2f54SBarry Smith if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1545a5ce6ee0Svictorle if (lda1>m || lda2>m) { 15460dbb7854Svictorle for (j=0; j<n; j++) { 1547a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1548a5ce6ee0Svictorle } 1549a5ce6ee0Svictorle } else { 1550d0f46423SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1551a5ce6ee0Svictorle } 1552273d9f13SBarry Smith PetscFunctionReturn(0); 1553273d9f13SBarry Smith } 1554273d9f13SBarry Smith 15554a2ae208SSatish Balay #undef __FUNCT__ 15564a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1557dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A) 1558273d9f13SBarry Smith { 1559dfbe8321SBarry Smith PetscErrorCode ierr; 1560273d9f13SBarry Smith 1561273d9f13SBarry Smith PetscFunctionBegin; 1562273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 15633a40ed3dSBarry Smith PetscFunctionReturn(0); 15644b0e389bSBarry Smith } 15654b0e389bSBarry Smith 1566284134d9SBarry Smith #undef __FUNCT__ 1567284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense" 1568284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N) 1569284134d9SBarry Smith { 1570284134d9SBarry Smith PetscFunctionBegin; 157121a2c019SBarry Smith /* this will not be called before lda, Mmax, and Nmax have been set */ 1572284134d9SBarry Smith m = PetscMax(m,M); 1573284134d9SBarry Smith n = PetscMax(n,N); 1574a868139aSShri Abhyankar 157586d161a7SShri 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); 157686d161a7SShri 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); 157786d161a7SShri Abhyankar */ 1578dc5cefdeSJed Brown A->rmap->n = A->rmap->N = m; 1579d0f46423SBarry Smith A->cmap->n = A->cmap->N = n; 1580284134d9SBarry Smith PetscFunctionReturn(0); 1581284134d9SBarry Smith } 1582170fe5c8SBarry Smith 1583ba337c44SJed Brown #undef __FUNCT__ 1584ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense" 1585ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A) 1586ba337c44SJed Brown { 1587ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1588ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1589ba337c44SJed Brown PetscScalar *aa = a->v; 1590ba337c44SJed Brown 1591ba337c44SJed Brown PetscFunctionBegin; 1592ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]); 1593ba337c44SJed Brown PetscFunctionReturn(0); 1594ba337c44SJed Brown } 1595ba337c44SJed Brown 1596ba337c44SJed Brown #undef __FUNCT__ 1597ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense" 1598ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A) 1599ba337c44SJed Brown { 1600ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1601ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1602ba337c44SJed Brown PetscScalar *aa = a->v; 1603ba337c44SJed Brown 1604ba337c44SJed Brown PetscFunctionBegin; 1605ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 1606ba337c44SJed Brown PetscFunctionReturn(0); 1607ba337c44SJed Brown } 1608ba337c44SJed Brown 1609ba337c44SJed Brown #undef __FUNCT__ 1610ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense" 1611ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A) 1612ba337c44SJed Brown { 1613ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1614ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1615ba337c44SJed Brown PetscScalar *aa = a->v; 1616ba337c44SJed Brown 1617ba337c44SJed Brown PetscFunctionBegin; 1618ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 1619ba337c44SJed Brown PetscFunctionReturn(0); 1620ba337c44SJed Brown } 1621284134d9SBarry Smith 1622a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/ 1623a9fe9ddaSSatish Balay #undef __FUNCT__ 1624a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1625a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1626a9fe9ddaSSatish Balay { 1627a9fe9ddaSSatish Balay PetscErrorCode ierr; 1628a9fe9ddaSSatish Balay 1629a9fe9ddaSSatish Balay PetscFunctionBegin; 1630a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1631a9fe9ddaSSatish Balay ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1632a9fe9ddaSSatish Balay } 1633a9fe9ddaSSatish Balay ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1634a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1635a9fe9ddaSSatish Balay } 1636a9fe9ddaSSatish Balay 1637a9fe9ddaSSatish Balay #undef __FUNCT__ 1638a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1639a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1640a9fe9ddaSSatish Balay { 1641ee16a9a1SHong Zhang PetscErrorCode ierr; 1642d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 1643ee16a9a1SHong Zhang Mat Cmat; 1644a9fe9ddaSSatish Balay 1645ee16a9a1SHong Zhang PetscFunctionBegin; 1646e32f2f54SBarry 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); 1647ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1648ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1649ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1650ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1651ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1652ee16a9a1SHong Zhang *C = Cmat; 1653ee16a9a1SHong Zhang PetscFunctionReturn(0); 1654ee16a9a1SHong Zhang } 1655a9fe9ddaSSatish Balay 165698a3b096SSatish Balay #undef __FUNCT__ 1657a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1658a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1659a9fe9ddaSSatish Balay { 1660a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1661a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1662a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 16630805154bSBarry Smith PetscBLASInt m,n,k; 1664a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1665a9fe9ddaSSatish Balay 1666a9fe9ddaSSatish Balay PetscFunctionBegin; 1667d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 1668d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1669d0f46423SBarry Smith k = PetscBLASIntCast(A->cmap->n); 1670a9fe9ddaSSatish Balay BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1671a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1672a9fe9ddaSSatish Balay } 1673a9fe9ddaSSatish Balay 1674a9fe9ddaSSatish Balay #undef __FUNCT__ 1675a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense" 1676a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1677a9fe9ddaSSatish Balay { 1678a9fe9ddaSSatish Balay PetscErrorCode ierr; 1679a9fe9ddaSSatish Balay 1680a9fe9ddaSSatish Balay PetscFunctionBegin; 1681a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1682a9fe9ddaSSatish Balay ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1683a9fe9ddaSSatish Balay } 1684a9fe9ddaSSatish Balay ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1685a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1686a9fe9ddaSSatish Balay } 1687a9fe9ddaSSatish Balay 1688a9fe9ddaSSatish Balay #undef __FUNCT__ 1689a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense" 1690a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1691a9fe9ddaSSatish Balay { 1692ee16a9a1SHong Zhang PetscErrorCode ierr; 1693d0f46423SBarry Smith PetscInt m=A->cmap->n,n=B->cmap->n; 1694ee16a9a1SHong Zhang Mat Cmat; 1695a9fe9ddaSSatish Balay 1696ee16a9a1SHong Zhang PetscFunctionBegin; 1697e32f2f54SBarry 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); 1698ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1699ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1700ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1701ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1702ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1703ee16a9a1SHong Zhang *C = Cmat; 1704ee16a9a1SHong Zhang PetscFunctionReturn(0); 1705ee16a9a1SHong Zhang } 1706a9fe9ddaSSatish Balay 1707a9fe9ddaSSatish Balay #undef __FUNCT__ 1708a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense" 1709a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1710a9fe9ddaSSatish Balay { 1711a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1712a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1713a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 17140805154bSBarry Smith PetscBLASInt m,n,k; 1715a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1716a9fe9ddaSSatish Balay 1717a9fe9ddaSSatish Balay PetscFunctionBegin; 1718d0f46423SBarry Smith m = PetscBLASIntCast(A->cmap->n); 1719d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1720d0f46423SBarry Smith k = PetscBLASIntCast(A->rmap->n); 17212fbe02b9SBarry Smith /* 17222fbe02b9SBarry Smith Note the m and n arguments below are the number rows and columns of A', not A! 17232fbe02b9SBarry Smith */ 1724a9fe9ddaSSatish Balay BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1725a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1726a9fe9ddaSSatish Balay } 1727985db425SBarry Smith 1728985db425SBarry Smith #undef __FUNCT__ 1729985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense" 1730985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 1731985db425SBarry Smith { 1732985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1733985db425SBarry Smith PetscErrorCode ierr; 1734d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1735985db425SBarry Smith PetscScalar *x; 1736985db425SBarry Smith MatScalar *aa = a->v; 1737985db425SBarry Smith 1738985db425SBarry Smith PetscFunctionBegin; 1739e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1740985db425SBarry Smith 1741985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1742985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1743985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1744e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1745985db425SBarry Smith for (i=0; i<m; i++) { 1746985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1747985db425SBarry Smith for (j=1; j<n; j++){ 1748985db425SBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1749985db425SBarry Smith } 1750985db425SBarry Smith } 1751985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1752985db425SBarry Smith PetscFunctionReturn(0); 1753985db425SBarry Smith } 1754985db425SBarry Smith 1755985db425SBarry Smith #undef __FUNCT__ 1756985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense" 1757985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 1758985db425SBarry Smith { 1759985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1760985db425SBarry Smith PetscErrorCode ierr; 1761d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1762985db425SBarry Smith PetscScalar *x; 1763985db425SBarry Smith PetscReal atmp; 1764985db425SBarry Smith MatScalar *aa = a->v; 1765985db425SBarry Smith 1766985db425SBarry Smith PetscFunctionBegin; 1767e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1768985db425SBarry Smith 1769985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1770985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1771985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1772e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1773985db425SBarry Smith for (i=0; i<m; i++) { 17749189402eSHong Zhang x[i] = PetscAbsScalar(aa[i]); 1775985db425SBarry Smith for (j=1; j<n; j++){ 1776985db425SBarry Smith atmp = PetscAbsScalar(aa[i+m*j]); 1777985db425SBarry Smith if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 1778985db425SBarry Smith } 1779985db425SBarry Smith } 1780985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1781985db425SBarry Smith PetscFunctionReturn(0); 1782985db425SBarry Smith } 1783985db425SBarry Smith 1784985db425SBarry Smith #undef __FUNCT__ 1785985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense" 1786985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 1787985db425SBarry Smith { 1788985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1789985db425SBarry Smith PetscErrorCode ierr; 1790d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1791985db425SBarry Smith PetscScalar *x; 1792985db425SBarry Smith MatScalar *aa = a->v; 1793985db425SBarry Smith 1794985db425SBarry Smith PetscFunctionBegin; 1795e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1796985db425SBarry Smith 1797985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1798985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1799985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1800e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1801985db425SBarry Smith for (i=0; i<m; i++) { 1802985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1803985db425SBarry Smith for (j=1; j<n; j++){ 1804985db425SBarry Smith if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1805985db425SBarry Smith } 1806985db425SBarry Smith } 1807985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1808985db425SBarry Smith PetscFunctionReturn(0); 1809985db425SBarry Smith } 1810985db425SBarry Smith 18118d0534beSBarry Smith #undef __FUNCT__ 18128d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense" 18138d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 18148d0534beSBarry Smith { 18158d0534beSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 18168d0534beSBarry Smith PetscErrorCode ierr; 18178d0534beSBarry Smith PetscScalar *x; 18188d0534beSBarry Smith 18198d0534beSBarry Smith PetscFunctionBegin; 1820e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 18218d0534beSBarry Smith 18228d0534beSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1823d0f46423SBarry Smith ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 18248d0534beSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 18258d0534beSBarry Smith PetscFunctionReturn(0); 18268d0534beSBarry Smith } 18278d0534beSBarry Smith 18280716a85fSBarry Smith 18290716a85fSBarry Smith #undef __FUNCT__ 18300716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqDense" 18310716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms) 18320716a85fSBarry Smith { 18330716a85fSBarry Smith PetscErrorCode ierr; 18340716a85fSBarry Smith PetscInt i,j,m,n; 18350716a85fSBarry Smith PetscScalar *a; 18360716a85fSBarry Smith 18370716a85fSBarry Smith PetscFunctionBegin; 18380716a85fSBarry Smith ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 18390716a85fSBarry Smith ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr); 18400716a85fSBarry Smith ierr = MatGetArray(A,&a);CHKERRQ(ierr); 18410716a85fSBarry Smith if (type == NORM_2) { 18420716a85fSBarry Smith for (i=0; i<n; i++ ){ 18430716a85fSBarry Smith for (j=0; j<m; j++) { 18440716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]*a[j]); 18450716a85fSBarry Smith } 18460716a85fSBarry Smith a += m; 18470716a85fSBarry Smith } 18480716a85fSBarry Smith } else if (type == NORM_1) { 18490716a85fSBarry Smith for (i=0; i<n; i++ ){ 18500716a85fSBarry Smith for (j=0; j<m; j++) { 18510716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]); 18520716a85fSBarry Smith } 18530716a85fSBarry Smith a += m; 18540716a85fSBarry Smith } 18550716a85fSBarry Smith } else if (type == NORM_INFINITY) { 18560716a85fSBarry Smith for (i=0; i<n; i++ ){ 18570716a85fSBarry Smith for (j=0; j<m; j++) { 18580716a85fSBarry Smith norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]); 18590716a85fSBarry Smith } 18600716a85fSBarry Smith a += m; 18610716a85fSBarry Smith } 18620716a85fSBarry Smith } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Unknown NormType"); 18630716a85fSBarry Smith if (type == NORM_2) { 18640716a85fSBarry Smith for (i=0; i<n; i++) norms[i] = sqrt(norms[i]); 18650716a85fSBarry Smith } 18660716a85fSBarry Smith PetscFunctionReturn(0); 18670716a85fSBarry Smith } 18680716a85fSBarry Smith 1869289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1870a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1871905e6a2fSBarry Smith MatGetRow_SeqDense, 1872905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1873905e6a2fSBarry Smith MatMult_SeqDense, 187497304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 18757c922b88SBarry Smith MatMultTranspose_SeqDense, 18767c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1877db4efbfdSBarry Smith 0, 1878db4efbfdSBarry Smith 0, 1879db4efbfdSBarry Smith 0, 1880db4efbfdSBarry Smith /*10*/ 0, 1881905e6a2fSBarry Smith MatLUFactor_SeqDense, 1882905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 188341f059aeSBarry Smith MatSOR_SeqDense, 1884ec8511deSBarry Smith MatTranspose_SeqDense, 188597304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1886905e6a2fSBarry Smith MatEqual_SeqDense, 1887905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1888905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1889905e6a2fSBarry Smith MatNorm_SeqDense, 1890c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense, 1891c0aa2d19SHong Zhang MatAssemblyEnd_SeqDense, 1892905e6a2fSBarry Smith MatSetOption_SeqDense, 1893905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1894d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense, 1895db4efbfdSBarry Smith 0, 1896db4efbfdSBarry Smith 0, 1897db4efbfdSBarry Smith 0, 1898db4efbfdSBarry Smith 0, 1899d519adbfSMatthew Knepley /*29*/ MatSetUpPreallocation_SeqDense, 1900273d9f13SBarry Smith 0, 1901905e6a2fSBarry Smith 0, 1902905e6a2fSBarry Smith MatGetArray_SeqDense, 1903905e6a2fSBarry Smith MatRestoreArray_SeqDense, 1904d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense, 1905a5ae1ecdSBarry Smith 0, 1906a5ae1ecdSBarry Smith 0, 1907a5ae1ecdSBarry Smith 0, 1908a5ae1ecdSBarry Smith 0, 1909d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense, 1910a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1911a5ae1ecdSBarry Smith 0, 19124b0e389bSBarry Smith MatGetValues_SeqDense, 1913a5ae1ecdSBarry Smith MatCopy_SeqDense, 1914d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense, 1915a5ae1ecdSBarry Smith MatScale_SeqDense, 1916a5ae1ecdSBarry Smith 0, 1917a5ae1ecdSBarry Smith 0, 1918a5ae1ecdSBarry Smith 0, 1919d519adbfSMatthew Knepley /*49*/ 0, 1920a5ae1ecdSBarry Smith 0, 1921a5ae1ecdSBarry Smith 0, 1922a5ae1ecdSBarry Smith 0, 1923a5ae1ecdSBarry Smith 0, 1924d519adbfSMatthew Knepley /*54*/ 0, 1925a5ae1ecdSBarry Smith 0, 1926a5ae1ecdSBarry Smith 0, 1927a5ae1ecdSBarry Smith 0, 1928a5ae1ecdSBarry Smith 0, 1929d519adbfSMatthew Knepley /*59*/ 0, 1930e03a110bSBarry Smith MatDestroy_SeqDense, 1931e03a110bSBarry Smith MatView_SeqDense, 1932357abbc8SBarry Smith 0, 193397304618SKris Buschelman 0, 1934d519adbfSMatthew Knepley /*64*/ 0, 193597304618SKris Buschelman 0, 193697304618SKris Buschelman 0, 193797304618SKris Buschelman 0, 193897304618SKris Buschelman 0, 1939d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense, 194097304618SKris Buschelman 0, 194197304618SKris Buschelman 0, 194297304618SKris Buschelman 0, 194397304618SKris Buschelman 0, 1944d519adbfSMatthew Knepley /*74*/ 0, 194597304618SKris Buschelman 0, 194697304618SKris Buschelman 0, 194797304618SKris Buschelman 0, 194897304618SKris Buschelman 0, 1949d519adbfSMatthew Knepley /*79*/ 0, 195097304618SKris Buschelman 0, 195197304618SKris Buschelman 0, 195297304618SKris Buschelman 0, 19535bba2384SShri Abhyankar /*83*/ MatLoad_SeqDense, 1954865e5f61SKris Buschelman 0, 19551cbb95d3SBarry Smith MatIsHermitian_SeqDense, 1956865e5f61SKris Buschelman 0, 1957865e5f61SKris Buschelman 0, 1958865e5f61SKris Buschelman 0, 1959d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense, 1960a9fe9ddaSSatish Balay MatMatMultSymbolic_SeqDense_SeqDense, 1961a9fe9ddaSSatish Balay MatMatMultNumeric_SeqDense_SeqDense, 1962865e5f61SKris Buschelman 0, 1963865e5f61SKris Buschelman 0, 1964d519adbfSMatthew Knepley /*94*/ 0, 1965a9fe9ddaSSatish Balay MatMatMultTranspose_SeqDense_SeqDense, 1966a9fe9ddaSSatish Balay MatMatMultTransposeSymbolic_SeqDense_SeqDense, 1967a9fe9ddaSSatish Balay MatMatMultTransposeNumeric_SeqDense_SeqDense, 1968284134d9SBarry Smith 0, 1969d519adbfSMatthew Knepley /*99*/ 0, 1970284134d9SBarry Smith 0, 1971284134d9SBarry Smith 0, 1972ba337c44SJed Brown MatConjugate_SeqDense, 1973985db425SBarry Smith MatSetSizes_SeqDense, 1974ba337c44SJed Brown /*104*/0, 1975ba337c44SJed Brown MatRealPart_SeqDense, 1976ba337c44SJed Brown MatImaginaryPart_SeqDense, 1977985db425SBarry Smith 0, 1978985db425SBarry Smith 0, 197985e2c93fSHong Zhang /*109*/MatMatSolve_SeqDense, 1980985db425SBarry Smith 0, 19818d0534beSBarry Smith MatGetRowMin_SeqDense, 1982aabbc4fbSShri Abhyankar MatGetColumnVector_SeqDense, 1983aabbc4fbSShri Abhyankar 0, 1984aabbc4fbSShri Abhyankar /*114*/0, 1985aabbc4fbSShri Abhyankar 0, 1986aabbc4fbSShri Abhyankar 0, 1987aabbc4fbSShri Abhyankar 0, 1988aabbc4fbSShri Abhyankar 0, 1989aabbc4fbSShri Abhyankar /*119*/0, 1990aabbc4fbSShri Abhyankar 0, 1991aabbc4fbSShri Abhyankar 0, 19920716a85fSBarry Smith 0, 19930716a85fSBarry Smith 0, 19940716a85fSBarry Smith /*124*/0, 19950716a85fSBarry Smith MatGetColumnNorms_SeqDense 1996985db425SBarry Smith }; 199790ace30eSBarry Smith 19984a2ae208SSatish Balay #undef __FUNCT__ 19994a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 20004b828684SBarry Smith /*@C 2001fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 2002d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 2003d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 2004289bc588SBarry Smith 2005db81eaa0SLois Curfman McInnes Collective on MPI_Comm 2006db81eaa0SLois Curfman McInnes 200720563c6bSBarry Smith Input Parameters: 2008db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 20090c775827SLois Curfman McInnes . m - number of rows 201018f449edSLois Curfman McInnes . n - number of columns 2011c0235b3cSMatthew Knepley - data - optional location of matrix data in column major order. Set data=PETSC_NULL for PETSc 2012dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 201320563c6bSBarry Smith 201420563c6bSBarry Smith Output Parameter: 201544cd7ae7SLois Curfman McInnes . A - the matrix 201620563c6bSBarry Smith 2017b259b22eSLois Curfman McInnes Notes: 201818f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 201918f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 2020b4fd4287SBarry Smith set data=PETSC_NULL. 202118f449edSLois Curfman McInnes 2022027ccd11SLois Curfman McInnes Level: intermediate 2023027ccd11SLois Curfman McInnes 2024dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 2025d65003e9SLois Curfman McInnes 2026db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 202720563c6bSBarry Smith @*/ 20287087cfbeSBarry Smith PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 2029289bc588SBarry Smith { 2030dfbe8321SBarry Smith PetscErrorCode ierr; 20313b2fbd54SBarry Smith 20323a40ed3dSBarry Smith PetscFunctionBegin; 2033f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 2034f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 2035273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 2036273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 2037273d9f13SBarry Smith PetscFunctionReturn(0); 2038273d9f13SBarry Smith } 2039273d9f13SBarry Smith 20404a2ae208SSatish Balay #undef __FUNCT__ 2041afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation" 2042273d9f13SBarry Smith /*@C 2043273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 2044273d9f13SBarry Smith 2045273d9f13SBarry Smith Collective on MPI_Comm 2046273d9f13SBarry Smith 2047273d9f13SBarry Smith Input Parameters: 2048273d9f13SBarry Smith + A - the matrix 2049273d9f13SBarry Smith - data - the array (or PETSC_NULL) 2050273d9f13SBarry Smith 2051273d9f13SBarry Smith Notes: 2052273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 2053273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 2054284134d9SBarry Smith need not call this routine. 2055273d9f13SBarry Smith 2056273d9f13SBarry Smith Level: intermediate 2057273d9f13SBarry Smith 2058273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 2059273d9f13SBarry Smith 2060867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA() 2061867c911aSBarry Smith 2062273d9f13SBarry Smith @*/ 20637087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 2064273d9f13SBarry Smith { 20654ac538c5SBarry Smith PetscErrorCode ierr; 2066a23d5eceSKris Buschelman 2067a23d5eceSKris Buschelman PetscFunctionBegin; 20684ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));CHKERRQ(ierr); 2069a23d5eceSKris Buschelman PetscFunctionReturn(0); 2070a23d5eceSKris Buschelman } 2071a23d5eceSKris Buschelman 2072a23d5eceSKris Buschelman EXTERN_C_BEGIN 2073a23d5eceSKris Buschelman #undef __FUNCT__ 2074afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense" 20757087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 2076a23d5eceSKris Buschelman { 2077273d9f13SBarry Smith Mat_SeqDense *b; 2078dfbe8321SBarry Smith PetscErrorCode ierr; 2079273d9f13SBarry Smith 2080273d9f13SBarry Smith PetscFunctionBegin; 2081273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 2082a868139aSShri Abhyankar 208334ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr); 208434ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr); 208534ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 208634ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 208734ef9618SShri Abhyankar 2088273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 208986d161a7SShri Abhyankar b->Mmax = B->rmap->n; 209086d161a7SShri Abhyankar b->Nmax = B->cmap->n; 209186d161a7SShri Abhyankar if(b->lda <= 0 || b->changelda) b->lda = B->rmap->n; 209286d161a7SShri Abhyankar 20939e8f95c4SLisandro Dalcin if (!data) { /* petsc-allocated storage */ 20949e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 20955afd5e0cSBarry Smith ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 2096284134d9SBarry Smith ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 2097284134d9SBarry Smith ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 20989e8f95c4SLisandro Dalcin b->user_alloc = PETSC_FALSE; 2099273d9f13SBarry Smith } else { /* user-allocated storage */ 21009e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 2101273d9f13SBarry Smith b->v = data; 2102273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 2103273d9f13SBarry Smith } 21040450473dSBarry Smith B->assembled = PETSC_TRUE; 2105273d9f13SBarry Smith PetscFunctionReturn(0); 2106273d9f13SBarry Smith } 2107a23d5eceSKris Buschelman EXTERN_C_END 2108273d9f13SBarry Smith 21091b807ce4Svictorle #undef __FUNCT__ 21101b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 21111b807ce4Svictorle /*@C 21121b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 21131b807ce4Svictorle 21141b807ce4Svictorle Input parameter: 21151b807ce4Svictorle + A - the matrix 21161b807ce4Svictorle - lda - the leading dimension 21171b807ce4Svictorle 21181b807ce4Svictorle Notes: 2119867c911aSBarry Smith This routine is to be used in conjunction with MatSeqDenseSetPreallocation(); 21201b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 21211b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 21221b807ce4Svictorle 21231b807ce4Svictorle Level: intermediate 21241b807ce4Svictorle 21251b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 21261b807ce4Svictorle 2127284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 2128867c911aSBarry Smith 21291b807ce4Svictorle @*/ 21307087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetLDA(Mat B,PetscInt lda) 21311b807ce4Svictorle { 21321b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 213321a2c019SBarry Smith 21341b807ce4Svictorle PetscFunctionBegin; 2135e32f2f54SBarry 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); 21361b807ce4Svictorle b->lda = lda; 213721a2c019SBarry Smith b->changelda = PETSC_FALSE; 213821a2c019SBarry Smith b->Mmax = PetscMax(b->Mmax,lda); 21391b807ce4Svictorle PetscFunctionReturn(0); 21401b807ce4Svictorle } 21411b807ce4Svictorle 21420bad9183SKris Buschelman /*MC 2143fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 21440bad9183SKris Buschelman 21450bad9183SKris Buschelman Options Database Keys: 21460bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 21470bad9183SKris Buschelman 21480bad9183SKris Buschelman Level: beginner 21490bad9183SKris Buschelman 215089665df3SBarry Smith .seealso: MatCreateSeqDense() 215189665df3SBarry Smith 21520bad9183SKris Buschelman M*/ 21530bad9183SKris Buschelman 2154273d9f13SBarry Smith EXTERN_C_BEGIN 21554a2ae208SSatish Balay #undef __FUNCT__ 21564a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 21577087cfbeSBarry Smith PetscErrorCode MatCreate_SeqDense(Mat B) 2158273d9f13SBarry Smith { 2159273d9f13SBarry Smith Mat_SeqDense *b; 2160dfbe8321SBarry Smith PetscErrorCode ierr; 21617c334f02SBarry Smith PetscMPIInt size; 2162273d9f13SBarry Smith 2163273d9f13SBarry Smith PetscFunctionBegin; 21647adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr); 2165e32f2f54SBarry Smith if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 216655659b69SBarry Smith 216738f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr); 2168549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 216944cd7ae7SLois Curfman McInnes B->data = (void*)b; 217018f449edSLois Curfman McInnes 217144cd7ae7SLois Curfman McInnes b->pivots = 0; 2172273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 2173273d9f13SBarry Smith b->v = 0; 217421a2c019SBarry Smith b->changelda = PETSC_FALSE; 21754e220ebcSLois Curfman McInnes 2176b24902e0SBarry Smith 2177ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C", 2178b24902e0SBarry Smith "MatGetFactor_seqdense_petsc", 2179b24902e0SBarry Smith MatGetFactor_seqdense_petsc);CHKERRQ(ierr); 2180a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 2181a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 2182a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 21834ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C", 21844ae313f4SHong Zhang "MatMatMult_SeqAIJ_SeqDense", 21854ae313f4SHong Zhang MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 21864ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C", 21874ae313f4SHong Zhang "MatMatMultSymbolic_SeqAIJ_SeqDense", 21884ae313f4SHong Zhang MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 21894ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C", 21904ae313f4SHong Zhang "MatMatMultNumeric_SeqAIJ_SeqDense", 21914ae313f4SHong Zhang MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 219217667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 21933a40ed3dSBarry Smith PetscFunctionReturn(0); 2194289bc588SBarry Smith } 2195273d9f13SBarry Smith EXTERN_C_END 2196