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
MatAXPY_ConstantDiagonal(Mat Y,PetscScalar a,Mat X,MatStructure str)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
MatEqual_ConstantDiagonal(Mat Y,Mat X,PetscBool * equal)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
MatGetRow_ConstantDiagonal(Mat A,PetscInt row,PetscInt * ncols,PetscInt * cols[],PetscScalar * vals[])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
MatRestoreRow_ConstantDiagonal(Mat A,PetscInt row,PetscInt * ncols,PetscInt * cols[],PetscScalar * vals[])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
MatMultAdd_ConstantDiagonal(Mat mat,Vec v1,Vec v2,Vec v3)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
MatMultHermitianTransposeAdd_ConstantDiagonal(Mat mat,Vec v1,Vec v2,Vec v3)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
MatNorm_ConstantDiagonal(Mat A,NormType type,PetscReal * nrm)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;
83*966bd95aSPierre Jolivet PetscCheck(type == NORM_FROBENIUS || type == NORM_2 || type == NORM_1 || type == NORM_INFINITY, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Unsupported norm");
84*966bd95aSPierre Jolivet *nrm = PetscAbsScalar(ctx->diag);
853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
86e5f3c4beSJose E. Roman }
87e5f3c4beSJose E. Roman
MatCreateSubMatrices_ConstantDiagonal(Mat A,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat * submat[])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 Mat B;
91d1ca2681SJose E. Roman
92d1ca2681SJose E. Roman PetscFunctionBegin;
939566063dSJacob Faibussowitsch PetscCall(MatConvert(A, MATAIJ, MAT_INITIAL_MATRIX, &B));
949566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B, n, irow, icol, scall, submat));
959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B));
963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
97d1ca2681SJose E. Roman }
98d1ca2681SJose E. Roman
MatDuplicate_ConstantDiagonal(Mat A,MatDuplicateOption op,Mat * B)99d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDuplicate_ConstantDiagonal(Mat A, MatDuplicateOption op, Mat *B)
100d71ae5a4SJacob Faibussowitsch {
101a8d4b458SStefano Zampini Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal *)A->data;
102a8d4b458SStefano Zampini
103a8d4b458SStefano Zampini PetscFunctionBegin;
1049566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), B));
1059566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*B, A->rmap->n, A->cmap->n, A->rmap->N, A->cmap->N));
1069566063dSJacob Faibussowitsch PetscCall(MatSetBlockSizesFromMats(*B, A, A));
1079566063dSJacob Faibussowitsch PetscCall(MatSetType(*B, MATCONSTANTDIAGONAL));
1089566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->rmap, &(*B)->rmap));
1099566063dSJacob Faibussowitsch PetscCall(PetscLayoutReference(A->cmap, &(*B)->cmap));
110a8d4b458SStefano Zampini if (op == MAT_COPY_VALUES) {
111a8d4b458SStefano Zampini Mat_ConstantDiagonal *bctx = (Mat_ConstantDiagonal *)(*B)->data;
112a8d4b458SStefano Zampini bctx->diag = actx->diag;
113a8d4b458SStefano Zampini }
1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
115a8d4b458SStefano Zampini }
116a8d4b458SStefano Zampini
MatDestroy_ConstantDiagonal(Mat mat)117d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatDestroy_ConstantDiagonal(Mat mat)
118d71ae5a4SJacob Faibussowitsch {
1193423f386SBarry Smith PetscFunctionBegin;
1209566063dSJacob Faibussowitsch PetscCall(PetscFree(mat->data));
1218a1df862SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConstantDiagonalGetConstant_C", NULL));
1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1233423f386SBarry Smith }
1243423f386SBarry Smith
MatView_ConstantDiagonal(Mat J,PetscViewer viewer)125d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_ConstantDiagonal(Mat J, PetscViewer viewer)
126d71ae5a4SJacob Faibussowitsch {
1273423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
1289f196a02SMartin Diehl PetscBool isascii;
1293423f386SBarry Smith
1303423f386SBarry Smith PetscFunctionBegin;
1319f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1329f196a02SMartin Diehl if (isascii) {
1333423f386SBarry Smith PetscViewerFormat format;
1343423f386SBarry Smith
1359566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format));
1363ba16761SJacob Faibussowitsch if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(PETSC_SUCCESS);
1378a1df862SToby Isaac if (PetscImaginaryPart(ctx->diag) == 0) {
1388a1df862SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g\n", (double)PetscRealPart(ctx->diag)));
1398a1df862SToby Isaac } else {
1409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g + i %g\n", (double)PetscRealPart(ctx->diag), (double)PetscImaginaryPart(ctx->diag)));
1418a1df862SToby Isaac }
1423423f386SBarry Smith }
1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1443423f386SBarry Smith }
1453423f386SBarry Smith
MatMult_ConstantDiagonal(Mat J,Vec x,Vec y)146d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_ConstantDiagonal(Mat J, Vec x, Vec y)
147d71ae5a4SJacob Faibussowitsch {
1483423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
1493423f386SBarry Smith
1503423f386SBarry Smith PetscFunctionBegin;
1519566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y, ctx->diag, 0.0, x));
1523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1533423f386SBarry Smith }
1543423f386SBarry Smith
MatMultHermitianTranspose_ConstantDiagonal(Mat J,Vec x,Vec y)155345a4b08SToby Isaac static PetscErrorCode MatMultHermitianTranspose_ConstantDiagonal(Mat J, Vec x, Vec y)
156345a4b08SToby Isaac {
157345a4b08SToby Isaac Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
158345a4b08SToby Isaac
159345a4b08SToby Isaac PetscFunctionBegin;
160345a4b08SToby Isaac PetscCall(VecAXPBY(y, PetscConj(ctx->diag), 0.0, x));
161345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
162345a4b08SToby Isaac }
163345a4b08SToby Isaac
MatGetDiagonal_ConstantDiagonal(Mat J,Vec x)164345a4b08SToby Isaac static PetscErrorCode MatGetDiagonal_ConstantDiagonal(Mat J, Vec x)
165d71ae5a4SJacob Faibussowitsch {
1663423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data;
1673423f386SBarry Smith
1683423f386SBarry Smith PetscFunctionBegin;
1699566063dSJacob Faibussowitsch PetscCall(VecSet(x, ctx->diag));
1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1713423f386SBarry Smith }
1723423f386SBarry Smith
MatShift_ConstantDiagonal(Mat Y,PetscScalar a)173d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatShift_ConstantDiagonal(Mat Y, PetscScalar a)
174d71ae5a4SJacob Faibussowitsch {
1753423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data;
1763423f386SBarry Smith
1773423f386SBarry Smith PetscFunctionBegin;
1783423f386SBarry Smith ctx->diag += a;
1793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1803423f386SBarry Smith }
1813423f386SBarry Smith
MatScale_ConstantDiagonal(Mat Y,PetscScalar a)182d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_ConstantDiagonal(Mat Y, PetscScalar a)
183d71ae5a4SJacob Faibussowitsch {
1843423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data;
1853423f386SBarry Smith
1863423f386SBarry Smith PetscFunctionBegin;
1873423f386SBarry Smith ctx->diag *= a;
1883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1893423f386SBarry Smith }
1903423f386SBarry Smith
MatZeroEntries_ConstantDiagonal(Mat Y)191d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatZeroEntries_ConstantDiagonal(Mat Y)
192d71ae5a4SJacob Faibussowitsch {
1933423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data;
1943423f386SBarry Smith
1953423f386SBarry Smith PetscFunctionBegin;
1963423f386SBarry Smith ctx->diag = 0.0;
1973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1983423f386SBarry Smith }
1993423f386SBarry Smith
MatConjugate_ConstantDiagonal(Mat Y)2008a1df862SToby Isaac static PetscErrorCode MatConjugate_ConstantDiagonal(Mat Y)
2018a1df862SToby Isaac {
2028a1df862SToby Isaac Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data;
2038a1df862SToby Isaac
2048a1df862SToby Isaac PetscFunctionBegin;
2058a1df862SToby Isaac ctx->diag = PetscConj(ctx->diag);
2068a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
2078a1df862SToby Isaac }
2088a1df862SToby Isaac
MatTranspose_ConstantDiagonal(Mat A,MatReuse reuse,Mat * matout)2098a1df862SToby Isaac static PetscErrorCode MatTranspose_ConstantDiagonal(Mat A, MatReuse reuse, Mat *matout)
2108a1df862SToby Isaac {
2118a1df862SToby Isaac Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data;
2128a1df862SToby Isaac
2138a1df862SToby Isaac PetscFunctionBegin;
2148a1df862SToby Isaac if (reuse == MAT_INPLACE_MATRIX) {
2158a1df862SToby Isaac PetscLayout tmplayout = A->rmap;
2168a1df862SToby Isaac
2178a1df862SToby Isaac A->rmap = A->cmap;
2188a1df862SToby Isaac A->cmap = tmplayout;
2198a1df862SToby Isaac } else {
2208a1df862SToby Isaac if (reuse == MAT_INITIAL_MATRIX) {
2218a1df862SToby Isaac PetscCall(MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A), A->cmap->n, A->rmap->n, A->cmap->N, A->rmap->N, ctx->diag, matout));
2228a1df862SToby Isaac } else {
2238a1df862SToby Isaac PetscCall(MatZeroEntries(*matout));
2248a1df862SToby Isaac PetscCall(MatShift(*matout, ctx->diag));
2258a1df862SToby Isaac }
2268a1df862SToby Isaac }
2278a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
2288a1df862SToby Isaac }
2298a1df862SToby Isaac
MatSetRandom_ConstantDiagonal(Mat A,PetscRandom rand)2308a1df862SToby Isaac static PetscErrorCode MatSetRandom_ConstantDiagonal(Mat A, PetscRandom rand)
2318a1df862SToby Isaac {
2328a1df862SToby Isaac PetscMPIInt rank;
2338a1df862SToby Isaac MPI_Comm comm;
2348a1df862SToby Isaac PetscScalar v = 0.0;
2358a1df862SToby Isaac Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data;
2368a1df862SToby Isaac
2378a1df862SToby Isaac PetscFunctionBegin;
2388a1df862SToby Isaac PetscCall(PetscObjectGetComm((PetscObject)A, &comm));
2398a1df862SToby Isaac PetscCallMPI(MPI_Comm_rank(comm, &rank));
2408a1df862SToby Isaac if (!rank) PetscCall(PetscRandomGetValue(rand, &v));
2418a1df862SToby Isaac PetscCallMPI(MPI_Bcast(&v, 1, MPIU_SCALAR, 0, comm));
2428a1df862SToby Isaac ctx->diag = v;
2438a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
2448a1df862SToby Isaac }
2458a1df862SToby Isaac
MatSolve_ConstantDiagonal(Mat matin,Vec b,Vec x)246345a4b08SToby Isaac static PetscErrorCode MatSolve_ConstantDiagonal(Mat matin, Vec b, Vec x)
247d71ae5a4SJacob Faibussowitsch {
2483423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)matin->data;
2493423f386SBarry Smith
2503423f386SBarry Smith PetscFunctionBegin;
2513423f386SBarry Smith if (ctx->diag == 0.0) matin->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
2523423f386SBarry Smith else matin->factorerrortype = MAT_FACTOR_NOERROR;
253345a4b08SToby Isaac PetscCall(VecAXPBY(x, 1.0 / ctx->diag, 0.0, b));
254345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
255345a4b08SToby Isaac }
256345a4b08SToby Isaac
MatSOR_ConstantDiagonal(Mat matin,Vec x,PetscReal omega,MatSORType flag,PetscReal fshift,PetscInt its,PetscInt lits,Vec y)25766976f2fSJacob Faibussowitsch static PetscErrorCode MatSOR_ConstantDiagonal(Mat matin, Vec x, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec y)
258345a4b08SToby Isaac {
259345a4b08SToby Isaac PetscFunctionBegin;
260345a4b08SToby Isaac PetscCall(MatSolve_ConstantDiagonal(matin, x, y));
2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2623423f386SBarry Smith }
2633423f386SBarry Smith
MatGetInfo_ConstantDiagonal(Mat A,MatInfoType flag,MatInfo * info)26466976f2fSJacob Faibussowitsch static PetscErrorCode MatGetInfo_ConstantDiagonal(Mat A, MatInfoType flag, MatInfo *info)
265d71ae5a4SJacob Faibussowitsch {
2663423f386SBarry Smith PetscFunctionBegin;
2673423f386SBarry Smith info->block_size = 1.0;
2683423f386SBarry Smith info->nz_allocated = 1.0;
2693423f386SBarry Smith info->nz_used = 1.0;
2703423f386SBarry Smith info->nz_unneeded = 0.0;
2713966268fSBarry Smith info->assemblies = A->num_ass;
2723423f386SBarry Smith info->mallocs = 0.0;
2734dfa11a4SJacob Faibussowitsch info->memory = 0; /* REVIEW ME */
2743423f386SBarry Smith if (A->factortype) {
2753423f386SBarry Smith info->fill_ratio_given = 1.0;
2763423f386SBarry Smith info->fill_ratio_needed = 1.0;
2773423f386SBarry Smith info->factor_mallocs = 0.0;
2783423f386SBarry Smith } else {
2793423f386SBarry Smith info->fill_ratio_given = 0;
2803423f386SBarry Smith info->fill_ratio_needed = 0;
2813423f386SBarry Smith info->factor_mallocs = 0;
2823423f386SBarry Smith }
2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2843423f386SBarry Smith }
2853423f386SBarry Smith
2863423f386SBarry Smith /*@
2873423f386SBarry Smith MatCreateConstantDiagonal - Creates a matrix with a uniform value along the diagonal
2883423f386SBarry Smith
289d083f849SBarry Smith Collective
2903423f386SBarry Smith
2913423f386SBarry Smith Input Parameters:
2923423f386SBarry Smith + comm - MPI communicator
2932ef1f0ffSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given)
2943423f386SBarry Smith This value should be the same as the local size used in creating the
2953423f386SBarry Smith y vector for the matrix-vector product y = Ax.
2963423f386SBarry Smith . n - This value should be the same as the local size used in creating the
2972ef1f0ffSBarry Smith x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have
2982ef1f0ffSBarry Smith calculated if `N` is given) For square matrices n is almost always `m`.
29911a5261eSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if m is given)
30011a5261eSBarry Smith . N - number of global columns (or `PETSC_DETERMINE` to have calculated if n is given)
3013423f386SBarry Smith - diag - the diagonal value
3023423f386SBarry Smith
3033423f386SBarry Smith Output Parameter:
3043423f386SBarry Smith . J - the diagonal matrix
3053423f386SBarry Smith
3063423f386SBarry Smith Level: advanced
3073423f386SBarry Smith
3083423f386SBarry Smith Notes:
3093423f386SBarry Smith Only supports square matrices with the same number of local rows and columns
3103423f386SBarry Smith
3111cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MATCONSTANTDIAGONAL`, `MatScale()`, `MatShift()`, `MatMult()`, `MatGetDiagonal()`, `MatGetFactor()`, `MatSolve()`
3123423f386SBarry Smith @*/
MatCreateConstantDiagonal(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscScalar diag,Mat * J)313d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateConstantDiagonal(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscScalar diag, Mat *J)
314d71ae5a4SJacob Faibussowitsch {
3153423f386SBarry Smith PetscFunctionBegin;
3169566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, J));
3179566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, m, n, M, N));
3189566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, MATCONSTANTDIAGONAL));
3199566063dSJacob Faibussowitsch PetscCall(MatShift(*J, diag));
3209566063dSJacob Faibussowitsch PetscCall(MatSetUp(*J));
3213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3223423f386SBarry Smith }
3233423f386SBarry Smith
3248a1df862SToby Isaac /*@
3258a1df862SToby Isaac MatConstantDiagonalGetConstant - Get the scalar constant of a constant diagonal matrix
3268a1df862SToby Isaac
3278a1df862SToby Isaac Not collective
3288a1df862SToby Isaac
3298a1df862SToby Isaac Input Parameter:
3308a1df862SToby Isaac . mat - a `MATCONSTANTDIAGONAL`
3318a1df862SToby Isaac
3328a1df862SToby Isaac Output Parameter:
3338a1df862SToby Isaac . value - the scalar value
3348a1df862SToby Isaac
3358a1df862SToby Isaac Level: developer
3368a1df862SToby Isaac
3378a1df862SToby Isaac .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MATCONSTANTDIAGONAL`
3388a1df862SToby Isaac @*/
MatConstantDiagonalGetConstant(Mat mat,PetscScalar * value)3398a1df862SToby Isaac PetscErrorCode MatConstantDiagonalGetConstant(Mat mat, PetscScalar *value)
3408a1df862SToby Isaac {
3418a1df862SToby Isaac PetscFunctionBegin;
3428a1df862SToby Isaac PetscUseMethod(mat, "MatConstantDiagonalGetConstant_C", (Mat, PetscScalar *), (mat, value));
3438a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
3448a1df862SToby Isaac }
3458a1df862SToby Isaac
MatConstantDiagonalGetConstant_ConstantDiagonal(Mat mat,PetscScalar * value)3468a1df862SToby Isaac static PetscErrorCode MatConstantDiagonalGetConstant_ConstantDiagonal(Mat mat, PetscScalar *value)
3478a1df862SToby Isaac {
3488a1df862SToby Isaac Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)mat->data;
3498a1df862SToby Isaac
3508a1df862SToby Isaac PetscFunctionBegin;
3518a1df862SToby Isaac *value = ctx->diag;
3528a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS);
3538a1df862SToby Isaac }
3548a1df862SToby Isaac
3554bf1a0e8SHansol Suh /*MC
3564bf1a0e8SHansol Suh MATCONSTANTDIAGONAL - "constant-diagonal" - A diagonal matrix type with a uniform value
3574bf1a0e8SHansol Suh along the diagonal.
3584bf1a0e8SHansol Suh
3594bf1a0e8SHansol Suh Level: advanced
3604bf1a0e8SHansol Suh
3614bf1a0e8SHansol Suh .seealso: [](ch_matrices), `Mat`, `MatCreateConstantDiagonal()`
3624bf1a0e8SHansol Suh M*/
MatCreate_ConstantDiagonal(Mat A)363d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_ConstantDiagonal(Mat A)
364d71ae5a4SJacob Faibussowitsch {
3653423f386SBarry Smith Mat_ConstantDiagonal *ctx;
3663423f386SBarry Smith
3673423f386SBarry Smith PetscFunctionBegin;
3689566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx));
3693423f386SBarry Smith ctx->diag = 0.0;
3703423f386SBarry Smith A->data = (void *)ctx;
3713423f386SBarry Smith
3723423f386SBarry Smith A->assembled = PETSC_TRUE;
3733423f386SBarry Smith A->preallocated = PETSC_TRUE;
374345a4b08SToby Isaac A->structurally_symmetric = PETSC_BOOL3_TRUE;
375345a4b08SToby Isaac A->structural_symmetry_eternal = PETSC_TRUE;
376345a4b08SToby Isaac A->symmetric = PETSC_BOOL3_TRUE;
377345a4b08SToby Isaac if (!PetscDefined(USE_COMPLEX)) A->hermitian = PETSC_BOOL3_TRUE;
378345a4b08SToby Isaac A->symmetry_eternal = PETSC_TRUE;
379a8d4b458SStefano Zampini
3803423f386SBarry Smith A->ops->mult = MatMult_ConstantDiagonal;
381a8d4b458SStefano Zampini A->ops->multadd = MatMultAdd_ConstantDiagonal;
382345a4b08SToby Isaac A->ops->multtranspose = MatMult_ConstantDiagonal;
383345a4b08SToby Isaac A->ops->multtransposeadd = MatMultAdd_ConstantDiagonal;
384345a4b08SToby Isaac A->ops->multhermitiantranspose = MatMultHermitianTranspose_ConstantDiagonal;
385345a4b08SToby Isaac A->ops->multhermitiantransposeadd = MatMultHermitianTransposeAdd_ConstantDiagonal;
386345a4b08SToby Isaac A->ops->solve = MatSolve_ConstantDiagonal;
387345a4b08SToby Isaac A->ops->solvetranspose = MatSolve_ConstantDiagonal;
388e5f3c4beSJose E. Roman A->ops->norm = MatNorm_ConstantDiagonal;
389d1ca2681SJose E. Roman A->ops->createsubmatrices = MatCreateSubMatrices_ConstantDiagonal;
390a8d4b458SStefano Zampini A->ops->duplicate = MatDuplicate_ConstantDiagonal;
391a8d4b458SStefano Zampini A->ops->getrow = MatGetRow_ConstantDiagonal;
392a8d4b458SStefano Zampini A->ops->restorerow = MatRestoreRow_ConstantDiagonal;
3933423f386SBarry Smith A->ops->sor = MatSOR_ConstantDiagonal;
3943423f386SBarry Smith A->ops->shift = MatShift_ConstantDiagonal;
3953423f386SBarry Smith A->ops->scale = MatScale_ConstantDiagonal;
3963423f386SBarry Smith A->ops->getdiagonal = MatGetDiagonal_ConstantDiagonal;
3973423f386SBarry Smith A->ops->view = MatView_ConstantDiagonal;
3983423f386SBarry Smith A->ops->zeroentries = MatZeroEntries_ConstantDiagonal;
3993423f386SBarry Smith A->ops->destroy = MatDestroy_ConstantDiagonal;
4003423f386SBarry Smith A->ops->getinfo = MatGetInfo_ConstantDiagonal;
401345a4b08SToby Isaac A->ops->equal = MatEqual_ConstantDiagonal;
402a01b8767SStefano Zampini A->ops->axpy = MatAXPY_ConstantDiagonal;
4038a1df862SToby Isaac A->ops->setrandom = MatSetRandom_ConstantDiagonal;
4048a1df862SToby Isaac A->ops->conjugate = MatConjugate_ConstantDiagonal;
4058a1df862SToby Isaac A->ops->transpose = MatTranspose_ConstantDiagonal;
406a8d4b458SStefano Zampini
4079566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATCONSTANTDIAGONAL));
4088a1df862SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConstantDiagonalGetConstant_C", MatConstantDiagonalGetConstant_ConstantDiagonal));
4093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4103423f386SBarry Smith }
4113423f386SBarry Smith
MatFactorNumeric_ConstantDiagonal(Mat fact,Mat A,const MatFactorInfo * info)412d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorNumeric_ConstantDiagonal(Mat fact, Mat A, const MatFactorInfo *info)
413d71ae5a4SJacob Faibussowitsch {
4143423f386SBarry Smith Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal *)A->data, *fctx = (Mat_ConstantDiagonal *)fact->data;
4153423f386SBarry Smith
4163423f386SBarry Smith PetscFunctionBegin;
4173423f386SBarry Smith if (actx->diag == 0.0) fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
4183423f386SBarry Smith else fact->factorerrortype = MAT_FACTOR_NOERROR;
4193423f386SBarry Smith fctx->diag = 1.0 / actx->diag;
4203423f386SBarry Smith fact->ops->solve = MatMult_ConstantDiagonal;
4213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4223423f386SBarry Smith }
4233423f386SBarry Smith
MatFactorSymbolic_LU_ConstantDiagonal(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo * info)424d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorSymbolic_LU_ConstantDiagonal(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info)
425d71ae5a4SJacob Faibussowitsch {
4263423f386SBarry Smith PetscFunctionBegin;
4273423f386SBarry Smith fact->ops->lufactornumeric = MatFactorNumeric_ConstantDiagonal;
4283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4293423f386SBarry Smith }
4303423f386SBarry Smith
MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact,Mat A,IS isrow,const MatFactorInfo * info)431d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact, Mat A, IS isrow, const MatFactorInfo *info)
432d71ae5a4SJacob Faibussowitsch {
4333423f386SBarry Smith PetscFunctionBegin;
434a01b8767SStefano Zampini fact->ops->choleskyfactornumeric = MatFactorNumeric_ConstantDiagonal;
4353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4363423f386SBarry Smith }
4373423f386SBarry Smith
MatGetFactor_constantdiagonal_petsc(Mat A,MatFactorType ftype,Mat * B)438d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_constantdiagonal_petsc(Mat A, MatFactorType ftype, Mat *B)
439d71ae5a4SJacob Faibussowitsch {
4403423f386SBarry Smith PetscInt n = A->rmap->n, N = A->rmap->N;
4413423f386SBarry Smith
4423423f386SBarry Smith PetscFunctionBegin;
4439566063dSJacob Faibussowitsch PetscCall(MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A), n, n, N, N, 0, B));
4443423f386SBarry Smith
4453423f386SBarry Smith (*B)->factortype = ftype;
4463423f386SBarry Smith (*B)->ops->ilufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal;
4473423f386SBarry Smith (*B)->ops->lufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal;
4483423f386SBarry Smith (*B)->ops->iccfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal;
4493423f386SBarry Smith (*B)->ops->choleskyfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal;
4503423f386SBarry Smith
4513423f386SBarry Smith (*B)->ops->shift = NULL;
4523423f386SBarry Smith (*B)->ops->scale = NULL;
4533423f386SBarry Smith (*B)->ops->mult = NULL;
4543423f386SBarry Smith (*B)->ops->sor = NULL;
4553423f386SBarry Smith (*B)->ops->zeroentries = NULL;
4563423f386SBarry Smith
4579566063dSJacob Faibussowitsch PetscCall(PetscFree((*B)->solvertype));
4589566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERPETSC, &(*B)->solvertype));
4593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4603423f386SBarry Smith }
461