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 8d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAXPY_ConstantDiagonal(Mat Y, PetscScalar a, Mat X, MatStructure str) 9d71ae5a4SJacob Faibussowitsch { 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; 153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16a01b8767SStefano Zampini } 17a01b8767SStefano Zampini 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[]) 19d71ae5a4SJacob Faibussowitsch { 20a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data; 213423f386SBarry Smith 22a8d4b458SStefano Zampini PetscFunctionBegin; 23a8d4b458SStefano Zampini if (ncols) *ncols = 1; 24a8d4b458SStefano Zampini if (cols) { 259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, cols)); 26a8d4b458SStefano Zampini (*cols)[0] = row; 27a8d4b458SStefano Zampini } 28a8d4b458SStefano Zampini if (vals) { 299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(1, vals)); 30a8d4b458SStefano Zampini (*vals)[0] = ctx->diag; 31a8d4b458SStefano Zampini } 323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33a8d4b458SStefano Zampini } 34a8d4b458SStefano Zampini 35d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[]) 36d71ae5a4SJacob Faibussowitsch { 37a8d4b458SStefano Zampini PetscFunctionBegin; 38a8d4b458SStefano Zampini if (ncols) *ncols = 0; 391baa6e33SBarry Smith if (cols) PetscCall(PetscFree(*cols)); 401baa6e33SBarry Smith if (vals) PetscCall(PetscFree(*vals)); 413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42a8d4b458SStefano Zampini } 43a8d4b458SStefano Zampini 44d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTranspose_ConstantDiagonal(Mat A, Vec x, Vec y) 45d71ae5a4SJacob Faibussowitsch { 46a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data; 47a8d4b458SStefano Zampini 48a8d4b458SStefano Zampini PetscFunctionBegin; 499566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y, ctx->diag, 0.0, x)); 503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51a8d4b458SStefano Zampini } 52a8d4b458SStefano Zampini 53d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_ConstantDiagonal(Mat mat, Vec v1, Vec v2, Vec v3) 54d71ae5a4SJacob Faibussowitsch { 55a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)mat->data; 56a8d4b458SStefano Zampini 57a8d4b458SStefano Zampini PetscFunctionBegin; 58a8d4b458SStefano Zampini if (v2 == v3) { 599566063dSJacob Faibussowitsch PetscCall(VecAXPBY(v3, ctx->diag, 1.0, v1)); 60a8d4b458SStefano Zampini } else { 619566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(v3, ctx->diag, 1.0, 0.0, v1, v2)); 62a8d4b458SStefano Zampini } 633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 64a8d4b458SStefano Zampini } 65a8d4b458SStefano Zampini 66d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultTransposeAdd_ConstantDiagonal(Mat mat, Vec v1, Vec v2, Vec v3) 67d71ae5a4SJacob Faibussowitsch { 68a8d4b458SStefano Zampini Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)mat->data; 69a8d4b458SStefano Zampini 70a8d4b458SStefano Zampini PetscFunctionBegin; 71a8d4b458SStefano Zampini if (v2 == v3) { 729566063dSJacob Faibussowitsch PetscCall(VecAXPBY(v3, ctx->diag, 1.0, v1)); 73a8d4b458SStefano Zampini } else { 749566063dSJacob Faibussowitsch PetscCall(VecAXPBYPCZ(v3, ctx->diag, 1.0, 0.0, v1, v2)); 75a8d4b458SStefano Zampini } 763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 77a8d4b458SStefano Zampini } 78a8d4b458SStefano Zampini 79d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNorm_ConstantDiagonal(Mat A, NormType type, PetscReal *nrm) 80d71ae5a4SJacob Faibussowitsch { 81e5f3c4beSJose E. Roman Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data; 82e5f3c4beSJose E. Roman 83e5f3c4beSJose E. Roman PetscFunctionBegin; 84e5f3c4beSJose E. Roman if (type == NORM_FROBENIUS || type == NORM_2 || type == NORM_1 || type == NORM_INFINITY) *nrm = PetscAbsScalar(ctx->diag); 85e5f3c4beSJose E. Roman else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Unsupported norm"); 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 87e5f3c4beSJose E. Roman } 88e5f3c4beSJose E. Roman 89d1ca2681SJose E. Roman static PetscErrorCode MatCreateSubMatrices_ConstantDiagonal(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[]) 90d1ca2681SJose E. Roman 91d1ca2681SJose E. Roman { 92d1ca2681SJose E. Roman Mat B; 93d1ca2681SJose E. Roman 94d1ca2681SJose E. Roman PetscFunctionBegin; 959566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B)); 969566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B, n, irow, icol, scall, submat)); 979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 99d1ca2681SJose E. Roman } 100d1ca2681SJose E. Roman 101d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDuplicate_ConstantDiagonal(Mat A, MatDuplicateOption op, Mat *B) 102d71ae5a4SJacob Faibussowitsch { 103a8d4b458SStefano Zampini Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal *)A->data; 104a8d4b458SStefano Zampini 105a8d4b458SStefano Zampini PetscFunctionBegin; 1069566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B)); 1079566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N)); 1089566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizesFromMats(*B, A, A)); 1099566063dSJacob Faibussowitsch PetscCall(MatSetType(*B, MATCONSTANTDIAGONAL)); 1109566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->rmap, &(*B)->rmap)); 1119566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap, &(*B)->cmap)); 112a8d4b458SStefano Zampini if (op == MAT_COPY_VALUES) { 113a8d4b458SStefano Zampini Mat_ConstantDiagonal *bctx = (Mat_ConstantDiagonal *)(*B)->data; 114a8d4b458SStefano Zampini bctx->diag = actx->diag; 115a8d4b458SStefano Zampini } 1163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 117a8d4b458SStefano Zampini } 118a8d4b458SStefano Zampini 119d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMissingDiagonal_ConstantDiagonal(Mat mat, PetscBool *missing, PetscInt *dd) 120d71ae5a4SJacob Faibussowitsch { 121a8d4b458SStefano Zampini PetscFunctionBegin; 122a8d4b458SStefano Zampini *missing = PETSC_FALSE; 1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124a8d4b458SStefano Zampini } 125a8d4b458SStefano Zampini 126d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_ConstantDiagonal(Mat mat) 127d71ae5a4SJacob Faibussowitsch { 1283423f386SBarry Smith PetscFunctionBegin; 1299566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data)); 1303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1313423f386SBarry Smith } 1323423f386SBarry Smith 133d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_ConstantDiagonal(Mat J, PetscViewer viewer) 134d71ae5a4SJacob Faibussowitsch { 1353423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data; 1363423f386SBarry Smith PetscBool iascii; 1373423f386SBarry Smith 1383423f386SBarry Smith PetscFunctionBegin; 1399566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1403423f386SBarry Smith if (iascii) { 1413423f386SBarry Smith PetscViewerFormat format; 1423423f386SBarry Smith 1439566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 1443ba16761SJacob Faibussowitsch if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(PETSC_SUCCESS); 145c0aa6a63SJacob Faibussowitsch #if defined(PETSC_USE_COMPLEX) 1469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g + i %g\n", (double)PetscRealPart(ctx->diag), (double)PetscImaginaryPart(ctx->diag))); 1473423f386SBarry Smith #else 1489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g\n", (double)(ctx->diag))); 1493423f386SBarry Smith #endif 1503423f386SBarry Smith } 1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1523423f386SBarry Smith } 1533423f386SBarry Smith 154d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_ConstantDiagonal(Mat J, MatAssemblyType mt) 155d71ae5a4SJacob Faibussowitsch { 1563423f386SBarry Smith PetscFunctionBegin; 1573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1583423f386SBarry Smith } 1593423f386SBarry Smith 160d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_ConstantDiagonal(Mat J, Vec x, Vec y) 161d71ae5a4SJacob Faibussowitsch { 1623423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data; 1633423f386SBarry Smith 1643423f386SBarry Smith PetscFunctionBegin; 1659566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y, ctx->diag, 0.0, x)); 1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1673423f386SBarry Smith } 1683423f386SBarry Smith 169d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetDiagonal_ConstantDiagonal(Mat J, Vec x) 170d71ae5a4SJacob Faibussowitsch { 1713423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data; 1723423f386SBarry Smith 1733423f386SBarry Smith PetscFunctionBegin; 1749566063dSJacob Faibussowitsch PetscCall(VecSet(x, ctx->diag)); 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1763423f386SBarry Smith } 1773423f386SBarry Smith 178d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatShift_ConstantDiagonal(Mat Y, PetscScalar a) 179d71ae5a4SJacob Faibussowitsch { 1803423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data; 1813423f386SBarry Smith 1823423f386SBarry Smith PetscFunctionBegin; 1833423f386SBarry Smith ctx->diag += a; 1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1853423f386SBarry Smith } 1863423f386SBarry Smith 187d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_ConstantDiagonal(Mat Y, PetscScalar a) 188d71ae5a4SJacob Faibussowitsch { 1893423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data; 1903423f386SBarry Smith 1913423f386SBarry Smith PetscFunctionBegin; 1923423f386SBarry Smith ctx->diag *= a; 1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1943423f386SBarry Smith } 1953423f386SBarry Smith 196d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatZeroEntries_ConstantDiagonal(Mat Y) 197d71ae5a4SJacob Faibussowitsch { 1983423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data; 1993423f386SBarry Smith 2003423f386SBarry Smith PetscFunctionBegin; 2013423f386SBarry Smith ctx->diag = 0.0; 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2033423f386SBarry Smith } 2043423f386SBarry Smith 205d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSOR_ConstantDiagonal(Mat matin, Vec x, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec y) 206d71ae5a4SJacob Faibussowitsch { 2073423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)matin->data; 2083423f386SBarry Smith 2093423f386SBarry Smith PetscFunctionBegin; 2103423f386SBarry Smith if (ctx->diag == 0.0) matin->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 2113423f386SBarry Smith else matin->factorerrortype = MAT_FACTOR_NOERROR; 2129566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y, 1.0 / ctx->diag, 0.0, x)); 2133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2143423f386SBarry Smith } 2153423f386SBarry Smith 216d71ae5a4SJacob Faibussowitsch PetscErrorCode MatGetInfo_ConstantDiagonal(Mat A, MatInfoType flag, MatInfo *info) 217d71ae5a4SJacob Faibussowitsch { 2183423f386SBarry Smith PetscFunctionBegin; 2193423f386SBarry Smith info->block_size = 1.0; 2203423f386SBarry Smith info->nz_allocated = 1.0; 2213423f386SBarry Smith info->nz_used = 1.0; 2223423f386SBarry Smith info->nz_unneeded = 0.0; 2233966268fSBarry Smith info->assemblies = A->num_ass; 2243423f386SBarry Smith info->mallocs = 0.0; 2254dfa11a4SJacob Faibussowitsch info->memory = 0; /* REVIEW ME */ 2263423f386SBarry Smith if (A->factortype) { 2273423f386SBarry Smith info->fill_ratio_given = 1.0; 2283423f386SBarry Smith info->fill_ratio_needed = 1.0; 2293423f386SBarry Smith info->factor_mallocs = 0.0; 2303423f386SBarry Smith } else { 2313423f386SBarry Smith info->fill_ratio_given = 0; 2323423f386SBarry Smith info->fill_ratio_needed = 0; 2333423f386SBarry Smith info->factor_mallocs = 0; 2343423f386SBarry Smith } 2353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2363423f386SBarry Smith } 2373423f386SBarry Smith 2383423f386SBarry Smith /*@ 2393423f386SBarry Smith MatCreateConstantDiagonal - Creates a matrix with a uniform value along the diagonal 2403423f386SBarry Smith 241d083f849SBarry Smith Collective 2423423f386SBarry Smith 2433423f386SBarry Smith Input Parameters: 2443423f386SBarry Smith + comm - MPI communicator 2452ef1f0ffSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given) 2463423f386SBarry Smith This value should be the same as the local size used in creating the 2473423f386SBarry Smith y vector for the matrix-vector product y = Ax. 2483423f386SBarry Smith . n - This value should be the same as the local size used in creating the 2492ef1f0ffSBarry Smith x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have 2502ef1f0ffSBarry Smith calculated if `N` is given) For square matrices n is almost always `m`. 25111a5261eSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if m is given) 25211a5261eSBarry Smith . N - number of global columns (or `PETSC_DETERMINE` to have calculated if n is given) 2533423f386SBarry Smith - diag - the diagonal value 2543423f386SBarry Smith 2553423f386SBarry Smith Output Parameter: 2563423f386SBarry Smith . J - the diagonal matrix 2573423f386SBarry Smith 2583423f386SBarry Smith Level: advanced 2593423f386SBarry Smith 2603423f386SBarry Smith Notes: 2613423f386SBarry Smith Only supports square matrices with the same number of local rows and columns 2623423f386SBarry Smith 263*1cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MATCONSTANTDIAGONAL`, `MatScale()`, `MatShift()`, `MatMult()`, `MatGetDiagonal()`, `MatGetFactor()`, `MatSolve()` 2643423f386SBarry Smith @*/ 265d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateConstantDiagonal(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscScalar diag, Mat *J) 266d71ae5a4SJacob Faibussowitsch { 2673423f386SBarry Smith PetscFunctionBegin; 2689566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, J)); 2699566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, m, n, M, N)); 2709566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, MATCONSTANTDIAGONAL)); 2719566063dSJacob Faibussowitsch PetscCall(MatShift(*J, diag)); 2729566063dSJacob Faibussowitsch PetscCall(MatSetUp(*J)); 2733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2743423f386SBarry Smith } 2753423f386SBarry Smith 276d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_ConstantDiagonal(Mat A) 277d71ae5a4SJacob Faibussowitsch { 2783423f386SBarry Smith Mat_ConstantDiagonal *ctx; 2793423f386SBarry Smith 2803423f386SBarry Smith PetscFunctionBegin; 2819566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 2823423f386SBarry Smith ctx->diag = 0.0; 2833423f386SBarry Smith A->data = (void *)ctx; 2843423f386SBarry Smith 2853423f386SBarry Smith A->assembled = PETSC_TRUE; 2863423f386SBarry Smith A->preallocated = PETSC_TRUE; 287a8d4b458SStefano Zampini 2883423f386SBarry Smith A->ops->mult = MatMult_ConstantDiagonal; 289a8d4b458SStefano Zampini A->ops->multadd = MatMultAdd_ConstantDiagonal; 290a8d4b458SStefano Zampini A->ops->multtranspose = MatMultTranspose_ConstantDiagonal; 291a8d4b458SStefano Zampini A->ops->multtransposeadd = MatMultTransposeAdd_ConstantDiagonal; 292e5f3c4beSJose E. Roman A->ops->norm = MatNorm_ConstantDiagonal; 293d1ca2681SJose E. Roman A->ops->createsubmatrices = MatCreateSubMatrices_ConstantDiagonal; 294a8d4b458SStefano Zampini A->ops->duplicate = MatDuplicate_ConstantDiagonal; 295a8d4b458SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_ConstantDiagonal; 296a8d4b458SStefano Zampini A->ops->getrow = MatGetRow_ConstantDiagonal; 297a8d4b458SStefano Zampini A->ops->restorerow = MatRestoreRow_ConstantDiagonal; 2983423f386SBarry Smith A->ops->sor = MatSOR_ConstantDiagonal; 2993423f386SBarry Smith A->ops->shift = MatShift_ConstantDiagonal; 3003423f386SBarry Smith A->ops->scale = MatScale_ConstantDiagonal; 3013423f386SBarry Smith A->ops->getdiagonal = MatGetDiagonal_ConstantDiagonal; 3023423f386SBarry Smith A->ops->view = MatView_ConstantDiagonal; 3033423f386SBarry Smith A->ops->zeroentries = MatZeroEntries_ConstantDiagonal; 3043423f386SBarry Smith A->ops->assemblyend = MatAssemblyEnd_ConstantDiagonal; 3053423f386SBarry Smith A->ops->destroy = MatDestroy_ConstantDiagonal; 3063423f386SBarry Smith A->ops->getinfo = MatGetInfo_ConstantDiagonal; 307a01b8767SStefano Zampini A->ops->axpy = MatAXPY_ConstantDiagonal; 308a8d4b458SStefano Zampini 3099566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATCONSTANTDIAGONAL)); 3103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3113423f386SBarry Smith } 3123423f386SBarry Smith 313d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorNumeric_ConstantDiagonal(Mat fact, Mat A, const MatFactorInfo *info) 314d71ae5a4SJacob Faibussowitsch { 3153423f386SBarry Smith Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal *)A->data, *fctx = (Mat_ConstantDiagonal *)fact->data; 3163423f386SBarry Smith 3173423f386SBarry Smith PetscFunctionBegin; 3183423f386SBarry Smith if (actx->diag == 0.0) fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 3193423f386SBarry Smith else fact->factorerrortype = MAT_FACTOR_NOERROR; 3203423f386SBarry Smith fctx->diag = 1.0 / actx->diag; 3213423f386SBarry Smith fact->ops->solve = MatMult_ConstantDiagonal; 3223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3233423f386SBarry Smith } 3243423f386SBarry Smith 325d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorSymbolic_LU_ConstantDiagonal(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info) 326d71ae5a4SJacob Faibussowitsch { 3273423f386SBarry Smith PetscFunctionBegin; 3283423f386SBarry Smith fact->ops->lufactornumeric = MatFactorNumeric_ConstantDiagonal; 3293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3303423f386SBarry Smith } 3313423f386SBarry Smith 332d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact, Mat A, IS isrow, const MatFactorInfo *info) 333d71ae5a4SJacob Faibussowitsch { 3343423f386SBarry Smith PetscFunctionBegin; 335a01b8767SStefano Zampini fact->ops->choleskyfactornumeric = MatFactorNumeric_ConstantDiagonal; 3363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3373423f386SBarry Smith } 3383423f386SBarry Smith 339d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_constantdiagonal_petsc(Mat A, MatFactorType ftype, Mat *B) 340d71ae5a4SJacob Faibussowitsch { 3413423f386SBarry Smith PetscInt n = A->rmap->n, N = A->rmap->N; 3423423f386SBarry Smith 3433423f386SBarry Smith PetscFunctionBegin; 3449566063dSJacob Faibussowitsch PetscCall(MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A), n, n, N, N, 0, B)); 3453423f386SBarry Smith 3463423f386SBarry Smith (*B)->factortype = ftype; 3473423f386SBarry Smith (*B)->ops->ilufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal; 3483423f386SBarry Smith (*B)->ops->lufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal; 3493423f386SBarry Smith (*B)->ops->iccfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal; 3503423f386SBarry Smith (*B)->ops->choleskyfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal; 3513423f386SBarry Smith 3523423f386SBarry Smith (*B)->ops->shift = NULL; 3533423f386SBarry Smith (*B)->ops->scale = NULL; 3543423f386SBarry Smith (*B)->ops->mult = NULL; 3553423f386SBarry Smith (*B)->ops->sor = NULL; 3563423f386SBarry Smith (*B)->ops->zeroentries = NULL; 3573423f386SBarry Smith 3589566063dSJacob Faibussowitsch PetscCall(PetscFree((*B)->solvertype)); 3599566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERPETSC, &(*B)->solvertype)); 3603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3613423f386SBarry Smith } 362