xref: /petsc/src/mat/impls/cdiagonal/cdiagonal.c (revision 4bf1a0e8ff80104a73c22b827149902ed229e956)
13423f386SBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
23423f386SBarry Smith 
33423f386SBarry Smith typedef struct {
43423f386SBarry Smith   PetscScalar diag;
53423f386SBarry Smith } Mat_ConstantDiagonal;
63423f386SBarry Smith 
7d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAXPY_ConstantDiagonal(Mat Y, PetscScalar a, Mat X, MatStructure str)
8d71ae5a4SJacob Faibussowitsch {
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;
143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15a01b8767SStefano Zampini }
16a01b8767SStefano Zampini 
17345a4b08SToby Isaac static PetscErrorCode MatEqual_ConstantDiagonal(Mat Y, Mat X, PetscBool *equal)
18345a4b08SToby Isaac {
19345a4b08SToby Isaac   Mat_ConstantDiagonal *yctx = (Mat_ConstantDiagonal *)Y->data;
20345a4b08SToby Isaac   Mat_ConstantDiagonal *xctx = (Mat_ConstantDiagonal *)X->data;
21345a4b08SToby Isaac 
22345a4b08SToby Isaac   PetscFunctionBegin;
23345a4b08SToby Isaac   *equal = (yctx->diag == xctx->diag) ? PETSC_TRUE : PETSC_FALSE;
24345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
25345a4b08SToby Isaac }
26345a4b08SToby Isaac 
27d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatGetRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[])
28d71ae5a4SJacob Faibussowitsch {
29a8d4b458SStefano Zampini   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data;
303423f386SBarry Smith 
31a8d4b458SStefano Zampini   PetscFunctionBegin;
32a8d4b458SStefano Zampini   if (ncols) *ncols = 1;
33a8d4b458SStefano Zampini   if (cols) {
349566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, cols));
35a8d4b458SStefano Zampini     (*cols)[0] = row;
36a8d4b458SStefano Zampini   }
37a8d4b458SStefano Zampini   if (vals) {
389566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(1, vals));
39a8d4b458SStefano Zampini     (*vals)[0] = ctx->diag;
40a8d4b458SStefano Zampini   }
413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42a8d4b458SStefano Zampini }
43a8d4b458SStefano Zampini 
44d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatRestoreRow_ConstantDiagonal(Mat A, PetscInt row, PetscInt *ncols, PetscInt *cols[], PetscScalar *vals[])
45d71ae5a4SJacob Faibussowitsch {
46a8d4b458SStefano Zampini   PetscFunctionBegin;
471baa6e33SBarry Smith   if (cols) PetscCall(PetscFree(*cols));
481baa6e33SBarry Smith   if (vals) PetscCall(PetscFree(*vals));
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50a8d4b458SStefano Zampini }
51a8d4b458SStefano Zampini 
52d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMultAdd_ConstantDiagonal(Mat mat, Vec v1, Vec v2, Vec v3)
53d71ae5a4SJacob Faibussowitsch {
54a8d4b458SStefano Zampini   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)mat->data;
55a8d4b458SStefano Zampini 
56a8d4b458SStefano Zampini   PetscFunctionBegin;
57a8d4b458SStefano Zampini   if (v2 == v3) {
589566063dSJacob Faibussowitsch     PetscCall(VecAXPBY(v3, ctx->diag, 1.0, v1));
59a8d4b458SStefano Zampini   } else {
609566063dSJacob Faibussowitsch     PetscCall(VecAXPBYPCZ(v3, ctx->diag, 1.0, 0.0, v1, v2));
61a8d4b458SStefano Zampini   }
623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
63a8d4b458SStefano Zampini }
64a8d4b458SStefano Zampini 
65345a4b08SToby Isaac static PetscErrorCode MatMultHermitianTransposeAdd_ConstantDiagonal(Mat mat, Vec v1, Vec v2, Vec v3)
66d71ae5a4SJacob Faibussowitsch {
67a8d4b458SStefano Zampini   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)mat->data;
68a8d4b458SStefano Zampini 
69a8d4b458SStefano Zampini   PetscFunctionBegin;
70a8d4b458SStefano Zampini   if (v2 == v3) {
71345a4b08SToby Isaac     PetscCall(VecAXPBY(v3, PetscConj(ctx->diag), 1.0, v1));
72a8d4b458SStefano Zampini   } else {
73345a4b08SToby Isaac     PetscCall(VecAXPBYPCZ(v3, PetscConj(ctx->diag), 1.0, 0.0, v1, v2));
74a8d4b458SStefano Zampini   }
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76a8d4b458SStefano Zampini }
77a8d4b458SStefano Zampini 
78d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatNorm_ConstantDiagonal(Mat A, NormType type, PetscReal *nrm)
79d71ae5a4SJacob Faibussowitsch {
80e5f3c4beSJose E. Roman   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data;
81e5f3c4beSJose E. Roman 
82e5f3c4beSJose E. Roman   PetscFunctionBegin;
83e5f3c4beSJose E. Roman   if (type == NORM_FROBENIUS || type == NORM_2 || type == NORM_1 || type == NORM_INFINITY) *nrm = PetscAbsScalar(ctx->diag);
84e5f3c4beSJose E. Roman   else SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Unsupported norm");
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86e5f3c4beSJose E. Roman }
87e5f3c4beSJose E. Roman 
88d1ca2681SJose E. Roman static PetscErrorCode MatCreateSubMatrices_ConstantDiagonal(Mat A, PetscInt n, const IS irow[], const IS icol[], MatReuse scall, Mat *submat[])
89d1ca2681SJose E. Roman 
90d1ca2681SJose E. Roman {
91d1ca2681SJose E. Roman   Mat B;
92d1ca2681SJose E. Roman 
93d1ca2681SJose E. Roman   PetscFunctionBegin;
949566063dSJacob Faibussowitsch   PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B));
959566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(B, n, irow, icol, scall, submat));
969566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
98d1ca2681SJose E. Roman }
99d1ca2681SJose E. Roman 
100d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDuplicate_ConstantDiagonal(Mat A, MatDuplicateOption op, Mat *B)
101d71ae5a4SJacob Faibussowitsch {
102a8d4b458SStefano Zampini   Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal *)A->data;
103a8d4b458SStefano Zampini 
104a8d4b458SStefano Zampini   PetscFunctionBegin;
1059566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
1069566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
1079566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSizesFromMats(*B, A, A));
1089566063dSJacob Faibussowitsch   PetscCall(MatSetType(*B, MATCONSTANTDIAGONAL));
1099566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->rmap, &(*B)->rmap));
1109566063dSJacob Faibussowitsch   PetscCall(PetscLayoutReference(A->cmap, &(*B)->cmap));
111a8d4b458SStefano Zampini   if (op == MAT_COPY_VALUES) {
112a8d4b458SStefano Zampini     Mat_ConstantDiagonal *bctx = (Mat_ConstantDiagonal *)(*B)->data;
113a8d4b458SStefano Zampini     bctx->diag                 = actx->diag;
114a8d4b458SStefano Zampini   }
1153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
116a8d4b458SStefano Zampini }
117a8d4b458SStefano Zampini 
118d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMissingDiagonal_ConstantDiagonal(Mat mat, PetscBool *missing, PetscInt *dd)
119d71ae5a4SJacob Faibussowitsch {
120a8d4b458SStefano Zampini   PetscFunctionBegin;
121a8d4b458SStefano Zampini   *missing = PETSC_FALSE;
1223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
123a8d4b458SStefano Zampini }
124a8d4b458SStefano Zampini 
125d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_ConstantDiagonal(Mat mat)
126d71ae5a4SJacob Faibussowitsch {
1273423f386SBarry Smith   PetscFunctionBegin;
1289566063dSJacob Faibussowitsch   PetscCall(PetscFree(mat->data));
129345a4b08SToby Isaac   mat->structural_symmetry_eternal = PETSC_FALSE;
130345a4b08SToby Isaac   mat->symmetry_eternal            = PETSC_FALSE;
1313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1323423f386SBarry Smith }
1333423f386SBarry Smith 
134d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_ConstantDiagonal(Mat J, PetscViewer viewer)
135d71ae5a4SJacob Faibussowitsch {
1363423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
1373423f386SBarry Smith   PetscBool             iascii;
1383423f386SBarry Smith 
1393423f386SBarry Smith   PetscFunctionBegin;
1409566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1413423f386SBarry Smith   if (iascii) {
1423423f386SBarry Smith     PetscViewerFormat format;
1433423f386SBarry Smith 
1449566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
1453ba16761SJacob Faibussowitsch     if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(PETSC_SUCCESS);
146c0aa6a63SJacob Faibussowitsch #if defined(PETSC_USE_COMPLEX)
1479566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g + i %g\n", (double)PetscRealPart(ctx->diag), (double)PetscImaginaryPart(ctx->diag)));
1483423f386SBarry Smith #else
149f4f49eeaSPierre Jolivet     PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g\n", (double)ctx->diag));
1503423f386SBarry Smith #endif
1513423f386SBarry Smith   }
1523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1533423f386SBarry Smith }
1543423f386SBarry Smith 
155d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatAssemblyEnd_ConstantDiagonal(Mat J, MatAssemblyType mt)
156d71ae5a4SJacob Faibussowitsch {
1573423f386SBarry Smith   PetscFunctionBegin;
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1593423f386SBarry Smith }
1603423f386SBarry Smith 
161d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_ConstantDiagonal(Mat J, Vec x, Vec y)
162d71ae5a4SJacob Faibussowitsch {
1633423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
1643423f386SBarry Smith 
1653423f386SBarry Smith   PetscFunctionBegin;
1669566063dSJacob Faibussowitsch   PetscCall(VecAXPBY(y, ctx->diag, 0.0, x));
1673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1683423f386SBarry Smith }
1693423f386SBarry Smith 
170345a4b08SToby Isaac static PetscErrorCode MatMultHermitianTranspose_ConstantDiagonal(Mat J, Vec x, Vec y)
171345a4b08SToby Isaac {
172345a4b08SToby Isaac   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
173345a4b08SToby Isaac 
174345a4b08SToby Isaac   PetscFunctionBegin;
175345a4b08SToby Isaac   PetscCall(VecAXPBY(y, PetscConj(ctx->diag), 0.0, x));
176345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
177345a4b08SToby Isaac }
178345a4b08SToby Isaac 
179345a4b08SToby Isaac static PetscErrorCode MatGetDiagonal_ConstantDiagonal(Mat J, Vec x)
180d71ae5a4SJacob Faibussowitsch {
1813423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
1823423f386SBarry Smith 
1833423f386SBarry Smith   PetscFunctionBegin;
1849566063dSJacob Faibussowitsch   PetscCall(VecSet(x, ctx->diag));
1853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1863423f386SBarry Smith }
1873423f386SBarry Smith 
188d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatShift_ConstantDiagonal(Mat Y, PetscScalar a)
189d71ae5a4SJacob Faibussowitsch {
1903423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data;
1913423f386SBarry Smith 
1923423f386SBarry Smith   PetscFunctionBegin;
1933423f386SBarry Smith   ctx->diag += a;
1943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1953423f386SBarry Smith }
1963423f386SBarry Smith 
197d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_ConstantDiagonal(Mat Y, PetscScalar a)
198d71ae5a4SJacob Faibussowitsch {
1993423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data;
2003423f386SBarry Smith 
2013423f386SBarry Smith   PetscFunctionBegin;
2023423f386SBarry Smith   ctx->diag *= a;
2033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2043423f386SBarry Smith }
2053423f386SBarry Smith 
206d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatZeroEntries_ConstantDiagonal(Mat Y)
207d71ae5a4SJacob Faibussowitsch {
2083423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data;
2093423f386SBarry Smith 
2103423f386SBarry Smith   PetscFunctionBegin;
2113423f386SBarry Smith   ctx->diag = 0.0;
2123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2133423f386SBarry Smith }
2143423f386SBarry Smith 
215345a4b08SToby Isaac static PetscErrorCode MatSolve_ConstantDiagonal(Mat matin, Vec b, Vec x)
216d71ae5a4SJacob Faibussowitsch {
2173423f386SBarry Smith   Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)matin->data;
2183423f386SBarry Smith 
2193423f386SBarry Smith   PetscFunctionBegin;
2203423f386SBarry Smith   if (ctx->diag == 0.0) matin->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
2213423f386SBarry Smith   else matin->factorerrortype = MAT_FACTOR_NOERROR;
222345a4b08SToby Isaac   PetscCall(VecAXPBY(x, 1.0 / ctx->diag, 0.0, b));
223345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
224345a4b08SToby Isaac }
225345a4b08SToby Isaac 
22666976f2fSJacob Faibussowitsch static PetscErrorCode MatSOR_ConstantDiagonal(Mat matin, Vec x, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec y)
227345a4b08SToby Isaac {
228345a4b08SToby Isaac   PetscFunctionBegin;
229345a4b08SToby Isaac   PetscCall(MatSolve_ConstantDiagonal(matin, x, y));
2303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2313423f386SBarry Smith }
2323423f386SBarry Smith 
23366976f2fSJacob Faibussowitsch static PetscErrorCode MatGetInfo_ConstantDiagonal(Mat A, MatInfoType flag, MatInfo *info)
234d71ae5a4SJacob Faibussowitsch {
2353423f386SBarry Smith   PetscFunctionBegin;
2363423f386SBarry Smith   info->block_size   = 1.0;
2373423f386SBarry Smith   info->nz_allocated = 1.0;
2383423f386SBarry Smith   info->nz_used      = 1.0;
2393423f386SBarry Smith   info->nz_unneeded  = 0.0;
2403966268fSBarry Smith   info->assemblies   = A->num_ass;
2413423f386SBarry Smith   info->mallocs      = 0.0;
2424dfa11a4SJacob Faibussowitsch   info->memory       = 0; /* REVIEW ME */
2433423f386SBarry Smith   if (A->factortype) {
2443423f386SBarry Smith     info->fill_ratio_given  = 1.0;
2453423f386SBarry Smith     info->fill_ratio_needed = 1.0;
2463423f386SBarry Smith     info->factor_mallocs    = 0.0;
2473423f386SBarry Smith   } else {
2483423f386SBarry Smith     info->fill_ratio_given  = 0;
2493423f386SBarry Smith     info->fill_ratio_needed = 0;
2503423f386SBarry Smith     info->factor_mallocs    = 0;
2513423f386SBarry Smith   }
2523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2533423f386SBarry Smith }
2543423f386SBarry Smith 
2553423f386SBarry Smith /*@
2563423f386SBarry Smith   MatCreateConstantDiagonal - Creates a matrix with a uniform value along the diagonal
2573423f386SBarry Smith 
258d083f849SBarry Smith   Collective
2593423f386SBarry Smith 
2603423f386SBarry Smith   Input Parameters:
2613423f386SBarry Smith + comm - MPI communicator
2622ef1f0ffSBarry Smith . m    - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
2633423f386SBarry Smith            This value should be the same as the local size used in creating the
2643423f386SBarry Smith            y vector for the matrix-vector product y = Ax.
2653423f386SBarry Smith . n    - This value should be the same as the local size used in creating the
2662ef1f0ffSBarry Smith        x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
2672ef1f0ffSBarry Smith        calculated if `N` is given) For square matrices n is almost always `m`.
26811a5261eSBarry Smith . M    - number of global rows (or `PETSC_DETERMINE` to have calculated if m is given)
26911a5261eSBarry Smith . N    - number of global columns (or `PETSC_DETERMINE` to have calculated if n is given)
2703423f386SBarry Smith - diag - the diagonal value
2713423f386SBarry Smith 
2723423f386SBarry Smith   Output Parameter:
2733423f386SBarry Smith . J - the diagonal matrix
2743423f386SBarry Smith 
2753423f386SBarry Smith   Level: advanced
2763423f386SBarry Smith 
2773423f386SBarry Smith   Notes:
2783423f386SBarry Smith   Only supports square matrices with the same number of local rows and columns
2793423f386SBarry Smith 
2801cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MATCONSTANTDIAGONAL`, `MatScale()`, `MatShift()`, `MatMult()`, `MatGetDiagonal()`, `MatGetFactor()`, `MatSolve()`
2813423f386SBarry Smith @*/
282d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateConstantDiagonal(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscScalar diag, Mat *J)
283d71ae5a4SJacob Faibussowitsch {
2843423f386SBarry Smith   PetscFunctionBegin;
2859566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, J));
2869566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, m, n, M, N));
2879566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, MATCONSTANTDIAGONAL));
2889566063dSJacob Faibussowitsch   PetscCall(MatShift(*J, diag));
2899566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*J));
2903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2913423f386SBarry Smith }
2923423f386SBarry Smith 
293*4bf1a0e8SHansol  Suh /*MC
294*4bf1a0e8SHansol  Suh    MATCONSTANTDIAGONAL - "constant-diagonal" - A diagonal matrix type with a uniform value
295*4bf1a0e8SHansol  Suh    along the diagonal.
296*4bf1a0e8SHansol  Suh 
297*4bf1a0e8SHansol  Suh   Level: advanced
298*4bf1a0e8SHansol  Suh 
299*4bf1a0e8SHansol  Suh .seealso: [](ch_matrices), `Mat`, `MatCreateConstantDiagonal()`
300*4bf1a0e8SHansol  Suh M*/
301d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_ConstantDiagonal(Mat A)
302d71ae5a4SJacob Faibussowitsch {
3033423f386SBarry Smith   Mat_ConstantDiagonal *ctx;
3043423f386SBarry Smith 
3053423f386SBarry Smith   PetscFunctionBegin;
3069566063dSJacob Faibussowitsch   PetscCall(PetscNew(&ctx));
3073423f386SBarry Smith   ctx->diag = 0.0;
3083423f386SBarry Smith   A->data   = (void *)ctx;
3093423f386SBarry Smith 
3103423f386SBarry Smith   A->assembled                   = PETSC_TRUE;
3113423f386SBarry Smith   A->preallocated                = PETSC_TRUE;
312345a4b08SToby Isaac   A->structurally_symmetric      = PETSC_BOOL3_TRUE;
313345a4b08SToby Isaac   A->structural_symmetry_eternal = PETSC_TRUE;
314345a4b08SToby Isaac   A->symmetric                   = PETSC_BOOL3_TRUE;
315345a4b08SToby Isaac   if (!PetscDefined(USE_COMPLEX)) A->hermitian = PETSC_BOOL3_TRUE;
316345a4b08SToby Isaac   A->symmetry_eternal = PETSC_TRUE;
317a8d4b458SStefano Zampini 
3183423f386SBarry Smith   A->ops->mult                      = MatMult_ConstantDiagonal;
319a8d4b458SStefano Zampini   A->ops->multadd                   = MatMultAdd_ConstantDiagonal;
320345a4b08SToby Isaac   A->ops->multtranspose             = MatMult_ConstantDiagonal;
321345a4b08SToby Isaac   A->ops->multtransposeadd          = MatMultAdd_ConstantDiagonal;
322345a4b08SToby Isaac   A->ops->multhermitiantranspose    = MatMultHermitianTranspose_ConstantDiagonal;
323345a4b08SToby Isaac   A->ops->multhermitiantransposeadd = MatMultHermitianTransposeAdd_ConstantDiagonal;
324345a4b08SToby Isaac   A->ops->solve                     = MatSolve_ConstantDiagonal;
325345a4b08SToby Isaac   A->ops->solvetranspose            = MatSolve_ConstantDiagonal;
326e5f3c4beSJose E. Roman   A->ops->norm                      = MatNorm_ConstantDiagonal;
327d1ca2681SJose E. Roman   A->ops->createsubmatrices         = MatCreateSubMatrices_ConstantDiagonal;
328a8d4b458SStefano Zampini   A->ops->duplicate                 = MatDuplicate_ConstantDiagonal;
329a8d4b458SStefano Zampini   A->ops->missingdiagonal           = MatMissingDiagonal_ConstantDiagonal;
330a8d4b458SStefano Zampini   A->ops->getrow                    = MatGetRow_ConstantDiagonal;
331a8d4b458SStefano Zampini   A->ops->restorerow                = MatRestoreRow_ConstantDiagonal;
3323423f386SBarry Smith   A->ops->sor                       = MatSOR_ConstantDiagonal;
3333423f386SBarry Smith   A->ops->shift                     = MatShift_ConstantDiagonal;
3343423f386SBarry Smith   A->ops->scale                     = MatScale_ConstantDiagonal;
3353423f386SBarry Smith   A->ops->getdiagonal               = MatGetDiagonal_ConstantDiagonal;
3363423f386SBarry Smith   A->ops->view                      = MatView_ConstantDiagonal;
3373423f386SBarry Smith   A->ops->zeroentries               = MatZeroEntries_ConstantDiagonal;
3383423f386SBarry Smith   A->ops->assemblyend               = MatAssemblyEnd_ConstantDiagonal;
3393423f386SBarry Smith   A->ops->destroy                   = MatDestroy_ConstantDiagonal;
3403423f386SBarry Smith   A->ops->getinfo                   = MatGetInfo_ConstantDiagonal;
341345a4b08SToby Isaac   A->ops->equal                     = MatEqual_ConstantDiagonal;
342a01b8767SStefano Zampini   A->ops->axpy                      = MatAXPY_ConstantDiagonal;
343a8d4b458SStefano Zampini 
3449566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATCONSTANTDIAGONAL));
3453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3463423f386SBarry Smith }
3473423f386SBarry Smith 
348d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorNumeric_ConstantDiagonal(Mat fact, Mat A, const MatFactorInfo *info)
349d71ae5a4SJacob Faibussowitsch {
3503423f386SBarry Smith   Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal *)A->data, *fctx = (Mat_ConstantDiagonal *)fact->data;
3513423f386SBarry Smith 
3523423f386SBarry Smith   PetscFunctionBegin;
3533423f386SBarry Smith   if (actx->diag == 0.0) fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
3543423f386SBarry Smith   else fact->factorerrortype = MAT_FACTOR_NOERROR;
3553423f386SBarry Smith   fctx->diag       = 1.0 / actx->diag;
3563423f386SBarry Smith   fact->ops->solve = MatMult_ConstantDiagonal;
3573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3583423f386SBarry Smith }
3593423f386SBarry Smith 
360d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorSymbolic_LU_ConstantDiagonal(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info)
361d71ae5a4SJacob Faibussowitsch {
3623423f386SBarry Smith   PetscFunctionBegin;
3633423f386SBarry Smith   fact->ops->lufactornumeric = MatFactorNumeric_ConstantDiagonal;
3643ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3653423f386SBarry Smith }
3663423f386SBarry Smith 
367d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact, Mat A, IS isrow, const MatFactorInfo *info)
368d71ae5a4SJacob Faibussowitsch {
3693423f386SBarry Smith   PetscFunctionBegin;
370a01b8767SStefano Zampini   fact->ops->choleskyfactornumeric = MatFactorNumeric_ConstantDiagonal;
3713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3723423f386SBarry Smith }
3733423f386SBarry Smith 
374d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_constantdiagonal_petsc(Mat A, MatFactorType ftype, Mat *B)
375d71ae5a4SJacob Faibussowitsch {
3763423f386SBarry Smith   PetscInt n = A->rmap->n, N = A->rmap->N;
3773423f386SBarry Smith 
3783423f386SBarry Smith   PetscFunctionBegin;
3799566063dSJacob Faibussowitsch   PetscCall(MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A), n, n, N, N, 0, B));
3803423f386SBarry Smith 
3813423f386SBarry Smith   (*B)->factortype                  = ftype;
3823423f386SBarry Smith   (*B)->ops->ilufactorsymbolic      = MatFactorSymbolic_LU_ConstantDiagonal;
3833423f386SBarry Smith   (*B)->ops->lufactorsymbolic       = MatFactorSymbolic_LU_ConstantDiagonal;
3843423f386SBarry Smith   (*B)->ops->iccfactorsymbolic      = MatFactorSymbolic_Cholesky_ConstantDiagonal;
3853423f386SBarry Smith   (*B)->ops->choleskyfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal;
3863423f386SBarry Smith 
3873423f386SBarry Smith   (*B)->ops->shift       = NULL;
3883423f386SBarry Smith   (*B)->ops->scale       = NULL;
3893423f386SBarry Smith   (*B)->ops->mult        = NULL;
3903423f386SBarry Smith   (*B)->ops->sor         = NULL;
3913423f386SBarry Smith   (*B)->ops->zeroentries = NULL;
3923423f386SBarry Smith 
3939566063dSJacob Faibussowitsch   PetscCall(PetscFree((*B)->solvertype));
3949566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERPETSC, &(*B)->solvertype));
3953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3963423f386SBarry Smith }
397