1be1d678aSKris Buschelman 267e560aaSBarry Smith /* 367e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 467e560aaSBarry Smith */ 5289bc588SBarry Smith 6c6db04a5SJed Brown #include <../src/mat/impls/dense/seq/dense.h> 7c6db04a5SJed Brown #include <petscblaslapack.h> 8289bc588SBarry Smith 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; 243bda8bf91SBarry Smith ierr = PetscTypeCompareAny((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"); 245bda8bf91SBarry Smith ierr = PetscTypeCompareAny((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); 25085e2c93fSHong Zhang ierr = MatGetArray(B,&b);CHKERRQ(ierr); 25185e2c93fSHong Zhang ierr = MatGetArray(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 27285e2c93fSHong Zhang ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 27385e2c93fSHong Zhang ierr = MatRestoreArray(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); 432db4efbfdSBarry Smith LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info); 433e32f2f54SBarry Smith if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization"); 434e32f2f54SBarry Smith if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 435db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 436db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 437db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 438db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 439d5f3da31SBarry Smith A->factortype = MAT_FACTOR_LU; 440db4efbfdSBarry Smith 441dc0b31edSSatish Balay ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr); 442db4efbfdSBarry Smith #endif 443db4efbfdSBarry Smith PetscFunctionReturn(0); 444db4efbfdSBarry Smith } 445db4efbfdSBarry Smith 446db4efbfdSBarry Smith #undef __FUNCT__ 447db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense" 4480481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo) 449db4efbfdSBarry Smith { 450db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 451db4efbfdSBarry Smith PetscFunctionBegin; 452e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 453db4efbfdSBarry Smith #else 454db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 455db4efbfdSBarry Smith PetscErrorCode ierr; 456db4efbfdSBarry Smith PetscBLASInt info,n = PetscBLASIntCast(A->cmap->n); 457db4efbfdSBarry Smith 458db4efbfdSBarry Smith PetscFunctionBegin; 459db4efbfdSBarry Smith ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 460db4efbfdSBarry Smith 461db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 462db4efbfdSBarry Smith LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info); 463e32f2f54SBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 464db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 465db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 466db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 467db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 468d5f3da31SBarry Smith A->factortype = MAT_FACTOR_CHOLESKY; 469dc0b31edSSatish Balay ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr); 470db4efbfdSBarry Smith #endif 471db4efbfdSBarry Smith PetscFunctionReturn(0); 472db4efbfdSBarry Smith } 473db4efbfdSBarry Smith 474db4efbfdSBarry Smith 475db4efbfdSBarry Smith #undef __FUNCT__ 476db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 4770481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 478db4efbfdSBarry Smith { 479db4efbfdSBarry Smith PetscErrorCode ierr; 480db4efbfdSBarry Smith MatFactorInfo info; 481db4efbfdSBarry Smith 482db4efbfdSBarry Smith PetscFunctionBegin; 483db4efbfdSBarry Smith info.fill = 1.0; 484c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 485719d5645SBarry Smith ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr); 486db4efbfdSBarry Smith PetscFunctionReturn(0); 487db4efbfdSBarry Smith } 488db4efbfdSBarry Smith 489db4efbfdSBarry Smith #undef __FUNCT__ 490db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 4910481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info) 492db4efbfdSBarry Smith { 493db4efbfdSBarry Smith PetscFunctionBegin; 494c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 495719d5645SBarry Smith fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense; 496db4efbfdSBarry Smith PetscFunctionReturn(0); 497db4efbfdSBarry Smith } 498db4efbfdSBarry Smith 499db4efbfdSBarry Smith #undef __FUNCT__ 500db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 5010481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 502db4efbfdSBarry Smith { 503db4efbfdSBarry Smith PetscFunctionBegin; 504c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 505719d5645SBarry Smith fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense; 506db4efbfdSBarry Smith PetscFunctionReturn(0); 507db4efbfdSBarry Smith } 508db4efbfdSBarry Smith 509bb5747d9SMatthew Knepley EXTERN_C_BEGIN 510db4efbfdSBarry Smith #undef __FUNCT__ 511db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc" 512db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact) 513db4efbfdSBarry Smith { 514db4efbfdSBarry Smith PetscErrorCode ierr; 515db4efbfdSBarry Smith 516db4efbfdSBarry Smith PetscFunctionBegin; 517db4efbfdSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr); 518db4efbfdSBarry Smith ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 519db4efbfdSBarry Smith ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr); 520db4efbfdSBarry Smith if (ftype == MAT_FACTOR_LU){ 521db4efbfdSBarry Smith (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense; 522db4efbfdSBarry Smith } else { 523db4efbfdSBarry Smith (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense; 524db4efbfdSBarry Smith } 525d5f3da31SBarry Smith (*fact)->factortype = ftype; 526db4efbfdSBarry Smith PetscFunctionReturn(0); 527db4efbfdSBarry Smith } 528bb5747d9SMatthew Knepley EXTERN_C_END 529db4efbfdSBarry Smith 530289bc588SBarry Smith /* ------------------------------------------------------------------*/ 5314a2ae208SSatish Balay #undef __FUNCT__ 53241f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense" 53341f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 534289bc588SBarry Smith { 535c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 53687828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 537dfbe8321SBarry Smith PetscErrorCode ierr; 538d0f46423SBarry Smith PetscInt m = A->rmap->n,i; 539aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 5400805154bSBarry Smith PetscBLASInt o = 1,bm = PetscBLASIntCast(m); 541bc1b551cSSatish Balay #endif 542289bc588SBarry Smith 5433a40ed3dSBarry Smith PetscFunctionBegin; 544289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 54571044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 5462dcb1b2aSMatthew Knepley ierr = VecSet(xx,zero);CHKERRQ(ierr); 547289bc588SBarry Smith } 5481ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5491ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 550b965ef7fSBarry Smith its = its*lits; 551e32f2f54SBarry 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); 552289bc588SBarry Smith while (its--) { 553fccaa45eSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 554289bc588SBarry Smith for (i=0; i<m; i++) { 555aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 556f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 557f1747703SBarry Smith not happy about returning a double complex */ 55813f74950SBarry Smith PetscInt _i; 55987828ca2SBarry Smith PetscScalar sum = b[i]; 560f1747703SBarry Smith for (_i=0; _i<m; _i++) { 5613f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 562f1747703SBarry Smith } 563f1747703SBarry Smith xt = sum; 564f1747703SBarry Smith #else 56571044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 566f1747703SBarry Smith #endif 56755a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 568289bc588SBarry Smith } 569289bc588SBarry Smith } 570fccaa45eSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 571289bc588SBarry Smith for (i=m-1; i>=0; i--) { 572aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 573f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 574f1747703SBarry Smith not happy about returning a double complex */ 57513f74950SBarry Smith PetscInt _i; 57687828ca2SBarry Smith PetscScalar sum = b[i]; 577f1747703SBarry Smith for (_i=0; _i<m; _i++) { 5783f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 579f1747703SBarry Smith } 580f1747703SBarry Smith xt = sum; 581f1747703SBarry Smith #else 58271044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 583f1747703SBarry Smith #endif 58455a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 585289bc588SBarry Smith } 586289bc588SBarry Smith } 587289bc588SBarry Smith } 5881ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 5891ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5903a40ed3dSBarry Smith PetscFunctionReturn(0); 591289bc588SBarry Smith } 592289bc588SBarry Smith 593289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5944a2ae208SSatish Balay #undef __FUNCT__ 5954a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 596dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 597289bc588SBarry Smith { 598c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 59987828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 600dfbe8321SBarry Smith PetscErrorCode ierr; 6010805154bSBarry Smith PetscBLASInt m, n,_One=1; 602ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 6033a40ed3dSBarry Smith 6043a40ed3dSBarry Smith PetscFunctionBegin; 605d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 606d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 607d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 6081ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6091ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 61071044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One); 6111ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6121ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 613dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 6143a40ed3dSBarry Smith PetscFunctionReturn(0); 615289bc588SBarry Smith } 616800995b7SMatthew Knepley 6174a2ae208SSatish Balay #undef __FUNCT__ 6184a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 619dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 620289bc588SBarry Smith { 621c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 62287828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 623dfbe8321SBarry Smith PetscErrorCode ierr; 6240805154bSBarry Smith PetscBLASInt m, n, _One=1; 6253a40ed3dSBarry Smith 6263a40ed3dSBarry Smith PetscFunctionBegin; 627d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 628d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 629d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 6301ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6311ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 63271044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 6331ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6341ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 635dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr); 6363a40ed3dSBarry Smith PetscFunctionReturn(0); 637289bc588SBarry Smith } 6386ee01492SSatish Balay 6394a2ae208SSatish Balay #undef __FUNCT__ 6404a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 641dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 642289bc588SBarry Smith { 643c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 64487828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 645dfbe8321SBarry Smith PetscErrorCode ierr; 6460805154bSBarry Smith PetscBLASInt m, n, _One=1; 6473a40ed3dSBarry Smith 6483a40ed3dSBarry Smith PetscFunctionBegin; 649d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 650d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 651d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 652600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6531ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6541ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 65571044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 6561ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6571ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 658dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6593a40ed3dSBarry Smith PetscFunctionReturn(0); 660289bc588SBarry Smith } 6616ee01492SSatish Balay 6624a2ae208SSatish Balay #undef __FUNCT__ 6634a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 664dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 665289bc588SBarry Smith { 666c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 66787828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 668dfbe8321SBarry Smith PetscErrorCode ierr; 6690805154bSBarry Smith PetscBLASInt m, n, _One=1; 67087828ca2SBarry Smith PetscScalar _DOne=1.0; 6713a40ed3dSBarry Smith 6723a40ed3dSBarry Smith PetscFunctionBegin; 673d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 674d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 675d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 676600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6771ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6781ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 67971044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 6801ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6811ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 682dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6833a40ed3dSBarry Smith PetscFunctionReturn(0); 684289bc588SBarry Smith } 685289bc588SBarry Smith 686289bc588SBarry Smith /* -----------------------------------------------------------------*/ 6874a2ae208SSatish Balay #undef __FUNCT__ 6884a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 68913f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 690289bc588SBarry Smith { 691c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 69287828ca2SBarry Smith PetscScalar *v; 6936849ba73SBarry Smith PetscErrorCode ierr; 69413f74950SBarry Smith PetscInt i; 69567e560aaSBarry Smith 6963a40ed3dSBarry Smith PetscFunctionBegin; 697d0f46423SBarry Smith *ncols = A->cmap->n; 698289bc588SBarry Smith if (cols) { 699d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr); 700d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) (*cols)[i] = i; 701289bc588SBarry Smith } 702289bc588SBarry Smith if (vals) { 703d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 704289bc588SBarry Smith v = mat->v + row; 705d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;} 706289bc588SBarry Smith } 7073a40ed3dSBarry Smith PetscFunctionReturn(0); 708289bc588SBarry Smith } 7096ee01492SSatish Balay 7104a2ae208SSatish Balay #undef __FUNCT__ 7114a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 71213f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 713289bc588SBarry Smith { 714dfbe8321SBarry Smith PetscErrorCode ierr; 715606d414cSSatish Balay PetscFunctionBegin; 716606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 717606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 7183a40ed3dSBarry Smith PetscFunctionReturn(0); 719289bc588SBarry Smith } 720289bc588SBarry Smith /* ----------------------------------------------------------------*/ 7214a2ae208SSatish Balay #undef __FUNCT__ 7224a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 72313f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 724289bc588SBarry Smith { 725c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 72613f74950SBarry Smith PetscInt i,j,idx=0; 727d6dfbf8fSBarry Smith 7283a40ed3dSBarry Smith PetscFunctionBegin; 72971fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 730289bc588SBarry Smith if (!mat->roworiented) { 731dbb450caSBarry Smith if (addv == INSERT_VALUES) { 732289bc588SBarry Smith for (j=0; j<n; j++) { 733cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 7342515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 735e32f2f54SBarry 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); 73658804f6dSBarry Smith #endif 737289bc588SBarry Smith for (i=0; i<m; i++) { 738cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 7392515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 740e32f2f54SBarry 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); 74158804f6dSBarry Smith #endif 742cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 743289bc588SBarry Smith } 744289bc588SBarry Smith } 7453a40ed3dSBarry Smith } else { 746289bc588SBarry Smith for (j=0; j<n; j++) { 747cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 7482515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 749e32f2f54SBarry 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); 75058804f6dSBarry Smith #endif 751289bc588SBarry Smith for (i=0; i<m; i++) { 752cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 7532515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 754e32f2f54SBarry 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); 75558804f6dSBarry Smith #endif 756cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 757289bc588SBarry Smith } 758289bc588SBarry Smith } 759289bc588SBarry Smith } 7603a40ed3dSBarry Smith } else { 761dbb450caSBarry Smith if (addv == INSERT_VALUES) { 762e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 763cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7642515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 765e32f2f54SBarry 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); 76658804f6dSBarry Smith #endif 767e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 768cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7692515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 770e32f2f54SBarry 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); 77158804f6dSBarry Smith #endif 772cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 773e8d4e0b9SBarry Smith } 774e8d4e0b9SBarry Smith } 7753a40ed3dSBarry Smith } else { 776289bc588SBarry Smith for (i=0; i<m; i++) { 777cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7782515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 779e32f2f54SBarry 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); 78058804f6dSBarry Smith #endif 781289bc588SBarry Smith for (j=0; j<n; j++) { 782cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7832515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 784e32f2f54SBarry 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); 78558804f6dSBarry Smith #endif 786cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 787289bc588SBarry Smith } 788289bc588SBarry Smith } 789289bc588SBarry Smith } 790e8d4e0b9SBarry Smith } 7913a40ed3dSBarry Smith PetscFunctionReturn(0); 792289bc588SBarry Smith } 793e8d4e0b9SBarry Smith 7944a2ae208SSatish Balay #undef __FUNCT__ 7954a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 79613f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 797ae80bb75SLois Curfman McInnes { 798ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 79913f74950SBarry Smith PetscInt i,j; 800ae80bb75SLois Curfman McInnes 8013a40ed3dSBarry Smith PetscFunctionBegin; 802ae80bb75SLois Curfman McInnes /* row-oriented output */ 803ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 80497e567efSBarry Smith if (indexm[i] < 0) {v += n;continue;} 805e32f2f54SBarry 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); 806ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 8076f31f424SBarry Smith if (indexn[j] < 0) {v++; continue;} 808e32f2f54SBarry 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); 80997e567efSBarry Smith *v++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 810ae80bb75SLois Curfman McInnes } 811ae80bb75SLois Curfman McInnes } 8123a40ed3dSBarry Smith PetscFunctionReturn(0); 813ae80bb75SLois Curfman McInnes } 814ae80bb75SLois Curfman McInnes 815289bc588SBarry Smith /* -----------------------------------------------------------------*/ 816289bc588SBarry Smith 8174a2ae208SSatish Balay #undef __FUNCT__ 8185bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_SeqDense" 819112444f4SShri Abhyankar PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer) 820aabbc4fbSShri Abhyankar { 821aabbc4fbSShri Abhyankar Mat_SeqDense *a; 822aabbc4fbSShri Abhyankar PetscErrorCode ierr; 823aabbc4fbSShri Abhyankar PetscInt *scols,i,j,nz,header[4]; 824aabbc4fbSShri Abhyankar int fd; 825aabbc4fbSShri Abhyankar PetscMPIInt size; 826aabbc4fbSShri Abhyankar PetscInt *rowlengths = 0,M,N,*cols,grows,gcols; 827aabbc4fbSShri Abhyankar PetscScalar *vals,*svals,*v,*w; 828aabbc4fbSShri Abhyankar MPI_Comm comm = ((PetscObject)viewer)->comm; 829aabbc4fbSShri Abhyankar 830aabbc4fbSShri Abhyankar PetscFunctionBegin; 831aabbc4fbSShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 832aabbc4fbSShri Abhyankar if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor"); 833aabbc4fbSShri Abhyankar ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 834aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 835aabbc4fbSShri Abhyankar if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 836aabbc4fbSShri Abhyankar M = header[1]; N = header[2]; nz = header[3]; 837aabbc4fbSShri Abhyankar 838aabbc4fbSShri Abhyankar /* set global size if not set already*/ 839aabbc4fbSShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 840aabbc4fbSShri Abhyankar ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr); 841aabbc4fbSShri Abhyankar } else { 842aabbc4fbSShri Abhyankar /* if sizes and type are already set, check if the vector global sizes are correct */ 843aabbc4fbSShri Abhyankar ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr); 844aabbc4fbSShri 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); 845aabbc4fbSShri Abhyankar } 846aabbc4fbSShri Abhyankar ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 847aabbc4fbSShri Abhyankar 848aabbc4fbSShri Abhyankar if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 849aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 850aabbc4fbSShri Abhyankar v = a->v; 851aabbc4fbSShri Abhyankar /* Allocate some temp space to read in the values and then flip them 852aabbc4fbSShri Abhyankar from row major to column major */ 853aabbc4fbSShri Abhyankar ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 854aabbc4fbSShri Abhyankar /* read in nonzero values */ 855aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 856aabbc4fbSShri Abhyankar /* now flip the values and store them in the matrix*/ 857aabbc4fbSShri Abhyankar for (j=0; j<N; j++) { 858aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 859aabbc4fbSShri Abhyankar *v++ =w[i*N+j]; 860aabbc4fbSShri Abhyankar } 861aabbc4fbSShri Abhyankar } 862aabbc4fbSShri Abhyankar ierr = PetscFree(w);CHKERRQ(ierr); 863aabbc4fbSShri Abhyankar } else { 864aabbc4fbSShri Abhyankar /* read row lengths */ 865aabbc4fbSShri Abhyankar ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 866aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 867aabbc4fbSShri Abhyankar 868aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 869aabbc4fbSShri Abhyankar v = a->v; 870aabbc4fbSShri Abhyankar 871aabbc4fbSShri Abhyankar /* read column indices and nonzeros */ 872aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr); 873aabbc4fbSShri Abhyankar cols = scols; 874aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 875aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 876aabbc4fbSShri Abhyankar vals = svals; 877aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 878aabbc4fbSShri Abhyankar 879aabbc4fbSShri Abhyankar /* insert into matrix */ 880aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 881aabbc4fbSShri Abhyankar for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 882aabbc4fbSShri Abhyankar svals += rowlengths[i]; scols += rowlengths[i]; 883aabbc4fbSShri Abhyankar } 884aabbc4fbSShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 885aabbc4fbSShri Abhyankar ierr = PetscFree(cols);CHKERRQ(ierr); 886aabbc4fbSShri Abhyankar ierr = PetscFree(rowlengths);CHKERRQ(ierr); 887aabbc4fbSShri Abhyankar } 888aabbc4fbSShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 889aabbc4fbSShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 890aabbc4fbSShri Abhyankar 891aabbc4fbSShri Abhyankar PetscFunctionReturn(0); 892aabbc4fbSShri Abhyankar } 893aabbc4fbSShri Abhyankar 894aabbc4fbSShri Abhyankar #undef __FUNCT__ 8954a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 8966849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 897289bc588SBarry Smith { 898932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 899dfbe8321SBarry Smith PetscErrorCode ierr; 90013f74950SBarry Smith PetscInt i,j; 9012dcb1b2aSMatthew Knepley const char *name; 90287828ca2SBarry Smith PetscScalar *v; 903f3ef73ceSBarry Smith PetscViewerFormat format; 9045f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX) 905ace3abfcSBarry Smith PetscBool allreal = PETSC_TRUE; 9065f481a85SSatish Balay #endif 907932b0c3eSLois Curfman McInnes 9083a40ed3dSBarry Smith PetscFunctionBegin; 909b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 910456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 9113a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 912fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 913d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 9147566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 915d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 91644cd7ae7SLois Curfman McInnes v = a->v + i; 91777431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 918d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 919aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 920329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 921a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 922329f5518SBarry Smith } else if (PetscRealPart(*v)) { 923a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr); 9246831982aSBarry Smith } 92580cd9d93SLois Curfman McInnes #else 9266831982aSBarry Smith if (*v) { 927a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr); 9286831982aSBarry Smith } 92980cd9d93SLois Curfman McInnes #endif 9301b807ce4Svictorle v += a->lda; 93180cd9d93SLois Curfman McInnes } 932b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 93380cd9d93SLois Curfman McInnes } 934d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 9353a40ed3dSBarry Smith } else { 936d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 937aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 93847989497SBarry Smith /* determine if matrix has all real values */ 93947989497SBarry Smith v = a->v; 940d0f46423SBarry Smith for (i=0; i<A->rmap->n*A->cmap->n; i++) { 941ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 94247989497SBarry Smith } 94347989497SBarry Smith #endif 944fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 9453a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 946d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr); 947d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 948fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 9497566de4bSShri Abhyankar } else { 9507566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 951ffac6cdbSBarry Smith } 952ffac6cdbSBarry Smith 953d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 954932b0c3eSLois Curfman McInnes v = a->v + i; 955d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 956aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 95747989497SBarry Smith if (allreal) { 958f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr); 95947989497SBarry Smith } else { 960f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 96147989497SBarry Smith } 962289bc588SBarry Smith #else 963f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr); 964289bc588SBarry Smith #endif 9651b807ce4Svictorle v += a->lda; 966289bc588SBarry Smith } 967b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 968289bc588SBarry Smith } 969fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 970b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 971ffac6cdbSBarry Smith } 972d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 973da3a660dSBarry Smith } 974b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 9753a40ed3dSBarry Smith PetscFunctionReturn(0); 976289bc588SBarry Smith } 977289bc588SBarry Smith 9784a2ae208SSatish Balay #undef __FUNCT__ 9794a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 9806849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 981932b0c3eSLois Curfman McInnes { 982932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9836849ba73SBarry Smith PetscErrorCode ierr; 98413f74950SBarry Smith int fd; 985d0f46423SBarry Smith PetscInt ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n; 986f4403165SShri Abhyankar PetscScalar *v,*anonz,*vals; 987f4403165SShri Abhyankar PetscViewerFormat format; 988932b0c3eSLois Curfman McInnes 9893a40ed3dSBarry Smith PetscFunctionBegin; 990b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 99190ace30eSBarry Smith 992f4403165SShri Abhyankar ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 993f4403165SShri Abhyankar if (format == PETSC_VIEWER_NATIVE) { 994f4403165SShri Abhyankar /* store the matrix as a dense matrix */ 995f4403165SShri Abhyankar ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 996f4403165SShri Abhyankar col_lens[0] = MAT_FILE_CLASSID; 997f4403165SShri Abhyankar col_lens[1] = m; 998f4403165SShri Abhyankar col_lens[2] = n; 999f4403165SShri Abhyankar col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 1000f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 1001f4403165SShri Abhyankar ierr = PetscFree(col_lens);CHKERRQ(ierr); 1002f4403165SShri Abhyankar 1003f4403165SShri Abhyankar /* write out matrix, by rows */ 1004f4403165SShri Abhyankar ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 1005f4403165SShri Abhyankar v = a->v; 1006f4403165SShri Abhyankar for (j=0; j<n; j++) { 1007f4403165SShri Abhyankar for (i=0; i<m; i++) { 1008f4403165SShri Abhyankar vals[j + i*n] = *v++; 1009f4403165SShri Abhyankar } 1010f4403165SShri Abhyankar } 1011f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 1012f4403165SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 1013f4403165SShri Abhyankar } else { 101413f74950SBarry Smith ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 10150700a824SBarry Smith col_lens[0] = MAT_FILE_CLASSID; 1016932b0c3eSLois Curfman McInnes col_lens[1] = m; 1017932b0c3eSLois Curfman McInnes col_lens[2] = n; 1018932b0c3eSLois Curfman McInnes col_lens[3] = nz; 1019932b0c3eSLois Curfman McInnes 1020932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 1021932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 10226f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 1023932b0c3eSLois Curfman McInnes 1024932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 1025932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 1026932b0c3eSLois Curfman McInnes ict = 0; 1027932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 1028932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 1029932b0c3eSLois Curfman McInnes } 10306f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 1031606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 1032932b0c3eSLois Curfman McInnes 1033932b0c3eSLois Curfman McInnes /* store nonzero values */ 103487828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 1035932b0c3eSLois Curfman McInnes ict = 0; 1036932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 1037932b0c3eSLois Curfman McInnes v = a->v + i; 1038932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 10391b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 1040932b0c3eSLois Curfman McInnes } 1041932b0c3eSLois Curfman McInnes } 10426f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 1043606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 1044f4403165SShri Abhyankar } 10453a40ed3dSBarry Smith PetscFunctionReturn(0); 1046932b0c3eSLois Curfman McInnes } 1047932b0c3eSLois Curfman McInnes 10484a2ae208SSatish Balay #undef __FUNCT__ 10494a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 1050dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 1051f1af5d2fSBarry Smith { 1052f1af5d2fSBarry Smith Mat A = (Mat) Aa; 1053f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 10546849ba73SBarry Smith PetscErrorCode ierr; 1055d0f46423SBarry Smith PetscInt m = A->rmap->n,n = A->cmap->n,color,i,j; 105687828ca2SBarry Smith PetscScalar *v = a->v; 1057b0a32e0cSBarry Smith PetscViewer viewer; 1058b0a32e0cSBarry Smith PetscDraw popup; 1059329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 1060f3ef73ceSBarry Smith PetscViewerFormat format; 1061f1af5d2fSBarry Smith 1062f1af5d2fSBarry Smith PetscFunctionBegin; 1063f1af5d2fSBarry Smith 1064f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 1065b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1066b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1067f1af5d2fSBarry Smith 1068f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 1069fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 1070f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 1071b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1072f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1073f1af5d2fSBarry Smith x_l = j; 1074f1af5d2fSBarry Smith x_r = x_l + 1.0; 1075f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1076f1af5d2fSBarry Smith y_l = m - i - 1.0; 1077f1af5d2fSBarry Smith y_r = y_l + 1.0; 1078f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 1079329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 1080b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1081329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 1082b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1083f1af5d2fSBarry Smith } else { 1084f1af5d2fSBarry Smith continue; 1085f1af5d2fSBarry Smith } 1086f1af5d2fSBarry Smith #else 1087f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 1088b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1089f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 1090b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1091f1af5d2fSBarry Smith } else { 1092f1af5d2fSBarry Smith continue; 1093f1af5d2fSBarry Smith } 1094f1af5d2fSBarry Smith #endif 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 } else { 1099f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 1100f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 1101f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 1102f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 1103f1af5d2fSBarry Smith } 1104b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 1105b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 1106b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 1107f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1108f1af5d2fSBarry Smith x_l = j; 1109f1af5d2fSBarry Smith x_r = x_l + 1.0; 1110f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1111f1af5d2fSBarry Smith y_l = m - i - 1.0; 1112f1af5d2fSBarry Smith y_r = y_l + 1.0; 1113b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 1114b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1115f1af5d2fSBarry Smith } 1116f1af5d2fSBarry Smith } 1117f1af5d2fSBarry Smith } 1118f1af5d2fSBarry Smith PetscFunctionReturn(0); 1119f1af5d2fSBarry Smith } 1120f1af5d2fSBarry Smith 11214a2ae208SSatish Balay #undef __FUNCT__ 11224a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 1123dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 1124f1af5d2fSBarry Smith { 1125b0a32e0cSBarry Smith PetscDraw draw; 1126ace3abfcSBarry Smith PetscBool isnull; 1127329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 1128dfbe8321SBarry Smith PetscErrorCode ierr; 1129f1af5d2fSBarry Smith 1130f1af5d2fSBarry Smith PetscFunctionBegin; 1131b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1132b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1133abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 1134f1af5d2fSBarry Smith 1135f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1136d0f46423SBarry Smith xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0; 1137f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 1138b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1139b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 1140f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 1141f1af5d2fSBarry Smith PetscFunctionReturn(0); 1142f1af5d2fSBarry Smith } 1143f1af5d2fSBarry Smith 11444a2ae208SSatish Balay #undef __FUNCT__ 11454a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 1146dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1147932b0c3eSLois Curfman McInnes { 1148dfbe8321SBarry Smith PetscErrorCode ierr; 1149ace3abfcSBarry Smith PetscBool iascii,isbinary,isdraw; 1150932b0c3eSLois Curfman McInnes 11513a40ed3dSBarry Smith PetscFunctionBegin; 11522692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 11532692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 11542692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 11550f5bd95cSBarry Smith 1156c45a1595SBarry Smith if (iascii) { 1157c45a1595SBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 11580f5bd95cSBarry Smith } else if (isbinary) { 11593a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1160f1af5d2fSBarry Smith } else if (isdraw) { 1161f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 11625cd90555SBarry Smith } else { 1163e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1164932b0c3eSLois Curfman McInnes } 11653a40ed3dSBarry Smith PetscFunctionReturn(0); 1166932b0c3eSLois Curfman McInnes } 1167289bc588SBarry Smith 11684a2ae208SSatish Balay #undef __FUNCT__ 11694a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1170dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1171289bc588SBarry Smith { 1172ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1173dfbe8321SBarry Smith PetscErrorCode ierr; 117490f02eecSBarry Smith 11753a40ed3dSBarry Smith PetscFunctionBegin; 1176aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1177d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n); 1178a5a9c739SBarry Smith #endif 117905b42c5fSBarry Smith ierr = PetscFree(l->pivots);CHKERRQ(ierr); 11806857c123SSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1181bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 1182dbd8c25aSHong Zhang 1183dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 1184901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 11854ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11864ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11874ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11883a40ed3dSBarry Smith PetscFunctionReturn(0); 1189289bc588SBarry Smith } 1190289bc588SBarry Smith 11914a2ae208SSatish Balay #undef __FUNCT__ 11924a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1193fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout) 1194289bc588SBarry Smith { 1195c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11966849ba73SBarry Smith PetscErrorCode ierr; 119713f74950SBarry Smith PetscInt k,j,m,n,M; 119887828ca2SBarry Smith PetscScalar *v,tmp; 119948b35521SBarry Smith 12003a40ed3dSBarry Smith PetscFunctionBegin; 1201d0f46423SBarry Smith v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n; 1202e9695a30SBarry Smith if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */ 1203e7e72b3dSBarry Smith if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1204e7e72b3dSBarry Smith else { 1205d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1206289bc588SBarry Smith for (k=0; k<j; k++) { 12071b807ce4Svictorle tmp = v[j + k*M]; 12081b807ce4Svictorle v[j + k*M] = v[k + j*M]; 12091b807ce4Svictorle v[k + j*M] = tmp; 1210289bc588SBarry Smith } 1211289bc588SBarry Smith } 1212d64ed03dSBarry Smith } 12133a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1214d3e5ee88SLois Curfman McInnes Mat tmat; 1215ec8511deSBarry Smith Mat_SeqDense *tmatd; 121687828ca2SBarry Smith PetscScalar *v2; 1217ea709b57SSatish Balay 1218fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 12197adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr); 1220d0f46423SBarry Smith ierr = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr); 12217adad957SLisandro Dalcin ierr = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 12225c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1223fc4dec0aSBarry Smith } else { 1224fc4dec0aSBarry Smith tmat = *matout; 1225fc4dec0aSBarry Smith } 1226ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 12270de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1228d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 12291b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1230d3e5ee88SLois Curfman McInnes } 12316d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12326d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1233d3e5ee88SLois Curfman McInnes *matout = tmat; 123448b35521SBarry Smith } 12353a40ed3dSBarry Smith PetscFunctionReturn(0); 1236289bc588SBarry Smith } 1237289bc588SBarry Smith 12384a2ae208SSatish Balay #undef __FUNCT__ 12394a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1240ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool *flg) 1241289bc588SBarry Smith { 1242c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1243c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 124413f74950SBarry Smith PetscInt i,j; 124587828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 12469ea5d5aeSSatish Balay 12473a40ed3dSBarry Smith PetscFunctionBegin; 1248d0f46423SBarry Smith if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1249d0f46423SBarry Smith if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1250d0f46423SBarry Smith for (i=0; i<A1->rmap->n; i++) { 12511b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1252d0f46423SBarry Smith for (j=0; j<A1->cmap->n; j++) { 12533a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 12541b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 12551b807ce4Svictorle } 1256289bc588SBarry Smith } 125777c4ece6SBarry Smith *flg = PETSC_TRUE; 12583a40ed3dSBarry Smith PetscFunctionReturn(0); 1259289bc588SBarry Smith } 1260289bc588SBarry Smith 12614a2ae208SSatish Balay #undef __FUNCT__ 12624a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1263dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1264289bc588SBarry Smith { 1265c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1266dfbe8321SBarry Smith PetscErrorCode ierr; 126713f74950SBarry Smith PetscInt i,n,len; 126887828ca2SBarry Smith PetscScalar *x,zero = 0.0; 126944cd7ae7SLois Curfman McInnes 12703a40ed3dSBarry Smith PetscFunctionBegin; 12712dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 12727a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 12731ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1274d0f46423SBarry Smith len = PetscMin(A->rmap->n,A->cmap->n); 1275e32f2f54SBarry Smith if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 127644cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 12771b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1278289bc588SBarry Smith } 12791ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 12803a40ed3dSBarry Smith PetscFunctionReturn(0); 1281289bc588SBarry Smith } 1282289bc588SBarry Smith 12834a2ae208SSatish Balay #undef __FUNCT__ 12844a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1285dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1286289bc588SBarry Smith { 1287c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 128887828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1289dfbe8321SBarry Smith PetscErrorCode ierr; 1290d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n; 129155659b69SBarry Smith 12923a40ed3dSBarry Smith PetscFunctionBegin; 129328988994SBarry Smith if (ll) { 12947a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 12951ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1296e32f2f54SBarry Smith if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1297da3a660dSBarry Smith for (i=0; i<m; i++) { 1298da3a660dSBarry Smith x = l[i]; 1299da3a660dSBarry Smith v = mat->v + i; 1300da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1301da3a660dSBarry Smith } 13021ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1303efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1304da3a660dSBarry Smith } 130528988994SBarry Smith if (rr) { 13067a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 13071ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1308e32f2f54SBarry Smith if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1309da3a660dSBarry Smith for (i=0; i<n; i++) { 1310da3a660dSBarry Smith x = r[i]; 1311da3a660dSBarry Smith v = mat->v + i*m; 1312da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1313da3a660dSBarry Smith } 13141ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1315efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1316da3a660dSBarry Smith } 13173a40ed3dSBarry Smith PetscFunctionReturn(0); 1318289bc588SBarry Smith } 1319289bc588SBarry Smith 13204a2ae208SSatish Balay #undef __FUNCT__ 13214a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1322dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1323289bc588SBarry Smith { 1324c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 132587828ca2SBarry Smith PetscScalar *v = mat->v; 1326329f5518SBarry Smith PetscReal sum = 0.0; 1327d0f46423SBarry Smith PetscInt lda=mat->lda,m=A->rmap->n,i,j; 1328efee365bSSatish Balay PetscErrorCode ierr; 132955659b69SBarry Smith 13303a40ed3dSBarry Smith PetscFunctionBegin; 1331289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1332a5ce6ee0Svictorle if (lda>m) { 1333d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1334a5ce6ee0Svictorle v = mat->v+j*lda; 1335a5ce6ee0Svictorle for (i=0; i<m; i++) { 1336a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1337a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1338a5ce6ee0Svictorle #else 1339a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1340a5ce6ee0Svictorle #endif 1341a5ce6ee0Svictorle } 1342a5ce6ee0Svictorle } 1343a5ce6ee0Svictorle } else { 1344d0f46423SBarry Smith for (i=0; i<A->cmap->n*A->rmap->n; i++) { 1345aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1346329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1347289bc588SBarry Smith #else 1348289bc588SBarry Smith sum += (*v)*(*v); v++; 1349289bc588SBarry Smith #endif 1350289bc588SBarry Smith } 1351a5ce6ee0Svictorle } 13528f1a2a5eSBarry Smith *nrm = PetscSqrtReal(sum); 1353dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13543a40ed3dSBarry Smith } else if (type == NORM_1) { 1355064f8208SBarry Smith *nrm = 0.0; 1356d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 13571b807ce4Svictorle v = mat->v + j*mat->lda; 1358289bc588SBarry Smith sum = 0.0; 1359d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 136033a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1361289bc588SBarry Smith } 1362064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1363289bc588SBarry Smith } 1364d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13653a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1366064f8208SBarry Smith *nrm = 0.0; 1367d0f46423SBarry Smith for (j=0; j<A->rmap->n; j++) { 1368289bc588SBarry Smith v = mat->v + j; 1369289bc588SBarry Smith sum = 0.0; 1370d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) { 13711b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1372289bc588SBarry Smith } 1373064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1374289bc588SBarry Smith } 1375d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 1376e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 13773a40ed3dSBarry Smith PetscFunctionReturn(0); 1378289bc588SBarry Smith } 1379289bc588SBarry Smith 13804a2ae208SSatish Balay #undef __FUNCT__ 13814a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 1382ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg) 1383289bc588SBarry Smith { 1384c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 138563ba0a88SBarry Smith PetscErrorCode ierr; 138667e560aaSBarry Smith 13873a40ed3dSBarry Smith PetscFunctionBegin; 1388b5a2b587SKris Buschelman switch (op) { 1389b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 13904e0d8c25SBarry Smith aij->roworiented = flg; 1391b5a2b587SKris Buschelman break; 1392512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 1393b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 13943971808eSMatthew Knepley case MAT_NEW_NONZERO_ALLOCATION_ERR: 13954e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 139613fa8e87SLisandro Dalcin case MAT_KEEP_NONZERO_PATTERN: 1397b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1398b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 139977e54ba9SKris Buschelman case MAT_SYMMETRIC: 140077e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 14019a4540c5SBarry Smith case MAT_HERMITIAN: 14029a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1403600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 1404290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 140577e54ba9SKris Buschelman break; 1406b5a2b587SKris Buschelman default: 1407e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 14083a40ed3dSBarry Smith } 14093a40ed3dSBarry Smith PetscFunctionReturn(0); 1410289bc588SBarry Smith } 1411289bc588SBarry Smith 14124a2ae208SSatish Balay #undef __FUNCT__ 14134a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1414dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 14156f0a148fSBarry Smith { 1416ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 14176849ba73SBarry Smith PetscErrorCode ierr; 1418d0f46423SBarry Smith PetscInt lda=l->lda,m=A->rmap->n,j; 14193a40ed3dSBarry Smith 14203a40ed3dSBarry Smith PetscFunctionBegin; 1421a5ce6ee0Svictorle if (lda>m) { 1422d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1423a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1424a5ce6ee0Svictorle } 1425a5ce6ee0Svictorle } else { 1426d0f46423SBarry Smith ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1427a5ce6ee0Svictorle } 14283a40ed3dSBarry Smith PetscFunctionReturn(0); 14296f0a148fSBarry Smith } 14306f0a148fSBarry Smith 14314a2ae208SSatish Balay #undef __FUNCT__ 14324a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 14332b40b63fSBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 14346f0a148fSBarry Smith { 143597b48c8fSBarry Smith PetscErrorCode ierr; 1436ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1437b9679d65SBarry Smith PetscInt m = l->lda, n = A->cmap->n, i,j; 143897b48c8fSBarry Smith PetscScalar *slot,*bb; 143997b48c8fSBarry Smith const PetscScalar *xx; 144055659b69SBarry Smith 14413a40ed3dSBarry Smith PetscFunctionBegin; 1442b9679d65SBarry Smith #if defined(PETSC_USE_DEBUG) 1443b9679d65SBarry Smith for (i=0; i<N; i++) { 1444b9679d65SBarry Smith if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed"); 1445b9679d65SBarry 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); 1446b9679d65SBarry Smith } 1447b9679d65SBarry Smith #endif 1448b9679d65SBarry Smith 144997b48c8fSBarry Smith /* fix right hand side if needed */ 145097b48c8fSBarry Smith if (x && b) { 145197b48c8fSBarry Smith ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 145297b48c8fSBarry Smith ierr = VecGetArray(b,&bb);CHKERRQ(ierr); 145397b48c8fSBarry Smith for (i=0; i<N; i++) { 145497b48c8fSBarry Smith bb[rows[i]] = diag*xx[rows[i]]; 145597b48c8fSBarry Smith } 145697b48c8fSBarry Smith ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 145797b48c8fSBarry Smith ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr); 145897b48c8fSBarry Smith } 145997b48c8fSBarry Smith 14606f0a148fSBarry Smith for (i=0; i<N; i++) { 14616f0a148fSBarry Smith slot = l->v + rows[i]; 1462b9679d65SBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += m;} 14636f0a148fSBarry Smith } 1464f4df32b1SMatthew Knepley if (diag != 0.0) { 1465b9679d65SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices"); 14666f0a148fSBarry Smith for (i=0; i<N; i++) { 1467b9679d65SBarry Smith slot = l->v + (m+1)*rows[i]; 1468f4df32b1SMatthew Knepley *slot = diag; 14696f0a148fSBarry Smith } 14706f0a148fSBarry Smith } 14713a40ed3dSBarry Smith PetscFunctionReturn(0); 14726f0a148fSBarry Smith } 1473557bce09SLois Curfman McInnes 14744a2ae208SSatish Balay #undef __FUNCT__ 14754a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 1476dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 147764e87e97SBarry Smith { 1478c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14793a40ed3dSBarry Smith 14803a40ed3dSBarry Smith PetscFunctionBegin; 1481e32f2f54SBarry 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"); 148264e87e97SBarry Smith *array = mat->v; 14833a40ed3dSBarry Smith PetscFunctionReturn(0); 148464e87e97SBarry Smith } 14850754003eSLois Curfman McInnes 14864a2ae208SSatish Balay #undef __FUNCT__ 14874a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 1488dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1489ff14e315SSatish Balay { 14903a40ed3dSBarry Smith PetscFunctionBegin; 149109b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 14923a40ed3dSBarry Smith PetscFunctionReturn(0); 1493ff14e315SSatish Balay } 14940754003eSLois Curfman McInnes 14954a2ae208SSatish Balay #undef __FUNCT__ 14964a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 149713f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 14980754003eSLois Curfman McInnes { 1499c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 15006849ba73SBarry Smith PetscErrorCode ierr; 15015d0c19d7SBarry Smith PetscInt i,j,nrows,ncols; 15025d0c19d7SBarry Smith const PetscInt *irow,*icol; 150387828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 15040754003eSLois Curfman McInnes Mat newmat; 15050754003eSLois Curfman McInnes 15063a40ed3dSBarry Smith PetscFunctionBegin; 150778b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 150878b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1509e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1510e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 15110754003eSLois Curfman McInnes 1512182d2002SSatish Balay /* Check submatrixcall */ 1513182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 151413f74950SBarry Smith PetscInt n_cols,n_rows; 1515182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 151621a2c019SBarry Smith if (n_rows != nrows || n_cols != ncols) { 1517f746d493SDmitry Karpeev /* resize the result matrix to match number of requested rows/columns */ 1518c61587bbSBarry Smith ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 151921a2c019SBarry Smith } 1520182d2002SSatish Balay newmat = *B; 1521182d2002SSatish Balay } else { 15220754003eSLois Curfman McInnes /* Create and fill new matrix */ 15237adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr); 1524f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 15257adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 15265c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1527182d2002SSatish Balay } 1528182d2002SSatish Balay 1529182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1530182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1531182d2002SSatish Balay 1532182d2002SSatish Balay for (i=0; i<ncols; i++) { 15336de62eeeSBarry Smith av = v + mat->lda*icol[i]; 1534182d2002SSatish Balay for (j=0; j<nrows; j++) { 1535182d2002SSatish Balay *bv++ = av[irow[j]]; 15360754003eSLois Curfman McInnes } 15370754003eSLois Curfman McInnes } 1538182d2002SSatish Balay 1539182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 15406d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15416d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15420754003eSLois Curfman McInnes 15430754003eSLois Curfman McInnes /* Free work space */ 154478b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 154578b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1546182d2002SSatish Balay *B = newmat; 15473a40ed3dSBarry Smith PetscFunctionReturn(0); 15480754003eSLois Curfman McInnes } 15490754003eSLois Curfman McInnes 15504a2ae208SSatish Balay #undef __FUNCT__ 15514a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 155213f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1553905e6a2fSBarry Smith { 15546849ba73SBarry Smith PetscErrorCode ierr; 155513f74950SBarry Smith PetscInt i; 1556905e6a2fSBarry Smith 15573a40ed3dSBarry Smith PetscFunctionBegin; 1558905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1559b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1560905e6a2fSBarry Smith } 1561905e6a2fSBarry Smith 1562905e6a2fSBarry Smith for (i=0; i<n; i++) { 15636a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1564905e6a2fSBarry Smith } 15653a40ed3dSBarry Smith PetscFunctionReturn(0); 1566905e6a2fSBarry Smith } 1567905e6a2fSBarry Smith 15684a2ae208SSatish Balay #undef __FUNCT__ 1569c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1570c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1571c0aa2d19SHong Zhang { 1572c0aa2d19SHong Zhang PetscFunctionBegin; 1573c0aa2d19SHong Zhang PetscFunctionReturn(0); 1574c0aa2d19SHong Zhang } 1575c0aa2d19SHong Zhang 1576c0aa2d19SHong Zhang #undef __FUNCT__ 1577c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1578c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1579c0aa2d19SHong Zhang { 1580c0aa2d19SHong Zhang PetscFunctionBegin; 1581c0aa2d19SHong Zhang PetscFunctionReturn(0); 1582c0aa2d19SHong Zhang } 1583c0aa2d19SHong Zhang 1584c0aa2d19SHong Zhang #undef __FUNCT__ 15854a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1586dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 15874b0e389bSBarry Smith { 15884b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 15896849ba73SBarry Smith PetscErrorCode ierr; 1590d0f46423SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j; 15913a40ed3dSBarry Smith 15923a40ed3dSBarry Smith PetscFunctionBegin; 159333f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 159433f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1595cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 15963a40ed3dSBarry Smith PetscFunctionReturn(0); 15973a40ed3dSBarry Smith } 1598e32f2f54SBarry Smith if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1599a5ce6ee0Svictorle if (lda1>m || lda2>m) { 16000dbb7854Svictorle for (j=0; j<n; j++) { 1601a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1602a5ce6ee0Svictorle } 1603a5ce6ee0Svictorle } else { 1604d0f46423SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1605a5ce6ee0Svictorle } 1606273d9f13SBarry Smith PetscFunctionReturn(0); 1607273d9f13SBarry Smith } 1608273d9f13SBarry Smith 16094a2ae208SSatish Balay #undef __FUNCT__ 1610*4994cf47SJed Brown #define __FUNCT__ "MatSetUp_SeqDense" 1611*4994cf47SJed Brown PetscErrorCode MatSetUp_SeqDense(Mat A) 1612273d9f13SBarry Smith { 1613dfbe8321SBarry Smith PetscErrorCode ierr; 1614273d9f13SBarry Smith 1615273d9f13SBarry Smith PetscFunctionBegin; 1616273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 16173a40ed3dSBarry Smith PetscFunctionReturn(0); 16184b0e389bSBarry Smith } 16194b0e389bSBarry Smith 1620284134d9SBarry Smith #undef __FUNCT__ 1621284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense" 1622284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N) 1623284134d9SBarry Smith { 1624284134d9SBarry Smith PetscFunctionBegin; 162521a2c019SBarry Smith /* this will not be called before lda, Mmax, and Nmax have been set */ 1626284134d9SBarry Smith m = PetscMax(m,M); 1627284134d9SBarry Smith n = PetscMax(n,N); 1628a868139aSShri Abhyankar 162986d161a7SShri Abhyankar /* if (m > a->Mmax) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot yet resize number rows of dense matrix larger then its initial size %d, requested %d",a->lda,(int)m); 163086d161a7SShri Abhyankar if (n > a->Nmax) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot yet resize number columns of dense matrix larger then its initial size %d, requested %d",a->Nmax,(int)n); 163186d161a7SShri Abhyankar */ 1632dc5cefdeSJed Brown A->rmap->n = A->rmap->N = m; 1633d0f46423SBarry Smith A->cmap->n = A->cmap->N = n; 1634284134d9SBarry Smith PetscFunctionReturn(0); 1635284134d9SBarry Smith } 1636170fe5c8SBarry Smith 1637ba337c44SJed Brown #undef __FUNCT__ 1638ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense" 1639ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A) 1640ba337c44SJed Brown { 1641ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1642ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1643ba337c44SJed Brown PetscScalar *aa = a->v; 1644ba337c44SJed Brown 1645ba337c44SJed Brown PetscFunctionBegin; 1646ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]); 1647ba337c44SJed Brown PetscFunctionReturn(0); 1648ba337c44SJed Brown } 1649ba337c44SJed Brown 1650ba337c44SJed Brown #undef __FUNCT__ 1651ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense" 1652ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A) 1653ba337c44SJed Brown { 1654ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1655ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1656ba337c44SJed Brown PetscScalar *aa = a->v; 1657ba337c44SJed Brown 1658ba337c44SJed Brown PetscFunctionBegin; 1659ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 1660ba337c44SJed Brown PetscFunctionReturn(0); 1661ba337c44SJed Brown } 1662ba337c44SJed Brown 1663ba337c44SJed Brown #undef __FUNCT__ 1664ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense" 1665ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A) 1666ba337c44SJed Brown { 1667ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1668ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1669ba337c44SJed Brown PetscScalar *aa = a->v; 1670ba337c44SJed Brown 1671ba337c44SJed Brown PetscFunctionBegin; 1672ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 1673ba337c44SJed Brown PetscFunctionReturn(0); 1674ba337c44SJed Brown } 1675284134d9SBarry Smith 1676a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/ 1677a9fe9ddaSSatish Balay #undef __FUNCT__ 1678a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1679a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1680a9fe9ddaSSatish Balay { 1681a9fe9ddaSSatish Balay PetscErrorCode ierr; 1682a9fe9ddaSSatish Balay 1683a9fe9ddaSSatish Balay PetscFunctionBegin; 1684a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1685a9fe9ddaSSatish Balay ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1686a9fe9ddaSSatish Balay } 1687a9fe9ddaSSatish Balay ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1688a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1689a9fe9ddaSSatish Balay } 1690a9fe9ddaSSatish Balay 1691a9fe9ddaSSatish Balay #undef __FUNCT__ 1692a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1693a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1694a9fe9ddaSSatish Balay { 1695ee16a9a1SHong Zhang PetscErrorCode ierr; 1696d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 1697ee16a9a1SHong Zhang Mat Cmat; 1698a9fe9ddaSSatish Balay 1699ee16a9a1SHong Zhang PetscFunctionBegin; 1700e32f2f54SBarry 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); 1701ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1702ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1703ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1704ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1705ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 17068cdbd757SHong Zhang Cmat->ops->matmult = MatMatMult_SeqDense_SeqDense; 1707ee16a9a1SHong Zhang *C = Cmat; 1708ee16a9a1SHong Zhang PetscFunctionReturn(0); 1709ee16a9a1SHong Zhang } 1710a9fe9ddaSSatish Balay 171198a3b096SSatish Balay #undef __FUNCT__ 1712a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1713a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1714a9fe9ddaSSatish Balay { 1715a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1716a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1717a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 17180805154bSBarry Smith PetscBLASInt m,n,k; 1719a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1720a9fe9ddaSSatish Balay 1721a9fe9ddaSSatish Balay PetscFunctionBegin; 1722d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 1723d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1724d0f46423SBarry Smith k = PetscBLASIntCast(A->cmap->n); 1725a9fe9ddaSSatish Balay BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1726a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1727a9fe9ddaSSatish Balay } 1728a9fe9ddaSSatish Balay 1729a9fe9ddaSSatish Balay #undef __FUNCT__ 173075648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqDense_SeqDense" 173175648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1732a9fe9ddaSSatish Balay { 1733a9fe9ddaSSatish Balay PetscErrorCode ierr; 1734a9fe9ddaSSatish Balay 1735a9fe9ddaSSatish Balay PetscFunctionBegin; 1736a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 173775648e8dSHong Zhang ierr = MatTransposeMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1738a9fe9ddaSSatish Balay } 173975648e8dSHong Zhang ierr = MatTransposeMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1740a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1741a9fe9ddaSSatish Balay } 1742a9fe9ddaSSatish Balay 1743a9fe9ddaSSatish Balay #undef __FUNCT__ 174475648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqDense_SeqDense" 174575648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1746a9fe9ddaSSatish Balay { 1747ee16a9a1SHong Zhang PetscErrorCode ierr; 1748d0f46423SBarry Smith PetscInt m=A->cmap->n,n=B->cmap->n; 1749ee16a9a1SHong Zhang Mat Cmat; 1750a9fe9ddaSSatish Balay 1751ee16a9a1SHong Zhang PetscFunctionBegin; 1752e32f2f54SBarry 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); 1753ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1754ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1755ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1756ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1757ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1758ee16a9a1SHong Zhang *C = Cmat; 1759ee16a9a1SHong Zhang PetscFunctionReturn(0); 1760ee16a9a1SHong Zhang } 1761a9fe9ddaSSatish Balay 1762a9fe9ddaSSatish Balay #undef __FUNCT__ 176375648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqDense_SeqDense" 176475648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1765a9fe9ddaSSatish Balay { 1766a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1767a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1768a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 17690805154bSBarry Smith PetscBLASInt m,n,k; 1770a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1771a9fe9ddaSSatish Balay 1772a9fe9ddaSSatish Balay PetscFunctionBegin; 1773d0f46423SBarry Smith m = PetscBLASIntCast(A->cmap->n); 1774d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1775d0f46423SBarry Smith k = PetscBLASIntCast(A->rmap->n); 17762fbe02b9SBarry Smith /* 17772fbe02b9SBarry Smith Note the m and n arguments below are the number rows and columns of A', not A! 17782fbe02b9SBarry Smith */ 1779a9fe9ddaSSatish Balay BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1780a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1781a9fe9ddaSSatish Balay } 1782985db425SBarry Smith 1783985db425SBarry Smith #undef __FUNCT__ 1784985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense" 1785985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 1786985db425SBarry Smith { 1787985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1788985db425SBarry Smith PetscErrorCode ierr; 1789d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1790985db425SBarry Smith PetscScalar *x; 1791985db425SBarry Smith MatScalar *aa = a->v; 1792985db425SBarry Smith 1793985db425SBarry Smith PetscFunctionBegin; 1794e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1795985db425SBarry Smith 1796985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1797985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1798985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1799e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1800985db425SBarry Smith for (i=0; i<m; i++) { 1801985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1802985db425SBarry Smith for (j=1; j<n; j++){ 1803985db425SBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1804985db425SBarry Smith } 1805985db425SBarry Smith } 1806985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1807985db425SBarry Smith PetscFunctionReturn(0); 1808985db425SBarry Smith } 1809985db425SBarry Smith 1810985db425SBarry Smith #undef __FUNCT__ 1811985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense" 1812985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 1813985db425SBarry Smith { 1814985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1815985db425SBarry Smith PetscErrorCode ierr; 1816d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1817985db425SBarry Smith PetscScalar *x; 1818985db425SBarry Smith PetscReal atmp; 1819985db425SBarry Smith MatScalar *aa = a->v; 1820985db425SBarry Smith 1821985db425SBarry Smith PetscFunctionBegin; 1822e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1823985db425SBarry Smith 1824985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1825985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1826985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1827e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1828985db425SBarry Smith for (i=0; i<m; i++) { 18299189402eSHong Zhang x[i] = PetscAbsScalar(aa[i]); 1830985db425SBarry Smith for (j=1; j<n; j++){ 1831985db425SBarry Smith atmp = PetscAbsScalar(aa[i+m*j]); 1832985db425SBarry Smith if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 1833985db425SBarry Smith } 1834985db425SBarry Smith } 1835985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1836985db425SBarry Smith PetscFunctionReturn(0); 1837985db425SBarry Smith } 1838985db425SBarry Smith 1839985db425SBarry Smith #undef __FUNCT__ 1840985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense" 1841985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 1842985db425SBarry Smith { 1843985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1844985db425SBarry Smith PetscErrorCode ierr; 1845d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1846985db425SBarry Smith PetscScalar *x; 1847985db425SBarry Smith MatScalar *aa = a->v; 1848985db425SBarry Smith 1849985db425SBarry Smith PetscFunctionBegin; 1850e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1851985db425SBarry Smith 1852985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1853985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1854985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1855e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1856985db425SBarry Smith for (i=0; i<m; i++) { 1857985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1858985db425SBarry Smith for (j=1; j<n; j++){ 1859985db425SBarry Smith if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1860985db425SBarry Smith } 1861985db425SBarry Smith } 1862985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1863985db425SBarry Smith PetscFunctionReturn(0); 1864985db425SBarry Smith } 1865985db425SBarry Smith 18668d0534beSBarry Smith #undef __FUNCT__ 18678d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense" 18688d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 18698d0534beSBarry Smith { 18708d0534beSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 18718d0534beSBarry Smith PetscErrorCode ierr; 18728d0534beSBarry Smith PetscScalar *x; 18738d0534beSBarry Smith 18748d0534beSBarry Smith PetscFunctionBegin; 1875e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 18768d0534beSBarry Smith 18778d0534beSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1878d0f46423SBarry Smith ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 18798d0534beSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 18808d0534beSBarry Smith PetscFunctionReturn(0); 18818d0534beSBarry Smith } 18828d0534beSBarry Smith 18830716a85fSBarry Smith 18840716a85fSBarry Smith #undef __FUNCT__ 18850716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqDense" 18860716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms) 18870716a85fSBarry Smith { 18880716a85fSBarry Smith PetscErrorCode ierr; 18890716a85fSBarry Smith PetscInt i,j,m,n; 18900716a85fSBarry Smith PetscScalar *a; 18910716a85fSBarry Smith 18920716a85fSBarry Smith PetscFunctionBegin; 18930716a85fSBarry Smith ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 18940716a85fSBarry Smith ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr); 18950716a85fSBarry Smith ierr = MatGetArray(A,&a);CHKERRQ(ierr); 18960716a85fSBarry Smith if (type == NORM_2) { 18970716a85fSBarry Smith for (i=0; i<n; i++ ){ 18980716a85fSBarry Smith for (j=0; j<m; j++) { 18990716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]*a[j]); 19000716a85fSBarry Smith } 19010716a85fSBarry Smith a += m; 19020716a85fSBarry Smith } 19030716a85fSBarry Smith } else if (type == NORM_1) { 19040716a85fSBarry Smith for (i=0; i<n; i++ ){ 19050716a85fSBarry Smith for (j=0; j<m; j++) { 19060716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]); 19070716a85fSBarry Smith } 19080716a85fSBarry Smith a += m; 19090716a85fSBarry Smith } 19100716a85fSBarry Smith } else if (type == NORM_INFINITY) { 19110716a85fSBarry Smith for (i=0; i<n; i++ ){ 19120716a85fSBarry Smith for (j=0; j<m; j++) { 19130716a85fSBarry Smith norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]); 19140716a85fSBarry Smith } 19150716a85fSBarry Smith a += m; 19160716a85fSBarry Smith } 19170716a85fSBarry Smith } else SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Unknown NormType"); 19180716a85fSBarry Smith if (type == NORM_2) { 19198f1a2a5eSBarry Smith for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]); 19200716a85fSBarry Smith } 19210716a85fSBarry Smith PetscFunctionReturn(0); 19220716a85fSBarry Smith } 19230716a85fSBarry Smith 1924289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1925a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1926905e6a2fSBarry Smith MatGetRow_SeqDense, 1927905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1928905e6a2fSBarry Smith MatMult_SeqDense, 192997304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 19307c922b88SBarry Smith MatMultTranspose_SeqDense, 19317c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1932db4efbfdSBarry Smith 0, 1933db4efbfdSBarry Smith 0, 1934db4efbfdSBarry Smith 0, 1935db4efbfdSBarry Smith /*10*/ 0, 1936905e6a2fSBarry Smith MatLUFactor_SeqDense, 1937905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 193841f059aeSBarry Smith MatSOR_SeqDense, 1939ec8511deSBarry Smith MatTranspose_SeqDense, 194097304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1941905e6a2fSBarry Smith MatEqual_SeqDense, 1942905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1943905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1944905e6a2fSBarry Smith MatNorm_SeqDense, 1945c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense, 1946c0aa2d19SHong Zhang MatAssemblyEnd_SeqDense, 1947905e6a2fSBarry Smith MatSetOption_SeqDense, 1948905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1949d519adbfSMatthew Knepley /*24*/ MatZeroRows_SeqDense, 1950db4efbfdSBarry Smith 0, 1951db4efbfdSBarry Smith 0, 1952db4efbfdSBarry Smith 0, 1953db4efbfdSBarry Smith 0, 1954*4994cf47SJed Brown /*29*/ MatSetUp_SeqDense, 1955273d9f13SBarry Smith 0, 1956905e6a2fSBarry Smith 0, 1957905e6a2fSBarry Smith MatGetArray_SeqDense, 1958905e6a2fSBarry Smith MatRestoreArray_SeqDense, 1959d519adbfSMatthew Knepley /*34*/ MatDuplicate_SeqDense, 1960a5ae1ecdSBarry Smith 0, 1961a5ae1ecdSBarry Smith 0, 1962a5ae1ecdSBarry Smith 0, 1963a5ae1ecdSBarry Smith 0, 1964d519adbfSMatthew Knepley /*39*/ MatAXPY_SeqDense, 1965a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1966a5ae1ecdSBarry Smith 0, 19674b0e389bSBarry Smith MatGetValues_SeqDense, 1968a5ae1ecdSBarry Smith MatCopy_SeqDense, 1969d519adbfSMatthew Knepley /*44*/ MatGetRowMax_SeqDense, 1970a5ae1ecdSBarry Smith MatScale_SeqDense, 1971a5ae1ecdSBarry Smith 0, 1972a5ae1ecdSBarry Smith 0, 1973a5ae1ecdSBarry Smith 0, 1974d519adbfSMatthew Knepley /*49*/ 0, 1975a5ae1ecdSBarry Smith 0, 1976a5ae1ecdSBarry Smith 0, 1977a5ae1ecdSBarry Smith 0, 1978a5ae1ecdSBarry Smith 0, 1979d519adbfSMatthew Knepley /*54*/ 0, 1980a5ae1ecdSBarry Smith 0, 1981a5ae1ecdSBarry Smith 0, 1982a5ae1ecdSBarry Smith 0, 1983a5ae1ecdSBarry Smith 0, 1984d519adbfSMatthew Knepley /*59*/ 0, 1985e03a110bSBarry Smith MatDestroy_SeqDense, 1986e03a110bSBarry Smith MatView_SeqDense, 1987357abbc8SBarry Smith 0, 198897304618SKris Buschelman 0, 1989d519adbfSMatthew Knepley /*64*/ 0, 199097304618SKris Buschelman 0, 199197304618SKris Buschelman 0, 199297304618SKris Buschelman 0, 199397304618SKris Buschelman 0, 1994d519adbfSMatthew Knepley /*69*/ MatGetRowMaxAbs_SeqDense, 199597304618SKris Buschelman 0, 199697304618SKris Buschelman 0, 199797304618SKris Buschelman 0, 199897304618SKris Buschelman 0, 1999d519adbfSMatthew Knepley /*74*/ 0, 200097304618SKris Buschelman 0, 200197304618SKris Buschelman 0, 200297304618SKris Buschelman 0, 200397304618SKris Buschelman 0, 2004d519adbfSMatthew Knepley /*79*/ 0, 200597304618SKris Buschelman 0, 200697304618SKris Buschelman 0, 200797304618SKris Buschelman 0, 20085bba2384SShri Abhyankar /*83*/ MatLoad_SeqDense, 2009865e5f61SKris Buschelman 0, 20101cbb95d3SBarry Smith MatIsHermitian_SeqDense, 2011865e5f61SKris Buschelman 0, 2012865e5f61SKris Buschelman 0, 2013865e5f61SKris Buschelman 0, 2014d519adbfSMatthew Knepley /*89*/ MatMatMult_SeqDense_SeqDense, 2015a9fe9ddaSSatish Balay MatMatMultSymbolic_SeqDense_SeqDense, 2016a9fe9ddaSSatish Balay MatMatMultNumeric_SeqDense_SeqDense, 2017865e5f61SKris Buschelman 0, 2018865e5f61SKris Buschelman 0, 2019d519adbfSMatthew Knepley /*94*/ 0, 20205df89d91SHong Zhang 0, 20215df89d91SHong Zhang 0, 20225df89d91SHong Zhang 0, 2023284134d9SBarry Smith 0, 2024d519adbfSMatthew Knepley /*99*/ 0, 2025284134d9SBarry Smith 0, 2026284134d9SBarry Smith 0, 2027ba337c44SJed Brown MatConjugate_SeqDense, 2028985db425SBarry Smith MatSetSizes_SeqDense, 2029ba337c44SJed Brown /*104*/0, 2030ba337c44SJed Brown MatRealPart_SeqDense, 2031ba337c44SJed Brown MatImaginaryPart_SeqDense, 2032985db425SBarry Smith 0, 2033985db425SBarry Smith 0, 203485e2c93fSHong Zhang /*109*/MatMatSolve_SeqDense, 2035985db425SBarry Smith 0, 20368d0534beSBarry Smith MatGetRowMin_SeqDense, 2037aabbc4fbSShri Abhyankar MatGetColumnVector_SeqDense, 2038aabbc4fbSShri Abhyankar 0, 2039aabbc4fbSShri Abhyankar /*114*/0, 2040aabbc4fbSShri Abhyankar 0, 2041aabbc4fbSShri Abhyankar 0, 2042aabbc4fbSShri Abhyankar 0, 2043aabbc4fbSShri Abhyankar 0, 2044aabbc4fbSShri Abhyankar /*119*/0, 2045aabbc4fbSShri Abhyankar 0, 2046aabbc4fbSShri Abhyankar 0, 20470716a85fSBarry Smith 0, 20480716a85fSBarry Smith 0, 20490716a85fSBarry Smith /*124*/0, 20505df89d91SHong Zhang MatGetColumnNorms_SeqDense, 20515df89d91SHong Zhang 0, 20525df89d91SHong Zhang 0, 20535df89d91SHong Zhang 0, 20545df89d91SHong Zhang /*129*/0, 205575648e8dSHong Zhang MatTransposeMatMult_SeqDense_SeqDense, 205675648e8dSHong Zhang MatTransposeMatMultSymbolic_SeqDense_SeqDense, 205775648e8dSHong Zhang MatTransposeMatMultNumeric_SeqDense_SeqDense, 2058985db425SBarry Smith }; 205990ace30eSBarry Smith 20604a2ae208SSatish Balay #undef __FUNCT__ 20614a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 20624b828684SBarry Smith /*@C 2063fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 2064d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 2065d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 2066289bc588SBarry Smith 2067db81eaa0SLois Curfman McInnes Collective on MPI_Comm 2068db81eaa0SLois Curfman McInnes 206920563c6bSBarry Smith Input Parameters: 2070db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 20710c775827SLois Curfman McInnes . m - number of rows 207218f449edSLois Curfman McInnes . n - number of columns 2073c0235b3cSMatthew Knepley - data - optional location of matrix data in column major order. Set data=PETSC_NULL for PETSc 2074dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 207520563c6bSBarry Smith 207620563c6bSBarry Smith Output Parameter: 207744cd7ae7SLois Curfman McInnes . A - the matrix 207820563c6bSBarry Smith 2079b259b22eSLois Curfman McInnes Notes: 208018f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 208118f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 2082b4fd4287SBarry Smith set data=PETSC_NULL. 208318f449edSLois Curfman McInnes 2084027ccd11SLois Curfman McInnes Level: intermediate 2085027ccd11SLois Curfman McInnes 2086dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 2087d65003e9SLois Curfman McInnes 2088db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 208920563c6bSBarry Smith @*/ 20907087cfbeSBarry Smith PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 2091289bc588SBarry Smith { 2092dfbe8321SBarry Smith PetscErrorCode ierr; 20933b2fbd54SBarry Smith 20943a40ed3dSBarry Smith PetscFunctionBegin; 2095f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 2096f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 2097273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 2098273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 2099273d9f13SBarry Smith PetscFunctionReturn(0); 2100273d9f13SBarry Smith } 2101273d9f13SBarry Smith 21024a2ae208SSatish Balay #undef __FUNCT__ 2103afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation" 2104273d9f13SBarry Smith /*@C 2105273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 2106273d9f13SBarry Smith 2107273d9f13SBarry Smith Collective on MPI_Comm 2108273d9f13SBarry Smith 2109273d9f13SBarry Smith Input Parameters: 2110273d9f13SBarry Smith + A - the matrix 2111273d9f13SBarry Smith - data - the array (or PETSC_NULL) 2112273d9f13SBarry Smith 2113273d9f13SBarry Smith Notes: 2114273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 2115273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 2116284134d9SBarry Smith need not call this routine. 2117273d9f13SBarry Smith 2118273d9f13SBarry Smith Level: intermediate 2119273d9f13SBarry Smith 2120273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 2121273d9f13SBarry Smith 2122867c911aSBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues(), MatSeqDenseSetLDA() 2123867c911aSBarry Smith 2124273d9f13SBarry Smith @*/ 21257087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 2126273d9f13SBarry Smith { 21274ac538c5SBarry Smith PetscErrorCode ierr; 2128a23d5eceSKris Buschelman 2129a23d5eceSKris Buschelman PetscFunctionBegin; 21304ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));CHKERRQ(ierr); 2131a23d5eceSKris Buschelman PetscFunctionReturn(0); 2132a23d5eceSKris Buschelman } 2133a23d5eceSKris Buschelman 2134a23d5eceSKris Buschelman EXTERN_C_BEGIN 2135a23d5eceSKris Buschelman #undef __FUNCT__ 2136afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense" 21377087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 2138a23d5eceSKris Buschelman { 2139273d9f13SBarry Smith Mat_SeqDense *b; 2140dfbe8321SBarry Smith PetscErrorCode ierr; 2141273d9f13SBarry Smith 2142273d9f13SBarry Smith PetscFunctionBegin; 2143273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 2144a868139aSShri Abhyankar 214534ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(B->rmap,1);CHKERRQ(ierr); 214634ef9618SShri Abhyankar ierr = PetscLayoutSetBlockSize(B->cmap,1);CHKERRQ(ierr); 214734ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 214834ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 214934ef9618SShri Abhyankar 2150273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 215186d161a7SShri Abhyankar b->Mmax = B->rmap->n; 215286d161a7SShri Abhyankar b->Nmax = B->cmap->n; 215386d161a7SShri Abhyankar if(b->lda <= 0 || b->changelda) b->lda = B->rmap->n; 215486d161a7SShri Abhyankar 21559e8f95c4SLisandro Dalcin if (!data) { /* petsc-allocated storage */ 21569e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 21575afd5e0cSBarry Smith ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 2158284134d9SBarry Smith ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 2159284134d9SBarry Smith ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 21609e8f95c4SLisandro Dalcin b->user_alloc = PETSC_FALSE; 2161273d9f13SBarry Smith } else { /* user-allocated storage */ 21629e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 2163273d9f13SBarry Smith b->v = data; 2164273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 2165273d9f13SBarry Smith } 21660450473dSBarry Smith B->assembled = PETSC_TRUE; 2167273d9f13SBarry Smith PetscFunctionReturn(0); 2168273d9f13SBarry Smith } 2169a23d5eceSKris Buschelman EXTERN_C_END 2170273d9f13SBarry Smith 21711b807ce4Svictorle #undef __FUNCT__ 21721b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 21731b807ce4Svictorle /*@C 21741b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 21751b807ce4Svictorle 21761b807ce4Svictorle Input parameter: 21771b807ce4Svictorle + A - the matrix 21781b807ce4Svictorle - lda - the leading dimension 21791b807ce4Svictorle 21801b807ce4Svictorle Notes: 2181867c911aSBarry Smith This routine is to be used in conjunction with MatSeqDenseSetPreallocation(); 21821b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 21831b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 21841b807ce4Svictorle 21851b807ce4Svictorle Level: intermediate 21861b807ce4Svictorle 21871b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 21881b807ce4Svictorle 2189284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 2190867c911aSBarry Smith 21911b807ce4Svictorle @*/ 21927087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetLDA(Mat B,PetscInt lda) 21931b807ce4Svictorle { 21941b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 219521a2c019SBarry Smith 21961b807ce4Svictorle PetscFunctionBegin; 2197e32f2f54SBarry 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); 21981b807ce4Svictorle b->lda = lda; 219921a2c019SBarry Smith b->changelda = PETSC_FALSE; 220021a2c019SBarry Smith b->Mmax = PetscMax(b->Mmax,lda); 22011b807ce4Svictorle PetscFunctionReturn(0); 22021b807ce4Svictorle } 22031b807ce4Svictorle 22040bad9183SKris Buschelman /*MC 2205fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 22060bad9183SKris Buschelman 22070bad9183SKris Buschelman Options Database Keys: 22080bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 22090bad9183SKris Buschelman 22100bad9183SKris Buschelman Level: beginner 22110bad9183SKris Buschelman 221289665df3SBarry Smith .seealso: MatCreateSeqDense() 221389665df3SBarry Smith 22140bad9183SKris Buschelman M*/ 22150bad9183SKris Buschelman 2216273d9f13SBarry Smith EXTERN_C_BEGIN 22174a2ae208SSatish Balay #undef __FUNCT__ 22184a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 22197087cfbeSBarry Smith PetscErrorCode MatCreate_SeqDense(Mat B) 2220273d9f13SBarry Smith { 2221273d9f13SBarry Smith Mat_SeqDense *b; 2222dfbe8321SBarry Smith PetscErrorCode ierr; 22237c334f02SBarry Smith PetscMPIInt size; 2224273d9f13SBarry Smith 2225273d9f13SBarry Smith PetscFunctionBegin; 22267adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr); 2227e32f2f54SBarry Smith if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 222855659b69SBarry Smith 222938f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr); 2230549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 223144cd7ae7SLois Curfman McInnes B->data = (void*)b; 223218f449edSLois Curfman McInnes 223344cd7ae7SLois Curfman McInnes b->pivots = 0; 2234273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 2235273d9f13SBarry Smith b->v = 0; 223621a2c019SBarry Smith b->changelda = PETSC_FALSE; 22374e220ebcSLois Curfman McInnes 2238b24902e0SBarry Smith 22396a63e612SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatConvert_seqdense_seqaij_C","MatConvert_SeqDense_SeqAIJ",MatConvert_SeqDense_SeqAIJ);CHKERRQ(ierr); 2240ec1065edSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_petsc_C", 2241b24902e0SBarry Smith "MatGetFactor_seqdense_petsc", 2242b24902e0SBarry Smith MatGetFactor_seqdense_petsc);CHKERRQ(ierr); 2243a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 2244a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 2245a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 22462356f84bSDmitry Karpeev 22474ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C", 22484ae313f4SHong Zhang "MatMatMult_SeqAIJ_SeqDense", 22494ae313f4SHong Zhang MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 22502356f84bSDmitry Karpeev 2251c0c8ee5eSDmitry Karpeev 22524ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C", 22534ae313f4SHong Zhang "MatMatMultSymbolic_SeqAIJ_SeqDense", 22544ae313f4SHong Zhang MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 22554ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C", 22564ae313f4SHong Zhang "MatMatMultNumeric_SeqAIJ_SeqDense", 22574ae313f4SHong Zhang MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 225817667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 22593a40ed3dSBarry Smith PetscFunctionReturn(0); 2260289bc588SBarry Smith } 2261273d9f13SBarry Smith EXTERN_C_END 2262