xref: /petsc/src/mat/impls/cdiagonal/cdiagonal.c (revision 66976f2f44dcc61d86a452a70219fb23b45d00f0)
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 
18345a4b08SToby Isaac static PetscErrorCode MatEqual_ConstantDiagonal(Mat Y, Mat X, PetscBool *equal)
19345a4b08SToby Isaac {
20345a4b08SToby Isaac   Mat_ConstantDiagonal *yctx = (Mat_ConstantDiagonal *)Y->data;
21345a4b08SToby Isaac   Mat_ConstantDiagonal *xctx = (Mat_ConstantDiagonal *)X->data;
22345a4b08SToby Isaac 
23345a4b08SToby Isaac   PetscFunctionBegin;
24345a4b08SToby Isaac   *equal = (yctx->diag == xctx->diag) ? PETSC_TRUE : PETSC_FALSE;
25345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
26345a4b08SToby Isaac }
27345a4b08SToby Isaac 
28d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[])
29d71ae5a4SJacob Faibussowitsch {
30a8d4b458SStefano Zampini   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data;
313423f386SBarry Smith 
32a8d4b458SStefano Zampini   PetscFunctionBegin;
33a8d4b458SStefano Zampini   if (ncols) *ncols = 1;
34a8d4b458SStefano Zampini   if (cols) {
359566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, cols));
36a8d4b458SStefano Zampini     (*cols)[0] = row;
37a8d4b458SStefano Zampini   }
38a8d4b458SStefano Zampini   if (vals) {
399566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, vals));
40a8d4b458SStefano Zampini     (*vals)[0] = ctx->diag;
41a8d4b458SStefano Zampini   }
423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43a8d4b458SStefano Zampini }
44a8d4b458SStefano Zampini 
45d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[])
46d71ae5a4SJacob Faibussowitsch {
47a8d4b458SStefano Zampini   PetscFunctionBegin;
48a8d4b458SStefano Zampini   if (ncols) *ncols = 0;
491baa6e33SBarry Smith   if (cols) PetscCall(PetscFree(*cols));
501baa6e33SBarry Smith   if (vals) PetscCall(PetscFree(*vals));
513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52a8d4b458SStefano Zampini }
53a8d4b458SStefano Zampini 
54d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_ConstantDiagonal(Mat mat, Vec v1, Vec v2, Vec v3)
55d71ae5a4SJacob Faibussowitsch {
56a8d4b458SStefano Zampini   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)mat->data;
57a8d4b458SStefano Zampini 
58a8d4b458SStefano Zampini   PetscFunctionBegin;
59a8d4b458SStefano Zampini   if (v2 == v3) {
609566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(v3, ctx->diag, 1.0, v1));
61a8d4b458SStefano Zampini   } else {
629566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(v3, ctx->diag, 1.0, 0.0, v1, v2));
63a8d4b458SStefano Zampini   }
643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
65a8d4b458SStefano Zampini }
66a8d4b458SStefano Zampini 
67345a4b08SToby Isaac static PetscErrorCode MatMultHermitianTransposeAdd_ConstantDiagonal(Mat mat, Vec v1, Vec v2, Vec v3)
68d71ae5a4SJacob Faibussowitsch {
69a8d4b458SStefano Zampini   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)mat->data;
70a8d4b458SStefano Zampini 
71a8d4b458SStefano Zampini   PetscFunctionBegin;
72a8d4b458SStefano Zampini   if (v2 == v3) {
73345a4b08SToby Isaac     PetscCall(VecAXPBY(v3, PetscConj(ctx->diag), 1.0, v1));
74a8d4b458SStefano Zampini   } else {
75345a4b08SToby Isaac     PetscCall(VecAXPBYPCZ(v3, PetscConj(ctx->diag), 1.0, 0.0, v1, v2));
76a8d4b458SStefano Zampini   }
773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
78a8d4b458SStefano Zampini }
79a8d4b458SStefano Zampini 
80d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNorm_ConstantDiagonal(Mat A, NormType type, PetscReal *nrm)
81d71ae5a4SJacob Faibussowitsch {
82e5f3c4beSJose E. Roman   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data;
83e5f3c4beSJose E. Roman 
84e5f3c4beSJose E. Roman   PetscFunctionBegin;
85e5f3c4beSJose E. Roman   if (type == NORM_FROBENIUS || type == NORM_2 || type == NORM_1 || type == NORM_INFINITY) *nrm = PetscAbsScalar(ctx->diag);
86e5f3c4beSJose E. Roman   else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Unsupported norm");
873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88e5f3c4beSJose E. Roman }
89e5f3c4beSJose E. Roman 
90d1ca2681SJose E. Roman static PetscErrorCode MatCreateSubMatrices_ConstantDiagonal(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
91d1ca2681SJose E. Roman 
92d1ca2681SJose E. Roman {
93d1ca2681SJose E. Roman   Mat B;
94d1ca2681SJose E. Roman 
95d1ca2681SJose E. Roman   PetscFunctionBegin;
969566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B));
979566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, n, irow, icol, scall, submat));
989566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
100d1ca2681SJose E. Roman }
101d1ca2681SJose E. Roman 
102d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDuplicate_ConstantDiagonal(Mat A, MatDuplicateOption op, Mat *B)
103d71ae5a4SJacob Faibussowitsch {
104a8d4b458SStefano Zampini   Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal *)A->data;
105a8d4b458SStefano Zampini 
106a8d4b458SStefano Zampini   PetscFunctionBegin;
1079566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
1089566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
1099566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizesFromMats(*B, A, A));
1109566063dSJacob Faibussowitsch   PetscCall(MatSetType(*B, MATCONSTANTDIAGONAL));
1119566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->rmap, &(*B)->rmap));
1129566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*B)->cmap));
113a8d4b458SStefano Zampini   if (op == MAT_COPY_VALUES) {
114a8d4b458SStefano Zampini     Mat_ConstantDiagonal *bctx = (Mat_ConstantDiagonal *)(*B)->data;
115a8d4b458SStefano Zampini     bctx->diag                 = actx->diag;
116a8d4b458SStefano Zampini   }
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118a8d4b458SStefano Zampini }
119a8d4b458SStefano Zampini 
120d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMissingDiagonal_ConstantDiagonal(Mat mat, PetscBool *missing, PetscInt *dd)
121d71ae5a4SJacob Faibussowitsch {
122a8d4b458SStefano Zampini   PetscFunctionBegin;
123a8d4b458SStefano Zampini   *missing = PETSC_FALSE;
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125a8d4b458SStefano Zampini }
126a8d4b458SStefano Zampini 
127d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_ConstantDiagonal(Mat mat)
128d71ae5a4SJacob Faibussowitsch {
1293423f386SBarry Smith   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
131345a4b08SToby Isaac   mat->structural_symmetry_eternal = PETSC_FALSE;
132345a4b08SToby Isaac   mat->symmetry_eternal            = PETSC_FALSE;
1333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1343423f386SBarry Smith }
1353423f386SBarry Smith 
136d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_ConstantDiagonal(Mat J, PetscViewer viewer)
137d71ae5a4SJacob Faibussowitsch {
1383423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
1393423f386SBarry Smith   PetscBool             iascii;
1403423f386SBarry Smith 
1413423f386SBarry Smith   PetscFunctionBegin;
1429566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1433423f386SBarry Smith   if (iascii) {
1443423f386SBarry Smith     PetscViewerFormat format;
1453423f386SBarry Smith 
1469566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
1473ba16761SJacob Faibussowitsch     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(PETSC_SUCCESS);
148c0aa6a63SJacob Faibussowitsch #if defined(PETSC_USE_COMPLEX)
1499566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g + i %g\n", (double)PetscRealPart(ctx->diag), (double)PetscImaginaryPart(ctx->diag)));
1503423f386SBarry Smith #else
1519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g\n", (double)(ctx->diag)));
1523423f386SBarry Smith #endif
1533423f386SBarry Smith   }
1543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1553423f386SBarry Smith }
1563423f386SBarry Smith 
157d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_ConstantDiagonal(Mat J, MatAssemblyType mt)
158d71ae5a4SJacob Faibussowitsch {
1593423f386SBarry Smith   PetscFunctionBegin;
1603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1613423f386SBarry Smith }
1623423f386SBarry Smith 
163d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_ConstantDiagonal(Mat J, Vec x, Vec y)
164d71ae5a4SJacob Faibussowitsch {
1653423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
1663423f386SBarry Smith 
1673423f386SBarry Smith   PetscFunctionBegin;
1689566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(y, ctx->diag, 0.0, x));
1693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1703423f386SBarry Smith }
1713423f386SBarry Smith 
172345a4b08SToby Isaac static PetscErrorCode MatMultHermitianTranspose_ConstantDiagonal(Mat J, Vec x, Vec y)
173345a4b08SToby Isaac {
174345a4b08SToby Isaac   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
175345a4b08SToby Isaac 
176345a4b08SToby Isaac   PetscFunctionBegin;
177345a4b08SToby Isaac   PetscCall(VecAXPBY(y, PetscConj(ctx->diag), 0.0, x));
178345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
179345a4b08SToby Isaac }
180345a4b08SToby Isaac 
181345a4b08SToby Isaac static PetscErrorCode MatGetDiagonal_ConstantDiagonal(Mat J, Vec x)
182d71ae5a4SJacob Faibussowitsch {
1833423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
1843423f386SBarry Smith 
1853423f386SBarry Smith   PetscFunctionBegin;
1869566063dSJacob Faibussowitsch   PetscCall(VecSet(x, ctx->diag));
1873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1883423f386SBarry Smith }
1893423f386SBarry Smith 
190d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatShift_ConstantDiagonal(Mat Y, PetscScalar a)
191d71ae5a4SJacob Faibussowitsch {
1923423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data;
1933423f386SBarry Smith 
1943423f386SBarry Smith   PetscFunctionBegin;
1953423f386SBarry Smith   ctx->diag += a;
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1973423f386SBarry Smith }
1983423f386SBarry Smith 
199d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_ConstantDiagonal(Mat Y, PetscScalar a)
200d71ae5a4SJacob Faibussowitsch {
2013423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data;
2023423f386SBarry Smith 
2033423f386SBarry Smith   PetscFunctionBegin;
2043423f386SBarry Smith   ctx->diag *= a;
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2063423f386SBarry Smith }
2073423f386SBarry Smith 
208d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatZeroEntries_ConstantDiagonal(Mat Y)
209d71ae5a4SJacob Faibussowitsch {
2103423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data;
2113423f386SBarry Smith 
2123423f386SBarry Smith   PetscFunctionBegin;
2133423f386SBarry Smith   ctx->diag = 0.0;
2143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2153423f386SBarry Smith }
2163423f386SBarry Smith 
217345a4b08SToby Isaac static PetscErrorCode MatSolve_ConstantDiagonal(Mat matin, Vec b, Vec x)
218d71ae5a4SJacob Faibussowitsch {
2193423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)matin->data;
2203423f386SBarry Smith 
2213423f386SBarry Smith   PetscFunctionBegin;
2223423f386SBarry Smith   if (ctx->diag == 0.0) matin->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
2233423f386SBarry Smith   else matin->factorerrortype = MAT_FACTOR_NOERROR;
224345a4b08SToby Isaac   PetscCall(VecAXPBY(x, 1.0 / ctx->diag, 0.0, b));
225345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
226345a4b08SToby Isaac }
227345a4b08SToby Isaac 
228*66976f2fSJacob Faibussowitsch static PetscErrorCode MatSOR_ConstantDiagonal(Mat matin, Vec x, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec y)
229345a4b08SToby Isaac {
230345a4b08SToby Isaac   PetscFunctionBegin;
231345a4b08SToby Isaac   PetscCall(MatSolve_ConstantDiagonal(matin, x, y));
2323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2333423f386SBarry Smith }
2343423f386SBarry Smith 
235*66976f2fSJacob Faibussowitsch static PetscErrorCode MatGetInfo_ConstantDiagonal(Mat A, MatInfoType flag, MatInfo *info)
236d71ae5a4SJacob Faibussowitsch {
2373423f386SBarry Smith   PetscFunctionBegin;
2383423f386SBarry Smith   info->block_size   = 1.0;
2393423f386SBarry Smith   info->nz_allocated = 1.0;
2403423f386SBarry Smith   info->nz_used      = 1.0;
2413423f386SBarry Smith   info->nz_unneeded  = 0.0;
2423966268fSBarry Smith   info->assemblies   = A->num_ass;
2433423f386SBarry Smith   info->mallocs      = 0.0;
2444dfa11a4SJacob Faibussowitsch   info->memory       = 0; /* REVIEW ME */
2453423f386SBarry Smith   if (A->factortype) {
2463423f386SBarry Smith     info->fill_ratio_given  = 1.0;
2473423f386SBarry Smith     info->fill_ratio_needed = 1.0;
2483423f386SBarry Smith     info->factor_mallocs    = 0.0;
2493423f386SBarry Smith   } else {
2503423f386SBarry Smith     info->fill_ratio_given  = 0;
2513423f386SBarry Smith     info->fill_ratio_needed = 0;
2523423f386SBarry Smith     info->factor_mallocs    = 0;
2533423f386SBarry Smith   }
2543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2553423f386SBarry Smith }
2563423f386SBarry Smith 
2573423f386SBarry Smith /*@
2583423f386SBarry Smith   MatCreateConstantDiagonal - Creates a matrix with a uniform value along the diagonal
2593423f386SBarry Smith 
260d083f849SBarry Smith   Collective
2613423f386SBarry Smith 
2623423f386SBarry Smith   Input Parameters:
2633423f386SBarry Smith + comm - MPI communicator
2642ef1f0ffSBarry Smith . m    - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
2653423f386SBarry Smith            This value should be the same as the local size used in creating the
2663423f386SBarry Smith            y vector for the matrix-vector product y = Ax.
2673423f386SBarry Smith . n    - This value should be the same as the local size used in creating the
2682ef1f0ffSBarry Smith        x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
2692ef1f0ffSBarry Smith        calculated if `N` is given) For square matrices n is almost always `m`.
27011a5261eSBarry Smith . M    - number of global rows (or `PETSC_DETERMINE` to have calculated if m is given)
27111a5261eSBarry Smith . N    - number of global columns (or `PETSC_DETERMINE` to have calculated if n is given)
2723423f386SBarry Smith - diag - the diagonal value
2733423f386SBarry Smith 
2743423f386SBarry Smith   Output Parameter:
2753423f386SBarry Smith . J - the diagonal matrix
2763423f386SBarry Smith 
2773423f386SBarry Smith   Level: advanced
2783423f386SBarry Smith 
2793423f386SBarry Smith   Notes:
2803423f386SBarry Smith   Only supports square matrices with the same number of local rows and columns
2813423f386SBarry Smith 
2821cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MATCONSTANTDIAGONAL`, `MatScale()`, `MatShift()`, `MatMult()`, `MatGetDiagonal()`, `MatGetFactor()`, `MatSolve()`
2833423f386SBarry Smith @*/
284d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateConstantDiagonal(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscScalar diag, Mat *J)
285d71ae5a4SJacob Faibussowitsch {
2863423f386SBarry Smith   PetscFunctionBegin;
2879566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, J));
2889566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, m, n, M, N));
2899566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATCONSTANTDIAGONAL));
2909566063dSJacob Faibussowitsch   PetscCall(MatShift(*J, diag));
2919566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*J));
2923ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2933423f386SBarry Smith }
2943423f386SBarry Smith 
295d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_ConstantDiagonal(Mat A)
296d71ae5a4SJacob Faibussowitsch {
2973423f386SBarry Smith   Mat_ConstantDiagonal *ctx;
2983423f386SBarry Smith 
2993423f386SBarry Smith   PetscFunctionBegin;
3009566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
3013423f386SBarry Smith   ctx->diag = 0.0;
3023423f386SBarry Smith   A->data   = (void *)ctx;
3033423f386SBarry Smith 
3043423f386SBarry Smith   A->assembled                   = PETSC_TRUE;
3053423f386SBarry Smith   A->preallocated                = PETSC_TRUE;
306345a4b08SToby Isaac   A->structurally_symmetric      = PETSC_BOOL3_TRUE;
307345a4b08SToby Isaac   A->structural_symmetry_eternal = PETSC_TRUE;
308345a4b08SToby Isaac   A->symmetric                   = PETSC_BOOL3_TRUE;
309345a4b08SToby Isaac   if (!PetscDefined(USE_COMPLEX)) A->hermitian = PETSC_BOOL3_TRUE;
310345a4b08SToby Isaac   A->symmetry_eternal = PETSC_TRUE;
311a8d4b458SStefano Zampini 
3123423f386SBarry Smith   A->ops->mult                      = MatMult_ConstantDiagonal;
313a8d4b458SStefano Zampini   A->ops->multadd                   = MatMultAdd_ConstantDiagonal;
314345a4b08SToby Isaac   A->ops->multtranspose             = MatMult_ConstantDiagonal;
315345a4b08SToby Isaac   A->ops->multtransposeadd          = MatMultAdd_ConstantDiagonal;
316345a4b08SToby Isaac   A->ops->multhermitiantranspose    = MatMultHermitianTranspose_ConstantDiagonal;
317345a4b08SToby Isaac   A->ops->multhermitiantransposeadd = MatMultHermitianTransposeAdd_ConstantDiagonal;
318345a4b08SToby Isaac   A->ops->solve                     = MatSolve_ConstantDiagonal;
319345a4b08SToby Isaac   A->ops->solvetranspose            = MatSolve_ConstantDiagonal;
320e5f3c4beSJose E. Roman   A->ops->norm                      = MatNorm_ConstantDiagonal;
321d1ca2681SJose E. Roman   A->ops->createsubmatrices         = MatCreateSubMatrices_ConstantDiagonal;
322a8d4b458SStefano Zampini   A->ops->duplicate                 = MatDuplicate_ConstantDiagonal;
323a8d4b458SStefano Zampini   A->ops->missingdiagonal           = MatMissingDiagonal_ConstantDiagonal;
324a8d4b458SStefano Zampini   A->ops->getrow                    = MatGetRow_ConstantDiagonal;
325a8d4b458SStefano Zampini   A->ops->restorerow                = MatRestoreRow_ConstantDiagonal;
3263423f386SBarry Smith   A->ops->sor                       = MatSOR_ConstantDiagonal;
3273423f386SBarry Smith   A->ops->shift                     = MatShift_ConstantDiagonal;
3283423f386SBarry Smith   A->ops->scale                     = MatScale_ConstantDiagonal;
3293423f386SBarry Smith   A->ops->getdiagonal               = MatGetDiagonal_ConstantDiagonal;
3303423f386SBarry Smith   A->ops->view                      = MatView_ConstantDiagonal;
3313423f386SBarry Smith   A->ops->zeroentries               = MatZeroEntries_ConstantDiagonal;
3323423f386SBarry Smith   A->ops->assemblyend               = MatAssemblyEnd_ConstantDiagonal;
3333423f386SBarry Smith   A->ops->destroy                   = MatDestroy_ConstantDiagonal;
3343423f386SBarry Smith   A->ops->getinfo                   = MatGetInfo_ConstantDiagonal;
335345a4b08SToby Isaac   A->ops->equal                     = MatEqual_ConstantDiagonal;
336a01b8767SStefano Zampini   A->ops->axpy                      = MatAXPY_ConstantDiagonal;
337a8d4b458SStefano Zampini 
3389566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATCONSTANTDIAGONAL));
3393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3403423f386SBarry Smith }
3413423f386SBarry Smith 
342d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorNumeric_ConstantDiagonal(Mat fact, Mat A, const MatFactorInfo *info)
343d71ae5a4SJacob Faibussowitsch {
3443423f386SBarry Smith   Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal *)A->data, *fctx = (Mat_ConstantDiagonal *)fact->data;
3453423f386SBarry Smith 
3463423f386SBarry Smith   PetscFunctionBegin;
3473423f386SBarry Smith   if (actx->diag == 0.0) fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3483423f386SBarry Smith   else fact->factorerrortype = MAT_FACTOR_NOERROR;
3493423f386SBarry Smith   fctx->diag       = 1.0 / actx->diag;
3503423f386SBarry Smith   fact->ops->solve = MatMult_ConstantDiagonal;
3513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3523423f386SBarry Smith }
3533423f386SBarry Smith 
354d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorSymbolic_LU_ConstantDiagonal(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info)
355d71ae5a4SJacob Faibussowitsch {
3563423f386SBarry Smith   PetscFunctionBegin;
3573423f386SBarry Smith   fact->ops->lufactornumeric = MatFactorNumeric_ConstantDiagonal;
3583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3593423f386SBarry Smith }
3603423f386SBarry Smith 
361d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact, Mat A, IS isrow, const MatFactorInfo *info)
362d71ae5a4SJacob Faibussowitsch {
3633423f386SBarry Smith   PetscFunctionBegin;
364a01b8767SStefano Zampini   fact->ops->choleskyfactornumeric = MatFactorNumeric_ConstantDiagonal;
3653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3663423f386SBarry Smith }
3673423f386SBarry Smith 
368d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_constantdiagonal_petsc(Mat A, MatFactorType ftype, Mat *B)
369d71ae5a4SJacob Faibussowitsch {
3703423f386SBarry Smith   PetscInt n = A->rmap->n, N = A->rmap->N;
3713423f386SBarry Smith 
3723423f386SBarry Smith   PetscFunctionBegin;
3739566063dSJacob Faibussowitsch   PetscCall(MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A), n, n, N, N, 0, B));
3743423f386SBarry Smith 
3753423f386SBarry Smith   (*B)->factortype                  = ftype;
3763423f386SBarry Smith   (*B)->ops->ilufactorsymbolic      = MatFactorSymbolic_LU_ConstantDiagonal;
3773423f386SBarry Smith   (*B)->ops->lufactorsymbolic       = MatFactorSymbolic_LU_ConstantDiagonal;
3783423f386SBarry Smith   (*B)->ops->iccfactorsymbolic      = MatFactorSymbolic_Cholesky_ConstantDiagonal;
3793423f386SBarry Smith   (*B)->ops->choleskyfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal;
3803423f386SBarry Smith 
3813423f386SBarry Smith   (*B)->ops->shift       = NULL;
3823423f386SBarry Smith   (*B)->ops->scale       = NULL;
3833423f386SBarry Smith   (*B)->ops->mult        = NULL;
3843423f386SBarry Smith   (*B)->ops->sor         = NULL;
3853423f386SBarry Smith   (*B)->ops->zeroentries = NULL;
3863423f386SBarry Smith 
3879566063dSJacob Faibussowitsch   PetscCall(PetscFree((*B)->solvertype));
3889566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERPETSC, &(*B)->solvertype));
3893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3903423f386SBarry Smith }
391