13423f386SBarry Smith 23423f386SBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/ 33423f386SBarry Smith 43423f386SBarry Smith typedef struct { 53423f386SBarry Smith PetscScalar diag; 63423f386SBarry Smith } Mat_ConstantDiagonal; 73423f386SBarry Smith 8*a8d4b458SStefano Zampini static PetscErrorCode MatGetRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[]) 9*a8d4b458SStefano Zampini { 10*a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)A->data; 11*a8d4b458SStefano Zampini PetscErrorCode ierr; 123423f386SBarry Smith 13*a8d4b458SStefano Zampini PetscFunctionBegin; 14*a8d4b458SStefano Zampini if (ncols) *ncols = 1; 15*a8d4b458SStefano Zampini if (cols) { 16*a8d4b458SStefano Zampini ierr = PetscMalloc1(1,cols);CHKERRQ(ierr); 17*a8d4b458SStefano Zampini (*cols)[0] = row; 18*a8d4b458SStefano Zampini } 19*a8d4b458SStefano Zampini if (vals) { 20*a8d4b458SStefano Zampini ierr = PetscMalloc1(1,vals);CHKERRQ(ierr); 21*a8d4b458SStefano Zampini (*vals)[0] = ctx->diag; 22*a8d4b458SStefano Zampini } 23*a8d4b458SStefano Zampini PetscFunctionReturn(0); 24*a8d4b458SStefano Zampini } 25*a8d4b458SStefano Zampini 26*a8d4b458SStefano Zampini static PetscErrorCode MatRestoreRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[]) 27*a8d4b458SStefano Zampini { 28*a8d4b458SStefano Zampini PetscErrorCode ierr; 29*a8d4b458SStefano Zampini 30*a8d4b458SStefano Zampini PetscFunctionBegin; 31*a8d4b458SStefano Zampini if (ncols) *ncols = 0; 32*a8d4b458SStefano Zampini if (cols) { 33*a8d4b458SStefano Zampini ierr = PetscFree(*cols);CHKERRQ(ierr); 34*a8d4b458SStefano Zampini } 35*a8d4b458SStefano Zampini if (vals) { 36*a8d4b458SStefano Zampini ierr = PetscFree(*vals);CHKERRQ(ierr); 37*a8d4b458SStefano Zampini } 38*a8d4b458SStefano Zampini PetscFunctionReturn(0); 39*a8d4b458SStefano Zampini } 40*a8d4b458SStefano Zampini 41*a8d4b458SStefano Zampini static PetscErrorCode MatMultTranspose_ConstantDiagonal(Mat A, Vec x, Vec y) 42*a8d4b458SStefano Zampini { 43*a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)A->data; 44*a8d4b458SStefano Zampini PetscErrorCode ierr; 45*a8d4b458SStefano Zampini 46*a8d4b458SStefano Zampini PetscFunctionBegin; 47*a8d4b458SStefano Zampini ierr = VecAXPBY(y,ctx->diag,0.0,x);CHKERRQ(ierr); 48*a8d4b458SStefano Zampini PetscFunctionReturn(0); 49*a8d4b458SStefano Zampini } 50*a8d4b458SStefano Zampini 51*a8d4b458SStefano Zampini static PetscErrorCode MatMultAdd_ConstantDiagonal(Mat mat,Vec v1,Vec v2,Vec v3) 52*a8d4b458SStefano Zampini { 53*a8d4b458SStefano Zampini PetscErrorCode ierr; 54*a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)mat->data; 55*a8d4b458SStefano Zampini 56*a8d4b458SStefano Zampini PetscFunctionBegin; 57*a8d4b458SStefano Zampini if (v2 == v3) { 58*a8d4b458SStefano Zampini ierr = VecAXPBY(v3,ctx->diag,1.0,v1);CHKERRQ(ierr); 59*a8d4b458SStefano Zampini } else { 60*a8d4b458SStefano Zampini ierr = VecAXPBYPCZ(v3,ctx->diag,1.0,0.0,v1,v2);CHKERRQ(ierr); 61*a8d4b458SStefano Zampini } 62*a8d4b458SStefano Zampini PetscFunctionReturn(0); 63*a8d4b458SStefano Zampini } 64*a8d4b458SStefano Zampini 65*a8d4b458SStefano Zampini static PetscErrorCode MatMultTransposeAdd_ConstantDiagonal(Mat mat,Vec v1,Vec v2,Vec v3) 66*a8d4b458SStefano Zampini { 67*a8d4b458SStefano Zampini PetscErrorCode ierr; 68*a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)mat->data; 69*a8d4b458SStefano Zampini 70*a8d4b458SStefano Zampini PetscFunctionBegin; 71*a8d4b458SStefano Zampini if (v2 == v3) { 72*a8d4b458SStefano Zampini ierr = VecAXPBY(v3,ctx->diag,1.0,v1);CHKERRQ(ierr); 73*a8d4b458SStefano Zampini } else { 74*a8d4b458SStefano Zampini ierr = VecAXPBYPCZ(v3,ctx->diag,1.0,0.0,v1,v2);CHKERRQ(ierr); 75*a8d4b458SStefano Zampini } 76*a8d4b458SStefano Zampini PetscFunctionReturn(0); 77*a8d4b458SStefano Zampini } 78*a8d4b458SStefano Zampini 79*a8d4b458SStefano Zampini static PetscErrorCode MatDuplicate_ConstantDiagonal(Mat A, MatDuplicateOption op, Mat *B) 80*a8d4b458SStefano Zampini { 81*a8d4b458SStefano Zampini PetscErrorCode ierr; 82*a8d4b458SStefano Zampini Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal*)A->data; 83*a8d4b458SStefano Zampini 84*a8d4b458SStefano Zampini PetscFunctionBegin; 85*a8d4b458SStefano Zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),B);CHKERRQ(ierr); 86*a8d4b458SStefano Zampini ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 87*a8d4b458SStefano Zampini ierr = MatSetBlockSizesFromMats(*B,A,A);CHKERRQ(ierr); 88*a8d4b458SStefano Zampini ierr = MatSetType(*B,MATCONSTANTDIAGONAL);CHKERRQ(ierr); 89*a8d4b458SStefano Zampini ierr = PetscLayoutReference(A->rmap,&(*B)->rmap);CHKERRQ(ierr); 90*a8d4b458SStefano Zampini ierr = PetscLayoutReference(A->cmap,&(*B)->cmap);CHKERRQ(ierr); 91*a8d4b458SStefano Zampini if (op == MAT_COPY_VALUES) { 92*a8d4b458SStefano Zampini Mat_ConstantDiagonal *bctx = (Mat_ConstantDiagonal*)(*B)->data; 93*a8d4b458SStefano Zampini bctx->diag = actx->diag; 94*a8d4b458SStefano Zampini } 95*a8d4b458SStefano Zampini PetscFunctionReturn(0); 96*a8d4b458SStefano Zampini } 97*a8d4b458SStefano Zampini 98*a8d4b458SStefano Zampini static PetscErrorCode MatMissingDiagonal_ConstantDiagonal(Mat mat,PetscBool *missing,PetscInt *dd) 99*a8d4b458SStefano Zampini { 100*a8d4b458SStefano Zampini PetscFunctionBegin; 101*a8d4b458SStefano Zampini *missing = PETSC_FALSE; 102*a8d4b458SStefano Zampini PetscFunctionReturn(0); 103*a8d4b458SStefano Zampini } 104*a8d4b458SStefano Zampini 1053423f386SBarry Smith static PetscErrorCode MatDestroy_ConstantDiagonal(Mat mat) 1063423f386SBarry Smith { 1073423f386SBarry Smith PetscErrorCode ierr; 1083423f386SBarry Smith 1093423f386SBarry Smith PetscFunctionBegin; 1103423f386SBarry Smith ierr = PetscFree(mat->data);CHKERRQ(ierr); 1113423f386SBarry Smith PetscFunctionReturn(0); 1123423f386SBarry Smith } 1133423f386SBarry Smith 1143423f386SBarry Smith static PetscErrorCode MatView_ConstantDiagonal(Mat J,PetscViewer viewer) 1153423f386SBarry Smith { 1163423f386SBarry Smith PetscErrorCode ierr; 1173423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data; 1183423f386SBarry Smith PetscBool iascii; 1193423f386SBarry Smith 1203423f386SBarry Smith PetscFunctionBegin; 1213423f386SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1223423f386SBarry Smith if (iascii) { 1233423f386SBarry Smith PetscViewerFormat format; 1243423f386SBarry Smith 1253423f386SBarry Smith ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 1263423f386SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(0); 1273423f386SBarry Smith #if !defined(PETSC_USE_COMPLEX) 1283423f386SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Diagonal value: %g\n",ctx->diag);CHKERRQ(ierr); 1293423f386SBarry Smith #else 1303423f386SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Diagonal value: %g + i %g\n",PetscRealPart(ctx->diag),PetscImaginaryPart(ctx->diag));CHKERRQ(ierr); 1313423f386SBarry Smith #endif 1323423f386SBarry Smith } 1333423f386SBarry Smith PetscFunctionReturn(0); 1343423f386SBarry Smith } 1353423f386SBarry Smith 1363423f386SBarry Smith static PetscErrorCode MatAssemblyEnd_ConstantDiagonal(Mat J,MatAssemblyType mt) 1373423f386SBarry Smith { 1383423f386SBarry Smith PetscFunctionBegin; 1393423f386SBarry Smith PetscFunctionReturn(0); 1403423f386SBarry Smith } 1413423f386SBarry Smith 1423423f386SBarry Smith static PetscErrorCode MatMult_ConstantDiagonal(Mat J,Vec x,Vec y) 1433423f386SBarry Smith { 1443423f386SBarry Smith PetscErrorCode ierr; 1453423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data; 1463423f386SBarry Smith 1473423f386SBarry Smith PetscFunctionBegin; 1483954eb43SBarry Smith ierr = VecAXPBY(y,ctx->diag,0.0,x);CHKERRQ(ierr); 1493423f386SBarry Smith PetscFunctionReturn(0); 1503423f386SBarry Smith } 1513423f386SBarry Smith 1523423f386SBarry Smith PetscErrorCode MatGetDiagonal_ConstantDiagonal(Mat J,Vec x) 1533423f386SBarry Smith { 1543423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data; 1553423f386SBarry Smith PetscErrorCode ierr; 1563423f386SBarry Smith 1573423f386SBarry Smith PetscFunctionBegin; 1583423f386SBarry Smith ierr = VecSet(x,ctx->diag);CHKERRQ(ierr); 1593423f386SBarry Smith PetscFunctionReturn(0); 1603423f386SBarry Smith } 1613423f386SBarry Smith 1623423f386SBarry Smith static PetscErrorCode MatShift_ConstantDiagonal(Mat Y,PetscScalar a) 1633423f386SBarry Smith { 1643954eb43SBarry Smith PetscErrorCode ierr; 1653423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)Y->data; 1663423f386SBarry Smith 1673423f386SBarry Smith PetscFunctionBegin; 1683954eb43SBarry Smith if (a != 0.) {ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);} 1693423f386SBarry Smith ctx->diag += a; 1703423f386SBarry Smith PetscFunctionReturn(0); 1713423f386SBarry Smith } 1723423f386SBarry Smith 1733423f386SBarry Smith static PetscErrorCode MatScale_ConstantDiagonal(Mat Y,PetscScalar a) 1743423f386SBarry Smith { 1754d43f095SSatish Balay PetscErrorCode ierr; 1763423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)Y->data; 1773423f386SBarry Smith 1783423f386SBarry Smith PetscFunctionBegin; 1794d43f095SSatish Balay if (a != 1.) {ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);} 1803423f386SBarry Smith ctx->diag *= a; 1813423f386SBarry Smith PetscFunctionReturn(0); 1823423f386SBarry Smith } 1833423f386SBarry Smith 1843423f386SBarry Smith static PetscErrorCode MatZeroEntries_ConstantDiagonal(Mat Y) 1853423f386SBarry Smith { 1864d43f095SSatish Balay PetscErrorCode ierr; 1873423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)Y->data; 1883423f386SBarry Smith 1893423f386SBarry Smith PetscFunctionBegin; 1904d43f095SSatish Balay if (ctx->diag != 0.0) {ierr = PetscObjectStateIncrease((PetscObject)Y);CHKERRQ(ierr);} 1913423f386SBarry Smith ctx->diag = 0.0; 1923423f386SBarry Smith PetscFunctionReturn(0); 1933423f386SBarry Smith } 1943423f386SBarry Smith 1953423f386SBarry Smith PetscErrorCode MatSOR_ConstantDiagonal(Mat matin,Vec x,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec y) 1963423f386SBarry Smith { 1973423f386SBarry Smith PetscErrorCode ierr; 1983423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)matin->data; 1993423f386SBarry Smith 2003423f386SBarry Smith PetscFunctionBegin; 2013423f386SBarry Smith if (ctx->diag == 0.0) matin->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 2023423f386SBarry Smith else matin->factorerrortype = MAT_FACTOR_NOERROR; 2033954eb43SBarry Smith ierr = VecAXPBY(y,1.0/ctx->diag,0.0,x);CHKERRQ(ierr); 2043423f386SBarry Smith PetscFunctionReturn(0); 2053423f386SBarry Smith } 2063423f386SBarry Smith 2073423f386SBarry Smith PetscErrorCode MatGetInfo_ConstantDiagonal(Mat A,MatInfoType flag,MatInfo *info) 2083423f386SBarry Smith { 2093423f386SBarry Smith PetscFunctionBegin; 2103423f386SBarry Smith info->block_size = 1.0; 2113423f386SBarry Smith info->nz_allocated = 1.0; 2123423f386SBarry Smith info->nz_used = 1.0; 2133423f386SBarry Smith info->nz_unneeded = 0.0; 2143966268fSBarry Smith info->assemblies = A->num_ass; 2153423f386SBarry Smith info->mallocs = 0.0; 2163423f386SBarry Smith info->memory = ((PetscObject)A)->mem; 2173423f386SBarry Smith if (A->factortype) { 2183423f386SBarry Smith info->fill_ratio_given = 1.0; 2193423f386SBarry Smith info->fill_ratio_needed = 1.0; 2203423f386SBarry Smith info->factor_mallocs = 0.0; 2213423f386SBarry Smith } else { 2223423f386SBarry Smith info->fill_ratio_given = 0; 2233423f386SBarry Smith info->fill_ratio_needed = 0; 2243423f386SBarry Smith info->factor_mallocs = 0; 2253423f386SBarry Smith } 2263423f386SBarry Smith PetscFunctionReturn(0); 2273423f386SBarry Smith } 2283423f386SBarry Smith 2293423f386SBarry Smith /*@ 2303423f386SBarry Smith MatCreateConstantDiagonal - Creates a matrix with a uniform value along the diagonal 2313423f386SBarry Smith 232d083f849SBarry Smith Collective 2333423f386SBarry Smith 2343423f386SBarry Smith Input Parameters: 2353423f386SBarry Smith + comm - MPI communicator 2363423f386SBarry Smith . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 2373423f386SBarry Smith This value should be the same as the local size used in creating the 2383423f386SBarry Smith y vector for the matrix-vector product y = Ax. 2393423f386SBarry Smith . n - This value should be the same as the local size used in creating the 2403423f386SBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 2413423f386SBarry Smith calculated if N is given) For square matrices n is almost always m. 2423423f386SBarry Smith . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 2433423f386SBarry Smith . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 2443423f386SBarry Smith - diag - the diagonal value 2453423f386SBarry Smith 2463423f386SBarry Smith Output Parameter: 2473423f386SBarry Smith . J - the diagonal matrix 2483423f386SBarry Smith 2493423f386SBarry Smith Level: advanced 2503423f386SBarry Smith 2513423f386SBarry Smith Notes: 2523423f386SBarry Smith Only supports square matrices with the same number of local rows and columns 2533423f386SBarry Smith 2543423f386SBarry Smith .seealso: MatDestroy(), MATCONSTANTDIAGONAL, MatScale(), MatShift(), MatMult(), MatGetDiagonal(), MatGetFactor(), MatSolve() 2553423f386SBarry Smith 2563423f386SBarry Smith @*/ 2573423f386SBarry Smith PetscErrorCode MatCreateConstantDiagonal(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar diag,Mat *J) 2583423f386SBarry Smith { 2593423f386SBarry Smith PetscErrorCode ierr; 2603423f386SBarry Smith 2613423f386SBarry Smith PetscFunctionBegin; 2623423f386SBarry Smith ierr = MatCreate(comm,J);CHKERRQ(ierr); 2633423f386SBarry Smith ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr); 2643423f386SBarry Smith ierr = MatSetType(*J,MATCONSTANTDIAGONAL);CHKERRQ(ierr); 2653423f386SBarry Smith ierr = MatShift(*J,diag);CHKERRQ(ierr); 2663423f386SBarry Smith ierr = MatSetUp(*J);CHKERRQ(ierr); 2673423f386SBarry Smith PetscFunctionReturn(0); 2683423f386SBarry Smith } 2693423f386SBarry Smith 2703423f386SBarry Smith PETSC_EXTERN PetscErrorCode MatCreate_ConstantDiagonal(Mat A) 2713423f386SBarry Smith { 2723423f386SBarry Smith PetscErrorCode ierr; 2733423f386SBarry Smith Mat_ConstantDiagonal *ctx; 2743423f386SBarry Smith 2753423f386SBarry Smith PetscFunctionBegin; 2763423f386SBarry Smith ierr = PetscNew(&ctx);CHKERRQ(ierr); 2773423f386SBarry Smith ctx->diag = 0.0; 2783423f386SBarry Smith A->data = (void*)ctx; 2793423f386SBarry Smith 2803423f386SBarry Smith A->assembled = PETSC_TRUE; 2813423f386SBarry Smith A->preallocated = PETSC_TRUE; 282*a8d4b458SStefano Zampini 2833423f386SBarry Smith A->ops->mult = MatMult_ConstantDiagonal; 284*a8d4b458SStefano Zampini A->ops->multadd = MatMultAdd_ConstantDiagonal; 285*a8d4b458SStefano Zampini A->ops->multtranspose = MatMultTranspose_ConstantDiagonal; 286*a8d4b458SStefano Zampini A->ops->multtransposeadd = MatMultTransposeAdd_ConstantDiagonal; 287*a8d4b458SStefano Zampini A->ops->duplicate = MatDuplicate_ConstantDiagonal; 288*a8d4b458SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_ConstantDiagonal; 289*a8d4b458SStefano Zampini A->ops->getrow = MatGetRow_ConstantDiagonal; 290*a8d4b458SStefano Zampini A->ops->restorerow = MatRestoreRow_ConstantDiagonal; 2913423f386SBarry Smith A->ops->sor = MatSOR_ConstantDiagonal; 2923423f386SBarry Smith A->ops->shift = MatShift_ConstantDiagonal; 2933423f386SBarry Smith A->ops->scale = MatScale_ConstantDiagonal; 2943423f386SBarry Smith A->ops->getdiagonal = MatGetDiagonal_ConstantDiagonal; 2953423f386SBarry Smith A->ops->view = MatView_ConstantDiagonal; 2963423f386SBarry Smith A->ops->zeroentries = MatZeroEntries_ConstantDiagonal; 2973423f386SBarry Smith A->ops->assemblyend = MatAssemblyEnd_ConstantDiagonal; 2983423f386SBarry Smith A->ops->destroy = MatDestroy_ConstantDiagonal; 2993423f386SBarry Smith A->ops->getinfo = MatGetInfo_ConstantDiagonal; 300*a8d4b458SStefano Zampini 3013423f386SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATCONSTANTDIAGONAL);CHKERRQ(ierr); 3023423f386SBarry Smith PetscFunctionReturn(0); 3033423f386SBarry Smith } 3043423f386SBarry Smith 3053423f386SBarry Smith static PetscErrorCode MatFactorNumeric_ConstantDiagonal(Mat fact,Mat A,const MatFactorInfo *info) 3063423f386SBarry Smith { 3074d43f095SSatish Balay PetscErrorCode ierr; 3083423f386SBarry Smith Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal*)A->data,*fctx = (Mat_ConstantDiagonal*)fact->data; 3093423f386SBarry Smith 3103423f386SBarry Smith PetscFunctionBegin; 3113423f386SBarry Smith if (actx->diag == 0.0) fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3123423f386SBarry Smith else fact->factorerrortype = MAT_FACTOR_NOERROR; 3133423f386SBarry Smith fctx->diag = 1.0/actx->diag; 3144d43f095SSatish Balay ierr = PetscObjectStateIncrease((PetscObject)fact);CHKERRQ(ierr); 3153423f386SBarry Smith fact->ops->solve = MatMult_ConstantDiagonal; 3163423f386SBarry Smith PetscFunctionReturn(0); 3173423f386SBarry Smith } 3183423f386SBarry Smith 3193423f386SBarry Smith static PetscErrorCode MatFactorSymbolic_LU_ConstantDiagonal(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 3203423f386SBarry Smith { 3213423f386SBarry Smith PetscFunctionBegin; 3223423f386SBarry Smith fact->ops->lufactornumeric = MatFactorNumeric_ConstantDiagonal; 3233423f386SBarry Smith PetscFunctionReturn(0); 3243423f386SBarry Smith } 3253423f386SBarry Smith 3263423f386SBarry Smith static PetscErrorCode MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact,Mat A,IS isrow,const MatFactorInfo *info) 3273423f386SBarry Smith { 3283423f386SBarry Smith PetscFunctionBegin; 3293423f386SBarry Smith fact->ops->choleskyfactornumeric = MatFactorNumeric_ConstantDiagonal; PetscFunctionReturn(0); 3303423f386SBarry Smith PetscFunctionReturn(0); 3313423f386SBarry Smith } 3323423f386SBarry Smith 3333423f386SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_constantdiagonal_petsc(Mat A,MatFactorType ftype,Mat *B) 3343423f386SBarry Smith { 3353423f386SBarry Smith PetscInt n = A->rmap->n, N = A->rmap->N; 3363423f386SBarry Smith PetscErrorCode ierr; 3373423f386SBarry Smith 3383423f386SBarry Smith PetscFunctionBegin; 3393423f386SBarry Smith ierr = MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A),n,n,N,N,0,B);CHKERRQ(ierr); 3403423f386SBarry Smith 3413423f386SBarry Smith (*B)->factortype = ftype; 3423423f386SBarry Smith (*B)->ops->ilufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal; 3433423f386SBarry Smith (*B)->ops->lufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal; 3443423f386SBarry Smith (*B)->ops->iccfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal; 3453423f386SBarry Smith (*B)->ops->choleskyfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal; 3463423f386SBarry Smith 3473423f386SBarry Smith (*B)->ops->shift = NULL; 3483423f386SBarry Smith (*B)->ops->scale = NULL; 3493423f386SBarry Smith (*B)->ops->mult = NULL; 3503423f386SBarry Smith (*B)->ops->sor = NULL; 3513423f386SBarry Smith (*B)->ops->zeroentries = NULL; 3523423f386SBarry Smith 3533423f386SBarry Smith ierr = PetscFree((*B)->solvertype);CHKERRQ(ierr); 3543423f386SBarry Smith ierr = PetscStrallocpy(MATSOLVERPETSC,&(*B)->solvertype);CHKERRQ(ierr); 3553423f386SBarry Smith PetscFunctionReturn(0); 3563423f386SBarry Smith } 357