1be1d678aSKris Buschelman 267e560aaSBarry Smith /* 367e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 467e560aaSBarry Smith */ 5289bc588SBarry Smith 6dec5eb66SMatthew G Knepley #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 7c6db04a5SJed Brown #include <petscblaslapack.h> 8289bc588SBarry Smith 96a63e612SBarry Smith #include <../src/mat/impls/aij/seq/aij.h> 106a63e612SBarry Smith EXTERN_C_BEGIN 116a63e612SBarry Smith #undef __FUNCT__ 126a63e612SBarry Smith #define __FUNCT__ "MatConvert_SeqDense_SeqAIJ" 136a63e612SBarry Smith PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 146a63e612SBarry Smith { 156a63e612SBarry Smith Mat B; 166a63e612SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 176a63e612SBarry Smith PetscErrorCode ierr; 186a63e612SBarry Smith PetscInt i; 196a63e612SBarry Smith PetscInt *rows; 206a63e612SBarry Smith MatScalar *aa = a->v; 216a63e612SBarry Smith 226a63e612SBarry Smith PetscFunctionBegin; 236a63e612SBarry Smith ierr = MatCreate(((PetscObject)A)->comm,&B);CHKERRQ(ierr); 246a63e612SBarry Smith ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 256a63e612SBarry Smith ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); 266a63e612SBarry Smith ierr = MatSeqAIJSetPreallocation(B,A->cmap->n,PETSC_NULL);CHKERRQ(ierr); 276a63e612SBarry Smith 286a63e612SBarry Smith ierr = PetscMalloc(A->rmap->n*sizeof(PetscInt),&rows);CHKERRQ(ierr); 296a63e612SBarry Smith for (i=0; i<A->rmap->n; i++) rows[i] = i; 306a63e612SBarry Smith 316a63e612SBarry Smith for (i=0; i<A->cmap->n; i++) { 326a63e612SBarry Smith ierr = MatSetValues(B,A->rmap->n,rows,1,&i,aa,INSERT_VALUES);CHKERRQ(ierr); 336a63e612SBarry Smith aa += a->lda; 346a63e612SBarry Smith } 356a63e612SBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 366a63e612SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 376a63e612SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 386a63e612SBarry Smith 396a63e612SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 406a63e612SBarry Smith ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 416a63e612SBarry Smith } else { 426a63e612SBarry Smith *newmat = B; 436a63e612SBarry Smith } 446a63e612SBarry Smith PetscFunctionReturn(0); 456a63e612SBarry Smith } 466a63e612SBarry Smith EXTERN_C_END 476a63e612SBarry Smith 484a2ae208SSatish Balay #undef __FUNCT__ 494a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense" 50f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 511987afe7SBarry Smith { 521987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 53f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 5413f74950SBarry Smith PetscInt j; 550805154bSBarry Smith PetscBLASInt N,m,ldax,lday,one = 1; 56efee365bSSatish Balay PetscErrorCode ierr; 573a40ed3dSBarry Smith 583a40ed3dSBarry Smith PetscFunctionBegin; 59d0f46423SBarry Smith N = PetscBLASIntCast(X->rmap->n*X->cmap->n); 60d0f46423SBarry Smith m = PetscBLASIntCast(X->rmap->n); 610805154bSBarry Smith ldax = PetscBLASIntCast(x->lda); 620805154bSBarry Smith lday = PetscBLASIntCast(y->lda); 63a5ce6ee0Svictorle if (ldax>m || lday>m) { 64d0f46423SBarry Smith for (j=0; j<X->cmap->n; j++) { 65f4df32b1SMatthew Knepley BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one); 66a5ce6ee0Svictorle } 67a5ce6ee0Svictorle } else { 68f4df32b1SMatthew Knepley BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one); 69a5ce6ee0Svictorle } 700450473dSBarry Smith ierr = PetscLogFlops(PetscMax(2*N-1,0));CHKERRQ(ierr); 713a40ed3dSBarry Smith PetscFunctionReturn(0); 721987afe7SBarry Smith } 731987afe7SBarry Smith 744a2ae208SSatish Balay #undef __FUNCT__ 754a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 76dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 77289bc588SBarry Smith { 78d0f46423SBarry Smith PetscInt N = A->rmap->n*A->cmap->n; 793a40ed3dSBarry Smith 803a40ed3dSBarry Smith PetscFunctionBegin; 814e220ebcSLois Curfman McInnes info->block_size = 1.0; 824e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 836de62eeeSBarry Smith info->nz_used = (double)N; 846de62eeeSBarry Smith info->nz_unneeded = (double)0; 854e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 864e220ebcSLois Curfman McInnes info->mallocs = 0; 877adad957SLisandro Dalcin info->memory = ((PetscObject)A)->mem; 884e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 894e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 904e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 913a40ed3dSBarry Smith PetscFunctionReturn(0); 92289bc588SBarry Smith } 93289bc588SBarry Smith 944a2ae208SSatish Balay #undef __FUNCT__ 954a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense" 96f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha) 9780cd9d93SLois Curfman McInnes { 98273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 99f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 100efee365bSSatish Balay PetscErrorCode ierr; 1010805154bSBarry Smith PetscBLASInt one = 1,j,nz,lda = PetscBLASIntCast(a->lda); 10280cd9d93SLois Curfman McInnes 1033a40ed3dSBarry Smith PetscFunctionBegin; 104d0f46423SBarry Smith if (lda>A->rmap->n) { 105d0f46423SBarry Smith nz = PetscBLASIntCast(A->rmap->n); 106d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 107f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v+j*lda,&one); 108a5ce6ee0Svictorle } 109a5ce6ee0Svictorle } else { 110d0f46423SBarry Smith nz = PetscBLASIntCast(A->rmap->n*A->cmap->n); 111f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v,&one); 112a5ce6ee0Svictorle } 113efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 1143a40ed3dSBarry Smith PetscFunctionReturn(0); 11580cd9d93SLois Curfman McInnes } 11680cd9d93SLois Curfman McInnes 1171cbb95d3SBarry Smith #undef __FUNCT__ 1181cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense" 119ace3abfcSBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool *fl) 1201cbb95d3SBarry Smith { 1211cbb95d3SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 122d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,N; 1231cbb95d3SBarry Smith PetscScalar *v = a->v; 1241cbb95d3SBarry Smith 1251cbb95d3SBarry Smith PetscFunctionBegin; 1261cbb95d3SBarry Smith *fl = PETSC_FALSE; 127d0f46423SBarry Smith if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0); 1281cbb95d3SBarry Smith N = a->lda; 1291cbb95d3SBarry Smith 1301cbb95d3SBarry Smith for (i=0; i<m; i++) { 1311cbb95d3SBarry Smith for (j=i+1; j<m; j++) { 1321cbb95d3SBarry Smith if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0); 1331cbb95d3SBarry Smith } 1341cbb95d3SBarry Smith } 1351cbb95d3SBarry Smith *fl = PETSC_TRUE; 1361cbb95d3SBarry Smith PetscFunctionReturn(0); 1371cbb95d3SBarry Smith } 1381cbb95d3SBarry Smith 139b24902e0SBarry Smith #undef __FUNCT__ 140b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense" 141719d5645SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues) 142b24902e0SBarry Smith { 143b24902e0SBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 144b24902e0SBarry Smith PetscErrorCode ierr; 145b24902e0SBarry Smith PetscInt lda = (PetscInt)mat->lda,j,m; 146b24902e0SBarry Smith 147b24902e0SBarry Smith PetscFunctionBegin; 148aa5ea44dSBarry Smith ierr = PetscLayoutReference(A->rmap,&newi->rmap);CHKERRQ(ierr); 149aa5ea44dSBarry Smith ierr = PetscLayoutReference(A->cmap,&newi->cmap);CHKERRQ(ierr); 150719d5645SBarry Smith ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr); 151b24902e0SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 152b24902e0SBarry Smith l = (Mat_SeqDense*)newi->data; 153d0f46423SBarry Smith if (lda>A->rmap->n) { 154d0f46423SBarry Smith m = A->rmap->n; 155d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 156b24902e0SBarry Smith ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 157b24902e0SBarry Smith } 158b24902e0SBarry Smith } else { 159d0f46423SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 160b24902e0SBarry Smith } 161b24902e0SBarry Smith } 162b24902e0SBarry Smith newi->assembled = PETSC_TRUE; 163b24902e0SBarry Smith PetscFunctionReturn(0); 164b24902e0SBarry Smith } 165b24902e0SBarry Smith 1664a2ae208SSatish Balay #undef __FUNCT__ 1674a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 168dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 16902cad45dSBarry Smith { 1706849ba73SBarry Smith PetscErrorCode ierr; 17102cad45dSBarry Smith 1723a40ed3dSBarry Smith PetscFunctionBegin; 1735c9eb25fSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr); 174d0f46423SBarry Smith ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 1755c9eb25fSBarry Smith ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 176719d5645SBarry Smith ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr); 177b24902e0SBarry Smith PetscFunctionReturn(0); 178b24902e0SBarry Smith } 179b24902e0SBarry Smith 1806ee01492SSatish Balay 1810481f469SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*); 182719d5645SBarry Smith 1834a2ae208SSatish Balay #undef __FUNCT__ 1844a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 1850481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 186289bc588SBarry Smith { 1874482741eSBarry Smith MatFactorInfo info; 188a093e273SMatthew Knepley PetscErrorCode ierr; 1893a40ed3dSBarry Smith 1903a40ed3dSBarry Smith PetscFunctionBegin; 191c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 192719d5645SBarry Smith ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr); 1933a40ed3dSBarry Smith PetscFunctionReturn(0); 194289bc588SBarry Smith } 1956ee01492SSatish Balay 1960b4b3355SBarry Smith #undef __FUNCT__ 1974a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 198dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 199289bc588SBarry Smith { 200c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2016849ba73SBarry Smith PetscErrorCode ierr; 20287828ca2SBarry Smith PetscScalar *x,*y; 203d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 20467e560aaSBarry Smith 2053a40ed3dSBarry Smith PetscFunctionBegin; 2061ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2071ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 208d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 209d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_LU) { 210ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 211e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 212ae7cfcebSSatish Balay #else 21371044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 214e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 215ae7cfcebSSatish Balay #endif 216d5f3da31SBarry Smith } else if (A->factortype == MAT_FACTOR_CHOLESKY){ 217ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 218e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 219ae7cfcebSSatish Balay #else 22071044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 221e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 222ae7cfcebSSatish Balay #endif 223289bc588SBarry Smith } 224e32f2f54SBarry Smith else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 2251ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2261ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 227dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 2283a40ed3dSBarry Smith PetscFunctionReturn(0); 229289bc588SBarry Smith } 2306ee01492SSatish Balay 2314a2ae208SSatish Balay #undef __FUNCT__ 23285e2c93fSHong Zhang #define __FUNCT__ "MatMatSolve_SeqDense" 23385e2c93fSHong Zhang PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X) 23485e2c93fSHong Zhang { 23585e2c93fSHong Zhang Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 23685e2c93fSHong Zhang PetscErrorCode ierr; 23785e2c93fSHong Zhang PetscScalar *b,*x; 238efb80c78SLisandro Dalcin PetscInt n; 23985e2c93fSHong Zhang PetscBLASInt nrhs,info,m=PetscBLASIntCast(A->rmap->n); 240bda8bf91SBarry Smith PetscBool flg; 24185e2c93fSHong Zhang 24285e2c93fSHong Zhang PetscFunctionBegin; 243251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr); 244bda8bf91SBarry Smith if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 245251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,PETSC_NULL);CHKERRQ(ierr); 246bda8bf91SBarry Smith if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 247bda8bf91SBarry Smith 248efb80c78SLisandro Dalcin ierr = MatGetSize(B,PETSC_NULL,&n);CHKERRQ(ierr); 249efb80c78SLisandro Dalcin nrhs = PetscBLASIntCast(n); 2508c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 2518c778c55SBarry Smith ierr = MatDenseGetArray(X,&x);CHKERRQ(ierr); 25285e2c93fSHong Zhang 253f272c775SHong Zhang ierr = PetscMemcpy(x,b,m*nrhs*sizeof(PetscScalar));CHKERRQ(ierr); 25485e2c93fSHong Zhang 25585e2c93fSHong Zhang if (A->factortype == MAT_FACTOR_LU) { 25685e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_GETRS) 25785e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 25885e2c93fSHong Zhang #else 25985e2c93fSHong Zhang LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info); 26085e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 26185e2c93fSHong Zhang #endif 26285e2c93fSHong Zhang } else if (A->factortype == MAT_FACTOR_CHOLESKY){ 26385e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_POTRS) 26485e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 26585e2c93fSHong Zhang #else 26685e2c93fSHong Zhang LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info); 26785e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 26885e2c93fSHong Zhang #endif 26985e2c93fSHong Zhang } 27085e2c93fSHong Zhang else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 27185e2c93fSHong Zhang 2728c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 2738c778c55SBarry Smith ierr = MatDenseRestoreArray(X,&x);CHKERRQ(ierr); 27485e2c93fSHong Zhang ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr); 27585e2c93fSHong Zhang PetscFunctionReturn(0); 27685e2c93fSHong Zhang } 27785e2c93fSHong Zhang 27885e2c93fSHong Zhang #undef __FUNCT__ 2794a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 280dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 281da3a660dSBarry Smith { 282c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 283dfbe8321SBarry Smith PetscErrorCode ierr; 28487828ca2SBarry Smith PetscScalar *x,*y; 285d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 28667e560aaSBarry Smith 2873a40ed3dSBarry Smith PetscFunctionBegin; 2881ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2891ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 290d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 291752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 292da3a660dSBarry Smith if (mat->pivots) { 293ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 294e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 295ae7cfcebSSatish Balay #else 29671044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 297e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 298ae7cfcebSSatish Balay #endif 2997a97a34bSBarry Smith } else { 300ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 301e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 302ae7cfcebSSatish Balay #else 30371044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 304e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 305ae7cfcebSSatish Balay #endif 306da3a660dSBarry Smith } 3071ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3081ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 309dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 3103a40ed3dSBarry Smith PetscFunctionReturn(0); 311da3a660dSBarry Smith } 3126ee01492SSatish Balay 3134a2ae208SSatish Balay #undef __FUNCT__ 3144a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 315dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 316da3a660dSBarry Smith { 317c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 318dfbe8321SBarry Smith PetscErrorCode ierr; 31987828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 320da3a660dSBarry Smith Vec tmp = 0; 321d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 32267e560aaSBarry Smith 3233a40ed3dSBarry Smith PetscFunctionBegin; 3241ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3251ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 326d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 327da3a660dSBarry Smith if (yy == zz) { 32878b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 32952e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 33078b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 331da3a660dSBarry Smith } 332d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 333752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 334da3a660dSBarry Smith if (mat->pivots) { 335ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 336e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 337ae7cfcebSSatish Balay #else 33871044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 339e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 340ae7cfcebSSatish Balay #endif 341a8c6a408SBarry Smith } else { 342ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 343e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 344ae7cfcebSSatish Balay #else 34571044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 346e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 347ae7cfcebSSatish Balay #endif 348da3a660dSBarry Smith } 3496bf464f9SBarry Smith if (tmp) { 3506bf464f9SBarry Smith ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 3516bf464f9SBarry Smith ierr = VecDestroy(&tmp);CHKERRQ(ierr); 3526bf464f9SBarry Smith } else { 3536bf464f9SBarry Smith ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 3546bf464f9SBarry Smith } 3551ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3561ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 357dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 3583a40ed3dSBarry Smith PetscFunctionReturn(0); 359da3a660dSBarry Smith } 36067e560aaSBarry Smith 3614a2ae208SSatish Balay #undef __FUNCT__ 3624a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 363dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 364da3a660dSBarry Smith { 365c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3666849ba73SBarry Smith PetscErrorCode ierr; 36787828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 368da3a660dSBarry Smith Vec tmp; 369d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 37067e560aaSBarry Smith 3713a40ed3dSBarry Smith PetscFunctionBegin; 372d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 3731ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3741ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 375da3a660dSBarry Smith if (yy == zz) { 37678b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 37752e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 37878b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 379da3a660dSBarry Smith } 380d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 381752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 382da3a660dSBarry Smith if (mat->pivots) { 383ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 384e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 385ae7cfcebSSatish Balay #else 38671044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 387e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 388ae7cfcebSSatish Balay #endif 3893a40ed3dSBarry Smith } else { 390ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 391e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 392ae7cfcebSSatish Balay #else 39371044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 394e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 395ae7cfcebSSatish Balay #endif 396da3a660dSBarry Smith } 39790f02eecSBarry Smith if (tmp) { 3982dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 3996bf464f9SBarry Smith ierr = VecDestroy(&tmp);CHKERRQ(ierr); 4003a40ed3dSBarry Smith } else { 4012dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 40290f02eecSBarry Smith } 4031ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4041ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 405dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 4063a40ed3dSBarry Smith PetscFunctionReturn(0); 407da3a660dSBarry Smith } 408db4efbfdSBarry Smith 409db4efbfdSBarry Smith /* ---------------------------------------------------------------*/ 410db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 411db4efbfdSBarry Smith rather than put it in the Mat->row slot.*/ 412db4efbfdSBarry Smith #undef __FUNCT__ 413db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense" 4140481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo) 415db4efbfdSBarry Smith { 416db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 417db4efbfdSBarry Smith PetscFunctionBegin; 418e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 419db4efbfdSBarry Smith #else 420db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 421db4efbfdSBarry Smith PetscErrorCode ierr; 422db4efbfdSBarry Smith PetscBLASInt n,m,info; 423db4efbfdSBarry Smith 424db4efbfdSBarry Smith PetscFunctionBegin; 425db4efbfdSBarry Smith n = PetscBLASIntCast(A->cmap->n); 426db4efbfdSBarry Smith m = PetscBLASIntCast(A->rmap->n); 427db4efbfdSBarry Smith if (!mat->pivots) { 428db4efbfdSBarry Smith ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr); 429db4efbfdSBarry Smith ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr); 430db4efbfdSBarry Smith } 431db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 4328e57ea43SSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 433db4efbfdSBarry Smith LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info); 4348e57ea43SSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 4358e57ea43SSatish Balay 436e32f2f54SBarry Smith if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization"); 437e32f2f54SBarry Smith if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 438db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 439db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 440db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 441db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 442d5f3da31SBarry Smith A->factortype = MAT_FACTOR_LU; 443db4efbfdSBarry Smith 444dc0b31edSSatish Balay ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr); 445db4efbfdSBarry Smith #endif 446db4efbfdSBarry Smith PetscFunctionReturn(0); 447db4efbfdSBarry Smith } 448db4efbfdSBarry Smith 449db4efbfdSBarry Smith #undef __FUNCT__ 450db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense" 4510481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo) 452db4efbfdSBarry Smith { 453db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 454db4efbfdSBarry Smith PetscFunctionBegin; 455e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 456db4efbfdSBarry Smith #else 457db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 458db4efbfdSBarry Smith PetscErrorCode ierr; 459db4efbfdSBarry Smith PetscBLASInt info,n = PetscBLASIntCast(A->cmap->n); 460db4efbfdSBarry Smith 461db4efbfdSBarry Smith PetscFunctionBegin; 462db4efbfdSBarry Smith ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 463db4efbfdSBarry Smith 464db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 465db4efbfdSBarry Smith LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info); 466e32f2f54SBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 467db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 468db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 469db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 470db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 471d5f3da31SBarry Smith A->factortype = MAT_FACTOR_CHOLESKY; 472dc0b31edSSatish Balay ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr); 473db4efbfdSBarry Smith #endif 474db4efbfdSBarry Smith PetscFunctionReturn(0); 475db4efbfdSBarry Smith } 476db4efbfdSBarry Smith 477db4efbfdSBarry Smith 478db4efbfdSBarry Smith #undef __FUNCT__ 479db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 4800481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 481db4efbfdSBarry Smith { 482db4efbfdSBarry Smith PetscErrorCode ierr; 483db4efbfdSBarry Smith MatFactorInfo info; 484db4efbfdSBarry Smith 485db4efbfdSBarry Smith PetscFunctionBegin; 486db4efbfdSBarry Smith info.fill = 1.0; 487c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 488719d5645SBarry Smith ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr); 489db4efbfdSBarry Smith PetscFunctionReturn(0); 490db4efbfdSBarry Smith } 491db4efbfdSBarry Smith 492db4efbfdSBarry Smith #undef __FUNCT__ 493db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 4940481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info) 495db4efbfdSBarry Smith { 496db4efbfdSBarry Smith PetscFunctionBegin; 497c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 4981bbcc794SSatish Balay fact->preallocated = PETSC_TRUE; 499719d5645SBarry Smith fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense; 500db4efbfdSBarry Smith PetscFunctionReturn(0); 501db4efbfdSBarry Smith } 502db4efbfdSBarry Smith 503db4efbfdSBarry Smith #undef __FUNCT__ 504db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 5050481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 506db4efbfdSBarry Smith { 507db4efbfdSBarry Smith PetscFunctionBegin; 508b66fe19dSMatthew G Knepley fact->preallocated = PETSC_TRUE; 509c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 510719d5645SBarry Smith fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense; 511db4efbfdSBarry Smith PetscFunctionReturn(0); 512db4efbfdSBarry Smith } 513db4efbfdSBarry Smith 514bb5747d9SMatthew Knepley EXTERN_C_BEGIN 515db4efbfdSBarry Smith #undef __FUNCT__ 516db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc" 517db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact) 518db4efbfdSBarry Smith { 519db4efbfdSBarry Smith PetscErrorCode ierr; 520db4efbfdSBarry Smith 521db4efbfdSBarry Smith PetscFunctionBegin; 522db4efbfdSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr); 523db4efbfdSBarry Smith ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 524db4efbfdSBarry Smith ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr); 525db4efbfdSBarry Smith if (ftype == MAT_FACTOR_LU){ 526db4efbfdSBarry Smith (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense; 527db4efbfdSBarry Smith } else { 528db4efbfdSBarry Smith (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense; 529db4efbfdSBarry Smith } 530d5f3da31SBarry Smith (*fact)->factortype = ftype; 531db4efbfdSBarry Smith PetscFunctionReturn(0); 532db4efbfdSBarry Smith } 533bb5747d9SMatthew Knepley EXTERN_C_END 534db4efbfdSBarry Smith 535289bc588SBarry Smith /* ------------------------------------------------------------------*/ 5364a2ae208SSatish Balay #undef __FUNCT__ 53741f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense" 53841f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 539289bc588SBarry Smith { 540c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 54187828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 542dfbe8321SBarry Smith PetscErrorCode ierr; 543d0f46423SBarry Smith PetscInt m = A->rmap->n,i; 5440805154bSBarry Smith PetscBLASInt o = 1,bm = PetscBLASIntCast(m); 545289bc588SBarry Smith 5463a40ed3dSBarry Smith PetscFunctionBegin; 547289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 54871044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 5492dcb1b2aSMatthew Knepley ierr = VecSet(xx,zero);CHKERRQ(ierr); 550289bc588SBarry Smith } 5511ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5521ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 553b965ef7fSBarry Smith its = its*lits; 554e32f2f54SBarry 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); 555289bc588SBarry Smith while (its--) { 556fccaa45eSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 557289bc588SBarry Smith for (i=0; i<m; i++) { 55871044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 55955a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 560289bc588SBarry Smith } 561289bc588SBarry Smith } 562fccaa45eSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 563289bc588SBarry Smith for (i=m-1; i>=0; i--) { 56471044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 56555a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 566289bc588SBarry Smith } 567289bc588SBarry Smith } 568289bc588SBarry Smith } 5691ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 5701ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5713a40ed3dSBarry Smith PetscFunctionReturn(0); 572289bc588SBarry Smith } 573289bc588SBarry Smith 574289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5754a2ae208SSatish Balay #undef __FUNCT__ 5764a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 577dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 578289bc588SBarry Smith { 579c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 58087828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 581dfbe8321SBarry Smith PetscErrorCode ierr; 5820805154bSBarry Smith PetscBLASInt m, n,_One=1; 583ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 5843a40ed3dSBarry Smith 5853a40ed3dSBarry Smith PetscFunctionBegin; 586d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 587d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 588d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5891ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5901ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 59171044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One); 5921ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5931ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 594dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 5953a40ed3dSBarry Smith PetscFunctionReturn(0); 596289bc588SBarry Smith } 597800995b7SMatthew Knepley 5984a2ae208SSatish Balay #undef __FUNCT__ 5994a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 600dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 601289bc588SBarry Smith { 602c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 60387828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 604dfbe8321SBarry Smith PetscErrorCode ierr; 6050805154bSBarry Smith PetscBLASInt m, n, _One=1; 6063a40ed3dSBarry Smith 6073a40ed3dSBarry Smith PetscFunctionBegin; 608d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 609d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 610d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 6111ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6121ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 61371044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 6141ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6151ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 616dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr); 6173a40ed3dSBarry Smith PetscFunctionReturn(0); 618289bc588SBarry Smith } 6196ee01492SSatish Balay 6204a2ae208SSatish Balay #undef __FUNCT__ 6214a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 622dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 623289bc588SBarry Smith { 624c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 62587828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 626dfbe8321SBarry Smith PetscErrorCode ierr; 6270805154bSBarry Smith PetscBLASInt m, n, _One=1; 6283a40ed3dSBarry Smith 6293a40ed3dSBarry Smith PetscFunctionBegin; 630d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 631d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 632d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 633600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6341ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6351ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 63671044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 6371ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6381ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 639dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6403a40ed3dSBarry Smith PetscFunctionReturn(0); 641289bc588SBarry Smith } 6426ee01492SSatish Balay 6434a2ae208SSatish Balay #undef __FUNCT__ 6444a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 645dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 646289bc588SBarry Smith { 647c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 64887828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 649dfbe8321SBarry Smith PetscErrorCode ierr; 6500805154bSBarry Smith PetscBLASInt m, n, _One=1; 65187828ca2SBarry Smith PetscScalar _DOne=1.0; 6523a40ed3dSBarry Smith 6533a40ed3dSBarry Smith PetscFunctionBegin; 654d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 655d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 656d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 657600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6581ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6591ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 66071044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 6611ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6621ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 663dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6643a40ed3dSBarry Smith PetscFunctionReturn(0); 665289bc588SBarry Smith } 666289bc588SBarry Smith 667289bc588SBarry Smith /* -----------------------------------------------------------------*/ 6684a2ae208SSatish Balay #undef __FUNCT__ 6694a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 67013f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 671289bc588SBarry Smith { 672c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 67387828ca2SBarry Smith PetscScalar *v; 6746849ba73SBarry Smith PetscErrorCode ierr; 67513f74950SBarry Smith PetscInt i; 67667e560aaSBarry Smith 6773a40ed3dSBarry Smith PetscFunctionBegin; 678d0f46423SBarry Smith *ncols = A->cmap->n; 679289bc588SBarry Smith if (cols) { 680d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr); 681d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) (*cols)[i] = i; 682289bc588SBarry Smith } 683289bc588SBarry Smith if (vals) { 684d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 685289bc588SBarry Smith v = mat->v + row; 686d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;} 687289bc588SBarry Smith } 6883a40ed3dSBarry Smith PetscFunctionReturn(0); 689289bc588SBarry Smith } 6906ee01492SSatish Balay 6914a2ae208SSatish Balay #undef __FUNCT__ 6924a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 69313f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 694289bc588SBarry Smith { 695dfbe8321SBarry Smith PetscErrorCode ierr; 696606d414cSSatish Balay PetscFunctionBegin; 697606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 698606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 6993a40ed3dSBarry Smith PetscFunctionReturn(0); 700289bc588SBarry Smith } 701289bc588SBarry Smith /* ----------------------------------------------------------------*/ 7024a2ae208SSatish Balay #undef __FUNCT__ 7034a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 70413f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 705289bc588SBarry Smith { 706c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 70713f74950SBarry Smith PetscInt i,j,idx=0; 708d6dfbf8fSBarry Smith 7093a40ed3dSBarry Smith PetscFunctionBegin; 71071fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 711289bc588SBarry Smith if (!mat->roworiented) { 712dbb450caSBarry Smith if (addv == INSERT_VALUES) { 713289bc588SBarry Smith for (j=0; j<n; j++) { 714cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 7152515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 716e32f2f54SBarry 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); 71758804f6dSBarry Smith #endif 718289bc588SBarry Smith for (i=0; i<m; i++) { 719cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 7202515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 721e32f2f54SBarry 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); 72258804f6dSBarry Smith #endif 723cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 724289bc588SBarry Smith } 725289bc588SBarry Smith } 7263a40ed3dSBarry Smith } else { 727289bc588SBarry Smith for (j=0; j<n; j++) { 728cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 7292515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 730e32f2f54SBarry 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); 73158804f6dSBarry Smith #endif 732289bc588SBarry Smith for (i=0; i<m; i++) { 733cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 7342515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 735e32f2f54SBarry 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); 73658804f6dSBarry Smith #endif 737cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 738289bc588SBarry Smith } 739289bc588SBarry Smith } 740289bc588SBarry Smith } 7413a40ed3dSBarry Smith } else { 742dbb450caSBarry Smith if (addv == INSERT_VALUES) { 743e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 744cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7452515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 746e32f2f54SBarry 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); 74758804f6dSBarry Smith #endif 748e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 749cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7502515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 751e32f2f54SBarry 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); 75258804f6dSBarry Smith #endif 753cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 754e8d4e0b9SBarry Smith } 755e8d4e0b9SBarry Smith } 7563a40ed3dSBarry Smith } else { 757289bc588SBarry Smith for (i=0; i<m; i++) { 758cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7592515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 760e32f2f54SBarry 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); 76158804f6dSBarry Smith #endif 762289bc588SBarry Smith for (j=0; j<n; j++) { 763cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7642515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 765e32f2f54SBarry 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); 76658804f6dSBarry Smith #endif 767cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 768289bc588SBarry Smith } 769289bc588SBarry Smith } 770289bc588SBarry Smith } 771e8d4e0b9SBarry Smith } 7723a40ed3dSBarry Smith PetscFunctionReturn(0); 773289bc588SBarry Smith } 774e8d4e0b9SBarry Smith 7754a2ae208SSatish Balay #undef __FUNCT__ 7764a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 77713f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 778ae80bb75SLois Curfman McInnes { 779ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 78013f74950SBarry Smith PetscInt i,j; 781ae80bb75SLois Curfman McInnes 7823a40ed3dSBarry Smith PetscFunctionBegin; 783ae80bb75SLois Curfman McInnes /* row-oriented output */ 784ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 78597e567efSBarry Smith if (indexm[i] < 0) {v += n;continue;} 786e32f2f54SBarry 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); 787ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 7886f31f424SBarry Smith if (indexn[j] < 0) {v++; continue;} 789e32f2f54SBarry 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); 79097e567efSBarry Smith *v++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 791ae80bb75SLois Curfman McInnes } 792ae80bb75SLois Curfman McInnes } 7933a40ed3dSBarry Smith PetscFunctionReturn(0); 794ae80bb75SLois Curfman McInnes } 795ae80bb75SLois Curfman McInnes 796289bc588SBarry Smith /* -----------------------------------------------------------------*/ 797289bc588SBarry Smith 7984a2ae208SSatish Balay #undef __FUNCT__ 7995bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_SeqDense" 800112444f4SShri Abhyankar PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer) 801aabbc4fbSShri Abhyankar { 802aabbc4fbSShri Abhyankar Mat_SeqDense *a; 803aabbc4fbSShri Abhyankar PetscErrorCode ierr; 804aabbc4fbSShri Abhyankar PetscInt *scols,i,j,nz,header[4]; 805aabbc4fbSShri Abhyankar int fd; 806aabbc4fbSShri Abhyankar PetscMPIInt size; 807aabbc4fbSShri Abhyankar PetscInt *rowlengths = 0,M,N,*cols,grows,gcols; 808aabbc4fbSShri Abhyankar PetscScalar *vals,*svals,*v,*w; 809aabbc4fbSShri Abhyankar MPI_Comm comm = ((PetscObject)viewer)->comm; 810aabbc4fbSShri Abhyankar 811aabbc4fbSShri Abhyankar PetscFunctionBegin; 812aabbc4fbSShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 813aabbc4fbSShri Abhyankar if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor"); 814aabbc4fbSShri Abhyankar ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 815aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 816aabbc4fbSShri Abhyankar if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 817aabbc4fbSShri Abhyankar M = header[1]; N = header[2]; nz = header[3]; 818aabbc4fbSShri Abhyankar 819aabbc4fbSShri Abhyankar /* set global size if not set already*/ 820aabbc4fbSShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 821aabbc4fbSShri Abhyankar ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr); 822aabbc4fbSShri Abhyankar } else { 823aabbc4fbSShri Abhyankar /* if sizes and type are already set, check if the vector global sizes are correct */ 824aabbc4fbSShri Abhyankar ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr); 825aabbc4fbSShri 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); 826aabbc4fbSShri Abhyankar } 827aabbc4fbSShri Abhyankar ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 828aabbc4fbSShri Abhyankar 829aabbc4fbSShri Abhyankar if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 830aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 831aabbc4fbSShri Abhyankar v = a->v; 832aabbc4fbSShri Abhyankar /* Allocate some temp space to read in the values and then flip them 833aabbc4fbSShri Abhyankar from row major to column major */ 834aabbc4fbSShri Abhyankar ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 835aabbc4fbSShri Abhyankar /* read in nonzero values */ 836aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 837aabbc4fbSShri Abhyankar /* now flip the values and store them in the matrix*/ 838aabbc4fbSShri Abhyankar for (j=0; j<N; j++) { 839aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 840aabbc4fbSShri Abhyankar *v++ =w[i*N+j]; 841aabbc4fbSShri Abhyankar } 842aabbc4fbSShri Abhyankar } 843aabbc4fbSShri Abhyankar ierr = PetscFree(w);CHKERRQ(ierr); 844aabbc4fbSShri Abhyankar } else { 845aabbc4fbSShri Abhyankar /* read row lengths */ 846aabbc4fbSShri Abhyankar ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 847aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 848aabbc4fbSShri Abhyankar 849aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 850aabbc4fbSShri Abhyankar v = a->v; 851aabbc4fbSShri Abhyankar 852aabbc4fbSShri Abhyankar /* read column indices and nonzeros */ 853aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr); 854aabbc4fbSShri Abhyankar cols = scols; 855aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 856aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 857aabbc4fbSShri Abhyankar vals = svals; 858aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 859aabbc4fbSShri Abhyankar 860aabbc4fbSShri Abhyankar /* insert into matrix */ 861aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 862aabbc4fbSShri Abhyankar for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 863aabbc4fbSShri Abhyankar svals += rowlengths[i]; scols += rowlengths[i]; 864aabbc4fbSShri Abhyankar } 865aabbc4fbSShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 866aabbc4fbSShri Abhyankar ierr = PetscFree(cols);CHKERRQ(ierr); 867aabbc4fbSShri Abhyankar ierr = PetscFree(rowlengths);CHKERRQ(ierr); 868aabbc4fbSShri Abhyankar } 869aabbc4fbSShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 870aabbc4fbSShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 871aabbc4fbSShri Abhyankar 872aabbc4fbSShri Abhyankar PetscFunctionReturn(0); 873aabbc4fbSShri Abhyankar } 874aabbc4fbSShri Abhyankar 875aabbc4fbSShri Abhyankar #undef __FUNCT__ 8764a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 8776849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 878289bc588SBarry Smith { 879932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 880dfbe8321SBarry Smith PetscErrorCode ierr; 88113f74950SBarry Smith PetscInt i,j; 8822dcb1b2aSMatthew Knepley const char *name; 88387828ca2SBarry Smith PetscScalar *v; 884f3ef73ceSBarry Smith PetscViewerFormat format; 8855f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX) 886ace3abfcSBarry Smith PetscBool allreal = PETSC_TRUE; 8875f481a85SSatish Balay #endif 888932b0c3eSLois Curfman McInnes 8893a40ed3dSBarry Smith PetscFunctionBegin; 890b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 891456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 8923a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 893fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 894d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 8957566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 896d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 89744cd7ae7SLois Curfman McInnes v = a->v + i; 89877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 899d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 900aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 901329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 902a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 903329f5518SBarry Smith } else if (PetscRealPart(*v)) { 904a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr); 9056831982aSBarry Smith } 90680cd9d93SLois Curfman McInnes #else 9076831982aSBarry Smith if (*v) { 908a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr); 9096831982aSBarry Smith } 91080cd9d93SLois Curfman McInnes #endif 9111b807ce4Svictorle v += a->lda; 91280cd9d93SLois Curfman McInnes } 913b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 91480cd9d93SLois Curfman McInnes } 915d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 9163a40ed3dSBarry Smith } else { 917d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 918aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 91947989497SBarry Smith /* determine if matrix has all real values */ 92047989497SBarry Smith v = a->v; 921d0f46423SBarry Smith for (i=0; i<A->rmap->n*A->cmap->n; i++) { 922ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 92347989497SBarry Smith } 92447989497SBarry Smith #endif 925fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 9263a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 927d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr); 928d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 929fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 9307566de4bSShri Abhyankar } else { 9317566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 932ffac6cdbSBarry Smith } 933ffac6cdbSBarry Smith 934d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 935932b0c3eSLois Curfman McInnes v = a->v + i; 936d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 937aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 93847989497SBarry Smith if (allreal) { 939c61cd2faSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));CHKERRQ(ierr); 94047989497SBarry Smith } else { 941c61cd2faSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));CHKERRQ(ierr); 94247989497SBarry Smith } 943289bc588SBarry Smith #else 944c61cd2faSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);CHKERRQ(ierr); 945289bc588SBarry Smith #endif 9461b807ce4Svictorle v += a->lda; 947289bc588SBarry Smith } 948b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 949289bc588SBarry Smith } 950fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 951b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 952ffac6cdbSBarry Smith } 953d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 954da3a660dSBarry Smith } 955b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 9563a40ed3dSBarry Smith PetscFunctionReturn(0); 957289bc588SBarry Smith } 958289bc588SBarry Smith 9594a2ae208SSatish Balay #undef __FUNCT__ 9604a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 9616849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 962932b0c3eSLois Curfman McInnes { 963932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9646849ba73SBarry Smith PetscErrorCode ierr; 96513f74950SBarry Smith int fd; 966d0f46423SBarry Smith PetscInt ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n; 967f4403165SShri Abhyankar PetscScalar *v,*anonz,*vals; 968f4403165SShri Abhyankar PetscViewerFormat format; 969932b0c3eSLois Curfman McInnes 9703a40ed3dSBarry Smith PetscFunctionBegin; 971b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 97290ace30eSBarry Smith 973f4403165SShri Abhyankar ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 974f4403165SShri Abhyankar if (format == PETSC_VIEWER_NATIVE) { 975f4403165SShri Abhyankar /* store the matrix as a dense matrix */ 976f4403165SShri Abhyankar ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 977f4403165SShri Abhyankar col_lens[0] = MAT_FILE_CLASSID; 978f4403165SShri Abhyankar col_lens[1] = m; 979f4403165SShri Abhyankar col_lens[2] = n; 980f4403165SShri Abhyankar col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 981f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 982f4403165SShri Abhyankar ierr = PetscFree(col_lens);CHKERRQ(ierr); 983f4403165SShri Abhyankar 984f4403165SShri Abhyankar /* write out matrix, by rows */ 985f4403165SShri Abhyankar ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 986f4403165SShri Abhyankar v = a->v; 987f4403165SShri Abhyankar for (j=0; j<n; j++) { 988f4403165SShri Abhyankar for (i=0; i<m; i++) { 989f4403165SShri Abhyankar vals[j + i*n] = *v++; 990f4403165SShri Abhyankar } 991f4403165SShri Abhyankar } 992f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 993f4403165SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 994f4403165SShri Abhyankar } else { 99513f74950SBarry Smith ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 9960700a824SBarry Smith col_lens[0] = MAT_FILE_CLASSID; 997932b0c3eSLois Curfman McInnes col_lens[1] = m; 998932b0c3eSLois Curfman McInnes col_lens[2] = n; 999932b0c3eSLois Curfman McInnes col_lens[3] = nz; 1000932b0c3eSLois Curfman McInnes 1001932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 1002932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 10036f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 1004932b0c3eSLois Curfman McInnes 1005932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 1006932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 1007932b0c3eSLois Curfman McInnes ict = 0; 1008932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 1009932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 1010932b0c3eSLois Curfman McInnes } 10116f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 1012606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 1013932b0c3eSLois Curfman McInnes 1014932b0c3eSLois Curfman McInnes /* store nonzero values */ 101587828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 1016932b0c3eSLois Curfman McInnes ict = 0; 1017932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 1018932b0c3eSLois Curfman McInnes v = a->v + i; 1019932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 10201b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 1021932b0c3eSLois Curfman McInnes } 1022932b0c3eSLois Curfman McInnes } 10236f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 1024606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 1025f4403165SShri Abhyankar } 10263a40ed3dSBarry Smith PetscFunctionReturn(0); 1027932b0c3eSLois Curfman McInnes } 1028932b0c3eSLois Curfman McInnes 10294a2ae208SSatish Balay #undef __FUNCT__ 10304a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 1031dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 1032f1af5d2fSBarry Smith { 1033f1af5d2fSBarry Smith Mat A = (Mat) Aa; 1034f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 10356849ba73SBarry Smith PetscErrorCode ierr; 1036d0f46423SBarry Smith PetscInt m = A->rmap->n,n = A->cmap->n,color,i,j; 103787828ca2SBarry Smith PetscScalar *v = a->v; 1038b0a32e0cSBarry Smith PetscViewer viewer; 1039b0a32e0cSBarry Smith PetscDraw popup; 1040329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 1041f3ef73ceSBarry Smith PetscViewerFormat format; 1042f1af5d2fSBarry Smith 1043f1af5d2fSBarry Smith PetscFunctionBegin; 1044f1af5d2fSBarry Smith 1045f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 1046b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1047b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1048f1af5d2fSBarry Smith 1049f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 1050fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 1051f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 1052b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1053f1af5d2fSBarry Smith for (j = 0; j < n; j++) { 1054f1af5d2fSBarry Smith x_l = j; 1055f1af5d2fSBarry Smith x_r = x_l + 1.0; 1056f1af5d2fSBarry Smith for (i = 0; i < m; i++) { 1057f1af5d2fSBarry Smith y_l = m - i - 1.0; 1058f1af5d2fSBarry Smith y_r = y_l + 1.0; 1059f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 1060329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 1061b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1062329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 1063b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1064f1af5d2fSBarry Smith } else { 1065f1af5d2fSBarry Smith continue; 1066f1af5d2fSBarry Smith } 1067f1af5d2fSBarry Smith #else 1068f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 1069b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1070f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 1071b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1072f1af5d2fSBarry Smith } else { 1073f1af5d2fSBarry Smith continue; 1074f1af5d2fSBarry Smith } 1075f1af5d2fSBarry Smith #endif 1076b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1077f1af5d2fSBarry Smith } 1078f1af5d2fSBarry Smith } 1079f1af5d2fSBarry Smith } else { 1080f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 1081f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 1082f1af5d2fSBarry Smith for (i = 0; i < m*n; i++) { 1083f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 1084f1af5d2fSBarry Smith } 1085b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 1086b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 1087b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 1088f1af5d2fSBarry Smith for (j = 0; j < n; j++) { 1089f1af5d2fSBarry Smith x_l = j; 1090f1af5d2fSBarry Smith x_r = x_l + 1.0; 1091f1af5d2fSBarry Smith for (i = 0; i < m; i++) { 1092f1af5d2fSBarry Smith y_l = m - i - 1.0; 1093f1af5d2fSBarry Smith y_r = y_l + 1.0; 1094b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 1095b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1096f1af5d2fSBarry Smith } 1097f1af5d2fSBarry Smith } 1098f1af5d2fSBarry Smith } 1099f1af5d2fSBarry Smith PetscFunctionReturn(0); 1100f1af5d2fSBarry Smith } 1101f1af5d2fSBarry Smith 11024a2ae208SSatish Balay #undef __FUNCT__ 11034a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 1104dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 1105f1af5d2fSBarry Smith { 1106b0a32e0cSBarry Smith PetscDraw draw; 1107ace3abfcSBarry Smith PetscBool isnull; 1108329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 1109dfbe8321SBarry Smith PetscErrorCode ierr; 1110f1af5d2fSBarry Smith 1111f1af5d2fSBarry Smith PetscFunctionBegin; 1112b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1113b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1114abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 1115f1af5d2fSBarry Smith 1116f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1117d0f46423SBarry Smith xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0; 1118f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 1119b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1120b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 1121f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 1122f1af5d2fSBarry Smith PetscFunctionReturn(0); 1123f1af5d2fSBarry Smith } 1124f1af5d2fSBarry Smith 11254a2ae208SSatish Balay #undef __FUNCT__ 11264a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 1127dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1128932b0c3eSLois Curfman McInnes { 1129dfbe8321SBarry Smith PetscErrorCode ierr; 1130ace3abfcSBarry Smith PetscBool iascii,isbinary,isdraw; 1131932b0c3eSLois Curfman McInnes 11323a40ed3dSBarry Smith PetscFunctionBegin; 1133251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1134251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 1135251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 11360f5bd95cSBarry Smith 1137c45a1595SBarry Smith if (iascii) { 1138c45a1595SBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 11390f5bd95cSBarry Smith } else if (isbinary) { 11403a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1141f1af5d2fSBarry Smith } else if (isdraw) { 1142f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 11435cd90555SBarry Smith } else { 1144e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1145932b0c3eSLois Curfman McInnes } 11463a40ed3dSBarry Smith PetscFunctionReturn(0); 1147932b0c3eSLois Curfman McInnes } 1148289bc588SBarry Smith 11494a2ae208SSatish Balay #undef __FUNCT__ 11504a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1151dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1152289bc588SBarry Smith { 1153ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1154dfbe8321SBarry Smith PetscErrorCode ierr; 115590f02eecSBarry Smith 11563a40ed3dSBarry Smith PetscFunctionBegin; 1157aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1158d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n); 1159a5a9c739SBarry Smith #endif 116005b42c5fSBarry Smith ierr = PetscFree(l->pivots);CHKERRQ(ierr); 11616857c123SSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1162bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 1163dbd8c25aSHong Zhang 1164dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 11658c778c55SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatDenseGetArray_C","",PETSC_NULL);CHKERRQ(ierr); 11668c778c55SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatDenseRestoreArray_C","",PETSC_NULL);CHKERRQ(ierr); 1167901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 11684ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11694ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11704ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11713a40ed3dSBarry Smith PetscFunctionReturn(0); 1172289bc588SBarry Smith } 1173289bc588SBarry Smith 11744a2ae208SSatish Balay #undef __FUNCT__ 11754a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1176fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout) 1177289bc588SBarry Smith { 1178c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11796849ba73SBarry Smith PetscErrorCode ierr; 118013f74950SBarry Smith PetscInt k,j,m,n,M; 118187828ca2SBarry Smith PetscScalar *v,tmp; 118248b35521SBarry Smith 11833a40ed3dSBarry Smith PetscFunctionBegin; 1184d0f46423SBarry Smith v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n; 1185e9695a30SBarry Smith if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */ 1186e7e72b3dSBarry Smith if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1187e7e72b3dSBarry Smith else { 1188d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1189289bc588SBarry Smith for (k=0; k<j; k++) { 11901b807ce4Svictorle tmp = v[j + k*M]; 11911b807ce4Svictorle v[j + k*M] = v[k + j*M]; 11921b807ce4Svictorle v[k + j*M] = tmp; 1193289bc588SBarry Smith } 1194289bc588SBarry Smith } 1195d64ed03dSBarry Smith } 11963a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1197d3e5ee88SLois Curfman McInnes Mat tmat; 1198ec8511deSBarry Smith Mat_SeqDense *tmatd; 119987828ca2SBarry Smith PetscScalar *v2; 1200ea709b57SSatish Balay 1201fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 12027adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr); 1203d0f46423SBarry Smith ierr = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr); 12047adad957SLisandro Dalcin ierr = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 12055c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1206fc4dec0aSBarry Smith } else { 1207fc4dec0aSBarry Smith tmat = *matout; 1208fc4dec0aSBarry Smith } 1209ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 12100de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1211d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 12121b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1213d3e5ee88SLois Curfman McInnes } 12146d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12156d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1216d3e5ee88SLois Curfman McInnes *matout = tmat; 121748b35521SBarry Smith } 12183a40ed3dSBarry Smith PetscFunctionReturn(0); 1219289bc588SBarry Smith } 1220289bc588SBarry Smith 12214a2ae208SSatish Balay #undef __FUNCT__ 12224a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1223ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool *flg) 1224289bc588SBarry Smith { 1225c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1226c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 122713f74950SBarry Smith PetscInt i,j; 1228*a2ea699eSBarry Smith PetscScalar *v1,*v2; 12299ea5d5aeSSatish Balay 12303a40ed3dSBarry Smith PetscFunctionBegin; 1231d0f46423SBarry Smith if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1232d0f46423SBarry Smith if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1233d0f46423SBarry Smith for (i=0; i<A1->rmap->n; i++) { 12341b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1235d0f46423SBarry Smith for (j=0; j<A1->cmap->n; j++) { 12363a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 12371b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 12381b807ce4Svictorle } 1239289bc588SBarry Smith } 124077c4ece6SBarry Smith *flg = PETSC_TRUE; 12413a40ed3dSBarry Smith PetscFunctionReturn(0); 1242289bc588SBarry Smith } 1243289bc588SBarry Smith 12444a2ae208SSatish Balay #undef __FUNCT__ 12454a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1246dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1247289bc588SBarry Smith { 1248c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1249dfbe8321SBarry Smith PetscErrorCode ierr; 125013f74950SBarry Smith PetscInt i,n,len; 125187828ca2SBarry Smith PetscScalar *x,zero = 0.0; 125244cd7ae7SLois Curfman McInnes 12533a40ed3dSBarry Smith PetscFunctionBegin; 12542dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 12557a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 12561ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1257d0f46423SBarry Smith len = PetscMin(A->rmap->n,A->cmap->n); 1258e32f2f54SBarry Smith if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 125944cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 12601b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1261289bc588SBarry Smith } 12621ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 12633a40ed3dSBarry Smith PetscFunctionReturn(0); 1264289bc588SBarry Smith } 1265289bc588SBarry Smith 12664a2ae208SSatish Balay #undef __FUNCT__ 12674a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1268dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1269289bc588SBarry Smith { 1270c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 127187828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1272dfbe8321SBarry Smith PetscErrorCode ierr; 1273d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n; 127455659b69SBarry Smith 12753a40ed3dSBarry Smith PetscFunctionBegin; 127628988994SBarry Smith if (ll) { 12777a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 12781ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1279e32f2f54SBarry Smith if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1280da3a660dSBarry Smith for (i=0; i<m; i++) { 1281da3a660dSBarry Smith x = l[i]; 1282da3a660dSBarry Smith v = mat->v + i; 1283da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1284da3a660dSBarry Smith } 12851ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1286efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1287da3a660dSBarry Smith } 128828988994SBarry Smith if (rr) { 12897a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 12901ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1291e32f2f54SBarry Smith if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1292da3a660dSBarry Smith for (i=0; i<n; i++) { 1293da3a660dSBarry Smith x = r[i]; 1294da3a660dSBarry Smith v = mat->v + i*m; 1295da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1296da3a660dSBarry Smith } 12971ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1298efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1299da3a660dSBarry Smith } 13003a40ed3dSBarry Smith PetscFunctionReturn(0); 1301289bc588SBarry Smith } 1302289bc588SBarry Smith 13034a2ae208SSatish Balay #undef __FUNCT__ 13044a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1305dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1306289bc588SBarry Smith { 1307c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 130887828ca2SBarry Smith PetscScalar *v = mat->v; 1309329f5518SBarry Smith PetscReal sum = 0.0; 1310d0f46423SBarry Smith PetscInt lda=mat->lda,m=A->rmap->n,i,j; 1311efee365bSSatish Balay PetscErrorCode ierr; 131255659b69SBarry Smith 13133a40ed3dSBarry Smith PetscFunctionBegin; 1314289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1315a5ce6ee0Svictorle if (lda>m) { 1316d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1317a5ce6ee0Svictorle v = mat->v+j*lda; 1318a5ce6ee0Svictorle for (i=0; i<m; i++) { 1319a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1320a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1321a5ce6ee0Svictorle #else 1322a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1323a5ce6ee0Svictorle #endif 1324a5ce6ee0Svictorle } 1325a5ce6ee0Svictorle } 1326a5ce6ee0Svictorle } else { 1327d0f46423SBarry Smith for (i=0; i<A->cmap->n*A->rmap->n; i++) { 1328aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1329329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1330289bc588SBarry Smith #else 1331289bc588SBarry Smith sum += (*v)*(*v); v++; 1332289bc588SBarry Smith #endif 1333289bc588SBarry Smith } 1334a5ce6ee0Svictorle } 13358f1a2a5eSBarry Smith *nrm = PetscSqrtReal(sum); 1336dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13373a40ed3dSBarry Smith } else if (type == NORM_1) { 1338064f8208SBarry Smith *nrm = 0.0; 1339d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 13401b807ce4Svictorle v = mat->v + j*mat->lda; 1341289bc588SBarry Smith sum = 0.0; 1342d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 134333a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1344289bc588SBarry Smith } 1345064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1346289bc588SBarry Smith } 1347d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13483a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1349064f8208SBarry Smith *nrm = 0.0; 1350d0f46423SBarry Smith for (j=0; j<A->rmap->n; j++) { 1351289bc588SBarry Smith v = mat->v + j; 1352289bc588SBarry Smith sum = 0.0; 1353d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) { 13541b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1355289bc588SBarry Smith } 1356064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1357289bc588SBarry Smith } 1358d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 1359e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 13603a40ed3dSBarry Smith PetscFunctionReturn(0); 1361289bc588SBarry Smith } 1362289bc588SBarry Smith 13634a2ae208SSatish Balay #undef __FUNCT__ 13644a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 1365ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg) 1366289bc588SBarry Smith { 1367c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 136863ba0a88SBarry Smith PetscErrorCode ierr; 136967e560aaSBarry Smith 13703a40ed3dSBarry Smith PetscFunctionBegin; 1371b5a2b587SKris Buschelman switch (op) { 1372b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 13734e0d8c25SBarry Smith aij->roworiented = flg; 1374b5a2b587SKris Buschelman break; 1375512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 1376b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 13773971808eSMatthew Knepley case MAT_NEW_NONZERO_ALLOCATION_ERR: 13784e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 137913fa8e87SLisandro Dalcin case MAT_KEEP_NONZERO_PATTERN: 1380b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1381b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 13825021d80fSJed Brown case MAT_IGNORE_LOWER_TRIANGULAR: 13835021d80fSJed Brown ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 13845021d80fSJed Brown break; 13855021d80fSJed Brown case MAT_SPD: 138677e54ba9SKris Buschelman case MAT_SYMMETRIC: 138777e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 13889a4540c5SBarry Smith case MAT_HERMITIAN: 13899a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 13905021d80fSJed Brown /* These options are handled directly by MatSetOption() */ 139177e54ba9SKris Buschelman break; 1392b5a2b587SKris Buschelman default: 1393e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 13943a40ed3dSBarry Smith } 13953a40ed3dSBarry Smith PetscFunctionReturn(0); 1396289bc588SBarry Smith } 1397289bc588SBarry Smith 13984a2ae208SSatish Balay #undef __FUNCT__ 13994a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1400dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 14016f0a148fSBarry Smith { 1402ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 14036849ba73SBarry Smith PetscErrorCode ierr; 1404d0f46423SBarry Smith PetscInt lda=l->lda,m=A->rmap->n,j; 14053a40ed3dSBarry Smith 14063a40ed3dSBarry Smith PetscFunctionBegin; 1407a5ce6ee0Svictorle if (lda>m) { 1408d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1409a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1410a5ce6ee0Svictorle } 1411a5ce6ee0Svictorle } else { 1412d0f46423SBarry Smith ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1413a5ce6ee0Svictorle } 14143a40ed3dSBarry Smith PetscFunctionReturn(0); 14156f0a148fSBarry Smith } 14166f0a148fSBarry Smith 14174a2ae208SSatish Balay #undef __FUNCT__ 14184a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 14192b40b63fSBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 14206f0a148fSBarry Smith { 142197b48c8fSBarry Smith PetscErrorCode ierr; 1422ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1423b9679d65SBarry Smith PetscInt m = l->lda, n = A->cmap->n, i,j; 142497b48c8fSBarry Smith PetscScalar *slot,*bb; 142597b48c8fSBarry Smith const PetscScalar *xx; 142655659b69SBarry Smith 14273a40ed3dSBarry Smith PetscFunctionBegin; 1428b9679d65SBarry Smith #if defined(PETSC_USE_DEBUG) 1429b9679d65SBarry Smith for (i=0; i<N; i++) { 1430b9679d65SBarry Smith if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed"); 1431b9679d65SBarry Smith if (rows[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested to be zeroed greater than or equal number of rows %D",rows[i],A->rmap->n); 1432b9679d65SBarry Smith } 1433b9679d65SBarry Smith #endif 1434b9679d65SBarry Smith 143597b48c8fSBarry Smith /* fix right hand side if needed */ 143697b48c8fSBarry Smith if (x && b) { 143797b48c8fSBarry Smith ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 143897b48c8fSBarry Smith ierr = VecGetArray(b,&bb);CHKERRQ(ierr); 143997b48c8fSBarry Smith for (i=0; i<N; i++) { 144097b48c8fSBarry Smith bb[rows[i]] = diag*xx[rows[i]]; 144197b48c8fSBarry Smith } 144297b48c8fSBarry Smith ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 144397b48c8fSBarry Smith ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr); 144497b48c8fSBarry Smith } 144597b48c8fSBarry Smith 14466f0a148fSBarry Smith for (i=0; i<N; i++) { 14476f0a148fSBarry Smith slot = l->v + rows[i]; 1448b9679d65SBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += m;} 14496f0a148fSBarry Smith } 1450f4df32b1SMatthew Knepley if (diag != 0.0) { 1451b9679d65SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices"); 14526f0a148fSBarry Smith for (i=0; i<N; i++) { 1453b9679d65SBarry Smith slot = l->v + (m+1)*rows[i]; 1454f4df32b1SMatthew Knepley *slot = diag; 14556f0a148fSBarry Smith } 14566f0a148fSBarry Smith } 14573a40ed3dSBarry Smith PetscFunctionReturn(0); 14586f0a148fSBarry Smith } 1459557bce09SLois Curfman McInnes 14604a2ae208SSatish Balay #undef __FUNCT__ 14618c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray_SeqDense" 14628c778c55SBarry Smith PetscErrorCode MatDenseGetArray_SeqDense(Mat A,PetscScalar *array[]) 146364e87e97SBarry Smith { 1464c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14653a40ed3dSBarry Smith 14663a40ed3dSBarry Smith PetscFunctionBegin; 1467e32f2f54SBarry 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"); 146864e87e97SBarry Smith *array = mat->v; 14693a40ed3dSBarry Smith PetscFunctionReturn(0); 147064e87e97SBarry Smith } 14710754003eSLois Curfman McInnes 14724a2ae208SSatish Balay #undef __FUNCT__ 14738c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray_SeqDense" 14748c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1475ff14e315SSatish Balay { 14763a40ed3dSBarry Smith PetscFunctionBegin; 147709b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 14783a40ed3dSBarry Smith PetscFunctionReturn(0); 1479ff14e315SSatish Balay } 14800754003eSLois Curfman McInnes 14814a2ae208SSatish Balay #undef __FUNCT__ 14828c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray" 1483dec5eb66SMatthew G Knepley /*@C 14848c778c55SBarry Smith MatDenseGetArray - gives access to the array where the data for a SeqDense matrix is stored 148573a71a0fSBarry Smith 148673a71a0fSBarry Smith Not Collective 148773a71a0fSBarry Smith 148873a71a0fSBarry Smith Input Parameter: 148973a71a0fSBarry Smith . mat - a MATSEQDENSE matrix 149073a71a0fSBarry Smith 149173a71a0fSBarry Smith Output Parameter: 149273a71a0fSBarry Smith . array - pointer to the data 149373a71a0fSBarry Smith 149473a71a0fSBarry Smith Level: intermediate 149573a71a0fSBarry Smith 14968c778c55SBarry Smith .seealso: MatDenseRestoreArray() 149773a71a0fSBarry Smith @*/ 14988c778c55SBarry Smith PetscErrorCode MatDenseGetArray(Mat A,PetscScalar **array) 149973a71a0fSBarry Smith { 150073a71a0fSBarry Smith PetscErrorCode ierr; 150173a71a0fSBarry Smith 150273a71a0fSBarry Smith PetscFunctionBegin; 15038c778c55SBarry Smith ierr = PetscUseMethod(A,"MatDenseGetArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr); 150473a71a0fSBarry Smith PetscFunctionReturn(0); 150573a71a0fSBarry Smith } 150673a71a0fSBarry Smith 150773a71a0fSBarry Smith #undef __FUNCT__ 15088c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray" 1509dec5eb66SMatthew G Knepley /*@C 15108c778c55SBarry Smith MatDenseRestoreArray - returns access to the array where the data for a SeqDense matrix is stored obtained by MatDenseGetArray() 151173a71a0fSBarry Smith 151273a71a0fSBarry Smith Not Collective 151373a71a0fSBarry Smith 151473a71a0fSBarry Smith Input Parameters: 151573a71a0fSBarry Smith . mat - a MATSEQDENSE matrix 151673a71a0fSBarry Smith . array - pointer to the data 151773a71a0fSBarry Smith 151873a71a0fSBarry Smith Level: intermediate 151973a71a0fSBarry Smith 15208c778c55SBarry Smith .seealso: MatDenseGetArray() 152173a71a0fSBarry Smith @*/ 15228c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray(Mat A,PetscScalar **array) 152373a71a0fSBarry Smith { 152473a71a0fSBarry Smith PetscErrorCode ierr; 152573a71a0fSBarry Smith 152673a71a0fSBarry Smith PetscFunctionBegin; 15278c778c55SBarry Smith ierr = PetscUseMethod(A,"MatDenseRestoreArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr); 152873a71a0fSBarry Smith PetscFunctionReturn(0); 152973a71a0fSBarry Smith } 153073a71a0fSBarry Smith 153173a71a0fSBarry Smith #undef __FUNCT__ 15324a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 153313f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 15340754003eSLois Curfman McInnes { 1535c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 15366849ba73SBarry Smith PetscErrorCode ierr; 15375d0c19d7SBarry Smith PetscInt i,j,nrows,ncols; 15385d0c19d7SBarry Smith const PetscInt *irow,*icol; 153987828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 15400754003eSLois Curfman McInnes Mat newmat; 15410754003eSLois Curfman McInnes 15423a40ed3dSBarry Smith PetscFunctionBegin; 154378b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 154478b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1545e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1546e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 15470754003eSLois Curfman McInnes 1548182d2002SSatish Balay /* Check submatrixcall */ 1549182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 155013f74950SBarry Smith PetscInt n_cols,n_rows; 1551182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 155221a2c019SBarry Smith if (n_rows != nrows || n_cols != ncols) { 1553f746d493SDmitry Karpeev /* resize the result matrix to match number of requested rows/columns */ 1554c61587bbSBarry Smith ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 155521a2c019SBarry Smith } 1556182d2002SSatish Balay newmat = *B; 1557182d2002SSatish Balay } else { 15580754003eSLois Curfman McInnes /* Create and fill new matrix */ 15597adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr); 1560f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 15617adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 15625c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1563182d2002SSatish Balay } 1564182d2002SSatish Balay 1565182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1566182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1567182d2002SSatish Balay 1568182d2002SSatish Balay for (i=0; i<ncols; i++) { 15696de62eeeSBarry Smith av = v + mat->lda*icol[i]; 1570182d2002SSatish Balay for (j=0; j<nrows; j++) { 1571182d2002SSatish Balay *bv++ = av[irow[j]]; 15720754003eSLois Curfman McInnes } 15730754003eSLois Curfman McInnes } 1574182d2002SSatish Balay 1575182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 15766d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15776d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15780754003eSLois Curfman McInnes 15790754003eSLois Curfman McInnes /* Free work space */ 158078b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 158178b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1582182d2002SSatish Balay *B = newmat; 15833a40ed3dSBarry Smith PetscFunctionReturn(0); 15840754003eSLois Curfman McInnes } 15850754003eSLois Curfman McInnes 15864a2ae208SSatish Balay #undef __FUNCT__ 15874a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 158813f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1589905e6a2fSBarry Smith { 15906849ba73SBarry Smith PetscErrorCode ierr; 159113f74950SBarry Smith PetscInt i; 1592905e6a2fSBarry Smith 15933a40ed3dSBarry Smith PetscFunctionBegin; 1594905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1595b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1596905e6a2fSBarry Smith } 1597905e6a2fSBarry Smith 1598905e6a2fSBarry Smith for (i=0; i<n; i++) { 15996a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1600905e6a2fSBarry Smith } 16013a40ed3dSBarry Smith PetscFunctionReturn(0); 1602905e6a2fSBarry Smith } 1603905e6a2fSBarry Smith 16044a2ae208SSatish Balay #undef __FUNCT__ 1605c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1606c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1607c0aa2d19SHong Zhang { 1608c0aa2d19SHong Zhang PetscFunctionBegin; 1609c0aa2d19SHong Zhang PetscFunctionReturn(0); 1610c0aa2d19SHong Zhang } 1611c0aa2d19SHong Zhang 1612c0aa2d19SHong Zhang #undef __FUNCT__ 1613c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1614c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1615c0aa2d19SHong Zhang { 1616c0aa2d19SHong Zhang PetscFunctionBegin; 1617c0aa2d19SHong Zhang PetscFunctionReturn(0); 1618c0aa2d19SHong Zhang } 1619c0aa2d19SHong Zhang 1620c0aa2d19SHong Zhang #undef __FUNCT__ 16214a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1622dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 16234b0e389bSBarry Smith { 16244b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 16256849ba73SBarry Smith PetscErrorCode ierr; 1626d0f46423SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j; 16273a40ed3dSBarry Smith 16283a40ed3dSBarry Smith PetscFunctionBegin; 162933f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 163033f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1631cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 16323a40ed3dSBarry Smith PetscFunctionReturn(0); 16333a40ed3dSBarry Smith } 1634e32f2f54SBarry Smith if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1635a5ce6ee0Svictorle if (lda1>m || lda2>m) { 16360dbb7854Svictorle for (j=0; j<n; j++) { 1637a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1638a5ce6ee0Svictorle } 1639a5ce6ee0Svictorle } else { 1640d0f46423SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1641a5ce6ee0Svictorle } 1642273d9f13SBarry Smith PetscFunctionReturn(0); 1643273d9f13SBarry Smith } 1644273d9f13SBarry Smith 16454a2ae208SSatish Balay #undef __FUNCT__ 16464994cf47SJed Brown #define __FUNCT__ "MatSetUp_SeqDense" 16474994cf47SJed Brown PetscErrorCode MatSetUp_SeqDense(Mat A) 1648273d9f13SBarry Smith { 1649dfbe8321SBarry Smith PetscErrorCode ierr; 1650273d9f13SBarry Smith 1651273d9f13SBarry Smith PetscFunctionBegin; 1652273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 16533a40ed3dSBarry Smith PetscFunctionReturn(0); 16544b0e389bSBarry Smith } 16554b0e389bSBarry Smith 1656284134d9SBarry Smith #undef __FUNCT__ 1657ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense" 1658ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A) 1659ba337c44SJed Brown { 1660ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1661ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1662ba337c44SJed Brown PetscScalar *aa = a->v; 1663ba337c44SJed Brown 1664ba337c44SJed Brown PetscFunctionBegin; 1665ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]); 1666ba337c44SJed Brown PetscFunctionReturn(0); 1667ba337c44SJed Brown } 1668ba337c44SJed Brown 1669ba337c44SJed Brown #undef __FUNCT__ 1670ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense" 1671ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A) 1672ba337c44SJed Brown { 1673ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1674ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1675ba337c44SJed Brown PetscScalar *aa = a->v; 1676ba337c44SJed Brown 1677ba337c44SJed Brown PetscFunctionBegin; 1678ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 1679ba337c44SJed Brown PetscFunctionReturn(0); 1680ba337c44SJed Brown } 1681ba337c44SJed Brown 1682ba337c44SJed Brown #undef __FUNCT__ 1683ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense" 1684ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A) 1685ba337c44SJed Brown { 1686ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1687ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1688ba337c44SJed Brown PetscScalar *aa = a->v; 1689ba337c44SJed Brown 1690ba337c44SJed Brown PetscFunctionBegin; 1691ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 1692ba337c44SJed Brown PetscFunctionReturn(0); 1693ba337c44SJed Brown } 1694284134d9SBarry Smith 1695a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/ 1696a9fe9ddaSSatish Balay #undef __FUNCT__ 1697a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1698a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1699a9fe9ddaSSatish Balay { 1700a9fe9ddaSSatish Balay PetscErrorCode ierr; 1701a9fe9ddaSSatish Balay 1702a9fe9ddaSSatish Balay PetscFunctionBegin; 1703a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1704a9fe9ddaSSatish Balay ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1705a9fe9ddaSSatish Balay } 1706a9fe9ddaSSatish Balay ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1707a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1708a9fe9ddaSSatish Balay } 1709a9fe9ddaSSatish Balay 1710a9fe9ddaSSatish Balay #undef __FUNCT__ 1711a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1712a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1713a9fe9ddaSSatish Balay { 1714ee16a9a1SHong Zhang PetscErrorCode ierr; 1715d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 1716ee16a9a1SHong Zhang Mat Cmat; 1717a9fe9ddaSSatish Balay 1718ee16a9a1SHong Zhang PetscFunctionBegin; 1719e32f2f54SBarry 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); 1720ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1721ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1722ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1723ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1724d73949e8SHong Zhang 1725ee16a9a1SHong Zhang *C = Cmat; 1726ee16a9a1SHong Zhang PetscFunctionReturn(0); 1727ee16a9a1SHong Zhang } 1728a9fe9ddaSSatish Balay 172998a3b096SSatish Balay #undef __FUNCT__ 1730a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1731a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1732a9fe9ddaSSatish Balay { 1733a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1734a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1735a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 17360805154bSBarry Smith PetscBLASInt m,n,k; 1737a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1738a9fe9ddaSSatish Balay 1739a9fe9ddaSSatish Balay PetscFunctionBegin; 1740d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 1741d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1742d0f46423SBarry Smith k = PetscBLASIntCast(A->cmap->n); 1743a9fe9ddaSSatish Balay BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1744a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1745a9fe9ddaSSatish Balay } 1746a9fe9ddaSSatish Balay 1747a9fe9ddaSSatish Balay #undef __FUNCT__ 174875648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqDense_SeqDense" 174975648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1750a9fe9ddaSSatish Balay { 1751a9fe9ddaSSatish Balay PetscErrorCode ierr; 1752a9fe9ddaSSatish Balay 1753a9fe9ddaSSatish Balay PetscFunctionBegin; 1754a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 175575648e8dSHong Zhang ierr = MatTransposeMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1756a9fe9ddaSSatish Balay } 175775648e8dSHong Zhang ierr = MatTransposeMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1758a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1759a9fe9ddaSSatish Balay } 1760a9fe9ddaSSatish Balay 1761a9fe9ddaSSatish Balay #undef __FUNCT__ 176275648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqDense_SeqDense" 176375648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1764a9fe9ddaSSatish Balay { 1765ee16a9a1SHong Zhang PetscErrorCode ierr; 1766d0f46423SBarry Smith PetscInt m=A->cmap->n,n=B->cmap->n; 1767ee16a9a1SHong Zhang Mat Cmat; 1768a9fe9ddaSSatish Balay 1769ee16a9a1SHong Zhang PetscFunctionBegin; 1770e32f2f54SBarry 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); 1771ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1772ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1773ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1774ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1775ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1776ee16a9a1SHong Zhang *C = Cmat; 1777ee16a9a1SHong Zhang PetscFunctionReturn(0); 1778ee16a9a1SHong Zhang } 1779a9fe9ddaSSatish Balay 1780a9fe9ddaSSatish Balay #undef __FUNCT__ 178175648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqDense_SeqDense" 178275648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1783a9fe9ddaSSatish Balay { 1784a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1785a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1786a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 17870805154bSBarry Smith PetscBLASInt m,n,k; 1788a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1789a9fe9ddaSSatish Balay 1790a9fe9ddaSSatish Balay PetscFunctionBegin; 1791d0f46423SBarry Smith m = PetscBLASIntCast(A->cmap->n); 1792d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1793d0f46423SBarry Smith k = PetscBLASIntCast(A->rmap->n); 17942fbe02b9SBarry Smith /* 17952fbe02b9SBarry Smith Note the m and n arguments below are the number rows and columns of A', not A! 17962fbe02b9SBarry Smith */ 1797a9fe9ddaSSatish Balay BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1798a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1799a9fe9ddaSSatish Balay } 1800985db425SBarry Smith 1801985db425SBarry Smith #undef __FUNCT__ 1802985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense" 1803985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 1804985db425SBarry Smith { 1805985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1806985db425SBarry Smith PetscErrorCode ierr; 1807d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1808985db425SBarry Smith PetscScalar *x; 1809985db425SBarry Smith MatScalar *aa = a->v; 1810985db425SBarry Smith 1811985db425SBarry Smith PetscFunctionBegin; 1812e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1813985db425SBarry Smith 1814985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1815985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1816985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1817e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1818985db425SBarry Smith for (i=0; i<m; i++) { 1819985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1820985db425SBarry Smith for (j=1; j<n; j++){ 1821985db425SBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1822985db425SBarry Smith } 1823985db425SBarry Smith } 1824985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1825985db425SBarry Smith PetscFunctionReturn(0); 1826985db425SBarry Smith } 1827985db425SBarry Smith 1828985db425SBarry Smith #undef __FUNCT__ 1829985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense" 1830985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 1831985db425SBarry Smith { 1832985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1833985db425SBarry Smith PetscErrorCode ierr; 1834d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1835985db425SBarry Smith PetscScalar *x; 1836985db425SBarry Smith PetscReal atmp; 1837985db425SBarry Smith MatScalar *aa = a->v; 1838985db425SBarry Smith 1839985db425SBarry Smith PetscFunctionBegin; 1840e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1841985db425SBarry Smith 1842985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1843985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1844985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1845e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1846985db425SBarry Smith for (i=0; i<m; i++) { 18479189402eSHong Zhang x[i] = PetscAbsScalar(aa[i]); 1848985db425SBarry Smith for (j=1; j<n; j++){ 1849985db425SBarry Smith atmp = PetscAbsScalar(aa[i+m*j]); 1850985db425SBarry Smith if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 1851985db425SBarry Smith } 1852985db425SBarry Smith } 1853985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1854985db425SBarry Smith PetscFunctionReturn(0); 1855985db425SBarry Smith } 1856985db425SBarry Smith 1857985db425SBarry Smith #undef __FUNCT__ 1858985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense" 1859985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 1860985db425SBarry Smith { 1861985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1862985db425SBarry Smith PetscErrorCode ierr; 1863d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1864985db425SBarry Smith PetscScalar *x; 1865985db425SBarry Smith MatScalar *aa = a->v; 1866985db425SBarry Smith 1867985db425SBarry Smith PetscFunctionBegin; 1868e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1869985db425SBarry Smith 1870985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1871985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1872985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1873e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1874985db425SBarry Smith for (i=0; i<m; i++) { 1875985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1876985db425SBarry Smith for (j=1; j<n; j++){ 1877985db425SBarry Smith if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1878985db425SBarry Smith } 1879985db425SBarry Smith } 1880985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1881985db425SBarry Smith PetscFunctionReturn(0); 1882985db425SBarry Smith } 1883985db425SBarry Smith 18848d0534beSBarry Smith #undef __FUNCT__ 18858d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense" 18868d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 18878d0534beSBarry Smith { 18888d0534beSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 18898d0534beSBarry Smith PetscErrorCode ierr; 18908d0534beSBarry Smith PetscScalar *x; 18918d0534beSBarry Smith 18928d0534beSBarry Smith PetscFunctionBegin; 1893e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 18948d0534beSBarry Smith 18958d0534beSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1896d0f46423SBarry Smith ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 18978d0534beSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 18988d0534beSBarry Smith PetscFunctionReturn(0); 18998d0534beSBarry Smith } 19008d0534beSBarry Smith 19010716a85fSBarry Smith 19020716a85fSBarry Smith #undef __FUNCT__ 19030716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqDense" 19040716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms) 19050716a85fSBarry Smith { 19060716a85fSBarry Smith PetscErrorCode ierr; 19070716a85fSBarry Smith PetscInt i,j,m,n; 19080716a85fSBarry Smith PetscScalar *a; 19090716a85fSBarry Smith 19100716a85fSBarry Smith PetscFunctionBegin; 19110716a85fSBarry Smith ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 19120716a85fSBarry Smith ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr); 19138c778c55SBarry Smith ierr = MatDenseGetArray(A,&a);CHKERRQ(ierr); 19140716a85fSBarry Smith if (type == NORM_2) { 19150716a85fSBarry Smith for (i=0; i<n; i++ ){ 19160716a85fSBarry Smith for (j=0; j<m; j++) { 19170716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]*a[j]); 19180716a85fSBarry Smith } 19190716a85fSBarry Smith a += m; 19200716a85fSBarry Smith } 19210716a85fSBarry Smith } else if (type == NORM_1) { 19220716a85fSBarry Smith for (i=0; i<n; i++ ){ 19230716a85fSBarry Smith for (j=0; j<m; j++) { 19240716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]); 19250716a85fSBarry Smith } 19260716a85fSBarry Smith a += m; 19270716a85fSBarry Smith } 19280716a85fSBarry Smith } else if (type == NORM_INFINITY) { 19290716a85fSBarry Smith for (i=0; i<n; i++ ){ 19300716a85fSBarry Smith for (j=0; j<m; j++) { 19310716a85fSBarry Smith norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]); 19320716a85fSBarry Smith } 19330716a85fSBarry Smith a += m; 19340716a85fSBarry Smith } 19350716a85fSBarry Smith } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Unknown NormType"); 19368c778c55SBarry Smith ierr = MatDenseRestoreArray(A,&a);CHKERRQ(ierr); 19370716a85fSBarry Smith if (type == NORM_2) { 19388f1a2a5eSBarry Smith for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]); 19390716a85fSBarry Smith } 19400716a85fSBarry Smith PetscFunctionReturn(0); 19410716a85fSBarry Smith } 19420716a85fSBarry Smith 194373a71a0fSBarry Smith #undef __FUNCT__ 194473a71a0fSBarry Smith #define __FUNCT__ "MatSetRandom_SeqDense" 194573a71a0fSBarry Smith static PetscErrorCode MatSetRandom_SeqDense(Mat x,PetscRandom rctx) 194673a71a0fSBarry Smith { 194773a71a0fSBarry Smith PetscErrorCode ierr; 194873a71a0fSBarry Smith PetscScalar *a; 194973a71a0fSBarry Smith PetscInt m,n,i; 195073a71a0fSBarry Smith 195173a71a0fSBarry Smith PetscFunctionBegin; 195273a71a0fSBarry Smith ierr = MatGetSize(x,&m,&n);CHKERRQ(ierr); 19538c778c55SBarry Smith ierr = MatDenseGetArray(x,&a);CHKERRQ(ierr); 195473a71a0fSBarry Smith for (i=0; i<m*n; i++) { 195573a71a0fSBarry Smith ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr); 195673a71a0fSBarry Smith } 19578c778c55SBarry Smith ierr = MatDenseRestoreArray(x,&a);CHKERRQ(ierr); 195873a71a0fSBarry Smith PetscFunctionReturn(0); 195973a71a0fSBarry Smith } 196073a71a0fSBarry Smith 196173a71a0fSBarry Smith 1962289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1963a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1964905e6a2fSBarry Smith MatGetRow_SeqDense, 1965905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1966905e6a2fSBarry Smith MatMult_SeqDense, 196797304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 19687c922b88SBarry Smith MatMultTranspose_SeqDense, 19697c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1970db4efbfdSBarry Smith 0, 1971db4efbfdSBarry Smith 0, 1972db4efbfdSBarry Smith 0, 1973db4efbfdSBarry Smith /*10*/ 0, 1974905e6a2fSBarry Smith MatLUFactor_SeqDense, 1975905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 197641f059aeSBarry Smith MatSOR_SeqDense, 1977ec8511deSBarry Smith MatTranspose_SeqDense, 197897304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1979905e6a2fSBarry Smith MatEqual_SeqDense, 1980905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1981905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1982905e6a2fSBarry Smith MatNorm_SeqDense, 1983c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense, 1984c0aa2d19SHong Zhang MatAssemblyEnd_SeqDense, 1985905e6a2fSBarry Smith MatSetOption_SeqDense, 1986905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1987d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense, 1988db4efbfdSBarry Smith 0, 1989db4efbfdSBarry Smith 0, 1990db4efbfdSBarry Smith 0, 1991db4efbfdSBarry Smith 0, 19924994cf47SJed Brown /*29*/ MatSetUp_SeqDense, 1993273d9f13SBarry Smith 0, 1994905e6a2fSBarry Smith 0, 199573a71a0fSBarry Smith 0, 199673a71a0fSBarry Smith 0, 1997d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense, 1998a5ae1ecdSBarry Smith 0, 1999a5ae1ecdSBarry Smith 0, 2000a5ae1ecdSBarry Smith 0, 2001a5ae1ecdSBarry Smith 0, 2002d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense, 2003a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 2004a5ae1ecdSBarry Smith 0, 20054b0e389bSBarry Smith MatGetValues_SeqDense, 2006a5ae1ecdSBarry Smith MatCopy_SeqDense, 2007d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense, 2008a5ae1ecdSBarry Smith MatScale_SeqDense, 2009a5ae1ecdSBarry Smith 0, 2010a5ae1ecdSBarry Smith 0, 2011a5ae1ecdSBarry Smith 0, 201273a71a0fSBarry Smith /*49*/ MatSetRandom_SeqDense, 2013a5ae1ecdSBarry Smith 0, 2014a5ae1ecdSBarry Smith 0, 2015a5ae1ecdSBarry Smith 0, 2016a5ae1ecdSBarry Smith 0, 2017d519adbfSMatthew Knepley /*54*/ 0, 2018a5ae1ecdSBarry Smith 0, 2019a5ae1ecdSBarry Smith 0, 2020a5ae1ecdSBarry Smith 0, 2021a5ae1ecdSBarry Smith 0, 2022d519adbfSMatthew Knepley /*59*/ 0, 2023e03a110bSBarry Smith MatDestroy_SeqDense, 2024e03a110bSBarry Smith MatView_SeqDense, 2025357abbc8SBarry Smith 0, 202697304618SKris Buschelman 0, 2027d519adbfSMatthew Knepley /*64*/ 0, 202897304618SKris Buschelman 0, 202997304618SKris Buschelman 0, 203097304618SKris Buschelman 0, 203197304618SKris Buschelman 0, 2032d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense, 203397304618SKris Buschelman 0, 203497304618SKris Buschelman 0, 203597304618SKris Buschelman 0, 203697304618SKris Buschelman 0, 2037d519adbfSMatthew Knepley /*74*/ 0, 203897304618SKris Buschelman 0, 203997304618SKris Buschelman 0, 204097304618SKris Buschelman 0, 204197304618SKris Buschelman 0, 2042d519adbfSMatthew Knepley /*79*/ 0, 204397304618SKris Buschelman 0, 204497304618SKris Buschelman 0, 204597304618SKris Buschelman 0, 20465bba2384SShri Abhyankar /*83*/ MatLoad_SeqDense, 2047865e5f61SKris Buschelman 0, 20481cbb95d3SBarry Smith MatIsHermitian_SeqDense, 2049865e5f61SKris Buschelman 0, 2050865e5f61SKris Buschelman 0, 2051865e5f61SKris Buschelman 0, 2052d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense, 2053a9fe9ddaSSatish Balay MatMatMultSymbolic_SeqDense_SeqDense, 2054a9fe9ddaSSatish Balay MatMatMultNumeric_SeqDense_SeqDense, 2055865e5f61SKris Buschelman 0, 2056865e5f61SKris Buschelman 0, 2057d519adbfSMatthew Knepley /*94*/ 0, 20585df89d91SHong Zhang 0, 20595df89d91SHong Zhang 0, 20605df89d91SHong Zhang 0, 2061284134d9SBarry Smith 0, 2062d519adbfSMatthew Knepley /*99*/ 0, 2063284134d9SBarry Smith 0, 2064284134d9SBarry Smith 0, 2065ba337c44SJed Brown MatConjugate_SeqDense, 2066f73d5cc4SBarry Smith 0, 2067ba337c44SJed Brown /*104*/0, 2068ba337c44SJed Brown MatRealPart_SeqDense, 2069ba337c44SJed Brown MatImaginaryPart_SeqDense, 2070985db425SBarry Smith 0, 2071985db425SBarry Smith 0, 207285e2c93fSHong Zhang /*109*/MatMatSolve_SeqDense, 2073985db425SBarry Smith 0, 20748d0534beSBarry Smith MatGetRowMin_SeqDense, 2075aabbc4fbSShri Abhyankar MatGetColumnVector_SeqDense, 2076aabbc4fbSShri Abhyankar 0, 2077aabbc4fbSShri Abhyankar /*114*/0, 2078aabbc4fbSShri Abhyankar 0, 2079aabbc4fbSShri Abhyankar 0, 2080aabbc4fbSShri Abhyankar 0, 2081aabbc4fbSShri Abhyankar 0, 2082aabbc4fbSShri Abhyankar /*119*/0, 2083aabbc4fbSShri Abhyankar 0, 2084aabbc4fbSShri Abhyankar 0, 20850716a85fSBarry Smith 0, 20860716a85fSBarry Smith 0, 20870716a85fSBarry Smith /*124*/0, 20885df89d91SHong Zhang MatGetColumnNorms_SeqDense, 20895df89d91SHong Zhang 0, 20905df89d91SHong Zhang 0, 20915df89d91SHong Zhang 0, 20925df89d91SHong Zhang /*129*/0, 209375648e8dSHong Zhang MatTransposeMatMult_SeqDense_SeqDense, 209475648e8dSHong Zhang MatTransposeMatMultSymbolic_SeqDense_SeqDense, 209575648e8dSHong Zhang MatTransposeMatMultNumeric_SeqDense_SeqDense, 2096985db425SBarry Smith }; 209790ace30eSBarry Smith 20984a2ae208SSatish Balay #undef __FUNCT__ 20994a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 21004b828684SBarry Smith /*@C 2101fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 2102d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 2103d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 2104289bc588SBarry Smith 2105db81eaa0SLois Curfman McInnes Collective on MPI_Comm 2106db81eaa0SLois Curfman McInnes 210720563c6bSBarry Smith Input Parameters: 2108db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 21090c775827SLois Curfman McInnes . m - number of rows 211018f449edSLois Curfman McInnes . n - number of columns 2111c0235b3cSMatthew Knepley - data - optional location of matrix data in column major order. Set data=PETSC_NULL for PETSc 2112dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 211320563c6bSBarry Smith 211420563c6bSBarry Smith Output Parameter: 211544cd7ae7SLois Curfman McInnes . A - the matrix 211620563c6bSBarry Smith 2117b259b22eSLois Curfman McInnes Notes: 211818f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 211918f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 2120b4fd4287SBarry Smith set data=PETSC_NULL. 212118f449edSLois Curfman McInnes 2122027ccd11SLois Curfman McInnes Level: intermediate 2123027ccd11SLois Curfman McInnes 2124dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 2125d65003e9SLois Curfman McInnes 212669b1f4b7SBarry Smith .seealso: MatCreate(), MatCreateDense(), MatSetValues() 212720563c6bSBarry Smith @*/ 21287087cfbeSBarry Smith PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 2129289bc588SBarry Smith { 2130dfbe8321SBarry Smith PetscErrorCode ierr; 21313b2fbd54SBarry Smith 21323a40ed3dSBarry Smith PetscFunctionBegin; 2133f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 2134f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 2135273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 2136273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 2137273d9f13SBarry Smith PetscFunctionReturn(0); 2138273d9f13SBarry Smith } 2139273d9f13SBarry Smith 21404a2ae208SSatish Balay #undef __FUNCT__ 2141afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation" 2142273d9f13SBarry Smith /*@C 2143273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 2144273d9f13SBarry Smith 2145273d9f13SBarry Smith Collective on MPI_Comm 2146273d9f13SBarry Smith 2147273d9f13SBarry Smith Input Parameters: 2148273d9f13SBarry Smith + A - the matrix 2149273d9f13SBarry Smith - data - the array (or PETSC_NULL) 2150273d9f13SBarry Smith 2151273d9f13SBarry Smith Notes: 2152273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 2153273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 2154284134d9SBarry Smith need not call this routine. 2155273d9f13SBarry Smith 2156273d9f13SBarry Smith Level: intermediate 2157273d9f13SBarry Smith 2158273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 2159273d9f13SBarry Smith 216069b1f4b7SBarry Smith .seealso: MatCreate(), MatCreateDense(), MatSetValues(), MatSeqDenseSetLDA() 2161867c911aSBarry Smith 2162273d9f13SBarry Smith @*/ 21637087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 2164273d9f13SBarry Smith { 21654ac538c5SBarry Smith PetscErrorCode ierr; 2166a23d5eceSKris Buschelman 2167a23d5eceSKris Buschelman PetscFunctionBegin; 21684ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));CHKERRQ(ierr); 2169a23d5eceSKris Buschelman PetscFunctionReturn(0); 2170a23d5eceSKris Buschelman } 2171a23d5eceSKris Buschelman 2172a23d5eceSKris Buschelman EXTERN_C_BEGIN 2173a23d5eceSKris Buschelman #undef __FUNCT__ 2174afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense" 21757087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 2176a23d5eceSKris Buschelman { 2177273d9f13SBarry Smith Mat_SeqDense *b; 2178dfbe8321SBarry Smith PetscErrorCode ierr; 2179273d9f13SBarry Smith 2180273d9f13SBarry Smith PetscFunctionBegin; 2181273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 2182a868139aSShri Abhyankar 218334ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 218434ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 218534ef9618SShri Abhyankar 2186273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 218786d161a7SShri Abhyankar b->Mmax = B->rmap->n; 218886d161a7SShri Abhyankar b->Nmax = B->cmap->n; 218986d161a7SShri Abhyankar if (b->lda <= 0 || b->changelda) b->lda = B->rmap->n; 219086d161a7SShri Abhyankar 21919e8f95c4SLisandro Dalcin if (!data) { /* petsc-allocated storage */ 21929e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 21935afd5e0cSBarry Smith ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 2194284134d9SBarry Smith ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 2195284134d9SBarry Smith ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 21969e8f95c4SLisandro Dalcin b->user_alloc = PETSC_FALSE; 2197273d9f13SBarry Smith } else { /* user-allocated storage */ 21989e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 2199273d9f13SBarry Smith b->v = data; 2200273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 2201273d9f13SBarry Smith } 22020450473dSBarry Smith B->assembled = PETSC_TRUE; 2203273d9f13SBarry Smith PetscFunctionReturn(0); 2204273d9f13SBarry Smith } 2205a23d5eceSKris Buschelman EXTERN_C_END 2206273d9f13SBarry Smith 22071b807ce4Svictorle #undef __FUNCT__ 22081b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 22091b807ce4Svictorle /*@C 22101b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 22111b807ce4Svictorle 22121b807ce4Svictorle Input parameter: 22131b807ce4Svictorle + A - the matrix 22141b807ce4Svictorle - lda - the leading dimension 22151b807ce4Svictorle 22161b807ce4Svictorle Notes: 2217867c911aSBarry Smith This routine is to be used in conjunction with MatSeqDenseSetPreallocation(); 22181b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 22191b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 22201b807ce4Svictorle 22211b807ce4Svictorle Level: intermediate 22221b807ce4Svictorle 22231b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 22241b807ce4Svictorle 2225284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 2226867c911aSBarry Smith 22271b807ce4Svictorle @*/ 22287087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetLDA(Mat B,PetscInt lda) 22291b807ce4Svictorle { 22301b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 223121a2c019SBarry Smith 22321b807ce4Svictorle PetscFunctionBegin; 2233e32f2f54SBarry 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); 22341b807ce4Svictorle b->lda = lda; 223521a2c019SBarry Smith b->changelda = PETSC_FALSE; 223621a2c019SBarry Smith b->Mmax = PetscMax(b->Mmax,lda); 22371b807ce4Svictorle PetscFunctionReturn(0); 22381b807ce4Svictorle } 22391b807ce4Svictorle 22400bad9183SKris Buschelman /*MC 2241fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 22420bad9183SKris Buschelman 22430bad9183SKris Buschelman Options Database Keys: 22440bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 22450bad9183SKris Buschelman 22460bad9183SKris Buschelman Level: beginner 22470bad9183SKris Buschelman 224889665df3SBarry Smith .seealso: MatCreateSeqDense() 224989665df3SBarry Smith 22500bad9183SKris Buschelman M*/ 22510bad9183SKris Buschelman 2252273d9f13SBarry Smith EXTERN_C_BEGIN 22534a2ae208SSatish Balay #undef __FUNCT__ 22544a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 22557087cfbeSBarry Smith PetscErrorCode MatCreate_SeqDense(Mat B) 2256273d9f13SBarry Smith { 2257273d9f13SBarry Smith Mat_SeqDense *b; 2258dfbe8321SBarry Smith PetscErrorCode ierr; 22597c334f02SBarry Smith PetscMPIInt size; 2260273d9f13SBarry Smith 2261273d9f13SBarry Smith PetscFunctionBegin; 22627adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr); 2263e32f2f54SBarry Smith if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 226455659b69SBarry Smith 226538f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr); 2266549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 226744cd7ae7SLois Curfman McInnes B->data = (void*)b; 226818f449edSLois Curfman McInnes 226944cd7ae7SLois Curfman McInnes b->pivots = 0; 2270273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 2271273d9f13SBarry Smith b->v = 0; 227221a2c019SBarry Smith b->changelda = PETSC_FALSE; 22734e220ebcSLois Curfman McInnes 22748c778c55SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatDenseGetArray_C","MatDenseGetArray_SeqDense",MatDenseGetArray_SeqDense);CHKERRQ(ierr); 22758c778c55SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatDenseRestoreArray_C","MatDenseRestoreArray_SeqDense",MatDenseRestoreArray_SeqDense);CHKERRQ(ierr); 2276b24902e0SBarry Smith 22776a63e612SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqdense_seqaij_C","MatConvert_SeqDense_SeqAIJ",MatConvert_SeqDense_SeqAIJ);CHKERRQ(ierr); 2278ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C", 2279b24902e0SBarry Smith "MatGetFactor_seqdense_petsc", 2280b24902e0SBarry Smith MatGetFactor_seqdense_petsc);CHKERRQ(ierr); 2281a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 2282a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 2283a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 22842356f84bSDmitry Karpeev 22854ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C", 22864ae313f4SHong Zhang "MatMatMult_SeqAIJ_SeqDense", 22874ae313f4SHong Zhang MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 22882356f84bSDmitry Karpeev 2289c0c8ee5eSDmitry Karpeev 22904ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C", 22914ae313f4SHong Zhang "MatMatMultSymbolic_SeqAIJ_SeqDense", 22924ae313f4SHong Zhang MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 22934ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C", 22944ae313f4SHong Zhang "MatMatMultNumeric_SeqAIJ_SeqDense", 22954ae313f4SHong Zhang MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 229617667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 22973a40ed3dSBarry Smith PetscFunctionReturn(0); 2298289bc588SBarry Smith } 2299273d9f13SBarry Smith EXTERN_C_END 2300