1be1d678aSKris Buschelman 267e560aaSBarry Smith /* 367e560aaSBarry Smith Defines the basic matrix operations for sequential dense. 467e560aaSBarry Smith */ 5289bc588SBarry Smith 6dec5eb66SMatthew G Knepley #include <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 7c6db04a5SJed Brown #include <petscblaslapack.h> 8289bc588SBarry Smith 96a63e612SBarry Smith #include <../src/mat/impls/aij/seq/aij.h> 10b2573a8aSBarry Smith 116a63e612SBarry Smith #undef __FUNCT__ 126a63e612SBarry Smith #define __FUNCT__ "MatConvert_SeqDense_SeqAIJ" 13*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatConvert_SeqDense_SeqAIJ(Mat A, MatType newtype,MatReuse reuse,Mat *newmat) 146a63e612SBarry Smith { 156a63e612SBarry Smith Mat B; 166a63e612SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 176a63e612SBarry Smith PetscErrorCode ierr; 186a63e612SBarry Smith PetscInt i; 196a63e612SBarry Smith PetscInt *rows; 206a63e612SBarry Smith MatScalar *aa = a->v; 216a63e612SBarry Smith 226a63e612SBarry Smith PetscFunctionBegin; 23ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&B);CHKERRQ(ierr); 246a63e612SBarry Smith ierr = MatSetSizes(B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 256a63e612SBarry Smith ierr = MatSetType(B,MATSEQAIJ);CHKERRQ(ierr); 260298fd71SBarry Smith ierr = MatSeqAIJSetPreallocation(B,A->cmap->n,NULL);CHKERRQ(ierr); 276a63e612SBarry Smith 286a63e612SBarry Smith ierr = PetscMalloc(A->rmap->n*sizeof(PetscInt),&rows);CHKERRQ(ierr); 296a63e612SBarry Smith for (i=0; i<A->rmap->n; i++) rows[i] = i; 306a63e612SBarry Smith 316a63e612SBarry Smith for (i=0; i<A->cmap->n; i++) { 326a63e612SBarry Smith ierr = MatSetValues(B,A->rmap->n,rows,1,&i,aa,INSERT_VALUES);CHKERRQ(ierr); 336a63e612SBarry Smith aa += a->lda; 346a63e612SBarry Smith } 356a63e612SBarry Smith ierr = PetscFree(rows);CHKERRQ(ierr); 366a63e612SBarry Smith ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 376a63e612SBarry Smith ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 386a63e612SBarry Smith 396a63e612SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 406a63e612SBarry Smith ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 416a63e612SBarry Smith } else { 426a63e612SBarry Smith *newmat = B; 436a63e612SBarry Smith } 446a63e612SBarry Smith PetscFunctionReturn(0); 456a63e612SBarry Smith } 466a63e612SBarry Smith 474a2ae208SSatish Balay #undef __FUNCT__ 484a2ae208SSatish Balay #define __FUNCT__ "MatAXPY_SeqDense" 49f4df32b1SMatthew Knepley PetscErrorCode MatAXPY_SeqDense(Mat Y,PetscScalar alpha,Mat X,MatStructure str) 501987afe7SBarry Smith { 511987afe7SBarry Smith Mat_SeqDense *x = (Mat_SeqDense*)X->data,*y = (Mat_SeqDense*)Y->data; 52f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 5313f74950SBarry Smith PetscInt j; 540805154bSBarry Smith PetscBLASInt N,m,ldax,lday,one = 1; 55efee365bSSatish Balay PetscErrorCode ierr; 563a40ed3dSBarry Smith 573a40ed3dSBarry Smith PetscFunctionBegin; 58c5df96a5SBarry Smith ierr = PetscBLASIntCast(X->rmap->n*X->cmap->n,&N);CHKERRQ(ierr); 59c5df96a5SBarry Smith ierr = PetscBLASIntCast(X->rmap->n,&m);CHKERRQ(ierr); 60c5df96a5SBarry Smith ierr = PetscBLASIntCast(x->lda,&ldax);CHKERRQ(ierr); 61c5df96a5SBarry Smith ierr = PetscBLASIntCast(y->lda,&lday);CHKERRQ(ierr); 62a5ce6ee0Svictorle if (ldax>m || lday>m) { 63d0f46423SBarry Smith for (j=0; j<X->cmap->n; j++) { 64a83cb05cSBarry Smith PetscStackCall("BLASaxpy",BLASaxpy_(&m,&oalpha,x->v+j*ldax,&one,y->v+j*lday,&one)); 65a5ce6ee0Svictorle } 66a5ce6ee0Svictorle } else { 67a83cb05cSBarry Smith PetscStackCall("BLASaxpy",BLASaxpy_(&N,&oalpha,x->v,&one,y->v,&one)); 68a5ce6ee0Svictorle } 690450473dSBarry Smith ierr = PetscLogFlops(PetscMax(2*N-1,0));CHKERRQ(ierr); 703a40ed3dSBarry Smith PetscFunctionReturn(0); 711987afe7SBarry Smith } 721987afe7SBarry Smith 734a2ae208SSatish Balay #undef __FUNCT__ 744a2ae208SSatish Balay #define __FUNCT__ "MatGetInfo_SeqDense" 75dfbe8321SBarry Smith PetscErrorCode MatGetInfo_SeqDense(Mat A,MatInfoType flag,MatInfo *info) 76289bc588SBarry Smith { 77d0f46423SBarry Smith PetscInt N = A->rmap->n*A->cmap->n; 783a40ed3dSBarry Smith 793a40ed3dSBarry Smith PetscFunctionBegin; 804e220ebcSLois Curfman McInnes info->block_size = 1.0; 814e220ebcSLois Curfman McInnes info->nz_allocated = (double)N; 826de62eeeSBarry Smith info->nz_used = (double)N; 836de62eeeSBarry Smith info->nz_unneeded = (double)0; 844e220ebcSLois Curfman McInnes info->assemblies = (double)A->num_ass; 854e220ebcSLois Curfman McInnes info->mallocs = 0; 867adad957SLisandro Dalcin info->memory = ((PetscObject)A)->mem; 874e220ebcSLois Curfman McInnes info->fill_ratio_given = 0; 884e220ebcSLois Curfman McInnes info->fill_ratio_needed = 0; 894e220ebcSLois Curfman McInnes info->factor_mallocs = 0; 903a40ed3dSBarry Smith PetscFunctionReturn(0); 91289bc588SBarry Smith } 92289bc588SBarry Smith 934a2ae208SSatish Balay #undef __FUNCT__ 944a2ae208SSatish Balay #define __FUNCT__ "MatScale_SeqDense" 95f4df32b1SMatthew Knepley PetscErrorCode MatScale_SeqDense(Mat A,PetscScalar alpha) 9680cd9d93SLois Curfman McInnes { 97273d9f13SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 98f4df32b1SMatthew Knepley PetscScalar oalpha = alpha; 99efee365bSSatish Balay PetscErrorCode ierr; 100c5df96a5SBarry Smith PetscBLASInt one = 1,j,nz,lda; 10180cd9d93SLois Curfman McInnes 1023a40ed3dSBarry Smith PetscFunctionBegin; 103c5df96a5SBarry Smith ierr = PetscBLASIntCast(a->lda,&lda);CHKERRQ(ierr); 104d0f46423SBarry Smith if (lda>A->rmap->n) { 105c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&nz);CHKERRQ(ierr); 106d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 107a83cb05cSBarry Smith PetscStackCall("BLASscal",BLASscal_(&nz,&oalpha,a->v+j*lda,&one)); 108a5ce6ee0Svictorle } 109a5ce6ee0Svictorle } else { 110c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n*A->cmap->n,&nz);CHKERRQ(ierr); 111a83cb05cSBarry Smith PetscStackCall("BLASscal",BLASscal_(&nz,&oalpha,a->v,&one)); 112a5ce6ee0Svictorle } 113efee365bSSatish Balay ierr = PetscLogFlops(nz);CHKERRQ(ierr); 1143a40ed3dSBarry Smith PetscFunctionReturn(0); 11580cd9d93SLois Curfman McInnes } 11680cd9d93SLois Curfman McInnes 1171cbb95d3SBarry Smith #undef __FUNCT__ 1181cbb95d3SBarry Smith #define __FUNCT__ "MatIsHermitian_SeqDense" 119ace3abfcSBarry Smith PetscErrorCode MatIsHermitian_SeqDense(Mat A,PetscReal rtol,PetscBool *fl) 1201cbb95d3SBarry Smith { 1211cbb95d3SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 122d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,N; 1231cbb95d3SBarry Smith PetscScalar *v = a->v; 1241cbb95d3SBarry Smith 1251cbb95d3SBarry Smith PetscFunctionBegin; 1261cbb95d3SBarry Smith *fl = PETSC_FALSE; 127d0f46423SBarry Smith if (A->rmap->n != A->cmap->n) PetscFunctionReturn(0); 1281cbb95d3SBarry Smith N = a->lda; 1291cbb95d3SBarry Smith 1301cbb95d3SBarry Smith for (i=0; i<m; i++) { 1311cbb95d3SBarry Smith for (j=i+1; j<m; j++) { 1321cbb95d3SBarry Smith if (PetscAbsScalar(v[i+j*N] - PetscConj(v[j+i*N])) > rtol) PetscFunctionReturn(0); 1331cbb95d3SBarry Smith } 1341cbb95d3SBarry Smith } 1351cbb95d3SBarry Smith *fl = PETSC_TRUE; 1361cbb95d3SBarry Smith PetscFunctionReturn(0); 1371cbb95d3SBarry Smith } 1381cbb95d3SBarry Smith 139b24902e0SBarry Smith #undef __FUNCT__ 140b24902e0SBarry Smith #define __FUNCT__ "MatDuplicateNoCreate_SeqDense" 141719d5645SBarry Smith PetscErrorCode MatDuplicateNoCreate_SeqDense(Mat newi,Mat A,MatDuplicateOption cpvalues) 142b24902e0SBarry Smith { 143b24902e0SBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data,*l; 144b24902e0SBarry Smith PetscErrorCode ierr; 145b24902e0SBarry Smith PetscInt lda = (PetscInt)mat->lda,j,m; 146b24902e0SBarry Smith 147b24902e0SBarry Smith PetscFunctionBegin; 148aa5ea44dSBarry Smith ierr = PetscLayoutReference(A->rmap,&newi->rmap);CHKERRQ(ierr); 149aa5ea44dSBarry Smith ierr = PetscLayoutReference(A->cmap,&newi->cmap);CHKERRQ(ierr); 1500298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(newi,NULL);CHKERRQ(ierr); 151b24902e0SBarry Smith if (cpvalues == MAT_COPY_VALUES) { 152b24902e0SBarry Smith l = (Mat_SeqDense*)newi->data; 153d0f46423SBarry Smith if (lda>A->rmap->n) { 154d0f46423SBarry Smith m = A->rmap->n; 155d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 156b24902e0SBarry Smith ierr = PetscMemcpy(l->v+j*m,mat->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 157b24902e0SBarry Smith } 158b24902e0SBarry Smith } else { 159d0f46423SBarry Smith ierr = PetscMemcpy(l->v,mat->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 160b24902e0SBarry Smith } 161b24902e0SBarry Smith } 162b24902e0SBarry Smith newi->assembled = PETSC_TRUE; 163b24902e0SBarry Smith PetscFunctionReturn(0); 164b24902e0SBarry Smith } 165b24902e0SBarry Smith 1664a2ae208SSatish Balay #undef __FUNCT__ 1674a2ae208SSatish Balay #define __FUNCT__ "MatDuplicate_SeqDense" 168dfbe8321SBarry Smith PetscErrorCode MatDuplicate_SeqDense(Mat A,MatDuplicateOption cpvalues,Mat *newmat) 16902cad45dSBarry Smith { 1706849ba73SBarry Smith PetscErrorCode ierr; 17102cad45dSBarry Smith 1723a40ed3dSBarry Smith PetscFunctionBegin; 173ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),newmat);CHKERRQ(ierr); 174d0f46423SBarry Smith ierr = MatSetSizes(*newmat,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 1755c9eb25fSBarry Smith ierr = MatSetType(*newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 176719d5645SBarry Smith ierr = MatDuplicateNoCreate_SeqDense(*newmat,A,cpvalues);CHKERRQ(ierr); 177b24902e0SBarry Smith PetscFunctionReturn(0); 178b24902e0SBarry Smith } 179b24902e0SBarry Smith 1806ee01492SSatish Balay 1810481f469SBarry Smith extern PetscErrorCode MatLUFactor_SeqDense(Mat,IS,IS,const MatFactorInfo*); 182719d5645SBarry Smith 1834a2ae208SSatish Balay #undef __FUNCT__ 1844a2ae208SSatish Balay #define __FUNCT__ "MatLUFactorNumeric_SeqDense" 1850481f469SBarry Smith PetscErrorCode MatLUFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 186289bc588SBarry Smith { 1874482741eSBarry Smith MatFactorInfo info; 188a093e273SMatthew Knepley PetscErrorCode ierr; 1893a40ed3dSBarry Smith 1903a40ed3dSBarry Smith PetscFunctionBegin; 191c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 192719d5645SBarry Smith ierr = MatLUFactor_SeqDense(fact,0,0,&info);CHKERRQ(ierr); 1933a40ed3dSBarry Smith PetscFunctionReturn(0); 194289bc588SBarry Smith } 1956ee01492SSatish Balay 1960b4b3355SBarry Smith #undef __FUNCT__ 1974a2ae208SSatish Balay #define __FUNCT__ "MatSolve_SeqDense" 198dfbe8321SBarry Smith PetscErrorCode MatSolve_SeqDense(Mat A,Vec xx,Vec yy) 199289bc588SBarry Smith { 200c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 2016849ba73SBarry Smith PetscErrorCode ierr; 20287828ca2SBarry Smith PetscScalar *x,*y; 203c5df96a5SBarry Smith PetscBLASInt one = 1,info,m; 20467e560aaSBarry Smith 2053a40ed3dSBarry Smith PetscFunctionBegin; 206c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 2071ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2081ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 209d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 210d5f3da31SBarry Smith if (A->factortype == MAT_FACTOR_LU) { 211ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 212e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 213ae7cfcebSSatish Balay #else 214f75e95b9SBarry Smith PetscStackCall("LAPACKgetrs",LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info)); 215e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 216ae7cfcebSSatish Balay #endif 217d5f3da31SBarry Smith } else if (A->factortype == MAT_FACTOR_CHOLESKY) { 218ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 219e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 220ae7cfcebSSatish Balay #else 221f75e95b9SBarry Smith PetscStackCall("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info)); 222e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 223ae7cfcebSSatish Balay #endif 2242205254eSKarl Rupp } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 2251ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2261ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 227dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 2283a40ed3dSBarry Smith PetscFunctionReturn(0); 229289bc588SBarry Smith } 2306ee01492SSatish Balay 2314a2ae208SSatish Balay #undef __FUNCT__ 23285e2c93fSHong Zhang #define __FUNCT__ "MatMatSolve_SeqDense" 23385e2c93fSHong Zhang PetscErrorCode MatMatSolve_SeqDense(Mat A,Mat B,Mat X) 23485e2c93fSHong Zhang { 23585e2c93fSHong Zhang Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 23685e2c93fSHong Zhang PetscErrorCode ierr; 23785e2c93fSHong Zhang PetscScalar *b,*x; 238efb80c78SLisandro Dalcin PetscInt n; 239783b601eSJed Brown PetscBLASInt nrhs,info,m; 240bda8bf91SBarry Smith PetscBool flg; 24185e2c93fSHong Zhang 24285e2c93fSHong Zhang PetscFunctionBegin; 243c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 2440298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)B,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 245ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix B must be MATDENSE matrix"); 2460298fd71SBarry Smith ierr = PetscObjectTypeCompareAny((PetscObject)X,&flg,MATSEQDENSE,MATMPIDENSE,NULL);CHKERRQ(ierr); 247ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix X must be MATDENSE matrix"); 248bda8bf91SBarry Smith 2490298fd71SBarry Smith ierr = MatGetSize(B,NULL,&n);CHKERRQ(ierr); 250c5df96a5SBarry Smith ierr = PetscBLASIntCast(n,&nrhs);CHKERRQ(ierr); 2518c778c55SBarry Smith ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 2528c778c55SBarry Smith ierr = MatDenseGetArray(X,&x);CHKERRQ(ierr); 25385e2c93fSHong Zhang 254f272c775SHong Zhang ierr = PetscMemcpy(x,b,m*nrhs*sizeof(PetscScalar));CHKERRQ(ierr); 25585e2c93fSHong Zhang 25685e2c93fSHong Zhang if (A->factortype == MAT_FACTOR_LU) { 25785e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_GETRS) 25885e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 25985e2c93fSHong Zhang #else 260f75e95b9SBarry Smith PetscStackCall("LAPACKgetrs",LAPACKgetrs_("N",&m,&nrhs,mat->v,&mat->lda,mat->pivots,x,&m,&info)); 26185e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"GETRS - Bad solve"); 26285e2c93fSHong Zhang #endif 26385e2c93fSHong Zhang } else if (A->factortype == MAT_FACTOR_CHOLESKY) { 26485e2c93fSHong Zhang #if defined(PETSC_MISSING_LAPACK_POTRS) 26585e2c93fSHong Zhang SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 26685e2c93fSHong Zhang #else 267f75e95b9SBarry Smith PetscStackCall("LAPACKpotrs",LAPACKpotrs_("L",&m,&nrhs,mat->v,&mat->lda,x,&m,&info)); 26885e2c93fSHong Zhang if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS Bad solve"); 26985e2c93fSHong Zhang #endif 2702205254eSKarl Rupp } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix must be factored to solve"); 27185e2c93fSHong Zhang 2728c778c55SBarry Smith ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 2738c778c55SBarry Smith ierr = MatDenseRestoreArray(X,&x);CHKERRQ(ierr); 27485e2c93fSHong Zhang ierr = PetscLogFlops(nrhs*(2.0*m*m - m));CHKERRQ(ierr); 27585e2c93fSHong Zhang PetscFunctionReturn(0); 27685e2c93fSHong Zhang } 27785e2c93fSHong Zhang 27885e2c93fSHong Zhang #undef __FUNCT__ 2794a2ae208SSatish Balay #define __FUNCT__ "MatSolveTranspose_SeqDense" 280dfbe8321SBarry Smith PetscErrorCode MatSolveTranspose_SeqDense(Mat A,Vec xx,Vec yy) 281da3a660dSBarry Smith { 282c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 283dfbe8321SBarry Smith PetscErrorCode ierr; 28487828ca2SBarry Smith PetscScalar *x,*y; 285c5df96a5SBarry Smith PetscBLASInt one = 1,info,m; 28667e560aaSBarry Smith 2873a40ed3dSBarry Smith PetscFunctionBegin; 288c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 2891ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2901ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 291d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 292752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 293da3a660dSBarry Smith if (mat->pivots) { 294ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 295e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 296ae7cfcebSSatish Balay #else 297f75e95b9SBarry Smith PetscStackCall("LAPACKgetrs",LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info)); 298e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 299ae7cfcebSSatish Balay #endif 3007a97a34bSBarry Smith } else { 301ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 302e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 303ae7cfcebSSatish Balay #else 304f75e95b9SBarry Smith PetscStackCall("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info)); 305e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"POTRS - Bad solve"); 306ae7cfcebSSatish Balay #endif 307da3a660dSBarry Smith } 3081ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3091ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 310dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 3113a40ed3dSBarry Smith PetscFunctionReturn(0); 312da3a660dSBarry Smith } 3136ee01492SSatish Balay 3144a2ae208SSatish Balay #undef __FUNCT__ 3154a2ae208SSatish Balay #define __FUNCT__ "MatSolveAdd_SeqDense" 316dfbe8321SBarry Smith PetscErrorCode MatSolveAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 317da3a660dSBarry Smith { 318c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 319dfbe8321SBarry Smith PetscErrorCode ierr; 32087828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 321da3a660dSBarry Smith Vec tmp = 0; 322c5df96a5SBarry Smith PetscBLASInt one = 1,info,m; 32367e560aaSBarry Smith 3243a40ed3dSBarry Smith PetscFunctionBegin; 325c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 3261ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3271ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 328d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 329da3a660dSBarry Smith if (yy == zz) { 33078b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 33152e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 33278b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 333da3a660dSBarry Smith } 334d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 335752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 336da3a660dSBarry Smith if (mat->pivots) { 337ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 338e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 339ae7cfcebSSatish Balay #else 340f75e95b9SBarry Smith PetscStackCall("LAPACKgetrs",LAPACKgetrs_("N",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info)); 341e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 342ae7cfcebSSatish Balay #endif 343a8c6a408SBarry Smith } else { 344ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 345e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 346ae7cfcebSSatish Balay #else 347f75e95b9SBarry Smith PetscStackCall("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info)); 348e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 349ae7cfcebSSatish Balay #endif 350da3a660dSBarry Smith } 3516bf464f9SBarry Smith if (tmp) { 3526bf464f9SBarry Smith ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 3536bf464f9SBarry Smith ierr = VecDestroy(&tmp);CHKERRQ(ierr); 3546bf464f9SBarry Smith } else { 3556bf464f9SBarry Smith ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 3566bf464f9SBarry Smith } 3571ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 3581ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 359dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 3603a40ed3dSBarry Smith PetscFunctionReturn(0); 361da3a660dSBarry Smith } 36267e560aaSBarry Smith 3634a2ae208SSatish Balay #undef __FUNCT__ 3644a2ae208SSatish Balay #define __FUNCT__ "MatSolveTransposeAdd_SeqDense" 365dfbe8321SBarry Smith PetscErrorCode MatSolveTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 366da3a660dSBarry Smith { 367c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 3686849ba73SBarry Smith PetscErrorCode ierr; 36987828ca2SBarry Smith PetscScalar *x,*y,sone = 1.0; 370da3a660dSBarry Smith Vec tmp; 371c5df96a5SBarry Smith PetscBLASInt one = 1,info,m; 37267e560aaSBarry Smith 3733a40ed3dSBarry Smith PetscFunctionBegin; 374c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 375d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 3761ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 3771ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 378da3a660dSBarry Smith if (yy == zz) { 37978b31e54SBarry Smith ierr = VecDuplicate(yy,&tmp);CHKERRQ(ierr); 38052e6d16bSBarry Smith ierr = PetscLogObjectParent(A,tmp);CHKERRQ(ierr); 38178b31e54SBarry Smith ierr = VecCopy(yy,tmp);CHKERRQ(ierr); 382da3a660dSBarry Smith } 383d0f46423SBarry Smith ierr = PetscMemcpy(y,x,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 384752f5567SLois Curfman McInnes /* assume if pivots exist then use LU; else Cholesky */ 385da3a660dSBarry Smith if (mat->pivots) { 386ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_GETRS) 387e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routine is unavailable."); 388ae7cfcebSSatish Balay #else 389f75e95b9SBarry Smith PetscStackCall("LAPACKgetrs",LAPACKgetrs_("T",&m,&one,mat->v,&mat->lda,mat->pivots,y,&m,&info)); 390e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 391ae7cfcebSSatish Balay #endif 3923a40ed3dSBarry Smith } else { 393ae7cfcebSSatish Balay #if defined(PETSC_MISSING_LAPACK_POTRS) 394e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRS - Lapack routine is unavailable."); 395ae7cfcebSSatish Balay #else 396f75e95b9SBarry Smith PetscStackCall("LAPACKpotrs",LAPACKpotrs_("L",&m,&one,mat->v,&mat->lda,y,&m,&info)); 397e32f2f54SBarry Smith if (info) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad solve"); 398ae7cfcebSSatish Balay #endif 399da3a660dSBarry Smith } 40090f02eecSBarry Smith if (tmp) { 4012dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,tmp);CHKERRQ(ierr); 4026bf464f9SBarry Smith ierr = VecDestroy(&tmp);CHKERRQ(ierr); 4033a40ed3dSBarry Smith } else { 4042dcb1b2aSMatthew Knepley ierr = VecAXPY(yy,sone,zz);CHKERRQ(ierr); 40590f02eecSBarry Smith } 4061ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 4071ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 408dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->cmap->n);CHKERRQ(ierr); 4093a40ed3dSBarry Smith PetscFunctionReturn(0); 410da3a660dSBarry Smith } 411db4efbfdSBarry Smith 412db4efbfdSBarry Smith /* ---------------------------------------------------------------*/ 413db4efbfdSBarry Smith /* COMMENT: I have chosen to hide row permutation in the pivots, 414db4efbfdSBarry Smith rather than put it in the Mat->row slot.*/ 415db4efbfdSBarry Smith #undef __FUNCT__ 416db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactor_SeqDense" 4170481f469SBarry Smith PetscErrorCode MatLUFactor_SeqDense(Mat A,IS row,IS col,const MatFactorInfo *minfo) 418db4efbfdSBarry Smith { 419db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_GETRF) 420db4efbfdSBarry Smith PetscFunctionBegin; 421e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable."); 422db4efbfdSBarry Smith #else 423db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 424db4efbfdSBarry Smith PetscErrorCode ierr; 425db4efbfdSBarry Smith PetscBLASInt n,m,info; 426db4efbfdSBarry Smith 427db4efbfdSBarry Smith PetscFunctionBegin; 428c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 429c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 430db4efbfdSBarry Smith if (!mat->pivots) { 431db4efbfdSBarry Smith ierr = PetscMalloc((A->rmap->n+1)*sizeof(PetscBLASInt),&mat->pivots);CHKERRQ(ierr); 432db4efbfdSBarry Smith ierr = PetscLogObjectMemory(A,A->rmap->n*sizeof(PetscBLASInt));CHKERRQ(ierr); 433db4efbfdSBarry Smith } 434db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 4358e57ea43SSatish Balay ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 436f75e95b9SBarry Smith PetscStackCall("LAPACKgetrf",LAPACKgetrf_(&m,&n,mat->v,&mat->lda,mat->pivots,&info)); 4378e57ea43SSatish Balay ierr = PetscFPTrapPop();CHKERRQ(ierr); 4388e57ea43SSatish Balay 439e32f2f54SBarry Smith if (info<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Bad argument to LU factorization"); 440e32f2f54SBarry Smith if (info>0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Bad LU factorization"); 441db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 442db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 443db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 444db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 445d5f3da31SBarry Smith A->factortype = MAT_FACTOR_LU; 446db4efbfdSBarry Smith 447dc0b31edSSatish Balay ierr = PetscLogFlops((2.0*A->cmap->n*A->cmap->n*A->cmap->n)/3);CHKERRQ(ierr); 448db4efbfdSBarry Smith #endif 449db4efbfdSBarry Smith PetscFunctionReturn(0); 450db4efbfdSBarry Smith } 451db4efbfdSBarry Smith 452db4efbfdSBarry Smith #undef __FUNCT__ 453db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactor_SeqDense" 4540481f469SBarry Smith PetscErrorCode MatCholeskyFactor_SeqDense(Mat A,IS perm,const MatFactorInfo *factinfo) 455db4efbfdSBarry Smith { 456db4efbfdSBarry Smith #if defined(PETSC_MISSING_LAPACK_POTRF) 457db4efbfdSBarry Smith PetscFunctionBegin; 458e32f2f54SBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"POTRF - Lapack routine is unavailable."); 459db4efbfdSBarry Smith #else 460db4efbfdSBarry Smith Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 461db4efbfdSBarry Smith PetscErrorCode ierr; 462c5df96a5SBarry Smith PetscBLASInt info,n; 463db4efbfdSBarry Smith 464db4efbfdSBarry Smith PetscFunctionBegin; 465c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 466db4efbfdSBarry Smith ierr = PetscFree(mat->pivots);CHKERRQ(ierr); 467db4efbfdSBarry Smith 468db4efbfdSBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 469f75e95b9SBarry Smith PetscStackCall("LAPACKpotrf",LAPACKpotrf_("L",&n,mat->v,&mat->lda,&info)); 470e32f2f54SBarry Smith if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_CH_ZRPVT,"Bad factorization: zero pivot in row %D",(PetscInt)info-1); 471db4efbfdSBarry Smith A->ops->solve = MatSolve_SeqDense; 472db4efbfdSBarry Smith A->ops->solvetranspose = MatSolveTranspose_SeqDense; 473db4efbfdSBarry Smith A->ops->solveadd = MatSolveAdd_SeqDense; 474db4efbfdSBarry Smith A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqDense; 475d5f3da31SBarry Smith A->factortype = MAT_FACTOR_CHOLESKY; 4762205254eSKarl Rupp 477dc0b31edSSatish Balay ierr = PetscLogFlops((A->cmap->n*A->cmap->n*A->cmap->n)/3.0);CHKERRQ(ierr); 478db4efbfdSBarry Smith #endif 479db4efbfdSBarry Smith PetscFunctionReturn(0); 480db4efbfdSBarry Smith } 481db4efbfdSBarry Smith 482db4efbfdSBarry Smith 483db4efbfdSBarry Smith #undef __FUNCT__ 484db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorNumeric_SeqDense" 4850481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqDense(Mat fact,Mat A,const MatFactorInfo *info_dummy) 486db4efbfdSBarry Smith { 487db4efbfdSBarry Smith PetscErrorCode ierr; 488db4efbfdSBarry Smith MatFactorInfo info; 489db4efbfdSBarry Smith 490db4efbfdSBarry Smith PetscFunctionBegin; 491db4efbfdSBarry Smith info.fill = 1.0; 4922205254eSKarl Rupp 493c3ef05f6SHong Zhang ierr = MatDuplicateNoCreate_SeqDense(fact,A,MAT_COPY_VALUES);CHKERRQ(ierr); 494719d5645SBarry Smith ierr = MatCholeskyFactor_SeqDense(fact,0,&info);CHKERRQ(ierr); 495db4efbfdSBarry Smith PetscFunctionReturn(0); 496db4efbfdSBarry Smith } 497db4efbfdSBarry Smith 498db4efbfdSBarry Smith #undef __FUNCT__ 499db4efbfdSBarry Smith #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqDense" 5000481f469SBarry Smith PetscErrorCode MatCholeskyFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,const MatFactorInfo *info) 501db4efbfdSBarry Smith { 502db4efbfdSBarry Smith PetscFunctionBegin; 503c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 5041bbcc794SSatish Balay fact->preallocated = PETSC_TRUE; 505719d5645SBarry Smith fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqDense; 506db4efbfdSBarry Smith PetscFunctionReturn(0); 507db4efbfdSBarry Smith } 508db4efbfdSBarry Smith 509db4efbfdSBarry Smith #undef __FUNCT__ 510db4efbfdSBarry Smith #define __FUNCT__ "MatLUFactorSymbolic_SeqDense" 5110481f469SBarry Smith PetscErrorCode MatLUFactorSymbolic_SeqDense(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 512db4efbfdSBarry Smith { 513db4efbfdSBarry Smith PetscFunctionBegin; 514b66fe19dSMatthew G Knepley fact->preallocated = PETSC_TRUE; 515c3ef05f6SHong Zhang fact->assembled = PETSC_TRUE; 516719d5645SBarry Smith fact->ops->lufactornumeric = MatLUFactorNumeric_SeqDense; 517db4efbfdSBarry Smith PetscFunctionReturn(0); 518db4efbfdSBarry Smith } 519db4efbfdSBarry Smith 520db4efbfdSBarry Smith #undef __FUNCT__ 521db4efbfdSBarry Smith #define __FUNCT__ "MatGetFactor_seqdense_petsc" 522*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatGetFactor_seqdense_petsc(Mat A,MatFactorType ftype,Mat *fact) 523db4efbfdSBarry Smith { 524db4efbfdSBarry Smith PetscErrorCode ierr; 525db4efbfdSBarry Smith 526db4efbfdSBarry Smith PetscFunctionBegin; 527ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),fact);CHKERRQ(ierr); 528db4efbfdSBarry Smith ierr = MatSetSizes(*fact,A->rmap->n,A->cmap->n,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 529db4efbfdSBarry Smith ierr = MatSetType(*fact,((PetscObject)A)->type_name);CHKERRQ(ierr); 530db4efbfdSBarry Smith if (ftype == MAT_FACTOR_LU) { 531db4efbfdSBarry Smith (*fact)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqDense; 532db4efbfdSBarry Smith } else { 533db4efbfdSBarry Smith (*fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqDense; 534db4efbfdSBarry Smith } 535d5f3da31SBarry Smith (*fact)->factortype = ftype; 536db4efbfdSBarry Smith PetscFunctionReturn(0); 537db4efbfdSBarry Smith } 538db4efbfdSBarry Smith 539289bc588SBarry Smith /* ------------------------------------------------------------------*/ 5404a2ae208SSatish Balay #undef __FUNCT__ 54141f059aeSBarry Smith #define __FUNCT__ "MatSOR_SeqDense" 54241f059aeSBarry Smith PetscErrorCode MatSOR_SeqDense(Mat A,Vec bb,PetscReal omega,MatSORType flag,PetscReal shift,PetscInt its,PetscInt lits,Vec xx) 543289bc588SBarry Smith { 544c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 54587828ca2SBarry Smith PetscScalar *x,*b,*v = mat->v,zero = 0.0,xt; 546dfbe8321SBarry Smith PetscErrorCode ierr; 547d0f46423SBarry Smith PetscInt m = A->rmap->n,i; 548c5df96a5SBarry Smith PetscBLASInt o = 1,bm; 549289bc588SBarry Smith 5503a40ed3dSBarry Smith PetscFunctionBegin; 551c5df96a5SBarry Smith ierr = PetscBLASIntCast(m,&bm);CHKERRQ(ierr); 552289bc588SBarry Smith if (flag & SOR_ZERO_INITIAL_GUESS) { 55371044d3cSBarry Smith /* this is a hack fix, should have another version without the second BLASdot */ 5542dcb1b2aSMatthew Knepley ierr = VecSet(xx,zero);CHKERRQ(ierr); 555289bc588SBarry Smith } 5561ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5571ebc52fbSHong Zhang ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 558b965ef7fSBarry Smith its = its*lits; 559e32f2f54SBarry Smith if (its <= 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Relaxation requires global its %D and local its %D both positive",its,lits); 560289bc588SBarry Smith while (its--) { 561fccaa45eSBarry Smith if (flag & SOR_FORWARD_SWEEP || flag & SOR_LOCAL_FORWARD_SWEEP) { 562289bc588SBarry Smith for (i=0; i<m; i++) { 563a83cb05cSBarry Smith PetscStackCall("BLASdot",xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o)); 56455a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 565289bc588SBarry Smith } 566289bc588SBarry Smith } 567fccaa45eSBarry Smith if (flag & SOR_BACKWARD_SWEEP || flag & SOR_LOCAL_BACKWARD_SWEEP) { 568289bc588SBarry Smith for (i=m-1; i>=0; i--) { 569a83cb05cSBarry Smith PetscStackCall("BLASdot",xt = b[i] - BLASdot_(&bm,v+i,&bm,x,&o)); 57055a1b374SBarry Smith x[i] = (1. - omega)*x[i] + omega*(xt+v[i + i*m]*x[i])/(v[i + i*m]+shift); 571289bc588SBarry Smith } 572289bc588SBarry Smith } 573289bc588SBarry Smith } 5741ebc52fbSHong Zhang ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 5751ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5763a40ed3dSBarry Smith PetscFunctionReturn(0); 577289bc588SBarry Smith } 578289bc588SBarry Smith 579289bc588SBarry Smith /* -----------------------------------------------------------------*/ 5804a2ae208SSatish Balay #undef __FUNCT__ 5814a2ae208SSatish Balay #define __FUNCT__ "MatMultTranspose_SeqDense" 582dfbe8321SBarry Smith PetscErrorCode MatMultTranspose_SeqDense(Mat A,Vec xx,Vec yy) 583289bc588SBarry Smith { 584c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 58587828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 586dfbe8321SBarry Smith PetscErrorCode ierr; 5870805154bSBarry Smith PetscBLASInt m, n,_One=1; 588ea709b57SSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 5893a40ed3dSBarry Smith 5903a40ed3dSBarry Smith PetscFunctionBegin; 591c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 592c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 593d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 5941ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 5951ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 596a83cb05cSBarry Smith PetscStackCall("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&mat->lda,x,&_One,&_DZero,y,&_One)); 5971ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 5981ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 599dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->cmap->n);CHKERRQ(ierr); 6003a40ed3dSBarry Smith PetscFunctionReturn(0); 601289bc588SBarry Smith } 602800995b7SMatthew Knepley 6034a2ae208SSatish Balay #undef __FUNCT__ 6044a2ae208SSatish Balay #define __FUNCT__ "MatMult_SeqDense" 605dfbe8321SBarry Smith PetscErrorCode MatMult_SeqDense(Mat A,Vec xx,Vec yy) 606289bc588SBarry Smith { 607c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 60887828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0,_DZero=0.0; 609dfbe8321SBarry Smith PetscErrorCode ierr; 6100805154bSBarry Smith PetscBLASInt m, n, _One=1; 6113a40ed3dSBarry Smith 6123a40ed3dSBarry Smith PetscFunctionBegin; 613c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 614c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 615d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 6161ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6171ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 618a83cb05cSBarry Smith PetscStackCall("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DZero,y,&_One)); 6191ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6201ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 621dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n - A->rmap->n);CHKERRQ(ierr); 6223a40ed3dSBarry Smith PetscFunctionReturn(0); 623289bc588SBarry Smith } 6246ee01492SSatish Balay 6254a2ae208SSatish Balay #undef __FUNCT__ 6264a2ae208SSatish Balay #define __FUNCT__ "MatMultAdd_SeqDense" 627dfbe8321SBarry Smith PetscErrorCode MatMultAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 628289bc588SBarry Smith { 629c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 63087828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y,_DOne=1.0; 631dfbe8321SBarry Smith PetscErrorCode ierr; 6320805154bSBarry Smith PetscBLASInt m, n, _One=1; 6333a40ed3dSBarry Smith 6343a40ed3dSBarry Smith PetscFunctionBegin; 635c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 636c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 637d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 638600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6391ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6401ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 641a83cb05cSBarry Smith PetscStackCall("BLASgemv",BLASgemv_("N",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One)); 6421ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6431ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 644dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6453a40ed3dSBarry Smith PetscFunctionReturn(0); 646289bc588SBarry Smith } 6476ee01492SSatish Balay 6484a2ae208SSatish Balay #undef __FUNCT__ 6494a2ae208SSatish Balay #define __FUNCT__ "MatMultTransposeAdd_SeqDense" 650dfbe8321SBarry Smith PetscErrorCode MatMultTransposeAdd_SeqDense(Mat A,Vec xx,Vec zz,Vec yy) 651289bc588SBarry Smith { 652c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 65387828ca2SBarry Smith PetscScalar *v = mat->v,*x,*y; 654dfbe8321SBarry Smith PetscErrorCode ierr; 6550805154bSBarry Smith PetscBLASInt m, n, _One=1; 65687828ca2SBarry Smith PetscScalar _DOne=1.0; 6573a40ed3dSBarry Smith 6583a40ed3dSBarry Smith PetscFunctionBegin; 659c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 660c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&n);CHKERRQ(ierr); 661d0f46423SBarry Smith if (!A->rmap->n || !A->cmap->n) PetscFunctionReturn(0); 662600d6d0bSBarry Smith if (zz != yy) {ierr = VecCopy(zz,yy);CHKERRQ(ierr);} 6631ebc52fbSHong Zhang ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 6641ebc52fbSHong Zhang ierr = VecGetArray(yy,&y);CHKERRQ(ierr); 665a83cb05cSBarry Smith PetscStackCall("BLASgemv",BLASgemv_("T",&m,&n,&_DOne,v,&(mat->lda),x,&_One,&_DOne,y,&_One)); 6661ebc52fbSHong Zhang ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 6671ebc52fbSHong Zhang ierr = VecRestoreArray(yy,&y);CHKERRQ(ierr); 668dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->rmap->n*A->cmap->n);CHKERRQ(ierr); 6693a40ed3dSBarry Smith PetscFunctionReturn(0); 670289bc588SBarry Smith } 671289bc588SBarry Smith 672289bc588SBarry Smith /* -----------------------------------------------------------------*/ 6734a2ae208SSatish Balay #undef __FUNCT__ 6744a2ae208SSatish Balay #define __FUNCT__ "MatGetRow_SeqDense" 67513f74950SBarry Smith PetscErrorCode MatGetRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 676289bc588SBarry Smith { 677c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 67887828ca2SBarry Smith PetscScalar *v; 6796849ba73SBarry Smith PetscErrorCode ierr; 68013f74950SBarry Smith PetscInt i; 68167e560aaSBarry Smith 6823a40ed3dSBarry Smith PetscFunctionBegin; 683d0f46423SBarry Smith *ncols = A->cmap->n; 684289bc588SBarry Smith if (cols) { 685d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscInt),cols);CHKERRQ(ierr); 686d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) (*cols)[i] = i; 687289bc588SBarry Smith } 688289bc588SBarry Smith if (vals) { 689d0f46423SBarry Smith ierr = PetscMalloc((A->cmap->n+1)*sizeof(PetscScalar),vals);CHKERRQ(ierr); 690289bc588SBarry Smith v = mat->v + row; 691d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) {(*vals)[i] = *v; v += mat->lda;} 692289bc588SBarry Smith } 6933a40ed3dSBarry Smith PetscFunctionReturn(0); 694289bc588SBarry Smith } 6956ee01492SSatish Balay 6964a2ae208SSatish Balay #undef __FUNCT__ 6974a2ae208SSatish Balay #define __FUNCT__ "MatRestoreRow_SeqDense" 69813f74950SBarry Smith PetscErrorCode MatRestoreRow_SeqDense(Mat A,PetscInt row,PetscInt *ncols,PetscInt **cols,PetscScalar **vals) 699289bc588SBarry Smith { 700dfbe8321SBarry Smith PetscErrorCode ierr; 7016e111a19SKarl Rupp 702606d414cSSatish Balay PetscFunctionBegin; 703606d414cSSatish Balay if (cols) {ierr = PetscFree(*cols);CHKERRQ(ierr);} 704606d414cSSatish Balay if (vals) {ierr = PetscFree(*vals);CHKERRQ(ierr); } 7053a40ed3dSBarry Smith PetscFunctionReturn(0); 706289bc588SBarry Smith } 707289bc588SBarry Smith /* ----------------------------------------------------------------*/ 7084a2ae208SSatish Balay #undef __FUNCT__ 7094a2ae208SSatish Balay #define __FUNCT__ "MatSetValues_SeqDense" 71013f74950SBarry Smith PetscErrorCode MatSetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],const PetscScalar v[],InsertMode addv) 711289bc588SBarry Smith { 712c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 71313f74950SBarry Smith PetscInt i,j,idx=0; 714d6dfbf8fSBarry Smith 7153a40ed3dSBarry Smith PetscFunctionBegin; 71671fd2e92SBarry Smith if (v) PetscValidScalarPointer(v,6); 717289bc588SBarry Smith if (!mat->roworiented) { 718dbb450caSBarry Smith if (addv == INSERT_VALUES) { 719289bc588SBarry Smith for (j=0; j<n; j++) { 720cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 7212515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 722e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 72358804f6dSBarry Smith #endif 724289bc588SBarry Smith for (i=0; i<m; i++) { 725cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 7262515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 727e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 72858804f6dSBarry Smith #endif 729cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 730289bc588SBarry Smith } 731289bc588SBarry Smith } 7323a40ed3dSBarry Smith } else { 733289bc588SBarry Smith for (j=0; j<n; j++) { 734cddbea37SSatish Balay if (indexn[j] < 0) {idx += m; continue;} 7352515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 736e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 73758804f6dSBarry Smith #endif 738289bc588SBarry Smith for (i=0; i<m; i++) { 739cddbea37SSatish Balay if (indexm[i] < 0) {idx++; continue;} 7402515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 741e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 74258804f6dSBarry Smith #endif 743cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 744289bc588SBarry Smith } 745289bc588SBarry Smith } 746289bc588SBarry Smith } 7473a40ed3dSBarry Smith } else { 748dbb450caSBarry Smith if (addv == INSERT_VALUES) { 749e8d4e0b9SBarry Smith for (i=0; i<m; i++) { 750cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7512515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 752e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 75358804f6dSBarry Smith #endif 754e8d4e0b9SBarry Smith for (j=0; j<n; j++) { 755cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7562515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 757e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 75858804f6dSBarry Smith #endif 759cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] = v[idx++]; 760e8d4e0b9SBarry Smith } 761e8d4e0b9SBarry Smith } 7623a40ed3dSBarry Smith } else { 763289bc588SBarry Smith for (i=0; i<m; i++) { 764cddbea37SSatish Balay if (indexm[i] < 0) { idx += n; continue;} 7652515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 766e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",indexm[i],A->rmap->n-1); 76758804f6dSBarry Smith #endif 768289bc588SBarry Smith for (j=0; j<n; j++) { 769cddbea37SSatish Balay if (indexn[j] < 0) { idx++; continue;} 7702515c552SBarry Smith #if defined(PETSC_USE_DEBUG) 771e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column too large: col %D max %D",indexn[j],A->cmap->n-1); 77258804f6dSBarry Smith #endif 773cddbea37SSatish Balay mat->v[indexn[j]*mat->lda + indexm[i]] += v[idx++]; 774289bc588SBarry Smith } 775289bc588SBarry Smith } 776289bc588SBarry Smith } 777e8d4e0b9SBarry Smith } 7783a40ed3dSBarry Smith PetscFunctionReturn(0); 779289bc588SBarry Smith } 780e8d4e0b9SBarry Smith 7814a2ae208SSatish Balay #undef __FUNCT__ 7824a2ae208SSatish Balay #define __FUNCT__ "MatGetValues_SeqDense" 78313f74950SBarry Smith PetscErrorCode MatGetValues_SeqDense(Mat A,PetscInt m,const PetscInt indexm[],PetscInt n,const PetscInt indexn[],PetscScalar v[]) 784ae80bb75SLois Curfman McInnes { 785ae80bb75SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 78613f74950SBarry Smith PetscInt i,j; 787ae80bb75SLois Curfman McInnes 7883a40ed3dSBarry Smith PetscFunctionBegin; 789ae80bb75SLois Curfman McInnes /* row-oriented output */ 790ae80bb75SLois Curfman McInnes for (i=0; i<m; i++) { 79197e567efSBarry Smith if (indexm[i] < 0) {v += n;continue;} 792e32f2f54SBarry Smith if (indexm[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested larger than number rows %D",indexm[i],A->rmap->n); 793ae80bb75SLois Curfman McInnes for (j=0; j<n; j++) { 7946f31f424SBarry Smith if (indexn[j] < 0) {v++; continue;} 795e32f2f54SBarry Smith if (indexn[j] >= A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column %D requested larger than number columns %D",indexn[j],A->cmap->n); 79697e567efSBarry Smith *v++ = mat->v[indexn[j]*mat->lda + indexm[i]]; 797ae80bb75SLois Curfman McInnes } 798ae80bb75SLois Curfman McInnes } 7993a40ed3dSBarry Smith PetscFunctionReturn(0); 800ae80bb75SLois Curfman McInnes } 801ae80bb75SLois Curfman McInnes 802289bc588SBarry Smith /* -----------------------------------------------------------------*/ 803289bc588SBarry Smith 8044a2ae208SSatish Balay #undef __FUNCT__ 8055bba2384SShri Abhyankar #define __FUNCT__ "MatLoad_SeqDense" 806112444f4SShri Abhyankar PetscErrorCode MatLoad_SeqDense(Mat newmat,PetscViewer viewer) 807aabbc4fbSShri Abhyankar { 808aabbc4fbSShri Abhyankar Mat_SeqDense *a; 809aabbc4fbSShri Abhyankar PetscErrorCode ierr; 810aabbc4fbSShri Abhyankar PetscInt *scols,i,j,nz,header[4]; 811aabbc4fbSShri Abhyankar int fd; 812aabbc4fbSShri Abhyankar PetscMPIInt size; 813aabbc4fbSShri Abhyankar PetscInt *rowlengths = 0,M,N,*cols,grows,gcols; 814aabbc4fbSShri Abhyankar PetscScalar *vals,*svals,*v,*w; 815ce94432eSBarry Smith MPI_Comm comm; 816aabbc4fbSShri Abhyankar 817aabbc4fbSShri Abhyankar PetscFunctionBegin; 818ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)viewer,&comm);CHKERRQ(ierr); 819aabbc4fbSShri Abhyankar ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 820aabbc4fbSShri Abhyankar if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"view must have one processor"); 821aabbc4fbSShri Abhyankar ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 822aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,header,4,PETSC_INT);CHKERRQ(ierr); 823aabbc4fbSShri Abhyankar if (header[0] != MAT_FILE_CLASSID) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED,"Not matrix object"); 824aabbc4fbSShri Abhyankar M = header[1]; N = header[2]; nz = header[3]; 825aabbc4fbSShri Abhyankar 826aabbc4fbSShri Abhyankar /* set global size if not set already*/ 827aabbc4fbSShri Abhyankar if (newmat->rmap->n < 0 && newmat->rmap->N < 0 && newmat->cmap->n < 0 && newmat->cmap->N < 0) { 828aabbc4fbSShri Abhyankar ierr = MatSetSizes(newmat,M,N,M,N);CHKERRQ(ierr); 829aabbc4fbSShri Abhyankar } else { 830aabbc4fbSShri Abhyankar /* if sizes and type are already set, check if the vector global sizes are correct */ 831aabbc4fbSShri Abhyankar ierr = MatGetSize(newmat,&grows,&gcols);CHKERRQ(ierr); 832aabbc4fbSShri Abhyankar if (M != grows || N != gcols) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Matrix in file of different length (%d, %d) than the input matrix (%d, %d)",M,N,grows,gcols); 833aabbc4fbSShri Abhyankar } 8340298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 835aabbc4fbSShri Abhyankar 836aabbc4fbSShri Abhyankar if (nz == MATRIX_BINARY_FORMAT_DENSE) { /* matrix in file is dense */ 837aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 838aabbc4fbSShri Abhyankar v = a->v; 839aabbc4fbSShri Abhyankar /* Allocate some temp space to read in the values and then flip them 840aabbc4fbSShri Abhyankar from row major to column major */ 841aabbc4fbSShri Abhyankar ierr = PetscMalloc((M*N > 0 ? M*N : 1)*sizeof(PetscScalar),&w);CHKERRQ(ierr); 842aabbc4fbSShri Abhyankar /* read in nonzero values */ 843aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,w,M*N,PETSC_SCALAR);CHKERRQ(ierr); 844aabbc4fbSShri Abhyankar /* now flip the values and store them in the matrix*/ 845aabbc4fbSShri Abhyankar for (j=0; j<N; j++) { 846aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 847aabbc4fbSShri Abhyankar *v++ =w[i*N+j]; 848aabbc4fbSShri Abhyankar } 849aabbc4fbSShri Abhyankar } 850aabbc4fbSShri Abhyankar ierr = PetscFree(w);CHKERRQ(ierr); 851aabbc4fbSShri Abhyankar } else { 852aabbc4fbSShri Abhyankar /* read row lengths */ 853aabbc4fbSShri Abhyankar ierr = PetscMalloc((M+1)*sizeof(PetscInt),&rowlengths);CHKERRQ(ierr); 854aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,rowlengths,M,PETSC_INT);CHKERRQ(ierr); 855aabbc4fbSShri Abhyankar 856aabbc4fbSShri Abhyankar a = (Mat_SeqDense*)newmat->data; 857aabbc4fbSShri Abhyankar v = a->v; 858aabbc4fbSShri Abhyankar 859aabbc4fbSShri Abhyankar /* read column indices and nonzeros */ 860aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&scols);CHKERRQ(ierr); 861aabbc4fbSShri Abhyankar cols = scols; 862aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,cols,nz,PETSC_INT);CHKERRQ(ierr); 863aabbc4fbSShri Abhyankar ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&svals);CHKERRQ(ierr); 864aabbc4fbSShri Abhyankar vals = svals; 865aabbc4fbSShri Abhyankar ierr = PetscBinaryRead(fd,vals,nz,PETSC_SCALAR);CHKERRQ(ierr); 866aabbc4fbSShri Abhyankar 867aabbc4fbSShri Abhyankar /* insert into matrix */ 868aabbc4fbSShri Abhyankar for (i=0; i<M; i++) { 869aabbc4fbSShri Abhyankar for (j=0; j<rowlengths[i]; j++) v[i+M*scols[j]] = svals[j]; 870aabbc4fbSShri Abhyankar svals += rowlengths[i]; scols += rowlengths[i]; 871aabbc4fbSShri Abhyankar } 872aabbc4fbSShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 873aabbc4fbSShri Abhyankar ierr = PetscFree(cols);CHKERRQ(ierr); 874aabbc4fbSShri Abhyankar ierr = PetscFree(rowlengths);CHKERRQ(ierr); 875aabbc4fbSShri Abhyankar } 876aabbc4fbSShri Abhyankar ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 877aabbc4fbSShri Abhyankar ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 878aabbc4fbSShri Abhyankar PetscFunctionReturn(0); 879aabbc4fbSShri Abhyankar } 880aabbc4fbSShri Abhyankar 881aabbc4fbSShri Abhyankar #undef __FUNCT__ 8824a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_ASCII" 8836849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_ASCII(Mat A,PetscViewer viewer) 884289bc588SBarry Smith { 885932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 886dfbe8321SBarry Smith PetscErrorCode ierr; 88713f74950SBarry Smith PetscInt i,j; 8882dcb1b2aSMatthew Knepley const char *name; 88987828ca2SBarry Smith PetscScalar *v; 890f3ef73ceSBarry Smith PetscViewerFormat format; 8915f481a85SSatish Balay #if defined(PETSC_USE_COMPLEX) 892ace3abfcSBarry Smith PetscBool allreal = PETSC_TRUE; 8935f481a85SSatish Balay #endif 894932b0c3eSLois Curfman McInnes 8953a40ed3dSBarry Smith PetscFunctionBegin; 896b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 897456192e2SBarry Smith if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 8983a40ed3dSBarry Smith PetscFunctionReturn(0); /* do nothing for now */ 899fb9695e5SSatish Balay } else if (format == PETSC_VIEWER_ASCII_COMMON) { 900d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 9017566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 902d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 90344cd7ae7SLois Curfman McInnes v = a->v + i; 90477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"row %D:",i);CHKERRQ(ierr); 905d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 906aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 907329f5518SBarry Smith if (PetscRealPart(*v) != 0.0 && PetscImaginaryPart(*v) != 0.0) { 908a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G + %G i) ",j,PetscRealPart(*v),PetscImaginaryPart(*v));CHKERRQ(ierr); 909329f5518SBarry Smith } else if (PetscRealPart(*v)) { 910a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,PetscRealPart(*v));CHKERRQ(ierr); 9116831982aSBarry Smith } 91280cd9d93SLois Curfman McInnes #else 9136831982aSBarry Smith if (*v) { 914a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," (%D, %G) ",j,*v);CHKERRQ(ierr); 9156831982aSBarry Smith } 91680cd9d93SLois Curfman McInnes #endif 9171b807ce4Svictorle v += a->lda; 91880cd9d93SLois Curfman McInnes } 919b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 92080cd9d93SLois Curfman McInnes } 921d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 9223a40ed3dSBarry Smith } else { 923d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 924aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 92547989497SBarry Smith /* determine if matrix has all real values */ 92647989497SBarry Smith v = a->v; 927d0f46423SBarry Smith for (i=0; i<A->rmap->n*A->cmap->n; i++) { 928ffac6cdbSBarry Smith if (PetscImaginaryPart(v[i])) { allreal = PETSC_FALSE; break;} 92947989497SBarry Smith } 93047989497SBarry Smith #endif 931fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 9323a7fca6bSBarry Smith ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 933d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%% Size = %D %D \n",A->rmap->n,A->cmap->n);CHKERRQ(ierr); 934d0f46423SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%s = zeros(%D,%D);\n",name,A->rmap->n,A->cmap->n);CHKERRQ(ierr); 935fb9695e5SSatish Balay ierr = PetscViewerASCIIPrintf(viewer,"%s = [\n",name);CHKERRQ(ierr); 9367566de4bSShri Abhyankar } else { 9377566de4bSShri Abhyankar ierr = PetscObjectPrintClassNamePrefixType((PetscObject)A,viewer,"Matrix Object");CHKERRQ(ierr); 938ffac6cdbSBarry Smith } 939ffac6cdbSBarry Smith 940d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 941932b0c3eSLois Curfman McInnes v = a->v + i; 942d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 943aa482453SBarry Smith #if defined(PETSC_USE_COMPLEX) 94447989497SBarry Smith if (allreal) { 945c61cd2faSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)PetscRealPart(*v));CHKERRQ(ierr); 94647989497SBarry Smith } else { 947c61cd2faSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"%18.16e + %18.16ei ",(double)PetscRealPart(*v),(double)PetscImaginaryPart(*v));CHKERRQ(ierr); 94847989497SBarry Smith } 949289bc588SBarry Smith #else 950c61cd2faSJed Brown ierr = PetscViewerASCIIPrintf(viewer,"%18.16e ",(double)*v);CHKERRQ(ierr); 951289bc588SBarry Smith #endif 9521b807ce4Svictorle v += a->lda; 953289bc588SBarry Smith } 954b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 955289bc588SBarry Smith } 956fb9695e5SSatish Balay if (format == PETSC_VIEWER_ASCII_MATLAB) { 957b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"];\n");CHKERRQ(ierr); 958ffac6cdbSBarry Smith } 959d00279f6SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 960da3a660dSBarry Smith } 961b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 9623a40ed3dSBarry Smith PetscFunctionReturn(0); 963289bc588SBarry Smith } 964289bc588SBarry Smith 9654a2ae208SSatish Balay #undef __FUNCT__ 9664a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Binary" 9676849ba73SBarry Smith static PetscErrorCode MatView_SeqDense_Binary(Mat A,PetscViewer viewer) 968932b0c3eSLois Curfman McInnes { 969932b0c3eSLois Curfman McInnes Mat_SeqDense *a = (Mat_SeqDense*)A->data; 9706849ba73SBarry Smith PetscErrorCode ierr; 97113f74950SBarry Smith int fd; 972d0f46423SBarry Smith PetscInt ict,j,n = A->cmap->n,m = A->rmap->n,i,*col_lens,nz = m*n; 973f4403165SShri Abhyankar PetscScalar *v,*anonz,*vals; 974f4403165SShri Abhyankar PetscViewerFormat format; 975932b0c3eSLois Curfman McInnes 9763a40ed3dSBarry Smith PetscFunctionBegin; 977b0a32e0cSBarry Smith ierr = PetscViewerBinaryGetDescriptor(viewer,&fd);CHKERRQ(ierr); 97890ace30eSBarry Smith 979f4403165SShri Abhyankar ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 980f4403165SShri Abhyankar if (format == PETSC_VIEWER_NATIVE) { 981f4403165SShri Abhyankar /* store the matrix as a dense matrix */ 982f4403165SShri Abhyankar ierr = PetscMalloc(4*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 9832205254eSKarl Rupp 984f4403165SShri Abhyankar col_lens[0] = MAT_FILE_CLASSID; 985f4403165SShri Abhyankar col_lens[1] = m; 986f4403165SShri Abhyankar col_lens[2] = n; 987f4403165SShri Abhyankar col_lens[3] = MATRIX_BINARY_FORMAT_DENSE; 9882205254eSKarl Rupp 989f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,col_lens,4,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 990f4403165SShri Abhyankar ierr = PetscFree(col_lens);CHKERRQ(ierr); 991f4403165SShri Abhyankar 992f4403165SShri Abhyankar /* write out matrix, by rows */ 993f4403165SShri Abhyankar ierr = PetscMalloc((m*n+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr); 994f4403165SShri Abhyankar v = a->v; 995f4403165SShri Abhyankar for (j=0; j<n; j++) { 996f4403165SShri Abhyankar for (i=0; i<m; i++) { 997f4403165SShri Abhyankar vals[j + i*n] = *v++; 998f4403165SShri Abhyankar } 999f4403165SShri Abhyankar } 1000f4403165SShri Abhyankar ierr = PetscBinaryWrite(fd,vals,n*m,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 1001f4403165SShri Abhyankar ierr = PetscFree(vals);CHKERRQ(ierr); 1002f4403165SShri Abhyankar } else { 100313f74950SBarry Smith ierr = PetscMalloc((4+nz)*sizeof(PetscInt),&col_lens);CHKERRQ(ierr); 10042205254eSKarl Rupp 10050700a824SBarry Smith col_lens[0] = MAT_FILE_CLASSID; 1006932b0c3eSLois Curfman McInnes col_lens[1] = m; 1007932b0c3eSLois Curfman McInnes col_lens[2] = n; 1008932b0c3eSLois Curfman McInnes col_lens[3] = nz; 1009932b0c3eSLois Curfman McInnes 1010932b0c3eSLois Curfman McInnes /* store lengths of each row and write (including header) to file */ 1011932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) col_lens[4+i] = n; 10126f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,4+m,PETSC_INT,PETSC_TRUE);CHKERRQ(ierr); 1013932b0c3eSLois Curfman McInnes 1014932b0c3eSLois Curfman McInnes /* Possibly should write in smaller increments, not whole matrix at once? */ 1015932b0c3eSLois Curfman McInnes /* store column indices (zero start index) */ 1016932b0c3eSLois Curfman McInnes ict = 0; 1017932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 1018932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) col_lens[ict++] = j; 1019932b0c3eSLois Curfman McInnes } 10206f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,col_lens,nz,PETSC_INT,PETSC_FALSE);CHKERRQ(ierr); 1021606d414cSSatish Balay ierr = PetscFree(col_lens);CHKERRQ(ierr); 1022932b0c3eSLois Curfman McInnes 1023932b0c3eSLois Curfman McInnes /* store nonzero values */ 102487828ca2SBarry Smith ierr = PetscMalloc((nz+1)*sizeof(PetscScalar),&anonz);CHKERRQ(ierr); 1025932b0c3eSLois Curfman McInnes ict = 0; 1026932b0c3eSLois Curfman McInnes for (i=0; i<m; i++) { 1027932b0c3eSLois Curfman McInnes v = a->v + i; 1028932b0c3eSLois Curfman McInnes for (j=0; j<n; j++) { 10291b807ce4Svictorle anonz[ict++] = *v; v += a->lda; 1030932b0c3eSLois Curfman McInnes } 1031932b0c3eSLois Curfman McInnes } 10326f69ff64SBarry Smith ierr = PetscBinaryWrite(fd,anonz,nz,PETSC_SCALAR,PETSC_FALSE);CHKERRQ(ierr); 1033606d414cSSatish Balay ierr = PetscFree(anonz);CHKERRQ(ierr); 1034f4403165SShri Abhyankar } 10353a40ed3dSBarry Smith PetscFunctionReturn(0); 1036932b0c3eSLois Curfman McInnes } 1037932b0c3eSLois Curfman McInnes 10389804daf3SBarry Smith #include <petscdraw.h> 10394a2ae208SSatish Balay #undef __FUNCT__ 10404a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw_Zoom" 1041dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw_Zoom(PetscDraw draw,void *Aa) 1042f1af5d2fSBarry Smith { 1043f1af5d2fSBarry Smith Mat A = (Mat) Aa; 1044f1af5d2fSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 10456849ba73SBarry Smith PetscErrorCode ierr; 1046d0f46423SBarry Smith PetscInt m = A->rmap->n,n = A->cmap->n,color,i,j; 104787828ca2SBarry Smith PetscScalar *v = a->v; 1048b0a32e0cSBarry Smith PetscViewer viewer; 1049b0a32e0cSBarry Smith PetscDraw popup; 1050329f5518SBarry Smith PetscReal xl,yl,xr,yr,x_l,x_r,y_l,y_r,scale,maxv = 0.0; 1051f3ef73ceSBarry Smith PetscViewerFormat format; 1052f1af5d2fSBarry Smith 1053f1af5d2fSBarry Smith PetscFunctionBegin; 1054f1af5d2fSBarry Smith ierr = PetscObjectQuery((PetscObject)A,"Zoomviewer",(PetscObject*)&viewer);CHKERRQ(ierr); 1055b0a32e0cSBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1056b0a32e0cSBarry Smith ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1057f1af5d2fSBarry Smith 1058f1af5d2fSBarry Smith /* Loop over matrix elements drawing boxes */ 1059fb9695e5SSatish Balay if (format != PETSC_VIEWER_DRAW_CONTOUR) { 1060f1af5d2fSBarry Smith /* Blue for negative and Red for positive */ 1061b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1062f1af5d2fSBarry Smith for (j = 0; j < n; j++) { 1063f1af5d2fSBarry Smith x_l = j; 1064f1af5d2fSBarry Smith x_r = x_l + 1.0; 1065f1af5d2fSBarry Smith for (i = 0; i < m; i++) { 1066f1af5d2fSBarry Smith y_l = m - i - 1.0; 1067f1af5d2fSBarry Smith y_r = y_l + 1.0; 1068329f5518SBarry Smith if (PetscRealPart(v[j*m+i]) > 0.) { 1069b0a32e0cSBarry Smith color = PETSC_DRAW_RED; 1070329f5518SBarry Smith } else if (PetscRealPart(v[j*m+i]) < 0.) { 1071b0a32e0cSBarry Smith color = PETSC_DRAW_BLUE; 1072f1af5d2fSBarry Smith } else { 1073f1af5d2fSBarry Smith continue; 1074f1af5d2fSBarry Smith } 1075b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1076f1af5d2fSBarry Smith } 1077f1af5d2fSBarry Smith } 1078f1af5d2fSBarry Smith } else { 1079f1af5d2fSBarry Smith /* use contour shading to indicate magnitude of values */ 1080f1af5d2fSBarry Smith /* first determine max of all nonzero values */ 1081f1af5d2fSBarry Smith for (i = 0; i < m*n; i++) { 1082f1af5d2fSBarry Smith if (PetscAbsScalar(v[i]) > maxv) maxv = PetscAbsScalar(v[i]); 1083f1af5d2fSBarry Smith } 1084b0a32e0cSBarry Smith scale = (245.0 - PETSC_DRAW_BASIC_COLORS)/maxv; 1085b0a32e0cSBarry Smith ierr = PetscDrawGetPopup(draw,&popup);CHKERRQ(ierr); 1086b0a32e0cSBarry Smith if (popup) {ierr = PetscDrawScalePopup(popup,0.0,maxv);CHKERRQ(ierr);} 1087f1af5d2fSBarry Smith for (j = 0; j < n; j++) { 1088f1af5d2fSBarry Smith x_l = j; 1089f1af5d2fSBarry Smith x_r = x_l + 1.0; 1090f1af5d2fSBarry Smith for (i = 0; i < m; i++) { 1091f1af5d2fSBarry Smith y_l = m - i - 1.0; 1092f1af5d2fSBarry Smith y_r = y_l + 1.0; 1093b0a32e0cSBarry Smith color = PETSC_DRAW_BASIC_COLORS + (int)(scale*PetscAbsScalar(v[j*m+i])); 1094b0a32e0cSBarry Smith ierr = PetscDrawRectangle(draw,x_l,y_l,x_r,y_r,color,color,color,color);CHKERRQ(ierr); 1095f1af5d2fSBarry Smith } 1096f1af5d2fSBarry Smith } 1097f1af5d2fSBarry Smith } 1098f1af5d2fSBarry Smith PetscFunctionReturn(0); 1099f1af5d2fSBarry Smith } 1100f1af5d2fSBarry Smith 11014a2ae208SSatish Balay #undef __FUNCT__ 11024a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense_Draw" 1103dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense_Draw(Mat A,PetscViewer viewer) 1104f1af5d2fSBarry Smith { 1105b0a32e0cSBarry Smith PetscDraw draw; 1106ace3abfcSBarry Smith PetscBool isnull; 1107329f5518SBarry Smith PetscReal xr,yr,xl,yl,h,w; 1108dfbe8321SBarry Smith PetscErrorCode ierr; 1109f1af5d2fSBarry Smith 1110f1af5d2fSBarry Smith PetscFunctionBegin; 1111b0a32e0cSBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 1112b0a32e0cSBarry Smith ierr = PetscDrawIsNull(draw,&isnull);CHKERRQ(ierr); 1113abc0a331SBarry Smith if (isnull) PetscFunctionReturn(0); 1114f1af5d2fSBarry Smith 1115f1af5d2fSBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",(PetscObject)viewer);CHKERRQ(ierr); 1116d0f46423SBarry Smith xr = A->cmap->n; yr = A->rmap->n; h = yr/10.0; w = xr/10.0; 1117f1af5d2fSBarry Smith xr += w; yr += h; xl = -w; yl = -h; 1118b0a32e0cSBarry Smith ierr = PetscDrawSetCoordinates(draw,xl,yl,xr,yr);CHKERRQ(ierr); 1119b0a32e0cSBarry Smith ierr = PetscDrawZoom(draw,MatView_SeqDense_Draw_Zoom,A);CHKERRQ(ierr); 11200298fd71SBarry Smith ierr = PetscObjectCompose((PetscObject)A,"Zoomviewer",NULL);CHKERRQ(ierr); 1121f1af5d2fSBarry Smith PetscFunctionReturn(0); 1122f1af5d2fSBarry Smith } 1123f1af5d2fSBarry Smith 11244a2ae208SSatish Balay #undef __FUNCT__ 11254a2ae208SSatish Balay #define __FUNCT__ "MatView_SeqDense" 1126dfbe8321SBarry Smith PetscErrorCode MatView_SeqDense(Mat A,PetscViewer viewer) 1127932b0c3eSLois Curfman McInnes { 1128dfbe8321SBarry Smith PetscErrorCode ierr; 1129ace3abfcSBarry Smith PetscBool iascii,isbinary,isdraw; 1130932b0c3eSLois Curfman McInnes 11313a40ed3dSBarry Smith PetscFunctionBegin; 1132251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1133251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 1134251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 11350f5bd95cSBarry Smith 1136c45a1595SBarry Smith if (iascii) { 1137c45a1595SBarry Smith ierr = MatView_SeqDense_ASCII(A,viewer);CHKERRQ(ierr); 11380f5bd95cSBarry Smith } else if (isbinary) { 11393a40ed3dSBarry Smith ierr = MatView_SeqDense_Binary(A,viewer);CHKERRQ(ierr); 1140f1af5d2fSBarry Smith } else if (isdraw) { 1141f1af5d2fSBarry Smith ierr = MatView_SeqDense_Draw(A,viewer);CHKERRQ(ierr); 1142932b0c3eSLois Curfman McInnes } 11433a40ed3dSBarry Smith PetscFunctionReturn(0); 1144932b0c3eSLois Curfman McInnes } 1145289bc588SBarry Smith 11464a2ae208SSatish Balay #undef __FUNCT__ 11474a2ae208SSatish Balay #define __FUNCT__ "MatDestroy_SeqDense" 1148dfbe8321SBarry Smith PetscErrorCode MatDestroy_SeqDense(Mat mat) 1149289bc588SBarry Smith { 1150ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)mat->data; 1151dfbe8321SBarry Smith PetscErrorCode ierr; 115290f02eecSBarry Smith 11533a40ed3dSBarry Smith PetscFunctionBegin; 1154aa482453SBarry Smith #if defined(PETSC_USE_LOG) 1155d0f46423SBarry Smith PetscLogObjectState((PetscObject)mat,"Rows %D Cols %D",mat->rmap->n,mat->cmap->n); 1156a5a9c739SBarry Smith #endif 115705b42c5fSBarry Smith ierr = PetscFree(l->pivots);CHKERRQ(ierr); 11586857c123SSatish Balay if (!l->user_alloc) {ierr = PetscFree(l->v);CHKERRQ(ierr);} 1159bf0cc555SLisandro Dalcin ierr = PetscFree(mat->data);CHKERRQ(ierr); 1160dbd8c25aSHong Zhang 1161dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 116200de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseGetArray_C","",NULL);CHKERRQ(ierr); 116300de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatDenseRestoreArray_C","",NULL);CHKERRQ(ierr); 116400de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSeqDenseSetPreallocation_C","",NULL);CHKERRQ(ierr); 116500de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMult_seqaij_seqdense_C","",NULL);CHKERRQ(ierr); 116600de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultSymbolic_seqaij_seqdense_C","",NULL);CHKERRQ(ierr); 116700de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMatMultNumeric_seqaij_seqdense_C","",NULL);CHKERRQ(ierr); 11683a40ed3dSBarry Smith PetscFunctionReturn(0); 1169289bc588SBarry Smith } 1170289bc588SBarry Smith 11714a2ae208SSatish Balay #undef __FUNCT__ 11724a2ae208SSatish Balay #define __FUNCT__ "MatTranspose_SeqDense" 1173fc4dec0aSBarry Smith PetscErrorCode MatTranspose_SeqDense(Mat A,MatReuse reuse,Mat *matout) 1174289bc588SBarry Smith { 1175c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 11766849ba73SBarry Smith PetscErrorCode ierr; 117713f74950SBarry Smith PetscInt k,j,m,n,M; 117887828ca2SBarry Smith PetscScalar *v,tmp; 117948b35521SBarry Smith 11803a40ed3dSBarry Smith PetscFunctionBegin; 1181d0f46423SBarry Smith v = mat->v; m = A->rmap->n; M = mat->lda; n = A->cmap->n; 1182e9695a30SBarry Smith if (reuse == MAT_REUSE_MATRIX && *matout == A) { /* in place transpose */ 1183e7e72b3dSBarry Smith if (m != n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Can not transpose non-square matrix in place"); 1184e7e72b3dSBarry Smith else { 1185d3e5ee88SLois Curfman McInnes for (j=0; j<m; j++) { 1186289bc588SBarry Smith for (k=0; k<j; k++) { 11871b807ce4Svictorle tmp = v[j + k*M]; 11881b807ce4Svictorle v[j + k*M] = v[k + j*M]; 11891b807ce4Svictorle v[k + j*M] = tmp; 1190289bc588SBarry Smith } 1191289bc588SBarry Smith } 1192d64ed03dSBarry Smith } 11933a40ed3dSBarry Smith } else { /* out-of-place transpose */ 1194d3e5ee88SLois Curfman McInnes Mat tmat; 1195ec8511deSBarry Smith Mat_SeqDense *tmatd; 119687828ca2SBarry Smith PetscScalar *v2; 1197ea709b57SSatish Balay 1198fc4dec0aSBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 1199ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&tmat);CHKERRQ(ierr); 1200d0f46423SBarry Smith ierr = MatSetSizes(tmat,A->cmap->n,A->rmap->n,A->cmap->n,A->rmap->n);CHKERRQ(ierr); 12017adad957SLisandro Dalcin ierr = MatSetType(tmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 12020298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(tmat,NULL);CHKERRQ(ierr); 1203fc4dec0aSBarry Smith } else { 1204fc4dec0aSBarry Smith tmat = *matout; 1205fc4dec0aSBarry Smith } 1206ec8511deSBarry Smith tmatd = (Mat_SeqDense*)tmat->data; 12070de55854SLois Curfman McInnes v = mat->v; v2 = tmatd->v; 1208d3e5ee88SLois Curfman McInnes for (j=0; j<n; j++) { 12091b807ce4Svictorle for (k=0; k<m; k++) v2[j + k*n] = v[k + j*M]; 1210d3e5ee88SLois Curfman McInnes } 12116d4a8577SBarry Smith ierr = MatAssemblyBegin(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 12126d4a8577SBarry Smith ierr = MatAssemblyEnd(tmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1213d3e5ee88SLois Curfman McInnes *matout = tmat; 121448b35521SBarry Smith } 12153a40ed3dSBarry Smith PetscFunctionReturn(0); 1216289bc588SBarry Smith } 1217289bc588SBarry Smith 12184a2ae208SSatish Balay #undef __FUNCT__ 12194a2ae208SSatish Balay #define __FUNCT__ "MatEqual_SeqDense" 1220ace3abfcSBarry Smith PetscErrorCode MatEqual_SeqDense(Mat A1,Mat A2,PetscBool *flg) 1221289bc588SBarry Smith { 1222c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat1 = (Mat_SeqDense*)A1->data; 1223c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat2 = (Mat_SeqDense*)A2->data; 122413f74950SBarry Smith PetscInt i,j; 1225a2ea699eSBarry Smith PetscScalar *v1,*v2; 12269ea5d5aeSSatish Balay 12273a40ed3dSBarry Smith PetscFunctionBegin; 1228d0f46423SBarry Smith if (A1->rmap->n != A2->rmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1229d0f46423SBarry Smith if (A1->cmap->n != A2->cmap->n) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 1230d0f46423SBarry Smith for (i=0; i<A1->rmap->n; i++) { 12311b807ce4Svictorle v1 = mat1->v+i; v2 = mat2->v+i; 1232d0f46423SBarry Smith for (j=0; j<A1->cmap->n; j++) { 12333a40ed3dSBarry Smith if (*v1 != *v2) {*flg = PETSC_FALSE; PetscFunctionReturn(0);} 12341b807ce4Svictorle v1 += mat1->lda; v2 += mat2->lda; 12351b807ce4Svictorle } 1236289bc588SBarry Smith } 123777c4ece6SBarry Smith *flg = PETSC_TRUE; 12383a40ed3dSBarry Smith PetscFunctionReturn(0); 1239289bc588SBarry Smith } 1240289bc588SBarry Smith 12414a2ae208SSatish Balay #undef __FUNCT__ 12424a2ae208SSatish Balay #define __FUNCT__ "MatGetDiagonal_SeqDense" 1243dfbe8321SBarry Smith PetscErrorCode MatGetDiagonal_SeqDense(Mat A,Vec v) 1244289bc588SBarry Smith { 1245c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 1246dfbe8321SBarry Smith PetscErrorCode ierr; 124713f74950SBarry Smith PetscInt i,n,len; 124887828ca2SBarry Smith PetscScalar *x,zero = 0.0; 124944cd7ae7SLois Curfman McInnes 12503a40ed3dSBarry Smith PetscFunctionBegin; 12512dcb1b2aSMatthew Knepley ierr = VecSet(v,zero);CHKERRQ(ierr); 12527a97a34bSBarry Smith ierr = VecGetSize(v,&n);CHKERRQ(ierr); 12531ebc52fbSHong Zhang ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1254d0f46423SBarry Smith len = PetscMin(A->rmap->n,A->cmap->n); 1255e32f2f54SBarry Smith if (n != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming mat and vec"); 125644cd7ae7SLois Curfman McInnes for (i=0; i<len; i++) { 12571b807ce4Svictorle x[i] = mat->v[i*mat->lda + i]; 1258289bc588SBarry Smith } 12591ebc52fbSHong Zhang ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 12603a40ed3dSBarry Smith PetscFunctionReturn(0); 1261289bc588SBarry Smith } 1262289bc588SBarry Smith 12634a2ae208SSatish Balay #undef __FUNCT__ 12644a2ae208SSatish Balay #define __FUNCT__ "MatDiagonalScale_SeqDense" 1265dfbe8321SBarry Smith PetscErrorCode MatDiagonalScale_SeqDense(Mat A,Vec ll,Vec rr) 1266289bc588SBarry Smith { 1267c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 126887828ca2SBarry Smith PetscScalar *l,*r,x,*v; 1269dfbe8321SBarry Smith PetscErrorCode ierr; 1270d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n; 127155659b69SBarry Smith 12723a40ed3dSBarry Smith PetscFunctionBegin; 127328988994SBarry Smith if (ll) { 12747a97a34bSBarry Smith ierr = VecGetSize(ll,&m);CHKERRQ(ierr); 12751ebc52fbSHong Zhang ierr = VecGetArray(ll,&l);CHKERRQ(ierr); 1276e32f2f54SBarry Smith if (m != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Left scaling vec wrong size"); 1277da3a660dSBarry Smith for (i=0; i<m; i++) { 1278da3a660dSBarry Smith x = l[i]; 1279da3a660dSBarry Smith v = mat->v + i; 1280da3a660dSBarry Smith for (j=0; j<n; j++) { (*v) *= x; v+= m;} 1281da3a660dSBarry Smith } 12821ebc52fbSHong Zhang ierr = VecRestoreArray(ll,&l);CHKERRQ(ierr); 1283efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1284da3a660dSBarry Smith } 128528988994SBarry Smith if (rr) { 12867a97a34bSBarry Smith ierr = VecGetSize(rr,&n);CHKERRQ(ierr); 12871ebc52fbSHong Zhang ierr = VecGetArray(rr,&r);CHKERRQ(ierr); 1288e32f2f54SBarry Smith if (n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Right scaling vec wrong size"); 1289da3a660dSBarry Smith for (i=0; i<n; i++) { 1290da3a660dSBarry Smith x = r[i]; 1291da3a660dSBarry Smith v = mat->v + i*m; 12922205254eSKarl Rupp for (j=0; j<m; j++) (*v++) *= x; 1293da3a660dSBarry Smith } 12941ebc52fbSHong Zhang ierr = VecRestoreArray(rr,&r);CHKERRQ(ierr); 1295efee365bSSatish Balay ierr = PetscLogFlops(n*m);CHKERRQ(ierr); 1296da3a660dSBarry Smith } 12973a40ed3dSBarry Smith PetscFunctionReturn(0); 1298289bc588SBarry Smith } 1299289bc588SBarry Smith 13004a2ae208SSatish Balay #undef __FUNCT__ 13014a2ae208SSatish Balay #define __FUNCT__ "MatNorm_SeqDense" 1302dfbe8321SBarry Smith PetscErrorCode MatNorm_SeqDense(Mat A,NormType type,PetscReal *nrm) 1303289bc588SBarry Smith { 1304c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 130587828ca2SBarry Smith PetscScalar *v = mat->v; 1306329f5518SBarry Smith PetscReal sum = 0.0; 1307d0f46423SBarry Smith PetscInt lda =mat->lda,m=A->rmap->n,i,j; 1308efee365bSSatish Balay PetscErrorCode ierr; 130955659b69SBarry Smith 13103a40ed3dSBarry Smith PetscFunctionBegin; 1311289bc588SBarry Smith if (type == NORM_FROBENIUS) { 1312a5ce6ee0Svictorle if (lda>m) { 1313d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1314a5ce6ee0Svictorle v = mat->v+j*lda; 1315a5ce6ee0Svictorle for (i=0; i<m; i++) { 1316a5ce6ee0Svictorle sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1317a5ce6ee0Svictorle } 1318a5ce6ee0Svictorle } 1319a5ce6ee0Svictorle } else { 1320d0f46423SBarry Smith for (i=0; i<A->cmap->n*A->rmap->n; i++) { 1321329f5518SBarry Smith sum += PetscRealPart(PetscConj(*v)*(*v)); v++; 1322289bc588SBarry Smith } 1323a5ce6ee0Svictorle } 13248f1a2a5eSBarry Smith *nrm = PetscSqrtReal(sum); 1325dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13263a40ed3dSBarry Smith } else if (type == NORM_1) { 1327064f8208SBarry Smith *nrm = 0.0; 1328d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 13291b807ce4Svictorle v = mat->v + j*mat->lda; 1330289bc588SBarry Smith sum = 0.0; 1331d0f46423SBarry Smith for (i=0; i<A->rmap->n; i++) { 133233a8263dSBarry Smith sum += PetscAbsScalar(*v); v++; 1333289bc588SBarry Smith } 1334064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1335289bc588SBarry Smith } 1336d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 13373a40ed3dSBarry Smith } else if (type == NORM_INFINITY) { 1338064f8208SBarry Smith *nrm = 0.0; 1339d0f46423SBarry Smith for (j=0; j<A->rmap->n; j++) { 1340289bc588SBarry Smith v = mat->v + j; 1341289bc588SBarry Smith sum = 0.0; 1342d0f46423SBarry Smith for (i=0; i<A->cmap->n; i++) { 13431b807ce4Svictorle sum += PetscAbsScalar(*v); v += mat->lda; 1344289bc588SBarry Smith } 1345064f8208SBarry Smith if (sum > *nrm) *nrm = sum; 1346289bc588SBarry Smith } 1347d0f46423SBarry Smith ierr = PetscLogFlops(A->cmap->n*A->rmap->n);CHKERRQ(ierr); 1348e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No two norm"); 13493a40ed3dSBarry Smith PetscFunctionReturn(0); 1350289bc588SBarry Smith } 1351289bc588SBarry Smith 13524a2ae208SSatish Balay #undef __FUNCT__ 13534a2ae208SSatish Balay #define __FUNCT__ "MatSetOption_SeqDense" 1354ace3abfcSBarry Smith PetscErrorCode MatSetOption_SeqDense(Mat A,MatOption op,PetscBool flg) 1355289bc588SBarry Smith { 1356c0bbcb79SLois Curfman McInnes Mat_SeqDense *aij = (Mat_SeqDense*)A->data; 135763ba0a88SBarry Smith PetscErrorCode ierr; 135867e560aaSBarry Smith 13593a40ed3dSBarry Smith PetscFunctionBegin; 1360b5a2b587SKris Buschelman switch (op) { 1361b5a2b587SKris Buschelman case MAT_ROW_ORIENTED: 13624e0d8c25SBarry Smith aij->roworiented = flg; 1363b5a2b587SKris Buschelman break; 1364512a5fc5SBarry Smith case MAT_NEW_NONZERO_LOCATIONS: 1365b5a2b587SKris Buschelman case MAT_NEW_NONZERO_LOCATION_ERR: 13663971808eSMatthew Knepley case MAT_NEW_NONZERO_ALLOCATION_ERR: 13674e0d8c25SBarry Smith case MAT_NEW_DIAGONALS: 136813fa8e87SLisandro Dalcin case MAT_KEEP_NONZERO_PATTERN: 1369b5a2b587SKris Buschelman case MAT_IGNORE_OFF_PROC_ENTRIES: 1370b5a2b587SKris Buschelman case MAT_USE_HASH_TABLE: 13715021d80fSJed Brown case MAT_IGNORE_LOWER_TRIANGULAR: 13725021d80fSJed Brown ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 13735021d80fSJed Brown break; 13745021d80fSJed Brown case MAT_SPD: 137577e54ba9SKris Buschelman case MAT_SYMMETRIC: 137677e54ba9SKris Buschelman case MAT_STRUCTURALLY_SYMMETRIC: 13779a4540c5SBarry Smith case MAT_HERMITIAN: 13789a4540c5SBarry Smith case MAT_SYMMETRY_ETERNAL: 13795021d80fSJed Brown /* These options are handled directly by MatSetOption() */ 138077e54ba9SKris Buschelman break; 1381b5a2b587SKris Buschelman default: 1382e32f2f54SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"unknown option %s",MatOptions[op]); 13833a40ed3dSBarry Smith } 13843a40ed3dSBarry Smith PetscFunctionReturn(0); 1385289bc588SBarry Smith } 1386289bc588SBarry Smith 13874a2ae208SSatish Balay #undef __FUNCT__ 13884a2ae208SSatish Balay #define __FUNCT__ "MatZeroEntries_SeqDense" 1389dfbe8321SBarry Smith PetscErrorCode MatZeroEntries_SeqDense(Mat A) 13906f0a148fSBarry Smith { 1391ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 13926849ba73SBarry Smith PetscErrorCode ierr; 1393d0f46423SBarry Smith PetscInt lda=l->lda,m=A->rmap->n,j; 13943a40ed3dSBarry Smith 13953a40ed3dSBarry Smith PetscFunctionBegin; 1396a5ce6ee0Svictorle if (lda>m) { 1397d0f46423SBarry Smith for (j=0; j<A->cmap->n; j++) { 1398a5ce6ee0Svictorle ierr = PetscMemzero(l->v+j*lda,m*sizeof(PetscScalar));CHKERRQ(ierr); 1399a5ce6ee0Svictorle } 1400a5ce6ee0Svictorle } else { 1401d0f46423SBarry Smith ierr = PetscMemzero(l->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1402a5ce6ee0Svictorle } 14033a40ed3dSBarry Smith PetscFunctionReturn(0); 14046f0a148fSBarry Smith } 14056f0a148fSBarry Smith 14064a2ae208SSatish Balay #undef __FUNCT__ 14074a2ae208SSatish Balay #define __FUNCT__ "MatZeroRows_SeqDense" 14082b40b63fSBarry Smith PetscErrorCode MatZeroRows_SeqDense(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 14096f0a148fSBarry Smith { 141097b48c8fSBarry Smith PetscErrorCode ierr; 1411ec8511deSBarry Smith Mat_SeqDense *l = (Mat_SeqDense*)A->data; 1412b9679d65SBarry Smith PetscInt m = l->lda, n = A->cmap->n, i,j; 141397b48c8fSBarry Smith PetscScalar *slot,*bb; 141497b48c8fSBarry Smith const PetscScalar *xx; 141555659b69SBarry Smith 14163a40ed3dSBarry Smith PetscFunctionBegin; 1417b9679d65SBarry Smith #if defined(PETSC_USE_DEBUG) 1418b9679d65SBarry Smith for (i=0; i<N; i++) { 1419b9679d65SBarry Smith if (rows[i] < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative row requested to be zeroed"); 1420b9679d65SBarry Smith if (rows[i] >= A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row %D requested to be zeroed greater than or equal number of rows %D",rows[i],A->rmap->n); 1421b9679d65SBarry Smith } 1422b9679d65SBarry Smith #endif 1423b9679d65SBarry Smith 142497b48c8fSBarry Smith /* fix right hand side if needed */ 142597b48c8fSBarry Smith if (x && b) { 142697b48c8fSBarry Smith ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr); 142797b48c8fSBarry Smith ierr = VecGetArray(b,&bb);CHKERRQ(ierr); 14282205254eSKarl Rupp for (i=0; i<N; i++) bb[rows[i]] = diag*xx[rows[i]]; 142997b48c8fSBarry Smith ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr); 143097b48c8fSBarry Smith ierr = VecRestoreArray(b,&bb);CHKERRQ(ierr); 143197b48c8fSBarry Smith } 143297b48c8fSBarry Smith 14336f0a148fSBarry Smith for (i=0; i<N; i++) { 14346f0a148fSBarry Smith slot = l->v + rows[i]; 1435b9679d65SBarry Smith for (j=0; j<n; j++) { *slot = 0.0; slot += m;} 14366f0a148fSBarry Smith } 1437f4df32b1SMatthew Knepley if (diag != 0.0) { 1438b9679d65SBarry Smith if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only coded for square matrices"); 14396f0a148fSBarry Smith for (i=0; i<N; i++) { 1440b9679d65SBarry Smith slot = l->v + (m+1)*rows[i]; 1441f4df32b1SMatthew Knepley *slot = diag; 14426f0a148fSBarry Smith } 14436f0a148fSBarry Smith } 14443a40ed3dSBarry Smith PetscFunctionReturn(0); 14456f0a148fSBarry Smith } 1446557bce09SLois Curfman McInnes 14474a2ae208SSatish Balay #undef __FUNCT__ 14488c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray_SeqDense" 14498c778c55SBarry Smith PetscErrorCode MatDenseGetArray_SeqDense(Mat A,PetscScalar *array[]) 145064e87e97SBarry Smith { 1451c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 14523a40ed3dSBarry Smith 14533a40ed3dSBarry Smith PetscFunctionBegin; 1454e32f2f54SBarry Smith if (mat->lda != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot get array for Dense matrices with LDA different from number of rows"); 145564e87e97SBarry Smith *array = mat->v; 14563a40ed3dSBarry Smith PetscFunctionReturn(0); 145764e87e97SBarry Smith } 14580754003eSLois Curfman McInnes 14594a2ae208SSatish Balay #undef __FUNCT__ 14608c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray_SeqDense" 14618c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray_SeqDense(Mat A,PetscScalar *array[]) 1462ff14e315SSatish Balay { 14633a40ed3dSBarry Smith PetscFunctionBegin; 146409b544d4SBarry Smith *array = 0; /* user cannot accidently use the array later */ 14653a40ed3dSBarry Smith PetscFunctionReturn(0); 1466ff14e315SSatish Balay } 14670754003eSLois Curfman McInnes 14684a2ae208SSatish Balay #undef __FUNCT__ 14698c778c55SBarry Smith #define __FUNCT__ "MatDenseGetArray" 1470dec5eb66SMatthew G Knepley /*@C 14718c778c55SBarry Smith MatDenseGetArray - gives access to the array where the data for a SeqDense matrix is stored 147273a71a0fSBarry Smith 147373a71a0fSBarry Smith Not Collective 147473a71a0fSBarry Smith 147573a71a0fSBarry Smith Input Parameter: 147673a71a0fSBarry Smith . mat - a MATSEQDENSE matrix 147773a71a0fSBarry Smith 147873a71a0fSBarry Smith Output Parameter: 147973a71a0fSBarry Smith . array - pointer to the data 148073a71a0fSBarry Smith 148173a71a0fSBarry Smith Level: intermediate 148273a71a0fSBarry Smith 14838c778c55SBarry Smith .seealso: MatDenseRestoreArray() 148473a71a0fSBarry Smith @*/ 14858c778c55SBarry Smith PetscErrorCode MatDenseGetArray(Mat A,PetscScalar **array) 148673a71a0fSBarry Smith { 148773a71a0fSBarry Smith PetscErrorCode ierr; 148873a71a0fSBarry Smith 148973a71a0fSBarry Smith PetscFunctionBegin; 14908c778c55SBarry Smith ierr = PetscUseMethod(A,"MatDenseGetArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr); 149173a71a0fSBarry Smith PetscFunctionReturn(0); 149273a71a0fSBarry Smith } 149373a71a0fSBarry Smith 149473a71a0fSBarry Smith #undef __FUNCT__ 14958c778c55SBarry Smith #define __FUNCT__ "MatDenseRestoreArray" 1496dec5eb66SMatthew G Knepley /*@C 14978c778c55SBarry Smith MatDenseRestoreArray - returns access to the array where the data for a SeqDense matrix is stored obtained by MatDenseGetArray() 149873a71a0fSBarry Smith 149973a71a0fSBarry Smith Not Collective 150073a71a0fSBarry Smith 150173a71a0fSBarry Smith Input Parameters: 150273a71a0fSBarry Smith . mat - a MATSEQDENSE matrix 150373a71a0fSBarry Smith . array - pointer to the data 150473a71a0fSBarry Smith 150573a71a0fSBarry Smith Level: intermediate 150673a71a0fSBarry Smith 15078c778c55SBarry Smith .seealso: MatDenseGetArray() 150873a71a0fSBarry Smith @*/ 15098c778c55SBarry Smith PetscErrorCode MatDenseRestoreArray(Mat A,PetscScalar **array) 151073a71a0fSBarry Smith { 151173a71a0fSBarry Smith PetscErrorCode ierr; 151273a71a0fSBarry Smith 151373a71a0fSBarry Smith PetscFunctionBegin; 15148c778c55SBarry Smith ierr = PetscUseMethod(A,"MatDenseRestoreArray_C",(Mat,PetscScalar**),(A,array));CHKERRQ(ierr); 151573a71a0fSBarry Smith PetscFunctionReturn(0); 151673a71a0fSBarry Smith } 151773a71a0fSBarry Smith 151873a71a0fSBarry Smith #undef __FUNCT__ 15194a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrix_SeqDense" 152013f74950SBarry Smith static PetscErrorCode MatGetSubMatrix_SeqDense(Mat A,IS isrow,IS iscol,PetscInt cs,MatReuse scall,Mat *B) 15210754003eSLois Curfman McInnes { 1522c0bbcb79SLois Curfman McInnes Mat_SeqDense *mat = (Mat_SeqDense*)A->data; 15236849ba73SBarry Smith PetscErrorCode ierr; 15245d0c19d7SBarry Smith PetscInt i,j,nrows,ncols; 15255d0c19d7SBarry Smith const PetscInt *irow,*icol; 152687828ca2SBarry Smith PetscScalar *av,*bv,*v = mat->v; 15270754003eSLois Curfman McInnes Mat newmat; 15280754003eSLois Curfman McInnes 15293a40ed3dSBarry Smith PetscFunctionBegin; 153078b31e54SBarry Smith ierr = ISGetIndices(isrow,&irow);CHKERRQ(ierr); 153178b31e54SBarry Smith ierr = ISGetIndices(iscol,&icol);CHKERRQ(ierr); 1532e03a110bSBarry Smith ierr = ISGetLocalSize(isrow,&nrows);CHKERRQ(ierr); 1533e03a110bSBarry Smith ierr = ISGetLocalSize(iscol,&ncols);CHKERRQ(ierr); 15340754003eSLois Curfman McInnes 1535182d2002SSatish Balay /* Check submatrixcall */ 1536182d2002SSatish Balay if (scall == MAT_REUSE_MATRIX) { 153713f74950SBarry Smith PetscInt n_cols,n_rows; 1538182d2002SSatish Balay ierr = MatGetSize(*B,&n_rows,&n_cols);CHKERRQ(ierr); 153921a2c019SBarry Smith if (n_rows != nrows || n_cols != ncols) { 1540f746d493SDmitry Karpeev /* resize the result matrix to match number of requested rows/columns */ 1541c61587bbSBarry Smith ierr = MatSetSizes(*B,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 154221a2c019SBarry Smith } 1543182d2002SSatish Balay newmat = *B; 1544182d2002SSatish Balay } else { 15450754003eSLois Curfman McInnes /* Create and fill new matrix */ 1546ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)A),&newmat);CHKERRQ(ierr); 1547f69a0ea3SMatthew Knepley ierr = MatSetSizes(newmat,nrows,ncols,nrows,ncols);CHKERRQ(ierr); 15487adad957SLisandro Dalcin ierr = MatSetType(newmat,((PetscObject)A)->type_name);CHKERRQ(ierr); 15490298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(newmat,NULL);CHKERRQ(ierr); 1550182d2002SSatish Balay } 1551182d2002SSatish Balay 1552182d2002SSatish Balay /* Now extract the data pointers and do the copy,column at a time */ 1553182d2002SSatish Balay bv = ((Mat_SeqDense*)newmat->data)->v; 1554182d2002SSatish Balay 1555182d2002SSatish Balay for (i=0; i<ncols; i++) { 15566de62eeeSBarry Smith av = v + mat->lda*icol[i]; 15572205254eSKarl Rupp for (j=0; j<nrows; j++) *bv++ = av[irow[j]]; 15580754003eSLois Curfman McInnes } 1559182d2002SSatish Balay 1560182d2002SSatish Balay /* Assemble the matrices so that the correct flags are set */ 15616d4a8577SBarry Smith ierr = MatAssemblyBegin(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15626d4a8577SBarry Smith ierr = MatAssemblyEnd(newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15630754003eSLois Curfman McInnes 15640754003eSLois Curfman McInnes /* Free work space */ 156578b31e54SBarry Smith ierr = ISRestoreIndices(isrow,&irow);CHKERRQ(ierr); 156678b31e54SBarry Smith ierr = ISRestoreIndices(iscol,&icol);CHKERRQ(ierr); 1567182d2002SSatish Balay *B = newmat; 15683a40ed3dSBarry Smith PetscFunctionReturn(0); 15690754003eSLois Curfman McInnes } 15700754003eSLois Curfman McInnes 15714a2ae208SSatish Balay #undef __FUNCT__ 15724a2ae208SSatish Balay #define __FUNCT__ "MatGetSubMatrices_SeqDense" 157313f74950SBarry Smith PetscErrorCode MatGetSubMatrices_SeqDense(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *B[]) 1574905e6a2fSBarry Smith { 15756849ba73SBarry Smith PetscErrorCode ierr; 157613f74950SBarry Smith PetscInt i; 1577905e6a2fSBarry Smith 15783a40ed3dSBarry Smith PetscFunctionBegin; 1579905e6a2fSBarry Smith if (scall == MAT_INITIAL_MATRIX) { 1580b0a32e0cSBarry Smith ierr = PetscMalloc((n+1)*sizeof(Mat),B);CHKERRQ(ierr); 1581905e6a2fSBarry Smith } 1582905e6a2fSBarry Smith 1583905e6a2fSBarry Smith for (i=0; i<n; i++) { 15846a6a5d1dSBarry Smith ierr = MatGetSubMatrix_SeqDense(A,irow[i],icol[i],PETSC_DECIDE,scall,&(*B)[i]);CHKERRQ(ierr); 1585905e6a2fSBarry Smith } 15863a40ed3dSBarry Smith PetscFunctionReturn(0); 1587905e6a2fSBarry Smith } 1588905e6a2fSBarry Smith 15894a2ae208SSatish Balay #undef __FUNCT__ 1590c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyBegin_SeqDense" 1591c0aa2d19SHong Zhang PetscErrorCode MatAssemblyBegin_SeqDense(Mat mat,MatAssemblyType mode) 1592c0aa2d19SHong Zhang { 1593c0aa2d19SHong Zhang PetscFunctionBegin; 1594c0aa2d19SHong Zhang PetscFunctionReturn(0); 1595c0aa2d19SHong Zhang } 1596c0aa2d19SHong Zhang 1597c0aa2d19SHong Zhang #undef __FUNCT__ 1598c0aa2d19SHong Zhang #define __FUNCT__ "MatAssemblyEnd_SeqDense" 1599c0aa2d19SHong Zhang PetscErrorCode MatAssemblyEnd_SeqDense(Mat mat,MatAssemblyType mode) 1600c0aa2d19SHong Zhang { 1601c0aa2d19SHong Zhang PetscFunctionBegin; 1602c0aa2d19SHong Zhang PetscFunctionReturn(0); 1603c0aa2d19SHong Zhang } 1604c0aa2d19SHong Zhang 1605c0aa2d19SHong Zhang #undef __FUNCT__ 16064a2ae208SSatish Balay #define __FUNCT__ "MatCopy_SeqDense" 1607dfbe8321SBarry Smith PetscErrorCode MatCopy_SeqDense(Mat A,Mat B,MatStructure str) 16084b0e389bSBarry Smith { 16094b0e389bSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data,*b = (Mat_SeqDense*)B->data; 16106849ba73SBarry Smith PetscErrorCode ierr; 1611d0f46423SBarry Smith PetscInt lda1=a->lda,lda2=b->lda, m=A->rmap->n,n=A->cmap->n, j; 16123a40ed3dSBarry Smith 16133a40ed3dSBarry Smith PetscFunctionBegin; 161433f4a19fSKris Buschelman /* If the two matrices don't have the same copy implementation, they aren't compatible for fast copy. */ 161533f4a19fSKris Buschelman if (A->ops->copy != B->ops->copy) { 1616cb5b572fSBarry Smith ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 16173a40ed3dSBarry Smith PetscFunctionReturn(0); 16183a40ed3dSBarry Smith } 1619e32f2f54SBarry Smith if (m != B->rmap->n || n != B->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"size(B) != size(A)"); 1620a5ce6ee0Svictorle if (lda1>m || lda2>m) { 16210dbb7854Svictorle for (j=0; j<n; j++) { 1622a5ce6ee0Svictorle ierr = PetscMemcpy(b->v+j*lda2,a->v+j*lda1,m*sizeof(PetscScalar));CHKERRQ(ierr); 1623a5ce6ee0Svictorle } 1624a5ce6ee0Svictorle } else { 1625d0f46423SBarry Smith ierr = PetscMemcpy(b->v,a->v,A->rmap->n*A->cmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 1626a5ce6ee0Svictorle } 1627273d9f13SBarry Smith PetscFunctionReturn(0); 1628273d9f13SBarry Smith } 1629273d9f13SBarry Smith 16304a2ae208SSatish Balay #undef __FUNCT__ 16314994cf47SJed Brown #define __FUNCT__ "MatSetUp_SeqDense" 16324994cf47SJed Brown PetscErrorCode MatSetUp_SeqDense(Mat A) 1633273d9f13SBarry Smith { 1634dfbe8321SBarry Smith PetscErrorCode ierr; 1635273d9f13SBarry Smith 1636273d9f13SBarry Smith PetscFunctionBegin; 1637273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(A,0);CHKERRQ(ierr); 16383a40ed3dSBarry Smith PetscFunctionReturn(0); 16394b0e389bSBarry Smith } 16404b0e389bSBarry Smith 1641284134d9SBarry Smith #undef __FUNCT__ 1642ba337c44SJed Brown #define __FUNCT__ "MatConjugate_SeqDense" 1643ba337c44SJed Brown static PetscErrorCode MatConjugate_SeqDense(Mat A) 1644ba337c44SJed Brown { 1645ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1646ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1647ba337c44SJed Brown PetscScalar *aa = a->v; 1648ba337c44SJed Brown 1649ba337c44SJed Brown PetscFunctionBegin; 1650ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscConj(aa[i]); 1651ba337c44SJed Brown PetscFunctionReturn(0); 1652ba337c44SJed Brown } 1653ba337c44SJed Brown 1654ba337c44SJed Brown #undef __FUNCT__ 1655ba337c44SJed Brown #define __FUNCT__ "MatRealPart_SeqDense" 1656ba337c44SJed Brown static PetscErrorCode MatRealPart_SeqDense(Mat A) 1657ba337c44SJed Brown { 1658ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1659ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1660ba337c44SJed Brown PetscScalar *aa = a->v; 1661ba337c44SJed Brown 1662ba337c44SJed Brown PetscFunctionBegin; 1663ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscRealPart(aa[i]); 1664ba337c44SJed Brown PetscFunctionReturn(0); 1665ba337c44SJed Brown } 1666ba337c44SJed Brown 1667ba337c44SJed Brown #undef __FUNCT__ 1668ba337c44SJed Brown #define __FUNCT__ "MatImaginaryPart_SeqDense" 1669ba337c44SJed Brown static PetscErrorCode MatImaginaryPart_SeqDense(Mat A) 1670ba337c44SJed Brown { 1671ba337c44SJed Brown Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1672ba337c44SJed Brown PetscInt i,nz = A->rmap->n*A->cmap->n; 1673ba337c44SJed Brown PetscScalar *aa = a->v; 1674ba337c44SJed Brown 1675ba337c44SJed Brown PetscFunctionBegin; 1676ba337c44SJed Brown for (i=0; i<nz; i++) aa[i] = PetscImaginaryPart(aa[i]); 1677ba337c44SJed Brown PetscFunctionReturn(0); 1678ba337c44SJed Brown } 1679284134d9SBarry Smith 1680a9fe9ddaSSatish Balay /* ----------------------------------------------------------------*/ 1681a9fe9ddaSSatish Balay #undef __FUNCT__ 1682a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMult_SeqDense_SeqDense" 1683a9fe9ddaSSatish Balay PetscErrorCode MatMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1684a9fe9ddaSSatish Balay { 1685a9fe9ddaSSatish Balay PetscErrorCode ierr; 1686a9fe9ddaSSatish Balay 1687a9fe9ddaSSatish Balay PetscFunctionBegin; 1688a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX) { 1689a9fe9ddaSSatish Balay ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1690a9fe9ddaSSatish Balay } 1691a9fe9ddaSSatish Balay ierr = MatMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1692a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1693a9fe9ddaSSatish Balay } 1694a9fe9ddaSSatish Balay 1695a9fe9ddaSSatish Balay #undef __FUNCT__ 1696a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultSymbolic_SeqDense_SeqDense" 1697a9fe9ddaSSatish Balay PetscErrorCode MatMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1698a9fe9ddaSSatish Balay { 1699ee16a9a1SHong Zhang PetscErrorCode ierr; 1700d0f46423SBarry Smith PetscInt m=A->rmap->n,n=B->cmap->n; 1701ee16a9a1SHong Zhang Mat Cmat; 1702a9fe9ddaSSatish Balay 1703ee16a9a1SHong Zhang PetscFunctionBegin; 1704e32f2f54SBarry Smith if (A->cmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->cmap->n %d != B->rmap->n %d\n",A->cmap->n,B->rmap->n); 1705ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1706ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1707ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 17080298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Cmat,NULL);CHKERRQ(ierr); 1709d73949e8SHong Zhang 1710ee16a9a1SHong Zhang *C = Cmat; 1711ee16a9a1SHong Zhang PetscFunctionReturn(0); 1712ee16a9a1SHong Zhang } 1713a9fe9ddaSSatish Balay 171498a3b096SSatish Balay #undef __FUNCT__ 1715a9fe9ddaSSatish Balay #define __FUNCT__ "MatMatMultNumeric_SeqDense_SeqDense" 1716a9fe9ddaSSatish Balay PetscErrorCode MatMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1717a9fe9ddaSSatish Balay { 1718a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1719a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1720a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 17210805154bSBarry Smith PetscBLASInt m,n,k; 1722a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1723c5df96a5SBarry Smith PetscErrorCode ierr; 1724a9fe9ddaSSatish Balay 1725a9fe9ddaSSatish Balay PetscFunctionBegin; 1726c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&m);CHKERRQ(ierr); 1727c5df96a5SBarry Smith ierr = PetscBLASIntCast(B->cmap->n,&n);CHKERRQ(ierr); 1728c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&k);CHKERRQ(ierr); 1729a83cb05cSBarry Smith PetscStackCall("BLASgemv",BLASgemm_("N","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda)); 1730a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1731a9fe9ddaSSatish Balay } 1732a9fe9ddaSSatish Balay 1733a9fe9ddaSSatish Balay #undef __FUNCT__ 173475648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMult_SeqDense_SeqDense" 173575648e8dSHong Zhang PetscErrorCode MatTransposeMatMult_SeqDense_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1736a9fe9ddaSSatish Balay { 1737a9fe9ddaSSatish Balay PetscErrorCode ierr; 1738a9fe9ddaSSatish Balay 1739a9fe9ddaSSatish Balay PetscFunctionBegin; 1740a9fe9ddaSSatish Balay if (scall == MAT_INITIAL_MATRIX) { 174175648e8dSHong Zhang ierr = MatTransposeMatMultSymbolic_SeqDense_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1742a9fe9ddaSSatish Balay } 174375648e8dSHong Zhang ierr = MatTransposeMatMultNumeric_SeqDense_SeqDense(A,B,*C);CHKERRQ(ierr); 1744a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1745a9fe9ddaSSatish Balay } 1746a9fe9ddaSSatish Balay 1747a9fe9ddaSSatish Balay #undef __FUNCT__ 174875648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqDense_SeqDense" 174975648e8dSHong Zhang PetscErrorCode MatTransposeMatMultSymbolic_SeqDense_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1750a9fe9ddaSSatish Balay { 1751ee16a9a1SHong Zhang PetscErrorCode ierr; 1752d0f46423SBarry Smith PetscInt m=A->cmap->n,n=B->cmap->n; 1753ee16a9a1SHong Zhang Mat Cmat; 1754a9fe9ddaSSatish Balay 1755ee16a9a1SHong Zhang PetscFunctionBegin; 1756e32f2f54SBarry Smith if (A->rmap->n != B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"A->rmap->n %d != B->rmap->n %d\n",A->rmap->n,B->rmap->n); 1757ee16a9a1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,&Cmat);CHKERRQ(ierr); 1758ee16a9a1SHong Zhang ierr = MatSetSizes(Cmat,m,n,m,n);CHKERRQ(ierr); 1759ee16a9a1SHong Zhang ierr = MatSetType(Cmat,MATSEQDENSE);CHKERRQ(ierr); 17600298fd71SBarry Smith ierr = MatSeqDenseSetPreallocation(Cmat,NULL);CHKERRQ(ierr); 17612205254eSKarl Rupp 1762ee16a9a1SHong Zhang Cmat->assembled = PETSC_TRUE; 17632205254eSKarl Rupp 1764ee16a9a1SHong Zhang *C = Cmat; 1765ee16a9a1SHong Zhang PetscFunctionReturn(0); 1766ee16a9a1SHong Zhang } 1767a9fe9ddaSSatish Balay 1768a9fe9ddaSSatish Balay #undef __FUNCT__ 176975648e8dSHong Zhang #define __FUNCT__ "MatTransposeMatMultNumeric_SeqDense_SeqDense" 177075648e8dSHong Zhang PetscErrorCode MatTransposeMatMultNumeric_SeqDense_SeqDense(Mat A,Mat B,Mat C) 1771a9fe9ddaSSatish Balay { 1772a9fe9ddaSSatish Balay Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1773a9fe9ddaSSatish Balay Mat_SeqDense *b = (Mat_SeqDense*)B->data; 1774a9fe9ddaSSatish Balay Mat_SeqDense *c = (Mat_SeqDense*)C->data; 17750805154bSBarry Smith PetscBLASInt m,n,k; 1776a9fe9ddaSSatish Balay PetscScalar _DOne=1.0,_DZero=0.0; 1777c5df96a5SBarry Smith PetscErrorCode ierr; 1778a9fe9ddaSSatish Balay 1779a9fe9ddaSSatish Balay PetscFunctionBegin; 1780c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->cmap->n,&m);CHKERRQ(ierr); 1781c5df96a5SBarry Smith ierr = PetscBLASIntCast(B->cmap->n,&n);CHKERRQ(ierr); 1782c5df96a5SBarry Smith ierr = PetscBLASIntCast(A->rmap->n,&k);CHKERRQ(ierr); 17832fbe02b9SBarry Smith /* 17842fbe02b9SBarry Smith Note the m and n arguments below are the number rows and columns of A', not A! 17852fbe02b9SBarry Smith */ 1786a83cb05cSBarry Smith PetscStackCall("BLASgemv",BLASgemm_("T","N",&m,&n,&k,&_DOne,a->v,&a->lda,b->v,&b->lda,&_DZero,c->v,&c->lda)); 1787a9fe9ddaSSatish Balay PetscFunctionReturn(0); 1788a9fe9ddaSSatish Balay } 1789985db425SBarry Smith 1790985db425SBarry Smith #undef __FUNCT__ 1791985db425SBarry Smith #define __FUNCT__ "MatGetRowMax_SeqDense" 1792985db425SBarry Smith PetscErrorCode MatGetRowMax_SeqDense(Mat A,Vec v,PetscInt idx[]) 1793985db425SBarry Smith { 1794985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1795985db425SBarry Smith PetscErrorCode ierr; 1796d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1797985db425SBarry Smith PetscScalar *x; 1798985db425SBarry Smith MatScalar *aa = a->v; 1799985db425SBarry Smith 1800985db425SBarry Smith PetscFunctionBegin; 1801e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1802985db425SBarry Smith 1803985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1804985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1805985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1806e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1807985db425SBarry Smith for (i=0; i<m; i++) { 1808985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1809985db425SBarry Smith for (j=1; j<n; j++) { 1810985db425SBarry Smith if (PetscRealPart(x[i]) < PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1811985db425SBarry Smith } 1812985db425SBarry Smith } 1813985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1814985db425SBarry Smith PetscFunctionReturn(0); 1815985db425SBarry Smith } 1816985db425SBarry Smith 1817985db425SBarry Smith #undef __FUNCT__ 1818985db425SBarry Smith #define __FUNCT__ "MatGetRowMaxAbs_SeqDense" 1819985db425SBarry Smith PetscErrorCode MatGetRowMaxAbs_SeqDense(Mat A,Vec v,PetscInt idx[]) 1820985db425SBarry Smith { 1821985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1822985db425SBarry Smith PetscErrorCode ierr; 1823d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1824985db425SBarry Smith PetscScalar *x; 1825985db425SBarry Smith PetscReal atmp; 1826985db425SBarry Smith MatScalar *aa = a->v; 1827985db425SBarry Smith 1828985db425SBarry Smith PetscFunctionBegin; 1829e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1830985db425SBarry Smith 1831985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1832985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1833985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1834e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1835985db425SBarry Smith for (i=0; i<m; i++) { 18369189402eSHong Zhang x[i] = PetscAbsScalar(aa[i]); 1837985db425SBarry Smith for (j=1; j<n; j++) { 1838985db425SBarry Smith atmp = PetscAbsScalar(aa[i+m*j]); 1839985db425SBarry Smith if (PetscAbsScalar(x[i]) < atmp) {x[i] = atmp; if (idx) idx[i] = j;} 1840985db425SBarry Smith } 1841985db425SBarry Smith } 1842985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1843985db425SBarry Smith PetscFunctionReturn(0); 1844985db425SBarry Smith } 1845985db425SBarry Smith 1846985db425SBarry Smith #undef __FUNCT__ 1847985db425SBarry Smith #define __FUNCT__ "MatGetRowMin_SeqDense" 1848985db425SBarry Smith PetscErrorCode MatGetRowMin_SeqDense(Mat A,Vec v,PetscInt idx[]) 1849985db425SBarry Smith { 1850985db425SBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 1851985db425SBarry Smith PetscErrorCode ierr; 1852d0f46423SBarry Smith PetscInt i,j,m = A->rmap->n,n = A->cmap->n,p; 1853985db425SBarry Smith PetscScalar *x; 1854985db425SBarry Smith MatScalar *aa = a->v; 1855985db425SBarry Smith 1856985db425SBarry Smith PetscFunctionBegin; 1857e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 1858985db425SBarry Smith 1859985db425SBarry Smith ierr = VecSet(v,0.0);CHKERRQ(ierr); 1860985db425SBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1861985db425SBarry Smith ierr = VecGetLocalSize(v,&p);CHKERRQ(ierr); 1862e32f2f54SBarry Smith if (p != A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Nonconforming matrix and vector"); 1863985db425SBarry Smith for (i=0; i<m; i++) { 1864985db425SBarry Smith x[i] = aa[i]; if (idx) idx[i] = 0; 1865985db425SBarry Smith for (j=1; j<n; j++) { 1866985db425SBarry Smith if (PetscRealPart(x[i]) > PetscRealPart(aa[i+m*j])) {x[i] = aa[i + m*j]; if (idx) idx[i] = j;} 1867985db425SBarry Smith } 1868985db425SBarry Smith } 1869985db425SBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 1870985db425SBarry Smith PetscFunctionReturn(0); 1871985db425SBarry Smith } 1872985db425SBarry Smith 18738d0534beSBarry Smith #undef __FUNCT__ 18748d0534beSBarry Smith #define __FUNCT__ "MatGetColumnVector_SeqDense" 18758d0534beSBarry Smith PetscErrorCode MatGetColumnVector_SeqDense(Mat A,Vec v,PetscInt col) 18768d0534beSBarry Smith { 18778d0534beSBarry Smith Mat_SeqDense *a = (Mat_SeqDense*)A->data; 18788d0534beSBarry Smith PetscErrorCode ierr; 18798d0534beSBarry Smith PetscScalar *x; 18808d0534beSBarry Smith 18818d0534beSBarry Smith PetscFunctionBegin; 1882e32f2f54SBarry Smith if (A->factortype) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 18838d0534beSBarry Smith 18848d0534beSBarry Smith ierr = VecGetArray(v,&x);CHKERRQ(ierr); 1885d0f46423SBarry Smith ierr = PetscMemcpy(x,a->v+col*a->lda,A->rmap->n*sizeof(PetscScalar));CHKERRQ(ierr); 18868d0534beSBarry Smith ierr = VecRestoreArray(v,&x);CHKERRQ(ierr); 18878d0534beSBarry Smith PetscFunctionReturn(0); 18888d0534beSBarry Smith } 18898d0534beSBarry Smith 18900716a85fSBarry Smith 18910716a85fSBarry Smith #undef __FUNCT__ 18920716a85fSBarry Smith #define __FUNCT__ "MatGetColumnNorms_SeqDense" 18930716a85fSBarry Smith PetscErrorCode MatGetColumnNorms_SeqDense(Mat A,NormType type,PetscReal *norms) 18940716a85fSBarry Smith { 18950716a85fSBarry Smith PetscErrorCode ierr; 18960716a85fSBarry Smith PetscInt i,j,m,n; 18970716a85fSBarry Smith PetscScalar *a; 18980716a85fSBarry Smith 18990716a85fSBarry Smith PetscFunctionBegin; 19000716a85fSBarry Smith ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr); 19010716a85fSBarry Smith ierr = PetscMemzero(norms,n*sizeof(PetscReal));CHKERRQ(ierr); 19028c778c55SBarry Smith ierr = MatDenseGetArray(A,&a);CHKERRQ(ierr); 19030716a85fSBarry Smith if (type == NORM_2) { 19040716a85fSBarry Smith for (i=0; i<n; i++) { 19050716a85fSBarry Smith for (j=0; j<m; j++) { 19060716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]*a[j]); 19070716a85fSBarry Smith } 19080716a85fSBarry Smith a += m; 19090716a85fSBarry Smith } 19100716a85fSBarry Smith } else if (type == NORM_1) { 19110716a85fSBarry Smith for (i=0; i<n; i++) { 19120716a85fSBarry Smith for (j=0; j<m; j++) { 19130716a85fSBarry Smith norms[i] += PetscAbsScalar(a[j]); 19140716a85fSBarry Smith } 19150716a85fSBarry Smith a += m; 19160716a85fSBarry Smith } 19170716a85fSBarry Smith } else if (type == NORM_INFINITY) { 19180716a85fSBarry Smith for (i=0; i<n; i++) { 19190716a85fSBarry Smith for (j=0; j<m; j++) { 19200716a85fSBarry Smith norms[i] = PetscMax(PetscAbsScalar(a[j]),norms[i]); 19210716a85fSBarry Smith } 19220716a85fSBarry Smith a += m; 19230716a85fSBarry Smith } 1924ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Unknown NormType"); 19258c778c55SBarry Smith ierr = MatDenseRestoreArray(A,&a);CHKERRQ(ierr); 19260716a85fSBarry Smith if (type == NORM_2) { 19278f1a2a5eSBarry Smith for (i=0; i<n; i++) norms[i] = PetscSqrtReal(norms[i]); 19280716a85fSBarry Smith } 19290716a85fSBarry Smith PetscFunctionReturn(0); 19300716a85fSBarry Smith } 19310716a85fSBarry Smith 193273a71a0fSBarry Smith #undef __FUNCT__ 193373a71a0fSBarry Smith #define __FUNCT__ "MatSetRandom_SeqDense" 193473a71a0fSBarry Smith static PetscErrorCode MatSetRandom_SeqDense(Mat x,PetscRandom rctx) 193573a71a0fSBarry Smith { 193673a71a0fSBarry Smith PetscErrorCode ierr; 193773a71a0fSBarry Smith PetscScalar *a; 193873a71a0fSBarry Smith PetscInt m,n,i; 193973a71a0fSBarry Smith 194073a71a0fSBarry Smith PetscFunctionBegin; 194173a71a0fSBarry Smith ierr = MatGetSize(x,&m,&n);CHKERRQ(ierr); 19428c778c55SBarry Smith ierr = MatDenseGetArray(x,&a);CHKERRQ(ierr); 194373a71a0fSBarry Smith for (i=0; i<m*n; i++) { 194473a71a0fSBarry Smith ierr = PetscRandomGetValue(rctx,a+i);CHKERRQ(ierr); 194573a71a0fSBarry Smith } 19468c778c55SBarry Smith ierr = MatDenseRestoreArray(x,&a);CHKERRQ(ierr); 194773a71a0fSBarry Smith PetscFunctionReturn(0); 194873a71a0fSBarry Smith } 194973a71a0fSBarry Smith 195073a71a0fSBarry Smith 1951289bc588SBarry Smith /* -------------------------------------------------------------------*/ 1952a5ae1ecdSBarry Smith static struct _MatOps MatOps_Values = { MatSetValues_SeqDense, 1953905e6a2fSBarry Smith MatGetRow_SeqDense, 1954905e6a2fSBarry Smith MatRestoreRow_SeqDense, 1955905e6a2fSBarry Smith MatMult_SeqDense, 195697304618SKris Buschelman /* 4*/ MatMultAdd_SeqDense, 19577c922b88SBarry Smith MatMultTranspose_SeqDense, 19587c922b88SBarry Smith MatMultTransposeAdd_SeqDense, 1959db4efbfdSBarry Smith 0, 1960db4efbfdSBarry Smith 0, 1961db4efbfdSBarry Smith 0, 1962db4efbfdSBarry Smith /* 10*/ 0, 1963905e6a2fSBarry Smith MatLUFactor_SeqDense, 1964905e6a2fSBarry Smith MatCholeskyFactor_SeqDense, 196541f059aeSBarry Smith MatSOR_SeqDense, 1966ec8511deSBarry Smith MatTranspose_SeqDense, 196797304618SKris Buschelman /* 15*/ MatGetInfo_SeqDense, 1968905e6a2fSBarry Smith MatEqual_SeqDense, 1969905e6a2fSBarry Smith MatGetDiagonal_SeqDense, 1970905e6a2fSBarry Smith MatDiagonalScale_SeqDense, 1971905e6a2fSBarry Smith MatNorm_SeqDense, 1972c0aa2d19SHong Zhang /* 20*/ MatAssemblyBegin_SeqDense, 1973c0aa2d19SHong Zhang MatAssemblyEnd_SeqDense, 1974905e6a2fSBarry Smith MatSetOption_SeqDense, 1975905e6a2fSBarry Smith MatZeroEntries_SeqDense, 1976d519adbfSMatthew Knepley /* 24*/ MatZeroRows_SeqDense, 1977db4efbfdSBarry Smith 0, 1978db4efbfdSBarry Smith 0, 1979db4efbfdSBarry Smith 0, 1980db4efbfdSBarry Smith 0, 19814994cf47SJed Brown /* 29*/ MatSetUp_SeqDense, 1982273d9f13SBarry Smith 0, 1983905e6a2fSBarry Smith 0, 198473a71a0fSBarry Smith 0, 198573a71a0fSBarry Smith 0, 1986d519adbfSMatthew Knepley /* 34*/ MatDuplicate_SeqDense, 1987a5ae1ecdSBarry Smith 0, 1988a5ae1ecdSBarry Smith 0, 1989a5ae1ecdSBarry Smith 0, 1990a5ae1ecdSBarry Smith 0, 1991d519adbfSMatthew Knepley /* 39*/ MatAXPY_SeqDense, 1992a5ae1ecdSBarry Smith MatGetSubMatrices_SeqDense, 1993a5ae1ecdSBarry Smith 0, 19944b0e389bSBarry Smith MatGetValues_SeqDense, 1995a5ae1ecdSBarry Smith MatCopy_SeqDense, 1996d519adbfSMatthew Knepley /* 44*/ MatGetRowMax_SeqDense, 1997a5ae1ecdSBarry Smith MatScale_SeqDense, 1998a5ae1ecdSBarry Smith 0, 1999a5ae1ecdSBarry Smith 0, 2000a5ae1ecdSBarry Smith 0, 200173a71a0fSBarry Smith /* 49*/ MatSetRandom_SeqDense, 2002a5ae1ecdSBarry Smith 0, 2003a5ae1ecdSBarry Smith 0, 2004a5ae1ecdSBarry Smith 0, 2005a5ae1ecdSBarry Smith 0, 2006d519adbfSMatthew Knepley /* 54*/ 0, 2007a5ae1ecdSBarry Smith 0, 2008a5ae1ecdSBarry Smith 0, 2009a5ae1ecdSBarry Smith 0, 2010a5ae1ecdSBarry Smith 0, 2011d519adbfSMatthew Knepley /* 59*/ 0, 2012e03a110bSBarry Smith MatDestroy_SeqDense, 2013e03a110bSBarry Smith MatView_SeqDense, 2014357abbc8SBarry Smith 0, 201597304618SKris Buschelman 0, 2016d519adbfSMatthew Knepley /* 64*/ 0, 201797304618SKris Buschelman 0, 201897304618SKris Buschelman 0, 201997304618SKris Buschelman 0, 202097304618SKris Buschelman 0, 2021d519adbfSMatthew Knepley /* 69*/ MatGetRowMaxAbs_SeqDense, 202297304618SKris Buschelman 0, 202397304618SKris Buschelman 0, 202497304618SKris Buschelman 0, 202597304618SKris Buschelman 0, 2026d519adbfSMatthew Knepley /* 74*/ 0, 202797304618SKris Buschelman 0, 202897304618SKris Buschelman 0, 202997304618SKris Buschelman 0, 203097304618SKris Buschelman 0, 2031d519adbfSMatthew Knepley /* 79*/ 0, 203297304618SKris Buschelman 0, 203397304618SKris Buschelman 0, 203497304618SKris Buschelman 0, 20355bba2384SShri Abhyankar /* 83*/ MatLoad_SeqDense, 2036865e5f61SKris Buschelman 0, 20371cbb95d3SBarry Smith MatIsHermitian_SeqDense, 2038865e5f61SKris Buschelman 0, 2039865e5f61SKris Buschelman 0, 2040865e5f61SKris Buschelman 0, 2041d519adbfSMatthew Knepley /* 89*/ MatMatMult_SeqDense_SeqDense, 2042a9fe9ddaSSatish Balay MatMatMultSymbolic_SeqDense_SeqDense, 2043a9fe9ddaSSatish Balay MatMatMultNumeric_SeqDense_SeqDense, 2044865e5f61SKris Buschelman 0, 2045865e5f61SKris Buschelman 0, 2046d519adbfSMatthew Knepley /* 94*/ 0, 20475df89d91SHong Zhang 0, 20485df89d91SHong Zhang 0, 20495df89d91SHong Zhang 0, 2050284134d9SBarry Smith 0, 2051d519adbfSMatthew Knepley /* 99*/ 0, 2052284134d9SBarry Smith 0, 2053284134d9SBarry Smith 0, 2054ba337c44SJed Brown MatConjugate_SeqDense, 2055f73d5cc4SBarry Smith 0, 2056ba337c44SJed Brown /*104*/ 0, 2057ba337c44SJed Brown MatRealPart_SeqDense, 2058ba337c44SJed Brown MatImaginaryPart_SeqDense, 2059985db425SBarry Smith 0, 2060985db425SBarry Smith 0, 206185e2c93fSHong Zhang /*109*/ MatMatSolve_SeqDense, 2062985db425SBarry Smith 0, 20638d0534beSBarry Smith MatGetRowMin_SeqDense, 2064aabbc4fbSShri Abhyankar MatGetColumnVector_SeqDense, 2065aabbc4fbSShri Abhyankar 0, 2066aabbc4fbSShri Abhyankar /*114*/ 0, 2067aabbc4fbSShri Abhyankar 0, 2068aabbc4fbSShri Abhyankar 0, 2069aabbc4fbSShri Abhyankar 0, 2070aabbc4fbSShri Abhyankar 0, 2071aabbc4fbSShri Abhyankar /*119*/ 0, 2072aabbc4fbSShri Abhyankar 0, 2073aabbc4fbSShri Abhyankar 0, 20740716a85fSBarry Smith 0, 20750716a85fSBarry Smith 0, 20760716a85fSBarry Smith /*124*/ 0, 20775df89d91SHong Zhang MatGetColumnNorms_SeqDense, 20785df89d91SHong Zhang 0, 20795df89d91SHong Zhang 0, 20805df89d91SHong Zhang 0, 20815df89d91SHong Zhang /*129*/ 0, 208275648e8dSHong Zhang MatTransposeMatMult_SeqDense_SeqDense, 208375648e8dSHong Zhang MatTransposeMatMultSymbolic_SeqDense_SeqDense, 208475648e8dSHong Zhang MatTransposeMatMultNumeric_SeqDense_SeqDense, 2085985db425SBarry Smith }; 208690ace30eSBarry Smith 20874a2ae208SSatish Balay #undef __FUNCT__ 20884a2ae208SSatish Balay #define __FUNCT__ "MatCreateSeqDense" 20894b828684SBarry Smith /*@C 2090fafbff53SBarry Smith MatCreateSeqDense - Creates a sequential dense matrix that 2091d65003e9SLois Curfman McInnes is stored in column major order (the usual Fortran 77 manner). Many 2092d65003e9SLois Curfman McInnes of the matrix operations use the BLAS and LAPACK routines. 2093289bc588SBarry Smith 2094db81eaa0SLois Curfman McInnes Collective on MPI_Comm 2095db81eaa0SLois Curfman McInnes 209620563c6bSBarry Smith Input Parameters: 2097db81eaa0SLois Curfman McInnes + comm - MPI communicator, set to PETSC_COMM_SELF 20980c775827SLois Curfman McInnes . m - number of rows 209918f449edSLois Curfman McInnes . n - number of columns 21000298fd71SBarry Smith - data - optional location of matrix data in column major order. Set data=NULL for PETSc 2101dfc5480cSLois Curfman McInnes to control all matrix memory allocation. 210220563c6bSBarry Smith 210320563c6bSBarry Smith Output Parameter: 210444cd7ae7SLois Curfman McInnes . A - the matrix 210520563c6bSBarry Smith 2106b259b22eSLois Curfman McInnes Notes: 210718f449edSLois Curfman McInnes The data input variable is intended primarily for Fortran programmers 210818f449edSLois Curfman McInnes who wish to allocate their own matrix memory space. Most users should 21090298fd71SBarry Smith set data=NULL. 211018f449edSLois Curfman McInnes 2111027ccd11SLois Curfman McInnes Level: intermediate 2112027ccd11SLois Curfman McInnes 2113dbd7a890SLois Curfman McInnes .keywords: dense, matrix, LAPACK, BLAS 2114d65003e9SLois Curfman McInnes 211569b1f4b7SBarry Smith .seealso: MatCreate(), MatCreateDense(), MatSetValues() 211620563c6bSBarry Smith @*/ 21177087cfbeSBarry Smith PetscErrorCode MatCreateSeqDense(MPI_Comm comm,PetscInt m,PetscInt n,PetscScalar *data,Mat *A) 2118289bc588SBarry Smith { 2119dfbe8321SBarry Smith PetscErrorCode ierr; 21203b2fbd54SBarry Smith 21213a40ed3dSBarry Smith PetscFunctionBegin; 2122f69a0ea3SMatthew Knepley ierr = MatCreate(comm,A);CHKERRQ(ierr); 2123f69a0ea3SMatthew Knepley ierr = MatSetSizes(*A,m,n,m,n);CHKERRQ(ierr); 2124273d9f13SBarry Smith ierr = MatSetType(*A,MATSEQDENSE);CHKERRQ(ierr); 2125273d9f13SBarry Smith ierr = MatSeqDenseSetPreallocation(*A,data);CHKERRQ(ierr); 2126273d9f13SBarry Smith PetscFunctionReturn(0); 2127273d9f13SBarry Smith } 2128273d9f13SBarry Smith 21294a2ae208SSatish Balay #undef __FUNCT__ 2130afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation" 2131273d9f13SBarry Smith /*@C 2132273d9f13SBarry Smith MatSeqDenseSetPreallocation - Sets the array used for storing the matrix elements 2133273d9f13SBarry Smith 2134273d9f13SBarry Smith Collective on MPI_Comm 2135273d9f13SBarry Smith 2136273d9f13SBarry Smith Input Parameters: 2137273d9f13SBarry Smith + A - the matrix 21380298fd71SBarry Smith - data - the array (or NULL) 2139273d9f13SBarry Smith 2140273d9f13SBarry Smith Notes: 2141273d9f13SBarry Smith The data input variable is intended primarily for Fortran programmers 2142273d9f13SBarry Smith who wish to allocate their own matrix memory space. Most users should 2143284134d9SBarry Smith need not call this routine. 2144273d9f13SBarry Smith 2145273d9f13SBarry Smith Level: intermediate 2146273d9f13SBarry Smith 2147273d9f13SBarry Smith .keywords: dense, matrix, LAPACK, BLAS 2148273d9f13SBarry Smith 214969b1f4b7SBarry Smith .seealso: MatCreate(), MatCreateDense(), MatSetValues(), MatSeqDenseSetLDA() 2150867c911aSBarry Smith 2151273d9f13SBarry Smith @*/ 21527087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation(Mat B,PetscScalar data[]) 2153273d9f13SBarry Smith { 21544ac538c5SBarry Smith PetscErrorCode ierr; 2155a23d5eceSKris Buschelman 2156a23d5eceSKris Buschelman PetscFunctionBegin; 21574ac538c5SBarry Smith ierr = PetscTryMethod(B,"MatSeqDenseSetPreallocation_C",(Mat,PetscScalar[]),(B,data));CHKERRQ(ierr); 2158a23d5eceSKris Buschelman PetscFunctionReturn(0); 2159a23d5eceSKris Buschelman } 2160a23d5eceSKris Buschelman 2161a23d5eceSKris Buschelman #undef __FUNCT__ 2162afc30d2aSLisandro Dalcin #define __FUNCT__ "MatSeqDenseSetPreallocation_SeqDense" 21637087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetPreallocation_SeqDense(Mat B,PetscScalar *data) 2164a23d5eceSKris Buschelman { 2165273d9f13SBarry Smith Mat_SeqDense *b; 2166dfbe8321SBarry Smith PetscErrorCode ierr; 2167273d9f13SBarry Smith 2168273d9f13SBarry Smith PetscFunctionBegin; 2169273d9f13SBarry Smith B->preallocated = PETSC_TRUE; 2170a868139aSShri Abhyankar 217134ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 217234ef9618SShri Abhyankar ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 217334ef9618SShri Abhyankar 2174273d9f13SBarry Smith b = (Mat_SeqDense*)B->data; 217586d161a7SShri Abhyankar b->Mmax = B->rmap->n; 217686d161a7SShri Abhyankar b->Nmax = B->cmap->n; 217786d161a7SShri Abhyankar if (b->lda <= 0 || b->changelda) b->lda = B->rmap->n; 217886d161a7SShri Abhyankar 21799e8f95c4SLisandro Dalcin if (!data) { /* petsc-allocated storage */ 21809e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 21815afd5e0cSBarry Smith ierr = PetscMalloc(b->lda*b->Nmax*sizeof(PetscScalar),&b->v);CHKERRQ(ierr); 2182284134d9SBarry Smith ierr = PetscMemzero(b->v,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 2183284134d9SBarry Smith ierr = PetscLogObjectMemory(B,b->lda*b->Nmax*sizeof(PetscScalar));CHKERRQ(ierr); 21842205254eSKarl Rupp 21859e8f95c4SLisandro Dalcin b->user_alloc = PETSC_FALSE; 2186273d9f13SBarry Smith } else { /* user-allocated storage */ 21879e8f95c4SLisandro Dalcin if (!b->user_alloc) { ierr = PetscFree(b->v);CHKERRQ(ierr); } 2188273d9f13SBarry Smith b->v = data; 2189273d9f13SBarry Smith b->user_alloc = PETSC_TRUE; 2190273d9f13SBarry Smith } 21910450473dSBarry Smith B->assembled = PETSC_TRUE; 2192273d9f13SBarry Smith PetscFunctionReturn(0); 2193273d9f13SBarry Smith } 2194273d9f13SBarry Smith 21951b807ce4Svictorle #undef __FUNCT__ 21961b807ce4Svictorle #define __FUNCT__ "MatSeqDenseSetLDA" 21971b807ce4Svictorle /*@C 21981b807ce4Svictorle MatSeqDenseSetLDA - Declare the leading dimension of the user-provided array 21991b807ce4Svictorle 22001b807ce4Svictorle Input parameter: 22011b807ce4Svictorle + A - the matrix 22021b807ce4Svictorle - lda - the leading dimension 22031b807ce4Svictorle 22041b807ce4Svictorle Notes: 2205867c911aSBarry Smith This routine is to be used in conjunction with MatSeqDenseSetPreallocation(); 22061b807ce4Svictorle it asserts that the preallocation has a leading dimension (the LDA parameter 22071b807ce4Svictorle of Blas and Lapack fame) larger than M, the first dimension of the matrix. 22081b807ce4Svictorle 22091b807ce4Svictorle Level: intermediate 22101b807ce4Svictorle 22111b807ce4Svictorle .keywords: dense, matrix, LAPACK, BLAS 22121b807ce4Svictorle 2213284134d9SBarry Smith .seealso: MatCreate(), MatCreateSeqDense(), MatSeqDenseSetPreallocation(), MatSetMaximumSize() 2214867c911aSBarry Smith 22151b807ce4Svictorle @*/ 22167087cfbeSBarry Smith PetscErrorCode MatSeqDenseSetLDA(Mat B,PetscInt lda) 22171b807ce4Svictorle { 22181b807ce4Svictorle Mat_SeqDense *b = (Mat_SeqDense*)B->data; 221921a2c019SBarry Smith 22201b807ce4Svictorle PetscFunctionBegin; 2221e32f2f54SBarry Smith if (lda < B->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"LDA %D must be at least matrix dimension %D",lda,B->rmap->n); 22221b807ce4Svictorle b->lda = lda; 222321a2c019SBarry Smith b->changelda = PETSC_FALSE; 222421a2c019SBarry Smith b->Mmax = PetscMax(b->Mmax,lda); 22251b807ce4Svictorle PetscFunctionReturn(0); 22261b807ce4Svictorle } 22271b807ce4Svictorle 22280bad9183SKris Buschelman /*MC 2229fafad747SKris Buschelman MATSEQDENSE - MATSEQDENSE = "seqdense" - A matrix type to be used for sequential dense matrices. 22300bad9183SKris Buschelman 22310bad9183SKris Buschelman Options Database Keys: 22320bad9183SKris Buschelman . -mat_type seqdense - sets the matrix type to "seqdense" during a call to MatSetFromOptions() 22330bad9183SKris Buschelman 22340bad9183SKris Buschelman Level: beginner 22350bad9183SKris Buschelman 223689665df3SBarry Smith .seealso: MatCreateSeqDense() 223789665df3SBarry Smith 22380bad9183SKris Buschelman M*/ 22390bad9183SKris Buschelman 22404a2ae208SSatish Balay #undef __FUNCT__ 22414a2ae208SSatish Balay #define __FUNCT__ "MatCreate_SeqDense" 2242*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_SeqDense(Mat B) 2243273d9f13SBarry Smith { 2244273d9f13SBarry Smith Mat_SeqDense *b; 2245dfbe8321SBarry Smith PetscErrorCode ierr; 22467c334f02SBarry Smith PetscMPIInt size; 2247273d9f13SBarry Smith 2248273d9f13SBarry Smith PetscFunctionBegin; 2249ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)B),&size);CHKERRQ(ierr); 2250e32f2f54SBarry Smith if (size > 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Comm must be of size 1"); 225155659b69SBarry Smith 225238f2d2fdSLisandro Dalcin ierr = PetscNewLog(B,Mat_SeqDense,&b);CHKERRQ(ierr); 2253549d3d68SSatish Balay ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 225444cd7ae7SLois Curfman McInnes B->data = (void*)b; 225518f449edSLois Curfman McInnes 225644cd7ae7SLois Curfman McInnes b->pivots = 0; 2257273d9f13SBarry Smith b->roworiented = PETSC_TRUE; 2258273d9f13SBarry Smith b->v = 0; 225921a2c019SBarry Smith b->changelda = PETSC_FALSE; 22604e220ebcSLois Curfman McInnes 226100de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseGetArray_C","MatDenseGetArray_SeqDense",MatDenseGetArray_SeqDense);CHKERRQ(ierr); 226200de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatDenseRestoreArray_C","MatDenseRestoreArray_SeqDense",MatDenseRestoreArray_SeqDense);CHKERRQ(ierr); 2263b24902e0SBarry Smith 226400de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_seqdense_seqaij_C","MatConvert_SeqDense_SeqAIJ",MatConvert_SeqDense_SeqAIJ);CHKERRQ(ierr); 226500de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatGetFactor_petsc_C","MatGetFactor_seqdense_petsc",MatGetFactor_seqdense_petsc);CHKERRQ(ierr); 226600de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatSeqDenseSetPreallocation_C","MatSeqDenseSetPreallocation_SeqDense",MatSeqDenseSetPreallocation_SeqDense);CHKERRQ(ierr); 226700de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMult_seqaij_seqdense_C","MatMatMult_SeqAIJ_SeqDense",MatMatMult_SeqAIJ_SeqDense);CHKERRQ(ierr); 226800de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultSymbolic_seqaij_seqdense_C","MatMatMultSymbolic_SeqAIJ_SeqDense",MatMatMultSymbolic_SeqAIJ_SeqDense);CHKERRQ(ierr); 226900de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)B,"MatMatMultNumeric_seqaij_seqdense_C","MatMatMultNumeric_SeqAIJ_SeqDense",MatMatMultNumeric_SeqAIJ_SeqDense);CHKERRQ(ierr); 227017667f90SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)B,MATSEQDENSE);CHKERRQ(ierr); 22713a40ed3dSBarry Smith PetscFunctionReturn(0); 2272289bc588SBarry Smith } 2273