1be1d678aSKris Buschelman #define PETSCMAT_DLL 2be1d678aSKris Buschelman 367e560aaSBarry Smith /* 467e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 567e560aaSBarry Smith */ 6289bc588SBarry Smith 770f55243SBarry Smith #include "src/mat/impls/dense/seq/dense.h" 8b4c862e0SSatish Balay #include "petscblaslapack.h" 9289bc588SBarry Smith 104a2ae208SSatish Balay #undef __FUNCT__ 114a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense" 12f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 131987afe7SBarry Smith { 141987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 15f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 1613f74950SBarry Smith PetscInt j; 170805154bSBarry Smith PetscBLASInt N,m,ldax,lday,one = 1; 18efee365bSSatish Balay PetscErrorCode ierr; 193a40ed3dSBarry Smith 203a40ed3dSBarry Smith PetscFunctionBegin; 21d0f46423SBarry Smith N = PetscBLASIntCast(X->rmap->n*X->cmap->n); 22d0f46423SBarry Smith m = PetscBLASIntCast(X->rmap->n); 230805154bSBarry Smith ldax = PetscBLASIntCast(x->lda); 240805154bSBarry Smith lday = PetscBLASIntCast(y->lda); 25a5ce6ee0Svictorle if (ldax>m || lday>m) { 26d0f46423SBarry Smith for (j=0; j<X->cmap->n; j++) { 27f4df32b1SMatthew Knepley BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one); 28a5ce6ee0Svictorle } 29a5ce6ee0Svictorle } else { 30f4df32b1SMatthew Knepley BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one); 31a5ce6ee0Svictorle } 32efee365bSSatish Balay ierr = PetscLogFlops(2*N-1);CHKERRQ(ierr); 333a40ed3dSBarry Smith PetscFunctionReturn(0); 341987afe7SBarry Smith } 351987afe7SBarry Smith 364a2ae208SSatish Balay #undef __FUNCT__ 374a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 38dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 39289bc588SBarry Smith { 40d0f46423SBarry Smith PetscInt N = A->rmap->n*A->cmap->n; 413a40ed3dSBarry Smith 423a40ed3dSBarry Smith PetscFunctionBegin; 43d0f46423SBarry Smith info->rows_global = (double)A->rmap->n; 44d0f46423SBarry Smith info->columns_global = (double)A->cmap->n; 45d0f46423SBarry Smith info->rows_local = (double)A->rmap->n; 46d0f46423SBarry Smith info->columns_local = (double)A->cmap->n; 474e220ebcSLois Curfman McInnes info->block_size = 1.0; 484e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 496de62eeeSBarry Smith info->nz_used = (double)N; 506de62eeeSBarry Smith info->nz_unneeded = (double)0; 514e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 524e220ebcSLois Curfman McInnes info->mallocs = 0; 537adad957SLisandro Dalcin info->memory = ((PetscObject)A)->mem; 544e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 554e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 564e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 573a40ed3dSBarry Smith PetscFunctionReturn(0); 58289bc588SBarry Smith } 59289bc588SBarry Smith 604a2ae208SSatish Balay #undef __FUNCT__ 614a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense" 62f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha) 6380cd9d93SLois Curfman McInnes { 64273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 65f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 66efee365bSSatish Balay PetscErrorCode ierr; 670805154bSBarry Smith PetscBLASInt one = 1,j,nz,lda = PetscBLASIntCast(a->lda); 6880cd9d93SLois Curfman McInnes 693a40ed3dSBarry Smith PetscFunctionBegin; 70d0f46423SBarry Smith if (lda>A->rmap->n) { 71d0f46423SBarry Smith nz = PetscBLASIntCast(A->rmap->n); 72d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 73f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v+j*lda,&one); 74a5ce6ee0Svictorle } 75a5ce6ee0Svictorle } else { 76d0f46423SBarry Smith nz = PetscBLASIntCast(A->rmap->n*A->cmap->n); 77f4df32b1SMatthew Knepley BLASscal_(&nz,&oalpha,a->v,&one); 78a5ce6ee0Svictorle } 79efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 803a40ed3dSBarry Smith PetscFunctionReturn(0); 8180cd9d93SLois Curfman McInnes } 8280cd9d93SLois Curfman McInnes 831cbb95d3SBarry Smith #undef __FUNCT__ 841cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense" 851cbb95d3SBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscTruth *fl) 861cbb95d3SBarry Smith { 871cbb95d3SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 88d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,N; 891cbb95d3SBarry Smith PetscScalar *v = a->v; 901cbb95d3SBarry Smith 911cbb95d3SBarry Smith PetscFunctionBegin; 921cbb95d3SBarry Smith *fl = PETSC_FALSE; 93d0f46423SBarry Smith if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0); 941cbb95d3SBarry Smith N = a->lda; 951cbb95d3SBarry Smith 961cbb95d3SBarry Smith for (i=0; i<m; i++) { 971cbb95d3SBarry Smith for (j=i+1; j<m; j++) { 981cbb95d3SBarry Smith if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0); 991cbb95d3SBarry Smith } 1001cbb95d3SBarry Smith } 1011cbb95d3SBarry Smith *fl = PETSC_TRUE; 1021cbb95d3SBarry Smith PetscFunctionReturn(0); 1031cbb95d3SBarry Smith } 1041cbb95d3SBarry Smith 1056ee01492SSatish Balay 106b24902e0SBarry Smith 107b24902e0SBarry Smith #undef __FUNCT__ 108b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense" 109*719d5645SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues) 110b24902e0SBarry Smith { 111b24902e0SBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 112b24902e0SBarry Smith PetscErrorCode ierr; 113b24902e0SBarry Smith PetscInt lda = (PetscInt)mat->lda,j,m; 114b24902e0SBarry Smith 115b24902e0SBarry Smith PetscFunctionBegin; 116*719d5645SBarry Smith ierr = MatSeqDenseSetPreallocation(newi,PETSC_NULL);CHKERRQ(ierr); 117b24902e0SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 118b24902e0SBarry Smith l = (Mat_SeqDense*)newi->data; 119d0f46423SBarry Smith if (lda>A->rmap->n) { 120d0f46423SBarry Smith m = A->rmap->n; 121d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 122b24902e0SBarry Smith ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 123b24902e0SBarry Smith } 124b24902e0SBarry Smith } else { 125d0f46423SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 126b24902e0SBarry Smith } 127b24902e0SBarry Smith } 128b24902e0SBarry Smith newi->assembled = PETSC_TRUE; 129b24902e0SBarry Smith PetscFunctionReturn(0); 130b24902e0SBarry Smith } 131b24902e0SBarry Smith 1324a2ae208SSatish Balay #undef __FUNCT__ 1334a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 134dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 13502cad45dSBarry Smith { 1366849ba73SBarry Smith PetscErrorCode ierr; 13702cad45dSBarry Smith 1383a40ed3dSBarry Smith PetscFunctionBegin; 1395c9eb25fSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,newmat);CHKERRQ(ierr); 140d0f46423SBarry Smith ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 1415c9eb25fSBarry Smith ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 142*719d5645SBarry Smith ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr); 143b24902e0SBarry Smith PetscFunctionReturn(0); 144b24902e0SBarry Smith } 145b24902e0SBarry Smith 1466ee01492SSatish Balay 147*719d5645SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,MatFactorInfo*); 148*719d5645SBarry Smith 1494a2ae208SSatish Balay #undef __FUNCT__ 1504a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 151*719d5645SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,MatFactorInfo *info_dummy) 152289bc588SBarry Smith { 153*719d5645SBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l = (Mat_SeqDense*)fact->data; 1546849ba73SBarry Smith PetscErrorCode ierr; 155d0f46423SBarry Smith PetscInt lda1=mat->lda,lda2=l->lda, m=A->rmap->n,n=A->cmap->n, j; 1564482741eSBarry Smith MatFactorInfo info; 1573a40ed3dSBarry Smith 1583a40ed3dSBarry Smith PetscFunctionBegin; 15902cad45dSBarry Smith /* copy the numerical values */ 1601b807ce4Svictorle if (lda1>m || lda2>m ) { 1611b807ce4Svictorle for (j=0; j<n; j++) { 1621b807ce4Svictorle ierr = PetscMemcpy(l->v+j*lda2,mat->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1631b807ce4Svictorle } 1641b807ce4Svictorle } else { 165d0f46423SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1661b807ce4Svictorle } 167*719d5645SBarry Smith ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr); 1683a40ed3dSBarry Smith PetscFunctionReturn(0); 169289bc588SBarry Smith } 1706ee01492SSatish Balay 1710b4b3355SBarry Smith 1720b4b3355SBarry Smith #undef __FUNCT__ 1734a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 174dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 175289bc588SBarry Smith { 176c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1776849ba73SBarry Smith PetscErrorCode ierr; 17887828ca2SBarry Smith PetscScalar *x,*y; 179d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 18067e560aaSBarry Smith 1813a40ed3dSBarry Smith PetscFunctionBegin; 1821ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1831ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 184d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1855c9eb25fSBarry Smith if (A->factor == MAT_FACTOR_LU) { 186ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 18780444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 188ae7cfcebSSatish Balay #else 18971044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 1904ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"GETRS - Bad solve"); 191ae7cfcebSSatish Balay #endif 1925c9eb25fSBarry Smith } else if (A->factor == MAT_FACTOR_CHOLESKY){ 193ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 19480444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 195ae7cfcebSSatish Balay #else 19671044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 1974ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS Bad solve"); 198ae7cfcebSSatish Balay #endif 199289bc588SBarry Smith } 20029bbc08cSBarry Smith else SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 2011ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2021ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 203d0f46423SBarry Smith ierr = PetscLogFlops(2*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 2043a40ed3dSBarry Smith PetscFunctionReturn(0); 205289bc588SBarry Smith } 2066ee01492SSatish Balay 2074a2ae208SSatish Balay #undef __FUNCT__ 2084a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 209dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 210da3a660dSBarry Smith { 211c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 212dfbe8321SBarry Smith PetscErrorCode ierr; 21387828ca2SBarry Smith PetscScalar *x,*y; 214d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 21567e560aaSBarry Smith 2163a40ed3dSBarry Smith PetscFunctionBegin; 2171ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2181ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 219d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 220752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 221da3a660dSBarry Smith if (mat->pivots) { 222ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 22380444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 224ae7cfcebSSatish Balay #else 22571044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 2264ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve"); 227ae7cfcebSSatish Balay #endif 2287a97a34bSBarry Smith } else { 229ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 23080444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 231ae7cfcebSSatish Balay #else 23271044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 2334ce68768SBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"POTRS - Bad solve"); 234ae7cfcebSSatish Balay #endif 235da3a660dSBarry Smith } 2361ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2371ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 238d0f46423SBarry Smith ierr = PetscLogFlops(2*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 2393a40ed3dSBarry Smith PetscFunctionReturn(0); 240da3a660dSBarry Smith } 2416ee01492SSatish Balay 2424a2ae208SSatish Balay #undef __FUNCT__ 2434a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 244dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 245da3a660dSBarry Smith { 246c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 247dfbe8321SBarry Smith PetscErrorCode ierr; 24887828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 249da3a660dSBarry Smith Vec tmp = 0; 250d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 25167e560aaSBarry Smith 2523a40ed3dSBarry Smith PetscFunctionBegin; 2531ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2541ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 255d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 256da3a660dSBarry Smith if (yy == zz) { 25778b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 25852e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 25978b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 260da3a660dSBarry Smith } 261d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 262752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 263da3a660dSBarry Smith if (mat->pivots) { 264ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 26580444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 266ae7cfcebSSatish Balay #else 26771044d3cSBarry Smith LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 2682ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 269ae7cfcebSSatish Balay #endif 270a8c6a408SBarry Smith } else { 271ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 27280444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 273ae7cfcebSSatish Balay #else 27471044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 2752ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 276ae7cfcebSSatish Balay #endif 277da3a660dSBarry Smith } 2782dcb1b2aSMatthew Knepley if (tmp) {ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); ierr = VecDestroy(tmp);CHKERRQ(ierr);} 2792dcb1b2aSMatthew Knepley else {ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr);} 2801ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2811ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 282d0f46423SBarry Smith ierr = PetscLogFlops(2*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 2833a40ed3dSBarry Smith PetscFunctionReturn(0); 284da3a660dSBarry Smith } 28567e560aaSBarry Smith 2864a2ae208SSatish Balay #undef __FUNCT__ 2874a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 288dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 289da3a660dSBarry Smith { 290c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2916849ba73SBarry Smith PetscErrorCode ierr; 29287828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 293da3a660dSBarry Smith Vec tmp; 294d0f46423SBarry Smith PetscBLASInt one = 1,info,m = PetscBLASIntCast(A->rmap->n); 29567e560aaSBarry Smith 2963a40ed3dSBarry Smith PetscFunctionBegin; 297d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 2981ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2991ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 300da3a660dSBarry Smith if (yy == zz) { 30178b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 30252e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 30378b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 304da3a660dSBarry Smith } 305d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 306752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 307da3a660dSBarry Smith if (mat->pivots) { 308ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 30980444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 310ae7cfcebSSatish Balay #else 31171044d3cSBarry Smith LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info); 3122ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 313ae7cfcebSSatish Balay #endif 3143a40ed3dSBarry Smith } else { 315ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 31680444007SBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 317ae7cfcebSSatish Balay #else 31871044d3cSBarry Smith LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info); 3192ffef20aSBarry Smith if (info) SETERRQ(PETSC_ERR_LIB,"Bad solve"); 320ae7cfcebSSatish Balay #endif 321da3a660dSBarry Smith } 32290f02eecSBarry Smith if (tmp) { 3232dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 32490f02eecSBarry Smith ierr = VecDestroy(tmp);CHKERRQ(ierr); 3253a40ed3dSBarry Smith } else { 3262dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 32790f02eecSBarry Smith } 3281ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3291ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 330d0f46423SBarry Smith ierr = PetscLogFlops(2*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 3313a40ed3dSBarry Smith PetscFunctionReturn(0); 332da3a660dSBarry Smith } 333db4efbfdSBarry Smith 334db4efbfdSBarry Smith /* ---------------------------------------------------------------*/ 335db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 336db4efbfdSBarry Smith rather than put it in the Mat->row slot.*/ 337db4efbfdSBarry Smith #undef __FUNCT__ 338db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense" 339db4efbfdSBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,MatFactorInfo *minfo) 340db4efbfdSBarry Smith { 341db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 342db4efbfdSBarry Smith PetscFunctionBegin; 343db4efbfdSBarry Smith SETERRQ(PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 344db4efbfdSBarry Smith #else 345db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 346db4efbfdSBarry Smith PetscErrorCode ierr; 347db4efbfdSBarry Smith PetscBLASInt n,m,info; 348db4efbfdSBarry Smith 349db4efbfdSBarry Smith PetscFunctionBegin; 350db4efbfdSBarry Smith n = PetscBLASIntCast(A->cmap->n); 351db4efbfdSBarry Smith m = PetscBLASIntCast(A->rmap->n); 352db4efbfdSBarry Smith if (!mat->pivots) { 353db4efbfdSBarry Smith ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr); 354db4efbfdSBarry Smith ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr); 355db4efbfdSBarry Smith } 356db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 357db4efbfdSBarry Smith LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info); 358db4efbfdSBarry Smith if (info<0) SETERRQ(PETSC_ERR_LIB,"Bad argument to LU factorization"); 359db4efbfdSBarry Smith if (info>0) SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 360db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 361db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 362db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 363db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 364db4efbfdSBarry Smith A->factor = MAT_FACTOR_LU; 365db4efbfdSBarry Smith 366db4efbfdSBarry Smith ierr = PetscLogFlops((2*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr); 367db4efbfdSBarry Smith #endif 368db4efbfdSBarry Smith PetscFunctionReturn(0); 369db4efbfdSBarry Smith } 370db4efbfdSBarry Smith 371db4efbfdSBarry Smith #undef __FUNCT__ 372db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense" 373db4efbfdSBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,MatFactorInfo *factinfo) 374db4efbfdSBarry Smith { 375db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 376db4efbfdSBarry Smith PetscFunctionBegin; 377db4efbfdSBarry Smith SETERRQ(PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 378db4efbfdSBarry Smith #else 379db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 380db4efbfdSBarry Smith PetscErrorCode ierr; 381db4efbfdSBarry Smith PetscBLASInt info,n = PetscBLASIntCast(A->cmap->n); 382db4efbfdSBarry Smith 383db4efbfdSBarry Smith PetscFunctionBegin; 384db4efbfdSBarry Smith ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 385db4efbfdSBarry Smith mat->pivots = 0; 386db4efbfdSBarry Smith 387db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 388db4efbfdSBarry Smith LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info); 389db4efbfdSBarry Smith if (info) SETERRQ1(PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 390db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 391db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 392db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 393db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 394db4efbfdSBarry Smith A->factor = MAT_FACTOR_CHOLESKY; 395db4efbfdSBarry Smith ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr); 396db4efbfdSBarry Smith #endif 397db4efbfdSBarry Smith PetscFunctionReturn(0); 398db4efbfdSBarry Smith } 399db4efbfdSBarry Smith 400db4efbfdSBarry Smith 401db4efbfdSBarry Smith #undef __FUNCT__ 402db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 403*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,MatFactorInfo *info_dummy) 404db4efbfdSBarry Smith { 405db4efbfdSBarry Smith PetscErrorCode ierr; 406db4efbfdSBarry Smith MatFactorInfo info; 407db4efbfdSBarry Smith 408db4efbfdSBarry Smith PetscFunctionBegin; 409db4efbfdSBarry Smith info.fill = 1.0; 410*719d5645SBarry Smith ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr); 411db4efbfdSBarry Smith PetscFunctionReturn(0); 412db4efbfdSBarry Smith } 413db4efbfdSBarry Smith 414db4efbfdSBarry Smith #undef __FUNCT__ 415db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 416*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,MatFactorInfo *info) 417db4efbfdSBarry Smith { 418db4efbfdSBarry Smith PetscErrorCode ierr; 419db4efbfdSBarry Smith 420db4efbfdSBarry Smith PetscFunctionBegin; 421*719d5645SBarry Smith ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 422*719d5645SBarry Smith fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense; 423db4efbfdSBarry Smith PetscFunctionReturn(0); 424db4efbfdSBarry Smith } 425db4efbfdSBarry Smith 426db4efbfdSBarry Smith #undef __FUNCT__ 427db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 428*719d5645SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,MatFactorInfo *info) 429db4efbfdSBarry Smith { 430db4efbfdSBarry Smith PetscErrorCode ierr; 431db4efbfdSBarry Smith 432db4efbfdSBarry Smith PetscFunctionBegin; 433*719d5645SBarry Smith ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_DO_NOT_COPY_VALUES);CHKERRQ(ierr); 434*719d5645SBarry Smith fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense; 435db4efbfdSBarry Smith PetscFunctionReturn(0); 436db4efbfdSBarry Smith } 437db4efbfdSBarry Smith 438db4efbfdSBarry Smith #undef __FUNCT__ 439db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc" 440db4efbfdSBarry Smith PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact) 441db4efbfdSBarry Smith { 442db4efbfdSBarry Smith PetscErrorCode ierr; 443db4efbfdSBarry Smith 444db4efbfdSBarry Smith PetscFunctionBegin; 445db4efbfdSBarry Smith ierr = MatCreate(((PetscObject)A)->comm,fact);CHKERRQ(ierr); 446db4efbfdSBarry Smith ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 447db4efbfdSBarry Smith ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr); 448db4efbfdSBarry Smith if (ftype == MAT_FACTOR_LU){ 449db4efbfdSBarry Smith (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense; 450db4efbfdSBarry Smith } else { 451db4efbfdSBarry Smith (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense; 452db4efbfdSBarry Smith } 453db4efbfdSBarry Smith (*fact)->factor = ftype; 454db4efbfdSBarry Smith PetscFunctionReturn(0); 455db4efbfdSBarry Smith } 456db4efbfdSBarry Smith 457289bc588SBarry Smith /* ------------------------------------------------------------------*/ 4584a2ae208SSatish Balay #undef __FUNCT__ 4594a2ae208SSatish Balay #define __FUNCT__ "MatRelax_SeqDense" 46013f74950SBarry Smith PetscErrorCode MatRelax_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 461289bc588SBarry Smith { 462c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 46387828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 464dfbe8321SBarry Smith PetscErrorCode ierr; 465d0f46423SBarry Smith PetscInt m = A->rmap->n,i; 466aa482453SBarry Smith #if !defined(PETSC_USE_COMPLEX) 4670805154bSBarry Smith PetscBLASInt o = 1,bm = PetscBLASIntCast(m); 468bc1b551cSSatish Balay #endif 469289bc588SBarry Smith 4703a40ed3dSBarry Smith PetscFunctionBegin; 471289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 47271044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 4732dcb1b2aSMatthew Knepley ierr = VecSet(xx,zero);CHKERRQ(ierr); 474289bc588SBarry Smith } 4751ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 4761ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 477b965ef7fSBarry Smith its = its*lits; 47877431f27SBarry Smith if (its <= 0) SETERRQ2(PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 479289bc588SBarry Smith while (its--) { 480fccaa45eSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP){ 481289bc588SBarry Smith for (i=0; i<m; i++) { 482aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 483f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 484f1747703SBarry Smith not happy about returning a double complex */ 48513f74950SBarry Smith PetscInt _i; 48687828ca2SBarry Smith PetscScalar sum = b[i]; 487f1747703SBarry Smith for (_i=0; _i<m; _i++) { 4883f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 489f1747703SBarry Smith } 490f1747703SBarry Smith xt = sum; 491f1747703SBarry Smith #else 49271044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 493f1747703SBarry Smith #endif 49455a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 495289bc588SBarry Smith } 496289bc588SBarry Smith } 497fccaa45eSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP){ 498289bc588SBarry Smith for (i=m-1; i>=0; i--) { 499aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 500f1747703SBarry Smith /* cannot use BLAS dot for complex because compiler/linker is 501f1747703SBarry Smith not happy about returning a double complex */ 50213f74950SBarry Smith PetscInt _i; 50387828ca2SBarry Smith PetscScalar sum = b[i]; 504f1747703SBarry Smith for (_i=0; _i<m; _i++) { 5053f6de6efSSatish Balay sum -= PetscConj(v[i+_i*m])*x[_i]; 506f1747703SBarry Smith } 507f1747703SBarry Smith xt = sum; 508f1747703SBarry Smith #else 50971044d3cSBarry Smith xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o); 510f1747703SBarry Smith #endif 51155a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 512289bc588SBarry Smith } 513289bc588SBarry Smith } 514289bc588SBarry Smith } 5151ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 5161ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5173a40ed3dSBarry Smith PetscFunctionReturn(0); 518289bc588SBarry Smith } 519289bc588SBarry Smith 520289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5214a2ae208SSatish Balay #undef __FUNCT__ 5224a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 523dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 524289bc588SBarry Smith { 525c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 52687828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 527dfbe8321SBarry Smith PetscErrorCode ierr; 5280805154bSBarry Smith PetscBLASInt m, n,_One=1; 529ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 5303a40ed3dSBarry Smith 5313a40ed3dSBarry Smith PetscFunctionBegin; 532d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 533d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 534d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5351ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5361ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 53771044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One); 5381ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5391ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 540d0f46423SBarry Smith ierr = PetscLogFlops(2*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 5413a40ed3dSBarry Smith PetscFunctionReturn(0); 542289bc588SBarry Smith } 5436ee01492SSatish Balay 5444a2ae208SSatish Balay #undef __FUNCT__ 5454a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 546dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 547289bc588SBarry Smith { 548c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 54987828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 550dfbe8321SBarry Smith PetscErrorCode ierr; 5510805154bSBarry Smith PetscBLASInt m, n, _One=1; 5523a40ed3dSBarry Smith 5533a40ed3dSBarry Smith PetscFunctionBegin; 554d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 555d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 556d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5571ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5581ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 55971044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One); 5601ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5611ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 562d0f46423SBarry Smith ierr = PetscLogFlops(2*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr); 5633a40ed3dSBarry Smith PetscFunctionReturn(0); 564289bc588SBarry Smith } 5656ee01492SSatish Balay 5664a2ae208SSatish Balay #undef __FUNCT__ 5674a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 568dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 569289bc588SBarry Smith { 570c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 57187828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 572dfbe8321SBarry Smith PetscErrorCode ierr; 5730805154bSBarry Smith PetscBLASInt m, n, _One=1; 5743a40ed3dSBarry Smith 5753a40ed3dSBarry Smith PetscFunctionBegin; 576d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 577d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 578d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 579600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 5801ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5811ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 58271044d3cSBarry Smith BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 5831ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5841ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 585d0f46423SBarry Smith ierr = PetscLogFlops(2*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 5863a40ed3dSBarry Smith PetscFunctionReturn(0); 587289bc588SBarry Smith } 5886ee01492SSatish Balay 5894a2ae208SSatish Balay #undef __FUNCT__ 5904a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 591dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 592289bc588SBarry Smith { 593c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 59487828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 595dfbe8321SBarry Smith PetscErrorCode ierr; 5960805154bSBarry Smith PetscBLASInt m, n, _One=1; 59787828ca2SBarry Smith PetscScalar _DOne=1.0; 5983a40ed3dSBarry Smith 5993a40ed3dSBarry Smith PetscFunctionBegin; 600d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 601d0f46423SBarry Smith n = PetscBLASIntCast(A->cmap->n); 602d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 603600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6041ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6051ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 60671044d3cSBarry Smith BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One); 6071ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6081ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 609d0f46423SBarry Smith ierr = PetscLogFlops(2*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6103a40ed3dSBarry Smith PetscFunctionReturn(0); 611289bc588SBarry Smith } 612289bc588SBarry Smith 613289bc588SBarry Smith /* -----------------------------------------------------------------*/ 6144a2ae208SSatish Balay #undef __FUNCT__ 6154a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 61613f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 617289bc588SBarry Smith { 618c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 61987828ca2SBarry Smith PetscScalar *v; 6206849ba73SBarry Smith PetscErrorCode ierr; 62113f74950SBarry Smith PetscInt i; 62267e560aaSBarry Smith 6233a40ed3dSBarry Smith PetscFunctionBegin; 624d0f46423SBarry Smith *ncols = A->cmap->n; 625289bc588SBarry Smith if (cols) { 626d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr); 627d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) (*cols)[i] = i; 628289bc588SBarry Smith } 629289bc588SBarry Smith if (vals) { 630d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 631289bc588SBarry Smith v = mat->v + row; 632d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;} 633289bc588SBarry Smith } 6343a40ed3dSBarry Smith PetscFunctionReturn(0); 635289bc588SBarry Smith } 6366ee01492SSatish Balay 6374a2ae208SSatish Balay #undef __FUNCT__ 6384a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 63913f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 640289bc588SBarry Smith { 641dfbe8321SBarry Smith PetscErrorCode ierr; 642606d414cSSatish Balay PetscFunctionBegin; 643606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 644606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 6453a40ed3dSBarry Smith PetscFunctionReturn(0); 646289bc588SBarry Smith } 647289bc588SBarry Smith /* ----------------------------------------------------------------*/ 6484a2ae208SSatish Balay #undef __FUNCT__ 6494a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 65013f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 651289bc588SBarry Smith { 652c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 65313f74950SBarry Smith PetscInt i,j,idx=0; 654d6dfbf8fSBarry Smith 6553a40ed3dSBarry Smith PetscFunctionBegin; 656289bc588SBarry Smith if (!mat->roworiented) { 657dbb450caSBarry Smith if (addv == INSERT_VALUES) { 658289bc588SBarry Smith for (j=0; j<n; j++) { 659cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6602515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 661d0f46423SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 66258804f6dSBarry Smith #endif 663289bc588SBarry Smith for (i=0; i<m; i++) { 664cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6652515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 666d0f46423SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 66758804f6dSBarry Smith #endif 668cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 669289bc588SBarry Smith } 670289bc588SBarry Smith } 6713a40ed3dSBarry Smith } else { 672289bc588SBarry Smith for (j=0; j<n; j++) { 673cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 6742515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 675d0f46423SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 67658804f6dSBarry Smith #endif 677289bc588SBarry Smith for (i=0; i<m; i++) { 678cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 6792515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 680d0f46423SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 68158804f6dSBarry Smith #endif 682cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 683289bc588SBarry Smith } 684289bc588SBarry Smith } 685289bc588SBarry Smith } 6863a40ed3dSBarry Smith } else { 687dbb450caSBarry Smith if (addv == INSERT_VALUES) { 688e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 689cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 6902515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 691d0f46423SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 69258804f6dSBarry Smith #endif 693e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 694cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 6952515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 696d0f46423SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 69758804f6dSBarry Smith #endif 698cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 699e8d4e0b9SBarry Smith } 700e8d4e0b9SBarry Smith } 7013a40ed3dSBarry Smith } else { 702289bc588SBarry Smith for (i=0; i<m; i++) { 703cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7042515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 705d0f46423SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 70658804f6dSBarry Smith #endif 707289bc588SBarry Smith for (j=0; j<n; j++) { 708cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7092515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 710d0f46423SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 71158804f6dSBarry Smith #endif 712cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 713289bc588SBarry Smith } 714289bc588SBarry Smith } 715289bc588SBarry Smith } 716e8d4e0b9SBarry Smith } 7173a40ed3dSBarry Smith PetscFunctionReturn(0); 718289bc588SBarry Smith } 719e8d4e0b9SBarry Smith 7204a2ae208SSatish Balay #undef __FUNCT__ 7214a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 72213f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 723ae80bb75SLois Curfman McInnes { 724ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 72513f74950SBarry Smith PetscInt i,j; 726ae80bb75SLois Curfman McInnes 7273a40ed3dSBarry Smith PetscFunctionBegin; 728ae80bb75SLois Curfman McInnes /* row-oriented output */ 729ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 73097e567efSBarry Smith if (indexm[i] < 0) {v += n;continue;} 731d0f46423SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap->n); 732ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 7336f31f424SBarry Smith if (indexn[j] < 0) {v++; continue;} 734d0f46423SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap->n); 73597e567efSBarry Smith *v++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 736ae80bb75SLois Curfman McInnes } 737ae80bb75SLois Curfman McInnes } 7383a40ed3dSBarry Smith PetscFunctionReturn(0); 739ae80bb75SLois Curfman McInnes } 740ae80bb75SLois Curfman McInnes 741289bc588SBarry Smith /* -----------------------------------------------------------------*/ 742289bc588SBarry Smith 743e090d566SSatish Balay #include "petscsys.h" 74488e32edaSLois Curfman McInnes 7454a2ae208SSatish Balay #undef __FUNCT__ 7464a2ae208SSatish Balay #define __FUNCT__ "MatLoad_SeqDense" 747a313700dSBarry Smith PetscErrorCode MatLoad_SeqDense(PetscViewer viewer, const MatType type,Mat *A) 74888e32edaSLois Curfman McInnes { 74988e32edaSLois Curfman McInnes Mat_SeqDense *a; 75088e32edaSLois Curfman McInnes Mat B; 7516849ba73SBarry Smith PetscErrorCode ierr; 75213f74950SBarry Smith PetscInt *scols,i,j,nz,header[4]; 75313f74950SBarry Smith int fd; 75413f74950SBarry Smith PetscMPIInt size; 75513f74950SBarry Smith PetscInt *rowlengths = 0,M,N,*cols; 75687828ca2SBarry Smith PetscScalar *vals,*svals,*v,*w; 75719bcc07fSBarry Smith MPI_Comm comm = ((PetscObject)viewer)->comm; 75888e32edaSLois Curfman McInnes 7593a40ed3dSBarry Smith PetscFunctionBegin; 760d132466eSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 76129bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"view must have one processor"); 762b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 7630752156aSBarry Smith ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 764552e946dSBarry Smith if (header[0] != MAT_FILE_COOKIE) SETERRQ(PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 76588e32edaSLois Curfman McInnes M = header[1]; N = header[2]; nz = header[3]; 76688e32edaSLois Curfman McInnes 76790ace30eSBarry Smith if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 768f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 769f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr); 770be5d1d56SKris Buschelman ierr = MatSetType(*A,type);CHKERRQ(ierr); 7715c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr); 77290ace30eSBarry Smith B = *A; 77390ace30eSBarry Smith a = (Mat_SeqDense*)B->data; 7744a41a4d3SSatish Balay v = a->v; 7754a41a4d3SSatish Balay /* Allocate some temp space to read in the values and then flip them 7764a41a4d3SSatish Balay from row major to column major */ 777b72907f3SBarry Smith ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 77890ace30eSBarry Smith /* read in nonzero values */ 7794a41a4d3SSatish Balay ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 7804a41a4d3SSatish Balay /* now flip the values and store them in the matrix*/ 7814a41a4d3SSatish Balay for (j=0; j<N; j++) { 7824a41a4d3SSatish Balay for (i=0; i<M; i++) { 7834a41a4d3SSatish Balay *v++ =w[i*N+j]; 7844a41a4d3SSatish Balay } 7854a41a4d3SSatish Balay } 786606d414cSSatish Balay ierr = PetscFree(w);CHKERRQ(ierr); 7876d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7886d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 78990ace30eSBarry Smith } else { 79088e32edaSLois Curfman McInnes /* read row lengths */ 79113f74950SBarry Smith ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 7920752156aSBarry Smith ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 79388e32edaSLois Curfman McInnes 79488e32edaSLois Curfman McInnes /* create our matrix */ 795f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 796f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,M,N,M,N);CHKERRQ(ierr); 797be5d1d56SKris Buschelman ierr = MatSetType(*A,type);CHKERRQ(ierr); 7985c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(*A,PETSC_NULL);CHKERRQ(ierr); 79988e32edaSLois Curfman McInnes B = *A; 80088e32edaSLois Curfman McInnes a = (Mat_SeqDense*)B->data; 80188e32edaSLois Curfman McInnes v = a->v; 80288e32edaSLois Curfman McInnes 80388e32edaSLois Curfman McInnes /* read column indices and nonzeros */ 80413f74950SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr); 805b0a32e0cSBarry Smith cols = scols; 8060752156aSBarry Smith ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 80787828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 808b0a32e0cSBarry Smith vals = svals; 8090752156aSBarry Smith ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 81088e32edaSLois Curfman McInnes 81188e32edaSLois Curfman McInnes /* insert into matrix */ 81288e32edaSLois Curfman McInnes for (i=0; i<M; i++) { 81388e32edaSLois Curfman McInnes for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 81488e32edaSLois Curfman McInnes svals += rowlengths[i]; scols += rowlengths[i]; 81588e32edaSLois Curfman McInnes } 816606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 817606d414cSSatish Balay ierr = PetscFree(cols);CHKERRQ(ierr); 818606d414cSSatish Balay ierr = PetscFree(rowlengths);CHKERRQ(ierr); 81988e32edaSLois Curfman McInnes 8206d4a8577SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8216d4a8577SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 82290ace30eSBarry Smith } 8233a40ed3dSBarry Smith PetscFunctionReturn(0); 82488e32edaSLois Curfman McInnes } 82588e32edaSLois Curfman McInnes 826e090d566SSatish Balay #include "petscsys.h" 827932b0c3eSLois Curfman McInnes 8284a2ae208SSatish Balay #undef __FUNCT__ 8294a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 8306849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 831289bc588SBarry Smith { 832932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 833dfbe8321SBarry Smith PetscErrorCode ierr; 83413f74950SBarry Smith PetscInt i,j; 8352dcb1b2aSMatthew Knepley const char *name; 83687828ca2SBarry Smith PetscScalar *v; 837f3ef73ceSBarry Smith PetscViewerFormat format; 8385f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX) 8395f481a85SSatish Balay PetscTruth allreal = PETSC_TRUE; 8405f481a85SSatish Balay #endif 841932b0c3eSLois Curfman McInnes 8423a40ed3dSBarry Smith PetscFunctionBegin; 843435da068SBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 844b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 845456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 8463a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 847fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 848b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 849d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 85044cd7ae7SLois Curfman McInnes v = a->v + i; 85177431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 852d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 853aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 854329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 855a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 856329f5518SBarry Smith } else if (PetscRealPart(*v)) { 857a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr); 8586831982aSBarry Smith } 85980cd9d93SLois Curfman McInnes #else 8606831982aSBarry Smith if (*v) { 861a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr); 8626831982aSBarry Smith } 86380cd9d93SLois Curfman McInnes #endif 8641b807ce4Svictorle v += a->lda; 86580cd9d93SLois Curfman McInnes } 866b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 86780cd9d93SLois Curfman McInnes } 868b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 8693a40ed3dSBarry Smith } else { 870b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_NO);CHKERRQ(ierr); 871aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 87247989497SBarry Smith /* determine if matrix has all real values */ 87347989497SBarry Smith v = a->v; 874d0f46423SBarry Smith for (i=0; i<A->rmap->n*A->cmap->n; i++) { 875ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break ;} 87647989497SBarry Smith } 87747989497SBarry Smith #endif 878fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 8793a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 880d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr); 881d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 882fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 883ffac6cdbSBarry Smith } 884ffac6cdbSBarry Smith 885d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 886932b0c3eSLois Curfman McInnes v = a->v + i; 887d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 888aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 88947989497SBarry Smith if (allreal) { 890f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",PetscRealPart(*v));CHKERRQ(ierr); 89147989497SBarry Smith } else { 892f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16e i ",PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 89347989497SBarry Smith } 894289bc588SBarry Smith #else 895f32d5d43SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",*v);CHKERRQ(ierr); 896289bc588SBarry Smith #endif 8971b807ce4Svictorle v += a->lda; 898289bc588SBarry Smith } 899b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 900289bc588SBarry Smith } 901fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 902b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 903ffac6cdbSBarry Smith } 904b0a32e0cSBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_YES);CHKERRQ(ierr); 905da3a660dSBarry Smith } 906b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 9073a40ed3dSBarry Smith PetscFunctionReturn(0); 908289bc588SBarry Smith } 909289bc588SBarry Smith 9104a2ae208SSatish Balay #undef __FUNCT__ 9114a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 9126849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 913932b0c3eSLois Curfman McInnes { 914932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9156849ba73SBarry Smith PetscErrorCode ierr; 91613f74950SBarry Smith int fd; 917d0f46423SBarry Smith PetscInt ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n; 91887828ca2SBarry Smith PetscScalar *v,*anonz,*vals; 919f3ef73ceSBarry Smith PetscViewerFormat format; 920932b0c3eSLois Curfman McInnes 9213a40ed3dSBarry Smith PetscFunctionBegin; 922b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 92390ace30eSBarry Smith 924b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 925fb9695e5SSatish Balay if (format == PETSC_VIEWER_BINARY_NATIVE) { 92690ace30eSBarry Smith /* store the matrix as a dense matrix */ 92713f74950SBarry Smith ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 928552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 92990ace30eSBarry Smith col_lens[1] = m; 93090ace30eSBarry Smith col_lens[2] = n; 93190ace30eSBarry Smith col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 9326f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 933606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 93490ace30eSBarry Smith 93590ace30eSBarry Smith /* write out matrix, by rows */ 93687828ca2SBarry Smith ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 93790ace30eSBarry Smith v = a->v; 93890ace30eSBarry Smith for (j=0; j<n; j++) { 939578230a0SSatish Balay for (i=0; i<m; i++) { 940578230a0SSatish Balay vals[j + i*n] = *v++; 94190ace30eSBarry Smith } 94290ace30eSBarry Smith } 9436f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 944606d414cSSatish Balay ierr = PetscFree(vals);CHKERRQ(ierr); 94590ace30eSBarry Smith } else { 94613f74950SBarry Smith ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 947552e946dSBarry Smith col_lens[0] = MAT_FILE_COOKIE; 948932b0c3eSLois Curfman McInnes col_lens[1] = m; 949932b0c3eSLois Curfman McInnes col_lens[2] = n; 950932b0c3eSLois Curfman McInnes col_lens[3] = nz; 951932b0c3eSLois Curfman McInnes 952932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 953932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 9546f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 955932b0c3eSLois Curfman McInnes 956932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 957932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 958932b0c3eSLois Curfman McInnes ict = 0; 959932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 960932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 961932b0c3eSLois Curfman McInnes } 9626f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 963606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 964932b0c3eSLois Curfman McInnes 965932b0c3eSLois Curfman McInnes /* store nonzero values */ 96687828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 967932b0c3eSLois Curfman McInnes ict = 0; 968932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 969932b0c3eSLois Curfman McInnes v = a->v + i; 970932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 9711b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 972932b0c3eSLois Curfman McInnes } 973932b0c3eSLois Curfman McInnes } 9746f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 975606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 97690ace30eSBarry Smith } 9773a40ed3dSBarry Smith PetscFunctionReturn(0); 978932b0c3eSLois Curfman McInnes } 979932b0c3eSLois Curfman McInnes 9804a2ae208SSatish Balay #undef __FUNCT__ 9814a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 982dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 983f1af5d2fSBarry Smith { 984f1af5d2fSBarry Smith Mat A = (Mat) Aa; 985f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9866849ba73SBarry Smith PetscErrorCode ierr; 987d0f46423SBarry Smith PetscInt m = A->rmap->n,n = A->cmap->n,color,i,j; 98887828ca2SBarry Smith PetscScalar *v = a->v; 989b0a32e0cSBarry Smith PetscViewer viewer; 990b0a32e0cSBarry Smith PetscDraw popup; 991329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 992f3ef73ceSBarry Smith PetscViewerFormat format; 993f1af5d2fSBarry Smith 994f1af5d2fSBarry Smith PetscFunctionBegin; 995f1af5d2fSBarry Smith 996f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 997b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 998b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 999f1af5d2fSBarry Smith 1000f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 1001fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 1002f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 1003b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1004f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1005f1af5d2fSBarry Smith x_l = j; 1006f1af5d2fSBarry Smith x_r = x_l + 1.0; 1007f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1008f1af5d2fSBarry Smith y_l = m - i - 1.0; 1009f1af5d2fSBarry Smith y_r = y_l + 1.0; 1010f1af5d2fSBarry Smith #if defined(PETSC_USE_COMPLEX) 1011329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 1012b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1013329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 1014b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1015f1af5d2fSBarry Smith } else { 1016f1af5d2fSBarry Smith continue; 1017f1af5d2fSBarry Smith } 1018f1af5d2fSBarry Smith #else 1019f1af5d2fSBarry Smith if (v[j*m+i] > 0.) { 1020b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1021f1af5d2fSBarry Smith } else if (v[j*m+i] < 0.) { 1022b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1023f1af5d2fSBarry Smith } else { 1024f1af5d2fSBarry Smith continue; 1025f1af5d2fSBarry Smith } 1026f1af5d2fSBarry Smith #endif 1027b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1028f1af5d2fSBarry Smith } 1029f1af5d2fSBarry Smith } 1030f1af5d2fSBarry Smith } else { 1031f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 1032f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 1033f1af5d2fSBarry Smith for(i = 0; i < m*n; i++) { 1034f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 1035f1af5d2fSBarry Smith } 1036b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 1037b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 1038b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 1039f1af5d2fSBarry Smith for(j = 0; j < n; j++) { 1040f1af5d2fSBarry Smith x_l = j; 1041f1af5d2fSBarry Smith x_r = x_l + 1.0; 1042f1af5d2fSBarry Smith for(i = 0; i < m; i++) { 1043f1af5d2fSBarry Smith y_l = m - i - 1.0; 1044f1af5d2fSBarry Smith y_r = y_l + 1.0; 1045b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 1046b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1047f1af5d2fSBarry Smith } 1048f1af5d2fSBarry Smith } 1049f1af5d2fSBarry Smith } 1050f1af5d2fSBarry Smith PetscFunctionReturn(0); 1051f1af5d2fSBarry Smith } 1052f1af5d2fSBarry Smith 10534a2ae208SSatish Balay #undef __FUNCT__ 10544a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 1055dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 1056f1af5d2fSBarry Smith { 1057b0a32e0cSBarry Smith PetscDraw draw; 1058f1af5d2fSBarry Smith PetscTruth isnull; 1059329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 1060dfbe8321SBarry Smith PetscErrorCode ierr; 1061f1af5d2fSBarry Smith 1062f1af5d2fSBarry Smith PetscFunctionBegin; 1063b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1064b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1065abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 1066f1af5d2fSBarry Smith 1067f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1068d0f46423SBarry Smith xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0; 1069f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 1070b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1071b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 1072f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",PETSC_NULL);CHKERRQ(ierr); 1073f1af5d2fSBarry Smith PetscFunctionReturn(0); 1074f1af5d2fSBarry Smith } 1075f1af5d2fSBarry Smith 10764a2ae208SSatish Balay #undef __FUNCT__ 10774a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 1078dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1079932b0c3eSLois Curfman McInnes { 1080dfbe8321SBarry Smith PetscErrorCode ierr; 10816805f65bSBarry Smith PetscTruth iascii,isbinary,isdraw; 1082932b0c3eSLois Curfman McInnes 10833a40ed3dSBarry Smith PetscFunctionBegin; 108432077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 1085fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_BINARY,&isbinary);CHKERRQ(ierr); 1086fb9695e5SSatish Balay ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);CHKERRQ(ierr); 10870f5bd95cSBarry Smith 1088c45a1595SBarry Smith if (iascii) { 1089c45a1595SBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 10900f5bd95cSBarry Smith } else if (isbinary) { 10913a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1092f1af5d2fSBarry Smith } else if (isdraw) { 1093f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 10945cd90555SBarry Smith } else { 1095958c9bccSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported by dense matrix",((PetscObject)viewer)->type_name); 1096932b0c3eSLois Curfman McInnes } 10973a40ed3dSBarry Smith PetscFunctionReturn(0); 1098932b0c3eSLois Curfman McInnes } 1099289bc588SBarry Smith 11004a2ae208SSatish Balay #undef __FUNCT__ 11014a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1102dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1103289bc588SBarry Smith { 1104ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1105dfbe8321SBarry Smith PetscErrorCode ierr; 110690f02eecSBarry Smith 11073a40ed3dSBarry Smith PetscFunctionBegin; 1108aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1109d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n); 1110a5a9c739SBarry Smith #endif 111105b42c5fSBarry Smith ierr = PetscFree(l->pivots);CHKERRQ(ierr); 11126857c123SSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1113606d414cSSatish Balay ierr = PetscFree(l);CHKERRQ(ierr); 1114dbd8c25aSHong Zhang 1115dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 1116901853e0SKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",PETSC_NULL);CHKERRQ(ierr); 11174ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11184ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11194ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",PETSC_NULL);CHKERRQ(ierr); 11203a40ed3dSBarry Smith PetscFunctionReturn(0); 1121289bc588SBarry Smith } 1122289bc588SBarry Smith 11234a2ae208SSatish Balay #undef __FUNCT__ 11244a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1125fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout) 1126289bc588SBarry Smith { 1127c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11286849ba73SBarry Smith PetscErrorCode ierr; 112913f74950SBarry Smith PetscInt k,j,m,n,M; 113087828ca2SBarry Smith PetscScalar *v,tmp; 113148b35521SBarry Smith 11323a40ed3dSBarry Smith PetscFunctionBegin; 1133d0f46423SBarry Smith v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n; 1134e9695a30SBarry Smith if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */ 1135a5ce6ee0Svictorle if (m != n) { 1136634064b4SBarry Smith SETERRQ(PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1137d64ed03dSBarry Smith } else { 1138d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1139289bc588SBarry Smith for (k=0; k<j; k++) { 11401b807ce4Svictorle tmp = v[j + k*M]; 11411b807ce4Svictorle v[j + k*M] = v[k + j*M]; 11421b807ce4Svictorle v[k + j*M] = tmp; 1143289bc588SBarry Smith } 1144289bc588SBarry Smith } 1145d64ed03dSBarry Smith } 11463a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1147d3e5ee88SLois Curfman McInnes Mat tmat; 1148ec8511deSBarry Smith Mat_SeqDense *tmatd; 114987828ca2SBarry Smith PetscScalar *v2; 1150ea709b57SSatish Balay 1151fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 11527adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&tmat);CHKERRQ(ierr); 1153d0f46423SBarry Smith ierr = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr); 11547adad957SLisandro Dalcin ierr = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 11555c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(tmat,PETSC_NULL);CHKERRQ(ierr); 1156fc4dec0aSBarry Smith } else { 1157fc4dec0aSBarry Smith tmat = *matout; 1158fc4dec0aSBarry Smith } 1159ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 11600de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1161d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 11621b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1163d3e5ee88SLois Curfman McInnes } 11646d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 11656d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1166d3e5ee88SLois Curfman McInnes *matout = tmat; 116748b35521SBarry Smith } 11683a40ed3dSBarry Smith PetscFunctionReturn(0); 1169289bc588SBarry Smith } 1170289bc588SBarry Smith 11714a2ae208SSatish Balay #undef __FUNCT__ 11724a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1173dfbe8321SBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscTruth *flg) 1174289bc588SBarry Smith { 1175c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1176c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 117713f74950SBarry Smith PetscInt i,j; 117887828ca2SBarry Smith PetscScalar *v1 = mat1->v,*v2 = mat2->v; 11799ea5d5aeSSatish Balay 11803a40ed3dSBarry Smith PetscFunctionBegin; 1181d0f46423SBarry Smith if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1182d0f46423SBarry Smith if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1183d0f46423SBarry Smith for (i=0; i<A1->rmap->n; i++) { 11841b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1185d0f46423SBarry Smith for (j=0; j<A1->cmap->n; j++) { 11863a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 11871b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 11881b807ce4Svictorle } 1189289bc588SBarry Smith } 119077c4ece6SBarry Smith *flg = PETSC_TRUE; 11913a40ed3dSBarry Smith PetscFunctionReturn(0); 1192289bc588SBarry Smith } 1193289bc588SBarry Smith 11944a2ae208SSatish Balay #undef __FUNCT__ 11954a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1196dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1197289bc588SBarry Smith { 1198c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1199dfbe8321SBarry Smith PetscErrorCode ierr; 120013f74950SBarry Smith PetscInt i,n,len; 120187828ca2SBarry Smith PetscScalar *x,zero = 0.0; 120244cd7ae7SLois Curfman McInnes 12033a40ed3dSBarry Smith PetscFunctionBegin; 12042dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 12057a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 12061ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1207d0f46423SBarry Smith len = PetscMin(A->rmap->n,A->cmap->n); 1208d0f46423SBarry Smith if (n != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 120944cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 12101b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1211289bc588SBarry Smith } 12121ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 12133a40ed3dSBarry Smith PetscFunctionReturn(0); 1214289bc588SBarry Smith } 1215289bc588SBarry Smith 12164a2ae208SSatish Balay #undef __FUNCT__ 12174a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1218dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1219289bc588SBarry Smith { 1220c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 122187828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1222dfbe8321SBarry Smith PetscErrorCode ierr; 1223d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n; 122455659b69SBarry Smith 12253a40ed3dSBarry Smith PetscFunctionBegin; 122628988994SBarry Smith if (ll) { 12277a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 12281ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1229d0f46423SBarry Smith if (m != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1230da3a660dSBarry Smith for (i=0; i<m; i++) { 1231da3a660dSBarry Smith x = l[i]; 1232da3a660dSBarry Smith v = mat->v + i; 1233da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1234da3a660dSBarry Smith } 12351ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1236efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1237da3a660dSBarry Smith } 123828988994SBarry Smith if (rr) { 12397a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 12401ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1241d0f46423SBarry Smith if (n != A->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1242da3a660dSBarry Smith for (i=0; i<n; i++) { 1243da3a660dSBarry Smith x = r[i]; 1244da3a660dSBarry Smith v = mat->v + i*m; 1245da3a660dSBarry Smith for (j=0; j<m; j++) { (*v++) *= x;} 1246da3a660dSBarry Smith } 12471ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1248efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1249da3a660dSBarry Smith } 12503a40ed3dSBarry Smith PetscFunctionReturn(0); 1251289bc588SBarry Smith } 1252289bc588SBarry Smith 12534a2ae208SSatish Balay #undef __FUNCT__ 12544a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1255dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1256289bc588SBarry Smith { 1257c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 125887828ca2SBarry Smith PetscScalar *v = mat->v; 1259329f5518SBarry Smith PetscReal sum = 0.0; 1260d0f46423SBarry Smith PetscInt lda=mat->lda,m=A->rmap->n,i,j; 1261efee365bSSatish Balay PetscErrorCode ierr; 126255659b69SBarry Smith 12633a40ed3dSBarry Smith PetscFunctionBegin; 1264289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1265a5ce6ee0Svictorle if (lda>m) { 1266d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1267a5ce6ee0Svictorle v = mat->v+j*lda; 1268a5ce6ee0Svictorle for (i=0; i<m; i++) { 1269a5ce6ee0Svictorle #if defined(PETSC_USE_COMPLEX) 1270a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1271a5ce6ee0Svictorle #else 1272a5ce6ee0Svictorle sum += (*v)*(*v); v++; 1273a5ce6ee0Svictorle #endif 1274a5ce6ee0Svictorle } 1275a5ce6ee0Svictorle } 1276a5ce6ee0Svictorle } else { 1277d0f46423SBarry Smith for (i=0; i<A->cmap->n*A->rmap->n; i++) { 1278aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 1279329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1280289bc588SBarry Smith #else 1281289bc588SBarry Smith sum += (*v)*(*v); v++; 1282289bc588SBarry Smith #endif 1283289bc588SBarry Smith } 1284a5ce6ee0Svictorle } 1285064f8208SBarry Smith *nrm = sqrt(sum); 1286d0f46423SBarry Smith ierr = PetscLogFlops(2*A->cmap->n*A->rmap->n);CHKERRQ(ierr); 12873a40ed3dSBarry Smith } else if (type == NORM_1) { 1288064f8208SBarry Smith *nrm = 0.0; 1289d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 12901b807ce4Svictorle v = mat->v + j*mat->lda; 1291289bc588SBarry Smith sum = 0.0; 1292d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 129333a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1294289bc588SBarry Smith } 1295064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1296289bc588SBarry Smith } 1297d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 12983a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1299064f8208SBarry Smith *nrm = 0.0; 1300d0f46423SBarry Smith for (j=0; j<A->rmap->n; j++) { 1301289bc588SBarry Smith v = mat->v + j; 1302289bc588SBarry Smith sum = 0.0; 1303d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) { 13041b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1305289bc588SBarry Smith } 1306064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1307289bc588SBarry Smith } 1308d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13093a40ed3dSBarry Smith } else { 131029bbc08cSBarry Smith SETERRQ(PETSC_ERR_SUP,"No two norm"); 1311289bc588SBarry Smith } 13123a40ed3dSBarry Smith PetscFunctionReturn(0); 1313289bc588SBarry Smith } 1314289bc588SBarry Smith 13154a2ae208SSatish Balay #undef __FUNCT__ 13164a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 13174e0d8c25SBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscTruth flg) 1318289bc588SBarry Smith { 1319c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 132063ba0a88SBarry Smith PetscErrorCode ierr; 132167e560aaSBarry Smith 13223a40ed3dSBarry Smith PetscFunctionBegin; 1323b5a2b587SKris Buschelman switch (op) { 1324b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 13254e0d8c25SBarry Smith aij->roworiented = flg; 1326b5a2b587SKris Buschelman break; 1327512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 1328b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 13293971808eSMatthew Knepley case MAT_NEW_NONZERO_ALLOCATION_ERR: 13304e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 1331b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1332b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 133377e54ba9SKris Buschelman case MAT_SYMMETRIC: 133477e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 13359a4540c5SBarry Smith case MAT_HERMITIAN: 13369a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 1337600fe468SBarry Smith case MAT_IGNORE_LOWER_TRIANGULAR: 1338290bbb0aSBarry Smith ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 133977e54ba9SKris Buschelman break; 1340b5a2b587SKris Buschelman default: 1341600fe468SBarry Smith SETERRQ1(PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 13423a40ed3dSBarry Smith } 13433a40ed3dSBarry Smith PetscFunctionReturn(0); 1344289bc588SBarry Smith } 1345289bc588SBarry Smith 13464a2ae208SSatish Balay #undef __FUNCT__ 13474a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1348dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 13496f0a148fSBarry Smith { 1350ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 13516849ba73SBarry Smith PetscErrorCode ierr; 1352d0f46423SBarry Smith PetscInt lda=l->lda,m=A->rmap->n,j; 13533a40ed3dSBarry Smith 13543a40ed3dSBarry Smith PetscFunctionBegin; 1355a5ce6ee0Svictorle if (lda>m) { 1356d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1357a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1358a5ce6ee0Svictorle } 1359a5ce6ee0Svictorle } else { 1360d0f46423SBarry Smith ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1361a5ce6ee0Svictorle } 13623a40ed3dSBarry Smith PetscFunctionReturn(0); 13636f0a148fSBarry Smith } 13646f0a148fSBarry Smith 13654a2ae208SSatish Balay #undef __FUNCT__ 13664a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 1367f4df32b1SMatthew Knepley PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag) 13686f0a148fSBarry Smith { 1369ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1370d0f46423SBarry Smith PetscInt n = A->cmap->n,i,j; 137187828ca2SBarry Smith PetscScalar *slot; 137255659b69SBarry Smith 13733a40ed3dSBarry Smith PetscFunctionBegin; 13746f0a148fSBarry Smith for (i=0; i<N; i++) { 13756f0a148fSBarry Smith slot = l->v + rows[i]; 13766f0a148fSBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += n;} 13776f0a148fSBarry Smith } 1378f4df32b1SMatthew Knepley if (diag != 0.0) { 13796f0a148fSBarry Smith for (i=0; i<N; i++) { 13806f0a148fSBarry Smith slot = l->v + (n+1)*rows[i]; 1381f4df32b1SMatthew Knepley *slot = diag; 13826f0a148fSBarry Smith } 13836f0a148fSBarry Smith } 13843a40ed3dSBarry Smith PetscFunctionReturn(0); 13856f0a148fSBarry Smith } 1386557bce09SLois Curfman McInnes 13874a2ae208SSatish Balay #undef __FUNCT__ 13884a2ae208SSatish Balay #define __FUNCT__ "MatGetArray_SeqDense" 1389dfbe8321SBarry Smith PetscErrorCode MatGetArray_SeqDense(Mat A,PetscScalar *array[]) 139064e87e97SBarry Smith { 1391c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 13923a40ed3dSBarry Smith 13933a40ed3dSBarry Smith PetscFunctionBegin; 1394d0f46423SBarry Smith if (mat->lda != A->rmap->n) SETERRQ(PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows"); 139564e87e97SBarry Smith *array = mat->v; 13963a40ed3dSBarry Smith PetscFunctionReturn(0); 139764e87e97SBarry Smith } 13980754003eSLois Curfman McInnes 13994a2ae208SSatish Balay #undef __FUNCT__ 14004a2ae208SSatish Balay #define __FUNCT__ "MatRestoreArray_SeqDense" 1401dfbe8321SBarry Smith PetscErrorCode MatRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1402ff14e315SSatish Balay { 14033a40ed3dSBarry Smith PetscFunctionBegin; 140409b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 14053a40ed3dSBarry Smith PetscFunctionReturn(0); 1406ff14e315SSatish Balay } 14070754003eSLois Curfman McInnes 14084a2ae208SSatish Balay #undef __FUNCT__ 14094a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 141013f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 14110754003eSLois Curfman McInnes { 1412c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14136849ba73SBarry Smith PetscErrorCode ierr; 141421a2c019SBarry Smith PetscInt i,j,*irow,*icol,nrows,ncols; 141587828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 14160754003eSLois Curfman McInnes Mat newmat; 14170754003eSLois Curfman McInnes 14183a40ed3dSBarry Smith PetscFunctionBegin; 141978b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 142078b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1421e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1422e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 14230754003eSLois Curfman McInnes 1424182d2002SSatish Balay /* Check submatrixcall */ 1425182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 142613f74950SBarry Smith PetscInt n_cols,n_rows; 1427182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 142821a2c019SBarry Smith if (n_rows != nrows || n_cols != ncols) { 142921a2c019SBarry Smith /* resize the result result matrix to match number of requested rows/columns */ 143021a2c019SBarry Smith ierr = MatSetSizes(*B,nrows,nrows,nrows,nrows);CHKERRQ(ierr); 143121a2c019SBarry Smith } 1432182d2002SSatish Balay newmat = *B; 1433182d2002SSatish Balay } else { 14340754003eSLois Curfman McInnes /* Create and fill new matrix */ 14357adad957SLisandro Dalcin ierr = MatCreate(((PetscObject)A)->comm,&newmat);CHKERRQ(ierr); 1436f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 14377adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 14385c5985e7SKris Buschelman ierr = MatSeqDenseSetPreallocation(newmat,PETSC_NULL);CHKERRQ(ierr); 1439182d2002SSatish Balay } 1440182d2002SSatish Balay 1441182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1442182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1443182d2002SSatish Balay 1444182d2002SSatish Balay for (i=0; i<ncols; i++) { 14456de62eeeSBarry Smith av = v + mat->lda*icol[i]; 1446182d2002SSatish Balay for (j=0; j<nrows; j++) { 1447182d2002SSatish Balay *bv++ = av[irow[j]]; 14480754003eSLois Curfman McInnes } 14490754003eSLois Curfman McInnes } 1450182d2002SSatish Balay 1451182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 14526d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14536d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14540754003eSLois Curfman McInnes 14550754003eSLois Curfman McInnes /* Free work space */ 145678b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 145778b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1458182d2002SSatish Balay *B = newmat; 14593a40ed3dSBarry Smith PetscFunctionReturn(0); 14600754003eSLois Curfman McInnes } 14610754003eSLois Curfman McInnes 14624a2ae208SSatish Balay #undef __FUNCT__ 14634a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 146413f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1465905e6a2fSBarry Smith { 14666849ba73SBarry Smith PetscErrorCode ierr; 146713f74950SBarry Smith PetscInt i; 1468905e6a2fSBarry Smith 14693a40ed3dSBarry Smith PetscFunctionBegin; 1470905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1471b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1472905e6a2fSBarry Smith } 1473905e6a2fSBarry Smith 1474905e6a2fSBarry Smith for (i=0; i<n; i++) { 14756a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1476905e6a2fSBarry Smith } 14773a40ed3dSBarry Smith PetscFunctionReturn(0); 1478905e6a2fSBarry Smith } 1479905e6a2fSBarry Smith 14804a2ae208SSatish Balay #undef __FUNCT__ 1481c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1482c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1483c0aa2d19SHong Zhang { 1484c0aa2d19SHong Zhang PetscFunctionBegin; 1485c0aa2d19SHong Zhang PetscFunctionReturn(0); 1486c0aa2d19SHong Zhang } 1487c0aa2d19SHong Zhang 1488c0aa2d19SHong Zhang #undef __FUNCT__ 1489c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1490c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1491c0aa2d19SHong Zhang { 1492c0aa2d19SHong Zhang PetscFunctionBegin; 1493c0aa2d19SHong Zhang PetscFunctionReturn(0); 1494c0aa2d19SHong Zhang } 1495c0aa2d19SHong Zhang 1496c0aa2d19SHong Zhang #undef __FUNCT__ 14974a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1498dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 14994b0e389bSBarry Smith { 15004b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense *)B->data; 15016849ba73SBarry Smith PetscErrorCode ierr; 1502d0f46423SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j; 15033a40ed3dSBarry Smith 15043a40ed3dSBarry Smith PetscFunctionBegin; 150533f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 150633f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1507cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 15083a40ed3dSBarry Smith PetscFunctionReturn(0); 15093a40ed3dSBarry Smith } 1510d0f46423SBarry Smith if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1511a5ce6ee0Svictorle if (lda1>m || lda2>m) { 15120dbb7854Svictorle for (j=0; j<n; j++) { 1513a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1514a5ce6ee0Svictorle } 1515a5ce6ee0Svictorle } else { 1516d0f46423SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1517a5ce6ee0Svictorle } 1518273d9f13SBarry Smith PetscFunctionReturn(0); 1519273d9f13SBarry Smith } 1520273d9f13SBarry Smith 15214a2ae208SSatish Balay #undef __FUNCT__ 15224a2ae208SSatish Balay #define __FUNCT__ "MatSetUpPreallocation_SeqDense" 1523dfbe8321SBarry Smith PetscErrorCode MatSetUpPreallocation_SeqDense(Mat A) 1524273d9f13SBarry Smith { 1525dfbe8321SBarry Smith PetscErrorCode ierr; 1526273d9f13SBarry Smith 1527273d9f13SBarry Smith PetscFunctionBegin; 1528273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 15293a40ed3dSBarry Smith PetscFunctionReturn(0); 15304b0e389bSBarry Smith } 15314b0e389bSBarry Smith 1532284134d9SBarry Smith #undef __FUNCT__ 1533284134d9SBarry Smith #define __FUNCT__ "MatSetSizes_SeqDense" 1534284134d9SBarry Smith PetscErrorCode MatSetSizes_SeqDense(Mat A,PetscInt m,PetscInt n,PetscInt M,PetscInt N) 1535284134d9SBarry Smith { 1536284134d9SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 153721a2c019SBarry Smith PetscErrorCode ierr; 1538284134d9SBarry Smith PetscFunctionBegin; 153921a2c019SBarry Smith /* this will not be called before lda, Mmax, and Nmax have been set */ 1540284134d9SBarry Smith m = PetscMax(m,M); 1541284134d9SBarry Smith n = PetscMax(n,N); 154221a2c019SBarry Smith if (m > a->Mmax) SETERRQ2(PETSC_ERR_SUP,"Cannot yet resize number rows of dense matrix larger then its initial size %d, requested %d",a->lda,(int)m); 1543284134d9SBarry Smith if (n > a->Nmax) SETERRQ2(PETSC_ERR_SUP,"Cannot yet resize number columns of dense matrix larger then its initial size %d, requested %d",a->Nmax,(int)n); 1544d0f46423SBarry Smith A->rmap->n = A->rmap->n = m; 1545d0f46423SBarry Smith A->cmap->n = A->cmap->N = n; 154621a2c019SBarry Smith if (a->changelda) a->lda = m; 154721a2c019SBarry Smith ierr = PetscMemzero(a->v,m*n*sizeof(PetscScalar));CHKERRQ(ierr); 1548284134d9SBarry Smith PetscFunctionReturn(0); 1549284134d9SBarry Smith } 1550170fe5c8SBarry Smith 1551284134d9SBarry Smith 1552a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/ 1553a9fe9ddaSSatish Balay #undef __FUNCT__ 1554a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1555a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1556a9fe9ddaSSatish Balay { 1557a9fe9ddaSSatish Balay PetscErrorCode ierr; 1558a9fe9ddaSSatish Balay 1559a9fe9ddaSSatish Balay PetscFunctionBegin; 1560a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1561a9fe9ddaSSatish Balay ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1562a9fe9ddaSSatish Balay } 1563a9fe9ddaSSatish Balay ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1564a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1565a9fe9ddaSSatish Balay } 1566a9fe9ddaSSatish Balay 1567a9fe9ddaSSatish Balay #undef __FUNCT__ 1568a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1569a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1570a9fe9ddaSSatish Balay { 1571ee16a9a1SHong Zhang PetscErrorCode ierr; 1572d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 1573ee16a9a1SHong Zhang Mat Cmat; 1574a9fe9ddaSSatish Balay 1575ee16a9a1SHong Zhang PetscFunctionBegin; 1576d0f46423SBarry Smith if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->rmap->n %d\n",A->cmap->n,B->rmap->n); 1577ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1578ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1579ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1580ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1581ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1582ee16a9a1SHong Zhang *C = Cmat; 1583ee16a9a1SHong Zhang PetscFunctionReturn(0); 1584ee16a9a1SHong Zhang } 1585a9fe9ddaSSatish Balay 158698a3b096SSatish Balay #undef __FUNCT__ 1587a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1588a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1589a9fe9ddaSSatish Balay { 1590a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1591a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1592a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 15930805154bSBarry Smith PetscBLASInt m,n,k; 1594a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1595a9fe9ddaSSatish Balay 1596a9fe9ddaSSatish Balay PetscFunctionBegin; 1597d0f46423SBarry Smith m = PetscBLASIntCast(A->rmap->n); 1598d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1599d0f46423SBarry Smith k = PetscBLASIntCast(A->cmap->n); 1600a9fe9ddaSSatish Balay BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1601a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1602a9fe9ddaSSatish Balay } 1603a9fe9ddaSSatish Balay 1604a9fe9ddaSSatish Balay #undef __FUNCT__ 1605a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTranspose_SeqDense_SeqDense" 1606a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTranspose_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1607a9fe9ddaSSatish Balay { 1608a9fe9ddaSSatish Balay PetscErrorCode ierr; 1609a9fe9ddaSSatish Balay 1610a9fe9ddaSSatish Balay PetscFunctionBegin; 1611a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX){ 1612a9fe9ddaSSatish Balay ierr = MatMatMultTransposeSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1613a9fe9ddaSSatish Balay } 1614a9fe9ddaSSatish Balay ierr = MatMatMultTransposeNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1615a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1616a9fe9ddaSSatish Balay } 1617a9fe9ddaSSatish Balay 1618a9fe9ddaSSatish Balay #undef __FUNCT__ 1619a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeSymbolic_SeqDense_SeqDense" 1620a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1621a9fe9ddaSSatish Balay { 1622ee16a9a1SHong Zhang PetscErrorCode ierr; 1623d0f46423SBarry Smith PetscInt m=A->cmap->n,n=B->cmap->n; 1624ee16a9a1SHong Zhang Mat Cmat; 1625a9fe9ddaSSatish Balay 1626ee16a9a1SHong Zhang PetscFunctionBegin; 1627d0f46423SBarry Smith if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"A->rmap->n %d != B->rmap->n %d\n",A->rmap->n,B->rmap->n); 1628ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1629ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1630ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 1631ee16a9a1SHong Zhang ierr = MatSeqDenseSetPreallocation(Cmat,PETSC_NULL);CHKERRQ(ierr); 1632ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 1633ee16a9a1SHong Zhang *C = Cmat; 1634ee16a9a1SHong Zhang PetscFunctionReturn(0); 1635ee16a9a1SHong Zhang } 1636a9fe9ddaSSatish Balay 1637a9fe9ddaSSatish Balay #undef __FUNCT__ 1638a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultTransposeNumeric_SeqDense_SeqDense" 1639a9fe9ddaSSatish Balay PetscErrorCode MatMatMultTransposeNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1640a9fe9ddaSSatish Balay { 1641a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1642a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1643a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 16440805154bSBarry Smith PetscBLASInt m,n,k; 1645a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1646a9fe9ddaSSatish Balay 1647a9fe9ddaSSatish Balay PetscFunctionBegin; 1648d0f46423SBarry Smith m = PetscBLASIntCast(A->cmap->n); 1649d0f46423SBarry Smith n = PetscBLASIntCast(B->cmap->n); 1650d0f46423SBarry Smith k = PetscBLASIntCast(A->rmap->n); 16512fbe02b9SBarry Smith /* 16522fbe02b9SBarry Smith Note the m and n arguments below are the number rows and columns of A', not A! 16532fbe02b9SBarry Smith */ 1654a9fe9ddaSSatish Balay BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda); 1655a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1656a9fe9ddaSSatish Balay } 1657985db425SBarry Smith 1658985db425SBarry Smith #undef __FUNCT__ 1659985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense" 1660985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 1661985db425SBarry Smith { 1662985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1663985db425SBarry Smith PetscErrorCode ierr; 1664d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1665985db425SBarry Smith PetscScalar *x; 1666985db425SBarry Smith MatScalar *aa = a->v; 1667985db425SBarry Smith 1668985db425SBarry Smith PetscFunctionBegin; 1669985db425SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1670985db425SBarry Smith 1671985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1672985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1673985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1674d0f46423SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1675985db425SBarry Smith for (i=0; i<m; i++) { 1676985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1677985db425SBarry Smith for (j=1; j<n; j++){ 1678985db425SBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1679985db425SBarry Smith } 1680985db425SBarry Smith } 1681985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1682985db425SBarry Smith PetscFunctionReturn(0); 1683985db425SBarry Smith } 1684985db425SBarry Smith 1685985db425SBarry Smith #undef __FUNCT__ 1686985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense" 1687985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 1688985db425SBarry Smith { 1689985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1690985db425SBarry Smith PetscErrorCode ierr; 1691d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1692985db425SBarry Smith PetscScalar *x; 1693985db425SBarry Smith PetscReal atmp; 1694985db425SBarry Smith MatScalar *aa = a->v; 1695985db425SBarry Smith 1696985db425SBarry Smith PetscFunctionBegin; 1697985db425SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1698985db425SBarry Smith 1699985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1700985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1701985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1702d0f46423SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1703985db425SBarry Smith for (i=0; i<m; i++) { 17049189402eSHong Zhang x[i] = PetscAbsScalar(aa[i]); 1705985db425SBarry Smith for (j=1; j<n; j++){ 1706985db425SBarry Smith atmp = PetscAbsScalar(aa[i+m*j]); 1707985db425SBarry Smith if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 1708985db425SBarry Smith } 1709985db425SBarry Smith } 1710985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1711985db425SBarry Smith PetscFunctionReturn(0); 1712985db425SBarry Smith } 1713985db425SBarry Smith 1714985db425SBarry Smith #undef __FUNCT__ 1715985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense" 1716985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 1717985db425SBarry Smith { 1718985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1719985db425SBarry Smith PetscErrorCode ierr; 1720d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1721985db425SBarry Smith PetscScalar *x; 1722985db425SBarry Smith MatScalar *aa = a->v; 1723985db425SBarry Smith 1724985db425SBarry Smith PetscFunctionBegin; 1725985db425SBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1726985db425SBarry Smith 1727985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1728985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1729985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1730d0f46423SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1731985db425SBarry Smith for (i=0; i<m; i++) { 1732985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1733985db425SBarry Smith for (j=1; j<n; j++){ 1734985db425SBarry Smith if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1735985db425SBarry Smith } 1736985db425SBarry Smith } 1737985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1738985db425SBarry Smith PetscFunctionReturn(0); 1739985db425SBarry Smith } 1740985db425SBarry Smith 17418d0534beSBarry Smith #undef __FUNCT__ 17428d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense" 17438d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 17448d0534beSBarry Smith { 17458d0534beSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 17468d0534beSBarry Smith PetscErrorCode ierr; 17478d0534beSBarry Smith PetscScalar *x; 17488d0534beSBarry Smith 17498d0534beSBarry Smith PetscFunctionBegin; 17508d0534beSBarry Smith if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 17518d0534beSBarry Smith 17528d0534beSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1753d0f46423SBarry Smith ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 17548d0534beSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 17558d0534beSBarry Smith PetscFunctionReturn(0); 17568d0534beSBarry Smith } 17578d0534beSBarry Smith 1758289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1759a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = {MatSetValues_SeqDense, 1760905e6a2fSBarry Smith MatGetRow_SeqDense, 1761905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1762905e6a2fSBarry Smith MatMult_SeqDense, 176397304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 17647c922b88SBarry Smith MatMultTranspose_SeqDense, 17657c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1766db4efbfdSBarry Smith 0, 1767db4efbfdSBarry Smith 0, 1768db4efbfdSBarry Smith 0, 1769db4efbfdSBarry Smith /*10*/ 0, 1770905e6a2fSBarry Smith MatLUFactor_SeqDense, 1771905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 1772ec8511deSBarry Smith MatRelax_SeqDense, 1773ec8511deSBarry Smith MatTranspose_SeqDense, 177497304618SKris Buschelman /*15*/ MatGetInfo_SeqDense, 1775905e6a2fSBarry Smith MatEqual_SeqDense, 1776905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1777905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1778905e6a2fSBarry Smith MatNorm_SeqDense, 1779c0aa2d19SHong Zhang /*20*/ MatAssemblyBegin_SeqDense, 1780c0aa2d19SHong Zhang MatAssemblyEnd_SeqDense, 1781905e6a2fSBarry Smith 0, 1782905e6a2fSBarry Smith MatSetOption_SeqDense, 1783905e6a2fSBarry Smith MatZeroEntries_SeqDense, 178497304618SKris Buschelman /*25*/ MatZeroRows_SeqDense, 1785db4efbfdSBarry Smith 0, 1786db4efbfdSBarry Smith 0, 1787db4efbfdSBarry Smith 0, 1788db4efbfdSBarry Smith 0, 178997304618SKris Buschelman /*30*/ MatSetUpPreallocation_SeqDense, 1790273d9f13SBarry Smith 0, 1791905e6a2fSBarry Smith 0, 1792905e6a2fSBarry Smith MatGetArray_SeqDense, 1793905e6a2fSBarry Smith MatRestoreArray_SeqDense, 179497304618SKris Buschelman /*35*/ MatDuplicate_SeqDense, 1795a5ae1ecdSBarry Smith 0, 1796a5ae1ecdSBarry Smith 0, 1797a5ae1ecdSBarry Smith 0, 1798a5ae1ecdSBarry Smith 0, 179997304618SKris Buschelman /*40*/ MatAXPY_SeqDense, 1800a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1801a5ae1ecdSBarry Smith 0, 18024b0e389bSBarry Smith MatGetValues_SeqDense, 1803a5ae1ecdSBarry Smith MatCopy_SeqDense, 1804985db425SBarry Smith /*45*/ MatGetRowMax_SeqDense, 1805a5ae1ecdSBarry Smith MatScale_SeqDense, 1806a5ae1ecdSBarry Smith 0, 1807a5ae1ecdSBarry Smith 0, 1808a5ae1ecdSBarry Smith 0, 1809521d7252SBarry Smith /*50*/ 0, 1810a5ae1ecdSBarry Smith 0, 1811a5ae1ecdSBarry Smith 0, 1812a5ae1ecdSBarry Smith 0, 1813a5ae1ecdSBarry Smith 0, 181497304618SKris Buschelman /*55*/ 0, 1815a5ae1ecdSBarry Smith 0, 1816a5ae1ecdSBarry Smith 0, 1817a5ae1ecdSBarry Smith 0, 1818a5ae1ecdSBarry Smith 0, 181997304618SKris Buschelman /*60*/ 0, 1820e03a110bSBarry Smith MatDestroy_SeqDense, 1821e03a110bSBarry Smith MatView_SeqDense, 1822357abbc8SBarry Smith 0, 182397304618SKris Buschelman 0, 182497304618SKris Buschelman /*65*/ 0, 182597304618SKris Buschelman 0, 182697304618SKris Buschelman 0, 182797304618SKris Buschelman 0, 182897304618SKris Buschelman 0, 1829985db425SBarry Smith /*70*/ MatGetRowMaxAbs_SeqDense, 183097304618SKris Buschelman 0, 183197304618SKris Buschelman 0, 183297304618SKris Buschelman 0, 183397304618SKris Buschelman 0, 183497304618SKris Buschelman /*75*/ 0, 183597304618SKris Buschelman 0, 183697304618SKris Buschelman 0, 183797304618SKris Buschelman 0, 183897304618SKris Buschelman 0, 183997304618SKris Buschelman /*80*/ 0, 184097304618SKris Buschelman 0, 184197304618SKris Buschelman 0, 184297304618SKris Buschelman 0, 1843865e5f61SKris Buschelman /*84*/ MatLoad_SeqDense, 1844865e5f61SKris Buschelman 0, 18451cbb95d3SBarry Smith MatIsHermitian_SeqDense, 1846865e5f61SKris Buschelman 0, 1847865e5f61SKris Buschelman 0, 1848865e5f61SKris Buschelman 0, 1849a9fe9ddaSSatish Balay /*90*/ MatMatMult_SeqDense_SeqDense, 1850a9fe9ddaSSatish Balay MatMatMultSymbolic_SeqDense_SeqDense, 1851a9fe9ddaSSatish Balay MatMatMultNumeric_SeqDense_SeqDense, 1852865e5f61SKris Buschelman 0, 1853865e5f61SKris Buschelman 0, 1854865e5f61SKris Buschelman /*95*/ 0, 1855a9fe9ddaSSatish Balay MatMatMultTranspose_SeqDense_SeqDense, 1856a9fe9ddaSSatish Balay MatMatMultTransposeSymbolic_SeqDense_SeqDense, 1857a9fe9ddaSSatish Balay MatMatMultTransposeNumeric_SeqDense_SeqDense, 1858284134d9SBarry Smith 0, 1859284134d9SBarry Smith /*100*/0, 1860284134d9SBarry Smith 0, 1861284134d9SBarry Smith 0, 1862284134d9SBarry Smith 0, 1863985db425SBarry Smith MatSetSizes_SeqDense, 1864985db425SBarry Smith 0, 1865985db425SBarry Smith 0, 1866985db425SBarry Smith 0, 1867985db425SBarry Smith 0, 1868985db425SBarry Smith 0, 1869985db425SBarry Smith /*110*/0, 1870985db425SBarry Smith 0, 18718d0534beSBarry Smith MatGetRowMin_SeqDense, 18728d0534beSBarry Smith MatGetColumnVector_SeqDense 1873985db425SBarry Smith }; 187490ace30eSBarry Smith 18754a2ae208SSatish Balay #undef __FUNCT__ 18764a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 18774b828684SBarry Smith /*@C 1878fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 1879d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 1880d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 1881289bc588SBarry Smith 1882db81eaa0SLois Curfman McInnes Collective on MPI_Comm 1883db81eaa0SLois Curfman McInnes 188420563c6bSBarry Smith Input Parameters: 1885db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 18860c775827SLois Curfman McInnes . m - number of rows 188718f449edSLois Curfman McInnes . n - number of columns 1888db81eaa0SLois Curfman McInnes - data - optional location of matrix data. Set data=PETSC_NULL for PETSc 1889dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 189020563c6bSBarry Smith 189120563c6bSBarry Smith Output Parameter: 189244cd7ae7SLois Curfman McInnes . A - the matrix 189320563c6bSBarry Smith 1894b259b22eSLois Curfman McInnes Notes: 189518f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 189618f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 1897b4fd4287SBarry Smith set data=PETSC_NULL. 189818f449edSLois Curfman McInnes 1899027ccd11SLois Curfman McInnes Level: intermediate 1900027ccd11SLois Curfman McInnes 1901dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 1902d65003e9SLois Curfman McInnes 1903db81eaa0SLois Curfman McInnes .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 190420563c6bSBarry Smith @*/ 1905be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 1906289bc588SBarry Smith { 1907dfbe8321SBarry Smith PetscErrorCode ierr; 19083b2fbd54SBarry Smith 19093a40ed3dSBarry Smith PetscFunctionBegin; 1910f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 1911f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 1912273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 1913273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 1914273d9f13SBarry Smith PetscFunctionReturn(0); 1915273d9f13SBarry Smith } 1916273d9f13SBarry Smith 19174a2ae208SSatish Balay #undef __FUNCT__ 1918afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation" 1919273d9f13SBarry Smith /*@C 1920273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 1921273d9f13SBarry Smith 1922273d9f13SBarry Smith Collective on MPI_Comm 1923273d9f13SBarry Smith 1924273d9f13SBarry Smith Input Parameters: 1925273d9f13SBarry Smith + A - the matrix 1926273d9f13SBarry Smith - data - the array (or PETSC_NULL) 1927273d9f13SBarry Smith 1928273d9f13SBarry Smith Notes: 1929273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 1930273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 1931284134d9SBarry Smith need not call this routine. 1932273d9f13SBarry Smith 1933273d9f13SBarry Smith Level: intermediate 1934273d9f13SBarry Smith 1935273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 1936273d9f13SBarry Smith 1937273d9f13SBarry Smith .seealso: MatCreate(), MatCreateMPIDense(), MatSetValues() 1938273d9f13SBarry Smith @*/ 1939be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 1940273d9f13SBarry Smith { 1941dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscScalar[]); 1942a23d5eceSKris Buschelman 1943a23d5eceSKris Buschelman PetscFunctionBegin; 1944a23d5eceSKris Buschelman ierr = PetscObjectQueryFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C",(void (**)(void))&f);CHKERRQ(ierr); 1945a23d5eceSKris Buschelman if (f) { 1946a23d5eceSKris Buschelman ierr = (*f)(B,data);CHKERRQ(ierr); 1947a23d5eceSKris Buschelman } 1948a23d5eceSKris Buschelman PetscFunctionReturn(0); 1949a23d5eceSKris Buschelman } 1950a23d5eceSKris Buschelman 1951a23d5eceSKris Buschelman EXTERN_C_BEGIN 1952a23d5eceSKris Buschelman #undef __FUNCT__ 1953afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense" 1954be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 1955a23d5eceSKris Buschelman { 1956273d9f13SBarry Smith Mat_SeqDense *b; 1957dfbe8321SBarry Smith PetscErrorCode ierr; 1958273d9f13SBarry Smith 1959273d9f13SBarry Smith PetscFunctionBegin; 1960273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 1961273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 1962d0f46423SBarry Smith if (b->lda <= 0) b->lda = B->rmap->n; 19639e8f95c4SLisandro Dalcin if (!data) { /* petsc-allocated storage */ 19649e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 19655afd5e0cSBarry Smith ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 1966284134d9SBarry Smith ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 1967284134d9SBarry Smith ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 19689e8f95c4SLisandro Dalcin b->user_alloc = PETSC_FALSE; 1969273d9f13SBarry Smith } else { /* user-allocated storage */ 19709e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 1971273d9f13SBarry Smith b->v = data; 1972273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 1973273d9f13SBarry Smith } 1974273d9f13SBarry Smith PetscFunctionReturn(0); 1975273d9f13SBarry Smith } 1976a23d5eceSKris Buschelman EXTERN_C_END 1977273d9f13SBarry Smith 19781b807ce4Svictorle #undef __FUNCT__ 19791b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 19801b807ce4Svictorle /*@C 19811b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 19821b807ce4Svictorle 19831b807ce4Svictorle Input parameter: 19841b807ce4Svictorle + A - the matrix 19851b807ce4Svictorle - lda - the leading dimension 19861b807ce4Svictorle 19871b807ce4Svictorle Notes: 19881b807ce4Svictorle This routine is to be used in conjunction with MatSeqDenseSetPreallocation; 19891b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 19901b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 19911b807ce4Svictorle 19921b807ce4Svictorle Level: intermediate 19931b807ce4Svictorle 19941b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 19951b807ce4Svictorle 1996284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 19971b807ce4Svictorle @*/ 1998be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatSeqDenseSetLDA(Mat B,PetscInt lda) 19991b807ce4Svictorle { 20001b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 200121a2c019SBarry Smith 20021b807ce4Svictorle PetscFunctionBegin; 2003d0f46423SBarry Smith if (lda < B->rmap->n) SETERRQ2(PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n); 20041b807ce4Svictorle b->lda = lda; 200521a2c019SBarry Smith b->changelda = PETSC_FALSE; 200621a2c019SBarry Smith b->Mmax = PetscMax(b->Mmax,lda); 20071b807ce4Svictorle PetscFunctionReturn(0); 20081b807ce4Svictorle } 20091b807ce4Svictorle 20100bad9183SKris Buschelman /*MC 2011fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 20120bad9183SKris Buschelman 20130bad9183SKris Buschelman Options Database Keys: 20140bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 20150bad9183SKris Buschelman 20160bad9183SKris Buschelman Level: beginner 20170bad9183SKris Buschelman 201889665df3SBarry Smith .seealso: MatCreateSeqDense() 201989665df3SBarry Smith 20200bad9183SKris Buschelman M*/ 20210bad9183SKris Buschelman 2022273d9f13SBarry Smith EXTERN_C_BEGIN 20234a2ae208SSatish Balay #undef __FUNCT__ 20244a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 2025be1d678aSKris Buschelman PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_SeqDense(Mat B) 2026273d9f13SBarry Smith { 2027273d9f13SBarry Smith Mat_SeqDense *b; 2028dfbe8321SBarry Smith PetscErrorCode ierr; 20297c334f02SBarry Smith PetscMPIInt size; 2030273d9f13SBarry Smith 2031273d9f13SBarry Smith PetscFunctionBegin; 20327adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)B)->comm,&size);CHKERRQ(ierr); 203329bbc08cSBarry Smith if (size > 1) SETERRQ(PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 203455659b69SBarry Smith 2035d0f46423SBarry Smith B->rmap->bs = B->cmap->bs = 1; 2036d0f46423SBarry Smith ierr = PetscMapSetUp(B->rmap);CHKERRQ(ierr); 2037d0f46423SBarry Smith ierr = PetscMapSetUp(B->cmap);CHKERRQ(ierr); 2038273d9f13SBarry Smith 203938f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr); 2040549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 204190f02eecSBarry Smith B->mapping = 0; 204244cd7ae7SLois Curfman McInnes B->data = (void*)b; 204318f449edSLois Curfman McInnes 2044a5ae1ecdSBarry Smith 204544cd7ae7SLois Curfman McInnes b->pivots = 0; 2046273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 2047273d9f13SBarry Smith b->v = 0; 2048d0f46423SBarry Smith b->lda = B->rmap->n; 204921a2c019SBarry Smith b->changelda = PETSC_FALSE; 2050d0f46423SBarry Smith b->Mmax = B->rmap->n; 2051d0f46423SBarry Smith b->Nmax = B->cmap->n; 20524e220ebcSLois Curfman McInnes 2053b24902e0SBarry Smith 2054b24902e0SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatGetFactor_seqdense_petsc_C", 2055b24902e0SBarry Smith "MatGetFactor_seqdense_petsc", 2056b24902e0SBarry Smith MatGetFactor_seqdense_petsc);CHKERRQ(ierr); 2057a23d5eceSKris Buschelman ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSeqDenseSetPreallocation_C", 2058a23d5eceSKris Buschelman "MatSeqDenseSetPreallocation_SeqDense", 2059a23d5eceSKris Buschelman MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 20604ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMult_seqaij_seqdense_C", 20614ae313f4SHong Zhang "MatMatMult_SeqAIJ_SeqDense", 20624ae313f4SHong Zhang MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 20634ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C", 20644ae313f4SHong Zhang "MatMatMultSymbolic_SeqAIJ_SeqDense", 20654ae313f4SHong Zhang MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 20664ae313f4SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C", 20674ae313f4SHong Zhang "MatMatMultNumeric_SeqAIJ_SeqDense", 20684ae313f4SHong Zhang MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 206917667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 20703a40ed3dSBarry Smith PetscFunctionReturn(0); 2071289bc588SBarry Smith } 2072c0aa2d19SHong Zhang 2073c0aa2d19SHong Zhang 2074273d9f13SBarry Smith EXTERN_C_END 2075