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 89371c9d4SSatish Balay static PetscErrorCode MatAXPY_ConstantDiagonal(Mat Y, PetscScalar a, Mat X, MatStructure str) { 9a01b8767SStefano Zampini Mat_ConstantDiagonal *yctx = (Mat_ConstantDiagonal *)Y->data; 10a01b8767SStefano Zampini Mat_ConstantDiagonal *xctx = (Mat_ConstantDiagonal *)X->data; 11a01b8767SStefano Zampini 12a01b8767SStefano Zampini PetscFunctionBegin; 13a01b8767SStefano Zampini yctx->diag += a * xctx->diag; 14a01b8767SStefano Zampini PetscFunctionReturn(0); 15a01b8767SStefano Zampini } 16a01b8767SStefano Zampini 179371c9d4SSatish Balay static PetscErrorCode MatGetRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[]) { 18a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data; 193423f386SBarry Smith 20a8d4b458SStefano Zampini PetscFunctionBegin; 21a8d4b458SStefano Zampini if (ncols) *ncols = 1; 22a8d4b458SStefano Zampini if (cols) { 239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, cols)); 24a8d4b458SStefano Zampini (*cols)[0] = row; 25a8d4b458SStefano Zampini } 26a8d4b458SStefano Zampini if (vals) { 279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, vals)); 28a8d4b458SStefano Zampini (*vals)[0] = ctx->diag; 29a8d4b458SStefano Zampini } 30a8d4b458SStefano Zampini PetscFunctionReturn(0); 31a8d4b458SStefano Zampini } 32a8d4b458SStefano Zampini 339371c9d4SSatish Balay static PetscErrorCode MatRestoreRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[]) { 34a8d4b458SStefano Zampini PetscFunctionBegin; 35a8d4b458SStefano Zampini if (ncols) *ncols = 0; 361baa6e33SBarry Smith if (cols) PetscCall(PetscFree(*cols)); 371baa6e33SBarry Smith if (vals) PetscCall(PetscFree(*vals)); 38a8d4b458SStefano Zampini PetscFunctionReturn(0); 39a8d4b458SStefano Zampini } 40a8d4b458SStefano Zampini 419371c9d4SSatish Balay static PetscErrorCode MatMultTranspose_ConstantDiagonal(Mat A, Vec x, Vec y) { 42a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data; 43a8d4b458SStefano Zampini 44a8d4b458SStefano Zampini PetscFunctionBegin; 459566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y, ctx->diag, 0.0, x)); 46a8d4b458SStefano Zampini PetscFunctionReturn(0); 47a8d4b458SStefano Zampini } 48a8d4b458SStefano Zampini 499371c9d4SSatish Balay static PetscErrorCode MatMultAdd_ConstantDiagonal(Mat mat, Vec v1, Vec v2, Vec v3) { 50a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)mat->data; 51a8d4b458SStefano Zampini 52a8d4b458SStefano Zampini PetscFunctionBegin; 53a8d4b458SStefano Zampini if (v2 == v3) { 549566063dSJacob Faibussowitsch PetscCall(VecAXPBY(v3, ctx->diag, 1.0, v1)); 55a8d4b458SStefano Zampini } else { 569566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(v3, ctx->diag, 1.0, 0.0, v1, v2)); 57a8d4b458SStefano Zampini } 58a8d4b458SStefano Zampini PetscFunctionReturn(0); 59a8d4b458SStefano Zampini } 60a8d4b458SStefano Zampini 619371c9d4SSatish Balay static PetscErrorCode MatMultTransposeAdd_ConstantDiagonal(Mat mat, Vec v1, Vec v2, Vec v3) { 62a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)mat->data; 63a8d4b458SStefano Zampini 64a8d4b458SStefano Zampini PetscFunctionBegin; 65a8d4b458SStefano Zampini if (v2 == v3) { 669566063dSJacob Faibussowitsch PetscCall(VecAXPBY(v3, ctx->diag, 1.0, v1)); 67a8d4b458SStefano Zampini } else { 689566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(v3, ctx->diag, 1.0, 0.0, v1, v2)); 69a8d4b458SStefano Zampini } 70a8d4b458SStefano Zampini PetscFunctionReturn(0); 71a8d4b458SStefano Zampini } 72a8d4b458SStefano Zampini 739371c9d4SSatish Balay static PetscErrorCode MatNorm_ConstantDiagonal(Mat A, NormType type, PetscReal *nrm) { 74e5f3c4beSJose E. Roman Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data; 75e5f3c4beSJose E. Roman 76e5f3c4beSJose E. Roman PetscFunctionBegin; 77e5f3c4beSJose E. Roman if (type == NORM_FROBENIUS || type == NORM_2 || type == NORM_1 || type == NORM_INFINITY) *nrm = PetscAbsScalar(ctx->diag); 78e5f3c4beSJose E. Roman else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Unsupported norm"); 79e5f3c4beSJose E. Roman PetscFunctionReturn(0); 80e5f3c4beSJose E. Roman } 81e5f3c4beSJose E. Roman 82d1ca2681SJose E. Roman static PetscErrorCode MatCreateSubMatrices_ConstantDiagonal(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 83d1ca2681SJose E. Roman 84d1ca2681SJose E. Roman { 85d1ca2681SJose E. Roman Mat B; 86d1ca2681SJose E. Roman 87d1ca2681SJose E. Roman PetscFunctionBegin; 889566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B)); 899566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B, n, irow, icol, scall, submat)); 909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 91d1ca2681SJose E. Roman PetscFunctionReturn(0); 92d1ca2681SJose E. Roman } 93d1ca2681SJose E. Roman 949371c9d4SSatish Balay static PetscErrorCode MatDuplicate_ConstantDiagonal(Mat A, MatDuplicateOption op, Mat *B) { 95a8d4b458SStefano Zampini Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal *)A->data; 96a8d4b458SStefano Zampini 97a8d4b458SStefano Zampini PetscFunctionBegin; 989566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B)); 999566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 1009566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizesFromMats(*B, A, A)); 1019566063dSJacob Faibussowitsch PetscCall(MatSetType(*B, MATCONSTANTDIAGONAL)); 1029566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->rmap, &(*B)->rmap)); 1039566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap, &(*B)->cmap)); 104a8d4b458SStefano Zampini if (op == MAT_COPY_VALUES) { 105a8d4b458SStefano Zampini Mat_ConstantDiagonal *bctx = (Mat_ConstantDiagonal *)(*B)->data; 106a8d4b458SStefano Zampini bctx->diag = actx->diag; 107a8d4b458SStefano Zampini } 108a8d4b458SStefano Zampini PetscFunctionReturn(0); 109a8d4b458SStefano Zampini } 110a8d4b458SStefano Zampini 1119371c9d4SSatish Balay static PetscErrorCode MatMissingDiagonal_ConstantDiagonal(Mat mat, PetscBool *missing, PetscInt *dd) { 112a8d4b458SStefano Zampini PetscFunctionBegin; 113a8d4b458SStefano Zampini *missing = PETSC_FALSE; 114a8d4b458SStefano Zampini PetscFunctionReturn(0); 115a8d4b458SStefano Zampini } 116a8d4b458SStefano Zampini 1179371c9d4SSatish Balay static PetscErrorCode MatDestroy_ConstantDiagonal(Mat mat) { 1183423f386SBarry Smith PetscFunctionBegin; 1199566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 1203423f386SBarry Smith PetscFunctionReturn(0); 1213423f386SBarry Smith } 1223423f386SBarry Smith 1239371c9d4SSatish Balay static PetscErrorCode MatView_ConstantDiagonal(Mat J, PetscViewer viewer) { 1243423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data; 1253423f386SBarry Smith PetscBool iascii; 1263423f386SBarry Smith 1273423f386SBarry Smith PetscFunctionBegin; 1289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1293423f386SBarry Smith if (iascii) { 1303423f386SBarry Smith PetscViewerFormat format; 1313423f386SBarry Smith 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 1333423f386SBarry Smith if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(0); 134c0aa6a63SJacob Faibussowitsch #if defined(PETSC_USE_COMPLEX) 1359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g + i %g\n", (double)PetscRealPart(ctx->diag), (double)PetscImaginaryPart(ctx->diag))); 1363423f386SBarry Smith #else 1379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g\n", (double)(ctx->diag))); 1383423f386SBarry Smith #endif 1393423f386SBarry Smith } 1403423f386SBarry Smith PetscFunctionReturn(0); 1413423f386SBarry Smith } 1423423f386SBarry Smith 1439371c9d4SSatish Balay static PetscErrorCode MatAssemblyEnd_ConstantDiagonal(Mat J, MatAssemblyType mt) { 1443423f386SBarry Smith PetscFunctionBegin; 1453423f386SBarry Smith PetscFunctionReturn(0); 1463423f386SBarry Smith } 1473423f386SBarry Smith 1489371c9d4SSatish Balay static PetscErrorCode MatMult_ConstantDiagonal(Mat J, Vec x, Vec y) { 1493423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data; 1503423f386SBarry Smith 1513423f386SBarry Smith PetscFunctionBegin; 1529566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y, ctx->diag, 0.0, x)); 1533423f386SBarry Smith PetscFunctionReturn(0); 1543423f386SBarry Smith } 1553423f386SBarry Smith 1569371c9d4SSatish Balay PetscErrorCode MatGetDiagonal_ConstantDiagonal(Mat J, Vec x) { 1573423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data; 1583423f386SBarry Smith 1593423f386SBarry Smith PetscFunctionBegin; 1609566063dSJacob Faibussowitsch PetscCall(VecSet(x, ctx->diag)); 1613423f386SBarry Smith PetscFunctionReturn(0); 1623423f386SBarry Smith } 1633423f386SBarry Smith 1649371c9d4SSatish Balay static PetscErrorCode MatShift_ConstantDiagonal(Mat Y, PetscScalar a) { 1653423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data; 1663423f386SBarry Smith 1673423f386SBarry Smith PetscFunctionBegin; 1683423f386SBarry Smith ctx->diag += a; 1693423f386SBarry Smith PetscFunctionReturn(0); 1703423f386SBarry Smith } 1713423f386SBarry Smith 1729371c9d4SSatish Balay static PetscErrorCode MatScale_ConstantDiagonal(Mat Y, PetscScalar a) { 1733423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data; 1743423f386SBarry Smith 1753423f386SBarry Smith PetscFunctionBegin; 1763423f386SBarry Smith ctx->diag *= a; 1773423f386SBarry Smith PetscFunctionReturn(0); 1783423f386SBarry Smith } 1793423f386SBarry Smith 1809371c9d4SSatish Balay static PetscErrorCode MatZeroEntries_ConstantDiagonal(Mat Y) { 1813423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data; 1823423f386SBarry Smith 1833423f386SBarry Smith PetscFunctionBegin; 1843423f386SBarry Smith ctx->diag = 0.0; 1853423f386SBarry Smith PetscFunctionReturn(0); 1863423f386SBarry Smith } 1873423f386SBarry Smith 1889371c9d4SSatish Balay PetscErrorCode MatSOR_ConstantDiagonal(Mat matin, Vec x, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec y) { 1893423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)matin->data; 1903423f386SBarry Smith 1913423f386SBarry Smith PetscFunctionBegin; 1923423f386SBarry Smith if (ctx->diag == 0.0) matin->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 1933423f386SBarry Smith else matin->factorerrortype = MAT_FACTOR_NOERROR; 1949566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y, 1.0 / ctx->diag, 0.0, x)); 1953423f386SBarry Smith PetscFunctionReturn(0); 1963423f386SBarry Smith } 1973423f386SBarry Smith 1989371c9d4SSatish Balay PetscErrorCode MatGetInfo_ConstantDiagonal(Mat A, MatInfoType flag, MatInfo *info) { 1993423f386SBarry Smith PetscFunctionBegin; 2003423f386SBarry Smith info->block_size = 1.0; 2013423f386SBarry Smith info->nz_allocated = 1.0; 2023423f386SBarry Smith info->nz_used = 1.0; 2033423f386SBarry Smith info->nz_unneeded = 0.0; 2043966268fSBarry Smith info->assemblies = A->num_ass; 2053423f386SBarry Smith info->mallocs = 0.0; 206*4dfa11a4SJacob Faibussowitsch info->memory = 0; /* REVIEW ME */ 2073423f386SBarry Smith if (A->factortype) { 2083423f386SBarry Smith info->fill_ratio_given = 1.0; 2093423f386SBarry Smith info->fill_ratio_needed = 1.0; 2103423f386SBarry Smith info->factor_mallocs = 0.0; 2113423f386SBarry Smith } else { 2123423f386SBarry Smith info->fill_ratio_given = 0; 2133423f386SBarry Smith info->fill_ratio_needed = 0; 2143423f386SBarry Smith info->factor_mallocs = 0; 2153423f386SBarry Smith } 2163423f386SBarry Smith PetscFunctionReturn(0); 2173423f386SBarry Smith } 2183423f386SBarry Smith 2193423f386SBarry Smith /*@ 2203423f386SBarry Smith MatCreateConstantDiagonal - Creates a matrix with a uniform value along the diagonal 2213423f386SBarry Smith 222d083f849SBarry Smith Collective 2233423f386SBarry Smith 2243423f386SBarry Smith Input Parameters: 2253423f386SBarry Smith + comm - MPI communicator 22611a5261eSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if M is given) 2273423f386SBarry Smith This value should be the same as the local size used in creating the 2283423f386SBarry Smith y vector for the matrix-vector product y = Ax. 2293423f386SBarry Smith . n - This value should be the same as the local size used in creating the 2303423f386SBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 2313423f386SBarry Smith calculated if N is given) For square matrices n is almost always m. 23211a5261eSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if m is given) 23311a5261eSBarry Smith . N - number of global columns (or `PETSC_DETERMINE` to have calculated if n is given) 2343423f386SBarry Smith - diag - the diagonal value 2353423f386SBarry Smith 2363423f386SBarry Smith Output Parameter: 2373423f386SBarry Smith . J - the diagonal matrix 2383423f386SBarry Smith 2393423f386SBarry Smith Level: advanced 2403423f386SBarry Smith 2413423f386SBarry Smith Notes: 2423423f386SBarry Smith Only supports square matrices with the same number of local rows and columns 2433423f386SBarry Smith 244db781477SPatrick Sanan .seealso: `MatDestroy()`, `MATCONSTANTDIAGONAL`, `MatScale()`, `MatShift()`, `MatMult()`, `MatGetDiagonal()`, `MatGetFactor()`, `MatSolve()` 2453423f386SBarry Smith @*/ 2469371c9d4SSatish Balay PetscErrorCode MatCreateConstantDiagonal(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscScalar diag, Mat *J) { 2473423f386SBarry Smith PetscFunctionBegin; 2489566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, J)); 2499566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, m, n, M, N)); 2509566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, MATCONSTANTDIAGONAL)); 2519566063dSJacob Faibussowitsch PetscCall(MatShift(*J, diag)); 2529566063dSJacob Faibussowitsch PetscCall(MatSetUp(*J)); 2533423f386SBarry Smith PetscFunctionReturn(0); 2543423f386SBarry Smith } 2553423f386SBarry Smith 2569371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode MatCreate_ConstantDiagonal(Mat A) { 2573423f386SBarry Smith Mat_ConstantDiagonal *ctx; 2583423f386SBarry Smith 2593423f386SBarry Smith PetscFunctionBegin; 2609566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 2613423f386SBarry Smith ctx->diag = 0.0; 2623423f386SBarry Smith A->data = (void *)ctx; 2633423f386SBarry Smith 2643423f386SBarry Smith A->assembled = PETSC_TRUE; 2653423f386SBarry Smith A->preallocated = PETSC_TRUE; 266a8d4b458SStefano Zampini 2673423f386SBarry Smith A->ops->mult = MatMult_ConstantDiagonal; 268a8d4b458SStefano Zampini A->ops->multadd = MatMultAdd_ConstantDiagonal; 269a8d4b458SStefano Zampini A->ops->multtranspose = MatMultTranspose_ConstantDiagonal; 270a8d4b458SStefano Zampini A->ops->multtransposeadd = MatMultTransposeAdd_ConstantDiagonal; 271e5f3c4beSJose E. Roman A->ops->norm = MatNorm_ConstantDiagonal; 272d1ca2681SJose E. Roman A->ops->createsubmatrices = MatCreateSubMatrices_ConstantDiagonal; 273a8d4b458SStefano Zampini A->ops->duplicate = MatDuplicate_ConstantDiagonal; 274a8d4b458SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_ConstantDiagonal; 275a8d4b458SStefano Zampini A->ops->getrow = MatGetRow_ConstantDiagonal; 276a8d4b458SStefano Zampini A->ops->restorerow = MatRestoreRow_ConstantDiagonal; 2773423f386SBarry Smith A->ops->sor = MatSOR_ConstantDiagonal; 2783423f386SBarry Smith A->ops->shift = MatShift_ConstantDiagonal; 2793423f386SBarry Smith A->ops->scale = MatScale_ConstantDiagonal; 2803423f386SBarry Smith A->ops->getdiagonal = MatGetDiagonal_ConstantDiagonal; 2813423f386SBarry Smith A->ops->view = MatView_ConstantDiagonal; 2823423f386SBarry Smith A->ops->zeroentries = MatZeroEntries_ConstantDiagonal; 2833423f386SBarry Smith A->ops->assemblyend = MatAssemblyEnd_ConstantDiagonal; 2843423f386SBarry Smith A->ops->destroy = MatDestroy_ConstantDiagonal; 2853423f386SBarry Smith A->ops->getinfo = MatGetInfo_ConstantDiagonal; 286a01b8767SStefano Zampini A->ops->axpy = MatAXPY_ConstantDiagonal; 287a8d4b458SStefano Zampini 2889566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATCONSTANTDIAGONAL)); 2893423f386SBarry Smith PetscFunctionReturn(0); 2903423f386SBarry Smith } 2913423f386SBarry Smith 2929371c9d4SSatish Balay static PetscErrorCode MatFactorNumeric_ConstantDiagonal(Mat fact, Mat A, const MatFactorInfo *info) { 2933423f386SBarry Smith Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal *)A->data, *fctx = (Mat_ConstantDiagonal *)fact->data; 2943423f386SBarry Smith 2953423f386SBarry Smith PetscFunctionBegin; 2963423f386SBarry Smith if (actx->diag == 0.0) fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 2973423f386SBarry Smith else fact->factorerrortype = MAT_FACTOR_NOERROR; 2983423f386SBarry Smith fctx->diag = 1.0 / actx->diag; 2993423f386SBarry Smith fact->ops->solve = MatMult_ConstantDiagonal; 3003423f386SBarry Smith PetscFunctionReturn(0); 3013423f386SBarry Smith } 3023423f386SBarry Smith 3039371c9d4SSatish Balay static PetscErrorCode MatFactorSymbolic_LU_ConstantDiagonal(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info) { 3043423f386SBarry Smith PetscFunctionBegin; 3053423f386SBarry Smith fact->ops->lufactornumeric = MatFactorNumeric_ConstantDiagonal; 3063423f386SBarry Smith PetscFunctionReturn(0); 3073423f386SBarry Smith } 3083423f386SBarry Smith 3099371c9d4SSatish Balay static PetscErrorCode MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact, Mat A, IS isrow, const MatFactorInfo *info) { 3103423f386SBarry Smith PetscFunctionBegin; 311a01b8767SStefano Zampini fact->ops->choleskyfactornumeric = MatFactorNumeric_ConstantDiagonal; 3123423f386SBarry Smith PetscFunctionReturn(0); 3133423f386SBarry Smith } 3143423f386SBarry Smith 3159371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatGetFactor_constantdiagonal_petsc(Mat A, MatFactorType ftype, Mat *B) { 3163423f386SBarry Smith PetscInt n = A->rmap->n, N = A->rmap->N; 3173423f386SBarry Smith 3183423f386SBarry Smith PetscFunctionBegin; 3199566063dSJacob Faibussowitsch PetscCall(MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A), n, n, N, N, 0, B)); 3203423f386SBarry Smith 3213423f386SBarry Smith (*B)->factortype = ftype; 3223423f386SBarry Smith (*B)->ops->ilufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal; 3233423f386SBarry Smith (*B)->ops->lufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal; 3243423f386SBarry Smith (*B)->ops->iccfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal; 3253423f386SBarry Smith (*B)->ops->choleskyfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal; 3263423f386SBarry Smith 3273423f386SBarry Smith (*B)->ops->shift = NULL; 3283423f386SBarry Smith (*B)->ops->scale = NULL; 3293423f386SBarry Smith (*B)->ops->mult = NULL; 3303423f386SBarry Smith (*B)->ops->sor = NULL; 3313423f386SBarry Smith (*B)->ops->zeroentries = NULL; 3323423f386SBarry Smith 3339566063dSJacob Faibussowitsch PetscCall(PetscFree((*B)->solvertype)); 3349566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERPETSC, &(*B)->solvertype)); 3353423f386SBarry Smith PetscFunctionReturn(0); 3363423f386SBarry Smith } 337