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 8a01b8767SStefano Zampini static PetscErrorCode MatAXPY_ConstantDiagonal(Mat Y, PetscScalar a, Mat X, MatStructure str) 9a01b8767SStefano Zampini { 10a01b8767SStefano Zampini Mat_ConstantDiagonal *yctx = (Mat_ConstantDiagonal*)Y->data; 11a01b8767SStefano Zampini Mat_ConstantDiagonal *xctx = (Mat_ConstantDiagonal*)X->data; 12a01b8767SStefano Zampini 13a01b8767SStefano Zampini PetscFunctionBegin; 14a01b8767SStefano Zampini yctx->diag += a*xctx->diag; 15a01b8767SStefano Zampini PetscFunctionReturn(0); 16a01b8767SStefano Zampini } 17a01b8767SStefano Zampini 18a8d4b458SStefano Zampini static PetscErrorCode MatGetRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[]) 19a8d4b458SStefano Zampini { 20a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)A->data; 213423f386SBarry Smith 22a8d4b458SStefano Zampini PetscFunctionBegin; 23a8d4b458SStefano Zampini if (ncols) *ncols = 1; 24a8d4b458SStefano Zampini if (cols) { 25*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1,cols)); 26a8d4b458SStefano Zampini (*cols)[0] = row; 27a8d4b458SStefano Zampini } 28a8d4b458SStefano Zampini if (vals) { 29*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1,vals)); 30a8d4b458SStefano Zampini (*vals)[0] = ctx->diag; 31a8d4b458SStefano Zampini } 32a8d4b458SStefano Zampini PetscFunctionReturn(0); 33a8d4b458SStefano Zampini } 34a8d4b458SStefano Zampini 35a8d4b458SStefano Zampini static PetscErrorCode MatRestoreRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[]) 36a8d4b458SStefano Zampini { 37a8d4b458SStefano Zampini PetscFunctionBegin; 38a8d4b458SStefano Zampini if (ncols) *ncols = 0; 39a8d4b458SStefano Zampini if (cols) { 40*9566063dSJacob Faibussowitsch PetscCall(PetscFree(*cols)); 41a8d4b458SStefano Zampini } 42a8d4b458SStefano Zampini if (vals) { 43*9566063dSJacob Faibussowitsch PetscCall(PetscFree(*vals)); 44a8d4b458SStefano Zampini } 45a8d4b458SStefano Zampini PetscFunctionReturn(0); 46a8d4b458SStefano Zampini } 47a8d4b458SStefano Zampini 48a8d4b458SStefano Zampini static PetscErrorCode MatMultTranspose_ConstantDiagonal(Mat A, Vec x, Vec y) 49a8d4b458SStefano Zampini { 50a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)A->data; 51a8d4b458SStefano Zampini 52a8d4b458SStefano Zampini PetscFunctionBegin; 53*9566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y,ctx->diag,0.0,x)); 54a8d4b458SStefano Zampini PetscFunctionReturn(0); 55a8d4b458SStefano Zampini } 56a8d4b458SStefano Zampini 57a8d4b458SStefano Zampini static PetscErrorCode MatMultAdd_ConstantDiagonal(Mat mat,Vec v1,Vec v2,Vec v3) 58a8d4b458SStefano Zampini { 59a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)mat->data; 60a8d4b458SStefano Zampini 61a8d4b458SStefano Zampini PetscFunctionBegin; 62a8d4b458SStefano Zampini if (v2 == v3) { 63*9566063dSJacob Faibussowitsch PetscCall(VecAXPBY(v3,ctx->diag,1.0,v1)); 64a8d4b458SStefano Zampini } else { 65*9566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(v3,ctx->diag,1.0,0.0,v1,v2)); 66a8d4b458SStefano Zampini } 67a8d4b458SStefano Zampini PetscFunctionReturn(0); 68a8d4b458SStefano Zampini } 69a8d4b458SStefano Zampini 70a8d4b458SStefano Zampini static PetscErrorCode MatMultTransposeAdd_ConstantDiagonal(Mat mat,Vec v1,Vec v2,Vec v3) 71a8d4b458SStefano Zampini { 72a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)mat->data; 73a8d4b458SStefano Zampini 74a8d4b458SStefano Zampini PetscFunctionBegin; 75a8d4b458SStefano Zampini if (v2 == v3) { 76*9566063dSJacob Faibussowitsch PetscCall(VecAXPBY(v3,ctx->diag,1.0,v1)); 77a8d4b458SStefano Zampini } else { 78*9566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(v3,ctx->diag,1.0,0.0,v1,v2)); 79a8d4b458SStefano Zampini } 80a8d4b458SStefano Zampini PetscFunctionReturn(0); 81a8d4b458SStefano Zampini } 82a8d4b458SStefano Zampini 83e5f3c4beSJose E. Roman static PetscErrorCode MatNorm_ConstantDiagonal(Mat A,NormType type,PetscReal *nrm) 84e5f3c4beSJose E. Roman { 85e5f3c4beSJose E. Roman Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)A->data; 86e5f3c4beSJose E. Roman 87e5f3c4beSJose E. Roman PetscFunctionBegin; 88e5f3c4beSJose E. Roman if (type == NORM_FROBENIUS || type == NORM_2 || type == NORM_1 || type == NORM_INFINITY) *nrm = PetscAbsScalar(ctx->diag); 89e5f3c4beSJose E. Roman else SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Unsupported norm"); 90e5f3c4beSJose E. Roman PetscFunctionReturn(0); 91e5f3c4beSJose E. Roman } 92e5f3c4beSJose E. Roman 93d1ca2681SJose E. Roman static PetscErrorCode MatCreateSubMatrices_ConstantDiagonal(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 94d1ca2681SJose E. Roman 95d1ca2681SJose E. Roman { 96d1ca2681SJose E. Roman Mat B; 97d1ca2681SJose E. Roman 98d1ca2681SJose E. Roman PetscFunctionBegin; 99*9566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B)); 100*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B,n,irow,icol,scall,submat)); 101*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 102d1ca2681SJose E. Roman PetscFunctionReturn(0); 103d1ca2681SJose E. Roman } 104d1ca2681SJose E. Roman 105a8d4b458SStefano Zampini static PetscErrorCode MatDuplicate_ConstantDiagonal(Mat A, MatDuplicateOption op, Mat *B) 106a8d4b458SStefano Zampini { 107a8d4b458SStefano Zampini Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal*)A->data; 108a8d4b458SStefano Zampini 109a8d4b458SStefano Zampini PetscFunctionBegin; 110*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A),B)); 111*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N)); 112*9566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizesFromMats(*B,A,A)); 113*9566063dSJacob Faibussowitsch PetscCall(MatSetType(*B,MATCONSTANTDIAGONAL)); 114*9566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->rmap,&(*B)->rmap)); 115*9566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap,&(*B)->cmap)); 116a8d4b458SStefano Zampini if (op == MAT_COPY_VALUES) { 117a8d4b458SStefano Zampini Mat_ConstantDiagonal *bctx = (Mat_ConstantDiagonal*)(*B)->data; 118a8d4b458SStefano Zampini bctx->diag = actx->diag; 119a8d4b458SStefano Zampini } 120a8d4b458SStefano Zampini PetscFunctionReturn(0); 121a8d4b458SStefano Zampini } 122a8d4b458SStefano Zampini 123a8d4b458SStefano Zampini static PetscErrorCode MatMissingDiagonal_ConstantDiagonal(Mat mat,PetscBool *missing,PetscInt *dd) 124a8d4b458SStefano Zampini { 125a8d4b458SStefano Zampini PetscFunctionBegin; 126a8d4b458SStefano Zampini *missing = PETSC_FALSE; 127a8d4b458SStefano Zampini PetscFunctionReturn(0); 128a8d4b458SStefano Zampini } 129a8d4b458SStefano Zampini 1303423f386SBarry Smith static PetscErrorCode MatDestroy_ConstantDiagonal(Mat mat) 1313423f386SBarry Smith { 1323423f386SBarry Smith PetscFunctionBegin; 133*9566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 1343423f386SBarry Smith PetscFunctionReturn(0); 1353423f386SBarry Smith } 1363423f386SBarry Smith 1373423f386SBarry Smith static PetscErrorCode MatView_ConstantDiagonal(Mat J,PetscViewer viewer) 1383423f386SBarry Smith { 1393423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data; 1403423f386SBarry Smith PetscBool iascii; 1413423f386SBarry Smith 1423423f386SBarry Smith PetscFunctionBegin; 143*9566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 1443423f386SBarry Smith if (iascii) { 1453423f386SBarry Smith PetscViewerFormat format; 1463423f386SBarry Smith 147*9566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 1483423f386SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(0); 149c0aa6a63SJacob Faibussowitsch #if defined(PETSC_USE_COMPLEX) 150*9566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Diagonal value: %g + i %g\n",(double)PetscRealPart(ctx->diag),(double)PetscImaginaryPart(ctx->diag))); 1513423f386SBarry Smith #else 152*9566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Diagonal value: %g\n",(double)(ctx->diag))); 1533423f386SBarry Smith #endif 1543423f386SBarry Smith } 1553423f386SBarry Smith PetscFunctionReturn(0); 1563423f386SBarry Smith } 1573423f386SBarry Smith 1583423f386SBarry Smith static PetscErrorCode MatAssemblyEnd_ConstantDiagonal(Mat J,MatAssemblyType mt) 1593423f386SBarry Smith { 1603423f386SBarry Smith PetscFunctionBegin; 1613423f386SBarry Smith PetscFunctionReturn(0); 1623423f386SBarry Smith } 1633423f386SBarry Smith 1643423f386SBarry Smith static PetscErrorCode MatMult_ConstantDiagonal(Mat J,Vec x,Vec y) 1653423f386SBarry Smith { 1663423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data; 1673423f386SBarry Smith 1683423f386SBarry Smith PetscFunctionBegin; 169*9566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y,ctx->diag,0.0,x)); 1703423f386SBarry Smith PetscFunctionReturn(0); 1713423f386SBarry Smith } 1723423f386SBarry Smith 1733423f386SBarry Smith PetscErrorCode MatGetDiagonal_ConstantDiagonal(Mat J,Vec x) 1743423f386SBarry Smith { 1753423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)J->data; 1763423f386SBarry Smith 1773423f386SBarry Smith PetscFunctionBegin; 178*9566063dSJacob Faibussowitsch PetscCall(VecSet(x,ctx->diag)); 1793423f386SBarry Smith PetscFunctionReturn(0); 1803423f386SBarry Smith } 1813423f386SBarry Smith 1823423f386SBarry Smith static PetscErrorCode MatShift_ConstantDiagonal(Mat Y,PetscScalar a) 1833423f386SBarry Smith { 1843423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)Y->data; 1853423f386SBarry Smith 1863423f386SBarry Smith PetscFunctionBegin; 1873423f386SBarry Smith ctx->diag += a; 1883423f386SBarry Smith PetscFunctionReturn(0); 1893423f386SBarry Smith } 1903423f386SBarry Smith 1913423f386SBarry Smith static PetscErrorCode MatScale_ConstantDiagonal(Mat Y,PetscScalar a) 1923423f386SBarry Smith { 1933423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)Y->data; 1943423f386SBarry Smith 1953423f386SBarry Smith PetscFunctionBegin; 1963423f386SBarry Smith ctx->diag *= a; 1973423f386SBarry Smith PetscFunctionReturn(0); 1983423f386SBarry Smith } 1993423f386SBarry Smith 2003423f386SBarry Smith static PetscErrorCode MatZeroEntries_ConstantDiagonal(Mat Y) 2013423f386SBarry Smith { 2023423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)Y->data; 2033423f386SBarry Smith 2043423f386SBarry Smith PetscFunctionBegin; 2053423f386SBarry Smith ctx->diag = 0.0; 2063423f386SBarry Smith PetscFunctionReturn(0); 2073423f386SBarry Smith } 2083423f386SBarry Smith 2093423f386SBarry Smith PetscErrorCode MatSOR_ConstantDiagonal(Mat matin,Vec x,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec y) 2103423f386SBarry Smith { 2113423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal*)matin->data; 2123423f386SBarry Smith 2133423f386SBarry Smith PetscFunctionBegin; 2143423f386SBarry Smith if (ctx->diag == 0.0) matin->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 2153423f386SBarry Smith else matin->factorerrortype = MAT_FACTOR_NOERROR; 216*9566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y,1.0/ctx->diag,0.0,x)); 2173423f386SBarry Smith PetscFunctionReturn(0); 2183423f386SBarry Smith } 2193423f386SBarry Smith 2203423f386SBarry Smith PetscErrorCode MatGetInfo_ConstantDiagonal(Mat A,MatInfoType flag,MatInfo *info) 2213423f386SBarry Smith { 2223423f386SBarry Smith PetscFunctionBegin; 2233423f386SBarry Smith info->block_size = 1.0; 2243423f386SBarry Smith info->nz_allocated = 1.0; 2253423f386SBarry Smith info->nz_used = 1.0; 2263423f386SBarry Smith info->nz_unneeded = 0.0; 2273966268fSBarry Smith info->assemblies = A->num_ass; 2283423f386SBarry Smith info->mallocs = 0.0; 2293423f386SBarry Smith info->memory = ((PetscObject)A)->mem; 2303423f386SBarry Smith if (A->factortype) { 2313423f386SBarry Smith info->fill_ratio_given = 1.0; 2323423f386SBarry Smith info->fill_ratio_needed = 1.0; 2333423f386SBarry Smith info->factor_mallocs = 0.0; 2343423f386SBarry Smith } else { 2353423f386SBarry Smith info->fill_ratio_given = 0; 2363423f386SBarry Smith info->fill_ratio_needed = 0; 2373423f386SBarry Smith info->factor_mallocs = 0; 2383423f386SBarry Smith } 2393423f386SBarry Smith PetscFunctionReturn(0); 2403423f386SBarry Smith } 2413423f386SBarry Smith 2423423f386SBarry Smith /*@ 2433423f386SBarry Smith MatCreateConstantDiagonal - Creates a matrix with a uniform value along the diagonal 2443423f386SBarry Smith 245d083f849SBarry Smith Collective 2463423f386SBarry Smith 2473423f386SBarry Smith Input Parameters: 2483423f386SBarry Smith + comm - MPI communicator 2493423f386SBarry Smith . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 2503423f386SBarry Smith This value should be the same as the local size used in creating the 2513423f386SBarry Smith y vector for the matrix-vector product y = Ax. 2523423f386SBarry Smith . n - This value should be the same as the local size used in creating the 2533423f386SBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 2543423f386SBarry Smith calculated if N is given) For square matrices n is almost always m. 2553423f386SBarry Smith . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 2563423f386SBarry Smith . N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 2573423f386SBarry Smith - diag - the diagonal value 2583423f386SBarry Smith 2593423f386SBarry Smith Output Parameter: 2603423f386SBarry Smith . J - the diagonal matrix 2613423f386SBarry Smith 2623423f386SBarry Smith Level: advanced 2633423f386SBarry Smith 2643423f386SBarry Smith Notes: 2653423f386SBarry Smith Only supports square matrices with the same number of local rows and columns 2663423f386SBarry Smith 2673423f386SBarry Smith .seealso: MatDestroy(), MATCONSTANTDIAGONAL, MatScale(), MatShift(), MatMult(), MatGetDiagonal(), MatGetFactor(), MatSolve() 2683423f386SBarry Smith 2693423f386SBarry Smith @*/ 2703423f386SBarry Smith PetscErrorCode MatCreateConstantDiagonal(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar diag,Mat *J) 2713423f386SBarry Smith { 2723423f386SBarry Smith PetscFunctionBegin; 273*9566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,J)); 274*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J,m,n,M,N)); 275*9566063dSJacob Faibussowitsch PetscCall(MatSetType(*J,MATCONSTANTDIAGONAL)); 276*9566063dSJacob Faibussowitsch PetscCall(MatShift(*J,diag)); 277*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(*J)); 2783423f386SBarry Smith PetscFunctionReturn(0); 2793423f386SBarry Smith } 2803423f386SBarry Smith 2813423f386SBarry Smith PETSC_EXTERN PetscErrorCode MatCreate_ConstantDiagonal(Mat A) 2823423f386SBarry Smith { 2833423f386SBarry Smith Mat_ConstantDiagonal *ctx; 2843423f386SBarry Smith 2853423f386SBarry Smith PetscFunctionBegin; 286*9566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 2873423f386SBarry Smith ctx->diag = 0.0; 2883423f386SBarry Smith A->data = (void*)ctx; 2893423f386SBarry Smith 2903423f386SBarry Smith A->assembled = PETSC_TRUE; 2913423f386SBarry Smith A->preallocated = PETSC_TRUE; 292a8d4b458SStefano Zampini 2933423f386SBarry Smith A->ops->mult = MatMult_ConstantDiagonal; 294a8d4b458SStefano Zampini A->ops->multadd = MatMultAdd_ConstantDiagonal; 295a8d4b458SStefano Zampini A->ops->multtranspose = MatMultTranspose_ConstantDiagonal; 296a8d4b458SStefano Zampini A->ops->multtransposeadd = MatMultTransposeAdd_ConstantDiagonal; 297e5f3c4beSJose E. Roman A->ops->norm = MatNorm_ConstantDiagonal; 298d1ca2681SJose E. Roman A->ops->createsubmatrices= MatCreateSubMatrices_ConstantDiagonal; 299a8d4b458SStefano Zampini A->ops->duplicate = MatDuplicate_ConstantDiagonal; 300a8d4b458SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_ConstantDiagonal; 301a8d4b458SStefano Zampini A->ops->getrow = MatGetRow_ConstantDiagonal; 302a8d4b458SStefano Zampini A->ops->restorerow = MatRestoreRow_ConstantDiagonal; 3033423f386SBarry Smith A->ops->sor = MatSOR_ConstantDiagonal; 3043423f386SBarry Smith A->ops->shift = MatShift_ConstantDiagonal; 3053423f386SBarry Smith A->ops->scale = MatScale_ConstantDiagonal; 3063423f386SBarry Smith A->ops->getdiagonal = MatGetDiagonal_ConstantDiagonal; 3073423f386SBarry Smith A->ops->view = MatView_ConstantDiagonal; 3083423f386SBarry Smith A->ops->zeroentries = MatZeroEntries_ConstantDiagonal; 3093423f386SBarry Smith A->ops->assemblyend = MatAssemblyEnd_ConstantDiagonal; 3103423f386SBarry Smith A->ops->destroy = MatDestroy_ConstantDiagonal; 3113423f386SBarry Smith A->ops->getinfo = MatGetInfo_ConstantDiagonal; 312a01b8767SStefano Zampini A->ops->axpy = MatAXPY_ConstantDiagonal; 313a8d4b458SStefano Zampini 314*9566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATCONSTANTDIAGONAL)); 3153423f386SBarry Smith PetscFunctionReturn(0); 3163423f386SBarry Smith } 3173423f386SBarry Smith 3183423f386SBarry Smith static PetscErrorCode MatFactorNumeric_ConstantDiagonal(Mat fact,Mat A,const MatFactorInfo *info) 3193423f386SBarry Smith { 3203423f386SBarry Smith Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal*)A->data,*fctx = (Mat_ConstantDiagonal*)fact->data; 3213423f386SBarry Smith 3223423f386SBarry Smith PetscFunctionBegin; 3233423f386SBarry Smith if (actx->diag == 0.0) fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3243423f386SBarry Smith else fact->factorerrortype = MAT_FACTOR_NOERROR; 3253423f386SBarry Smith fctx->diag = 1.0/actx->diag; 3263423f386SBarry Smith fact->ops->solve = MatMult_ConstantDiagonal; 3273423f386SBarry Smith PetscFunctionReturn(0); 3283423f386SBarry Smith } 3293423f386SBarry Smith 3303423f386SBarry Smith static PetscErrorCode MatFactorSymbolic_LU_ConstantDiagonal(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 3313423f386SBarry Smith { 3323423f386SBarry Smith PetscFunctionBegin; 3333423f386SBarry Smith fact->ops->lufactornumeric = MatFactorNumeric_ConstantDiagonal; 3343423f386SBarry Smith PetscFunctionReturn(0); 3353423f386SBarry Smith } 3363423f386SBarry Smith 3373423f386SBarry Smith static PetscErrorCode MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact,Mat A,IS isrow,const MatFactorInfo *info) 3383423f386SBarry Smith { 3393423f386SBarry Smith PetscFunctionBegin; 340a01b8767SStefano Zampini fact->ops->choleskyfactornumeric = MatFactorNumeric_ConstantDiagonal; 3413423f386SBarry Smith PetscFunctionReturn(0); 3423423f386SBarry Smith } 3433423f386SBarry Smith 3443423f386SBarry Smith PETSC_INTERN PetscErrorCode MatGetFactor_constantdiagonal_petsc(Mat A,MatFactorType ftype,Mat *B) 3453423f386SBarry Smith { 3463423f386SBarry Smith PetscInt n = A->rmap->n, N = A->rmap->N; 3473423f386SBarry Smith 3483423f386SBarry Smith PetscFunctionBegin; 349*9566063dSJacob Faibussowitsch PetscCall(MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A),n,n,N,N,0,B)); 3503423f386SBarry Smith 3513423f386SBarry Smith (*B)->factortype = ftype; 3523423f386SBarry Smith (*B)->ops->ilufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal; 3533423f386SBarry Smith (*B)->ops->lufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal; 3543423f386SBarry Smith (*B)->ops->iccfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal; 3553423f386SBarry Smith (*B)->ops->choleskyfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal; 3563423f386SBarry Smith 3573423f386SBarry Smith (*B)->ops->shift = NULL; 3583423f386SBarry Smith (*B)->ops->scale = NULL; 3593423f386SBarry Smith (*B)->ops->mult = NULL; 3603423f386SBarry Smith (*B)->ops->sor = NULL; 3613423f386SBarry Smith (*B)->ops->zeroentries = NULL; 3623423f386SBarry Smith 363*9566063dSJacob Faibussowitsch PetscCall(PetscFree((*B)->solvertype)); 364*9566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERPETSC,&(*B)->solvertype)); 3653423f386SBarry Smith PetscFunctionReturn(0); 3663423f386SBarry Smith } 367