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; 1318a1df862SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)mat, "MatConstantDiagonalGetConstant_C", NULL)); 1323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1333423f386SBarry Smith } 1343423f386SBarry Smith 135d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_ConstantDiagonal(Mat J, PetscViewer viewer) 136d71ae5a4SJacob Faibussowitsch { 1373423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data; 138*9f196a02SMartin Diehl PetscBool isascii; 1393423f386SBarry Smith 1403423f386SBarry Smith PetscFunctionBegin; 141*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 142*9f196a02SMartin Diehl if (isascii) { 1433423f386SBarry Smith PetscViewerFormat format; 1443423f386SBarry Smith 1459566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 1463ba16761SJacob Faibussowitsch if (format == PETSC_VIEWER_ASCII_FACTOR_INFO || format == PETSC_VIEWER_ASCII_INFO) PetscFunctionReturn(PETSC_SUCCESS); 1478a1df862SToby Isaac if (PetscImaginaryPart(ctx->diag) == 0) { 1488a1df862SToby Isaac PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g\n", (double)PetscRealPart(ctx->diag))); 1498a1df862SToby Isaac } else { 1509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Diagonal value: %g + i %g\n", (double)PetscRealPart(ctx->diag), (double)PetscImaginaryPart(ctx->diag))); 1518a1df862SToby Isaac } 1523423f386SBarry Smith } 1533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1543423f386SBarry Smith } 1553423f386SBarry Smith 156d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_ConstantDiagonal(Mat J, Vec x, Vec y) 157d71ae5a4SJacob Faibussowitsch { 1583423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data; 1593423f386SBarry Smith 1603423f386SBarry Smith PetscFunctionBegin; 1619566063dSJacob Faibussowitsch PetscCall(VecAXPBY(y, ctx->diag, 0.0, x)); 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1633423f386SBarry Smith } 1643423f386SBarry Smith 165345a4b08SToby Isaac static PetscErrorCode MatMultHermitianTranspose_ConstantDiagonal(Mat J, Vec x, Vec y) 166345a4b08SToby Isaac { 167345a4b08SToby Isaac Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data; 168345a4b08SToby Isaac 169345a4b08SToby Isaac PetscFunctionBegin; 170345a4b08SToby Isaac PetscCall(VecAXPBY(y, PetscConj(ctx->diag), 0.0, x)); 171345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 172345a4b08SToby Isaac } 173345a4b08SToby Isaac 174345a4b08SToby Isaac static PetscErrorCode MatGetDiagonal_ConstantDiagonal(Mat J, Vec x) 175d71ae5a4SJacob Faibussowitsch { 1763423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)J->data; 1773423f386SBarry Smith 1783423f386SBarry Smith PetscFunctionBegin; 1799566063dSJacob Faibussowitsch PetscCall(VecSet(x, ctx->diag)); 1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1813423f386SBarry Smith } 1823423f386SBarry Smith 183d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatShift_ConstantDiagonal(Mat Y, PetscScalar a) 184d71ae5a4SJacob Faibussowitsch { 1853423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data; 1863423f386SBarry Smith 1873423f386SBarry Smith PetscFunctionBegin; 1883423f386SBarry Smith ctx->diag += a; 1893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1903423f386SBarry Smith } 1913423f386SBarry Smith 192d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatScale_ConstantDiagonal(Mat Y, PetscScalar a) 193d71ae5a4SJacob Faibussowitsch { 1943423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data; 1953423f386SBarry Smith 1963423f386SBarry Smith PetscFunctionBegin; 1973423f386SBarry Smith ctx->diag *= a; 1983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1993423f386SBarry Smith } 2003423f386SBarry Smith 201d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatZeroEntries_ConstantDiagonal(Mat Y) 202d71ae5a4SJacob Faibussowitsch { 2033423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data; 2043423f386SBarry Smith 2053423f386SBarry Smith PetscFunctionBegin; 2063423f386SBarry Smith ctx->diag = 0.0; 2073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2083423f386SBarry Smith } 2093423f386SBarry Smith 2108a1df862SToby Isaac static PetscErrorCode MatConjugate_ConstantDiagonal(Mat Y) 2118a1df862SToby Isaac { 2128a1df862SToby Isaac Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)Y->data; 2138a1df862SToby Isaac 2148a1df862SToby Isaac PetscFunctionBegin; 2158a1df862SToby Isaac ctx->diag = PetscConj(ctx->diag); 2168a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 2178a1df862SToby Isaac } 2188a1df862SToby Isaac 2198a1df862SToby Isaac static PetscErrorCode MatTranspose_ConstantDiagonal(Mat A, MatReuse reuse, Mat *matout) 2208a1df862SToby Isaac { 2218a1df862SToby Isaac Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data; 2228a1df862SToby Isaac 2238a1df862SToby Isaac PetscFunctionBegin; 2248a1df862SToby Isaac if (reuse == MAT_INPLACE_MATRIX) { 2258a1df862SToby Isaac PetscLayout tmplayout = A->rmap; 2268a1df862SToby Isaac 2278a1df862SToby Isaac A->rmap = A->cmap; 2288a1df862SToby Isaac A->cmap = tmplayout; 2298a1df862SToby Isaac } else { 2308a1df862SToby Isaac if (reuse == MAT_INITIAL_MATRIX) { 2318a1df862SToby Isaac PetscCall(MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A), A->cmap->n, A->rmap->n, A->cmap->N, A->rmap->N, ctx->diag, matout)); 2328a1df862SToby Isaac } else { 2338a1df862SToby Isaac PetscCall(MatZeroEntries(*matout)); 2348a1df862SToby Isaac PetscCall(MatShift(*matout, ctx->diag)); 2358a1df862SToby Isaac } 2368a1df862SToby Isaac } 2378a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 2388a1df862SToby Isaac } 2398a1df862SToby Isaac 2408a1df862SToby Isaac static PetscErrorCode MatSetRandom_ConstantDiagonal(Mat A, PetscRandom rand) 2418a1df862SToby Isaac { 2428a1df862SToby Isaac PetscMPIInt rank; 2438a1df862SToby Isaac MPI_Comm comm; 2448a1df862SToby Isaac PetscScalar v = 0.0; 2458a1df862SToby Isaac Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)A->data; 2468a1df862SToby Isaac 2478a1df862SToby Isaac PetscFunctionBegin; 2488a1df862SToby Isaac PetscCall(PetscObjectGetComm((PetscObject)A, &comm)); 2498a1df862SToby Isaac PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2508a1df862SToby Isaac if (!rank) PetscCall(PetscRandomGetValue(rand, &v)); 2518a1df862SToby Isaac PetscCallMPI(MPI_Bcast(&v, 1, MPIU_SCALAR, 0, comm)); 2528a1df862SToby Isaac ctx->diag = v; 2538a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 2548a1df862SToby Isaac } 2558a1df862SToby Isaac 256345a4b08SToby Isaac static PetscErrorCode MatSolve_ConstantDiagonal(Mat matin, Vec b, Vec x) 257d71ae5a4SJacob Faibussowitsch { 2583423f386SBarry Smith Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)matin->data; 2593423f386SBarry Smith 2603423f386SBarry Smith PetscFunctionBegin; 2613423f386SBarry Smith if (ctx->diag == 0.0) matin->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 2623423f386SBarry Smith else matin->factorerrortype = MAT_FACTOR_NOERROR; 263345a4b08SToby Isaac PetscCall(VecAXPBY(x, 1.0 / ctx->diag, 0.0, b)); 264345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 265345a4b08SToby Isaac } 266345a4b08SToby Isaac 26766976f2fSJacob Faibussowitsch static PetscErrorCode MatSOR_ConstantDiagonal(Mat matin, Vec x, PetscReal omega, MatSORType flag, PetscReal fshift, PetscInt its, PetscInt lits, Vec y) 268345a4b08SToby Isaac { 269345a4b08SToby Isaac PetscFunctionBegin; 270345a4b08SToby Isaac PetscCall(MatSolve_ConstantDiagonal(matin, x, y)); 2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2723423f386SBarry Smith } 2733423f386SBarry Smith 27466976f2fSJacob Faibussowitsch static PetscErrorCode MatGetInfo_ConstantDiagonal(Mat A, MatInfoType flag, MatInfo *info) 275d71ae5a4SJacob Faibussowitsch { 2763423f386SBarry Smith PetscFunctionBegin; 2773423f386SBarry Smith info->block_size = 1.0; 2783423f386SBarry Smith info->nz_allocated = 1.0; 2793423f386SBarry Smith info->nz_used = 1.0; 2803423f386SBarry Smith info->nz_unneeded = 0.0; 2813966268fSBarry Smith info->assemblies = A->num_ass; 2823423f386SBarry Smith info->mallocs = 0.0; 2834dfa11a4SJacob Faibussowitsch info->memory = 0; /* REVIEW ME */ 2843423f386SBarry Smith if (A->factortype) { 2853423f386SBarry Smith info->fill_ratio_given = 1.0; 2863423f386SBarry Smith info->fill_ratio_needed = 1.0; 2873423f386SBarry Smith info->factor_mallocs = 0.0; 2883423f386SBarry Smith } else { 2893423f386SBarry Smith info->fill_ratio_given = 0; 2903423f386SBarry Smith info->fill_ratio_needed = 0; 2913423f386SBarry Smith info->factor_mallocs = 0; 2923423f386SBarry Smith } 2933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2943423f386SBarry Smith } 2953423f386SBarry Smith 2963423f386SBarry Smith /*@ 2973423f386SBarry Smith MatCreateConstantDiagonal - Creates a matrix with a uniform value along the diagonal 2983423f386SBarry Smith 299d083f849SBarry Smith Collective 3003423f386SBarry Smith 3013423f386SBarry Smith Input Parameters: 3023423f386SBarry Smith + comm - MPI communicator 3032ef1f0ffSBarry Smith . m - number of local rows (or `PETSC_DECIDE` to have calculated if `M` is given) 3043423f386SBarry Smith This value should be the same as the local size used in creating the 3053423f386SBarry Smith y vector for the matrix-vector product y = Ax. 3063423f386SBarry Smith . n - This value should be the same as the local size used in creating the 3072ef1f0ffSBarry Smith x vector for the matrix-vector product y = Ax. (or `PETSC_DECIDE` to have 3082ef1f0ffSBarry Smith calculated if `N` is given) For square matrices n is almost always `m`. 30911a5261eSBarry Smith . M - number of global rows (or `PETSC_DETERMINE` to have calculated if m is given) 31011a5261eSBarry Smith . N - number of global columns (or `PETSC_DETERMINE` to have calculated if n is given) 3113423f386SBarry Smith - diag - the diagonal value 3123423f386SBarry Smith 3133423f386SBarry Smith Output Parameter: 3143423f386SBarry Smith . J - the diagonal matrix 3153423f386SBarry Smith 3163423f386SBarry Smith Level: advanced 3173423f386SBarry Smith 3183423f386SBarry Smith Notes: 3193423f386SBarry Smith Only supports square matrices with the same number of local rows and columns 3203423f386SBarry Smith 3211cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MATCONSTANTDIAGONAL`, `MatScale()`, `MatShift()`, `MatMult()`, `MatGetDiagonal()`, `MatGetFactor()`, `MatSolve()` 3223423f386SBarry Smith @*/ 323d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCreateConstantDiagonal(MPI_Comm comm, PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscScalar diag, Mat *J) 324d71ae5a4SJacob Faibussowitsch { 3253423f386SBarry Smith PetscFunctionBegin; 3269566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, J)); 3279566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, m, n, M, N)); 3289566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, MATCONSTANTDIAGONAL)); 3299566063dSJacob Faibussowitsch PetscCall(MatShift(*J, diag)); 3309566063dSJacob Faibussowitsch PetscCall(MatSetUp(*J)); 3313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3323423f386SBarry Smith } 3333423f386SBarry Smith 3348a1df862SToby Isaac /*@ 3358a1df862SToby Isaac MatConstantDiagonalGetConstant - Get the scalar constant of a constant diagonal matrix 3368a1df862SToby Isaac 3378a1df862SToby Isaac Not collective 3388a1df862SToby Isaac 3398a1df862SToby Isaac Input Parameter: 3408a1df862SToby Isaac . mat - a `MATCONSTANTDIAGONAL` 3418a1df862SToby Isaac 3428a1df862SToby Isaac Output Parameter: 3438a1df862SToby Isaac . value - the scalar value 3448a1df862SToby Isaac 3458a1df862SToby Isaac Level: developer 3468a1df862SToby Isaac 3478a1df862SToby Isaac .seealso: [](ch_matrices), `Mat`, `MatDestroy()`, `MATCONSTANTDIAGONAL` 3488a1df862SToby Isaac @*/ 3498a1df862SToby Isaac PetscErrorCode MatConstantDiagonalGetConstant(Mat mat, PetscScalar *value) 3508a1df862SToby Isaac { 3518a1df862SToby Isaac PetscFunctionBegin; 3528a1df862SToby Isaac PetscUseMethod(mat, "MatConstantDiagonalGetConstant_C", (Mat, PetscScalar *), (mat, value)); 3538a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 3548a1df862SToby Isaac } 3558a1df862SToby Isaac 3568a1df862SToby Isaac static PetscErrorCode MatConstantDiagonalGetConstant_ConstantDiagonal(Mat mat, PetscScalar *value) 3578a1df862SToby Isaac { 3588a1df862SToby Isaac Mat_ConstantDiagonal *ctx = (Mat_ConstantDiagonal *)mat->data; 3598a1df862SToby Isaac 3608a1df862SToby Isaac PetscFunctionBegin; 3618a1df862SToby Isaac *value = ctx->diag; 3628a1df862SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 3638a1df862SToby Isaac } 3648a1df862SToby Isaac 3654bf1a0e8SHansol Suh /*MC 3664bf1a0e8SHansol Suh MATCONSTANTDIAGONAL - "constant-diagonal" - A diagonal matrix type with a uniform value 3674bf1a0e8SHansol Suh along the diagonal. 3684bf1a0e8SHansol Suh 3694bf1a0e8SHansol Suh Level: advanced 3704bf1a0e8SHansol Suh 3714bf1a0e8SHansol Suh .seealso: [](ch_matrices), `Mat`, `MatCreateConstantDiagonal()` 3724bf1a0e8SHansol Suh M*/ 373d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreate_ConstantDiagonal(Mat A) 374d71ae5a4SJacob Faibussowitsch { 3753423f386SBarry Smith Mat_ConstantDiagonal *ctx; 3763423f386SBarry Smith 3773423f386SBarry Smith PetscFunctionBegin; 3789566063dSJacob Faibussowitsch PetscCall(PetscNew(&ctx)); 3793423f386SBarry Smith ctx->diag = 0.0; 3803423f386SBarry Smith A->data = (void *)ctx; 3813423f386SBarry Smith 3823423f386SBarry Smith A->assembled = PETSC_TRUE; 3833423f386SBarry Smith A->preallocated = PETSC_TRUE; 384345a4b08SToby Isaac A->structurally_symmetric = PETSC_BOOL3_TRUE; 385345a4b08SToby Isaac A->structural_symmetry_eternal = PETSC_TRUE; 386345a4b08SToby Isaac A->symmetric = PETSC_BOOL3_TRUE; 387345a4b08SToby Isaac if (!PetscDefined(USE_COMPLEX)) A->hermitian = PETSC_BOOL3_TRUE; 388345a4b08SToby Isaac A->symmetry_eternal = PETSC_TRUE; 389a8d4b458SStefano Zampini 3903423f386SBarry Smith A->ops->mult = MatMult_ConstantDiagonal; 391a8d4b458SStefano Zampini A->ops->multadd = MatMultAdd_ConstantDiagonal; 392345a4b08SToby Isaac A->ops->multtranspose = MatMult_ConstantDiagonal; 393345a4b08SToby Isaac A->ops->multtransposeadd = MatMultAdd_ConstantDiagonal; 394345a4b08SToby Isaac A->ops->multhermitiantranspose = MatMultHermitianTranspose_ConstantDiagonal; 395345a4b08SToby Isaac A->ops->multhermitiantransposeadd = MatMultHermitianTransposeAdd_ConstantDiagonal; 396345a4b08SToby Isaac A->ops->solve = MatSolve_ConstantDiagonal; 397345a4b08SToby Isaac A->ops->solvetranspose = MatSolve_ConstantDiagonal; 398e5f3c4beSJose E. Roman A->ops->norm = MatNorm_ConstantDiagonal; 399d1ca2681SJose E. Roman A->ops->createsubmatrices = MatCreateSubMatrices_ConstantDiagonal; 400a8d4b458SStefano Zampini A->ops->duplicate = MatDuplicate_ConstantDiagonal; 401a8d4b458SStefano Zampini A->ops->missingdiagonal = MatMissingDiagonal_ConstantDiagonal; 402a8d4b458SStefano Zampini A->ops->getrow = MatGetRow_ConstantDiagonal; 403a8d4b458SStefano Zampini A->ops->restorerow = MatRestoreRow_ConstantDiagonal; 4043423f386SBarry Smith A->ops->sor = MatSOR_ConstantDiagonal; 4053423f386SBarry Smith A->ops->shift = MatShift_ConstantDiagonal; 4063423f386SBarry Smith A->ops->scale = MatScale_ConstantDiagonal; 4073423f386SBarry Smith A->ops->getdiagonal = MatGetDiagonal_ConstantDiagonal; 4083423f386SBarry Smith A->ops->view = MatView_ConstantDiagonal; 4093423f386SBarry Smith A->ops->zeroentries = MatZeroEntries_ConstantDiagonal; 4103423f386SBarry Smith A->ops->destroy = MatDestroy_ConstantDiagonal; 4113423f386SBarry Smith A->ops->getinfo = MatGetInfo_ConstantDiagonal; 412345a4b08SToby Isaac A->ops->equal = MatEqual_ConstantDiagonal; 413a01b8767SStefano Zampini A->ops->axpy = MatAXPY_ConstantDiagonal; 4148a1df862SToby Isaac A->ops->setrandom = MatSetRandom_ConstantDiagonal; 4158a1df862SToby Isaac A->ops->conjugate = MatConjugate_ConstantDiagonal; 4168a1df862SToby Isaac A->ops->transpose = MatTranspose_ConstantDiagonal; 417a8d4b458SStefano Zampini 4189566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A, MATCONSTANTDIAGONAL)); 4198a1df862SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)A, "MatConstantDiagonalGetConstant_C", MatConstantDiagonalGetConstant_ConstantDiagonal)); 4203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4213423f386SBarry Smith } 4223423f386SBarry Smith 423d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorNumeric_ConstantDiagonal(Mat fact, Mat A, const MatFactorInfo *info) 424d71ae5a4SJacob Faibussowitsch { 4253423f386SBarry Smith Mat_ConstantDiagonal *actx = (Mat_ConstantDiagonal *)A->data, *fctx = (Mat_ConstantDiagonal *)fact->data; 4263423f386SBarry Smith 4273423f386SBarry Smith PetscFunctionBegin; 4283423f386SBarry Smith if (actx->diag == 0.0) fact->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 4293423f386SBarry Smith else fact->factorerrortype = MAT_FACTOR_NOERROR; 4303423f386SBarry Smith fctx->diag = 1.0 / actx->diag; 4313423f386SBarry Smith fact->ops->solve = MatMult_ConstantDiagonal; 4323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4333423f386SBarry Smith } 4343423f386SBarry Smith 435d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorSymbolic_LU_ConstantDiagonal(Mat fact, Mat A, IS isrow, IS iscol, const MatFactorInfo *info) 436d71ae5a4SJacob Faibussowitsch { 4373423f386SBarry Smith PetscFunctionBegin; 4383423f386SBarry Smith fact->ops->lufactornumeric = MatFactorNumeric_ConstantDiagonal; 4393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4403423f386SBarry Smith } 4413423f386SBarry Smith 442d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorSymbolic_Cholesky_ConstantDiagonal(Mat fact, Mat A, IS isrow, const MatFactorInfo *info) 443d71ae5a4SJacob Faibussowitsch { 4443423f386SBarry Smith PetscFunctionBegin; 445a01b8767SStefano Zampini fact->ops->choleskyfactornumeric = MatFactorNumeric_ConstantDiagonal; 4463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4473423f386SBarry Smith } 4483423f386SBarry Smith 449d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_constantdiagonal_petsc(Mat A, MatFactorType ftype, Mat *B) 450d71ae5a4SJacob Faibussowitsch { 4513423f386SBarry Smith PetscInt n = A->rmap->n, N = A->rmap->N; 4523423f386SBarry Smith 4533423f386SBarry Smith PetscFunctionBegin; 4549566063dSJacob Faibussowitsch PetscCall(MatCreateConstantDiagonal(PetscObjectComm((PetscObject)A), n, n, N, N, 0, B)); 4553423f386SBarry Smith 4563423f386SBarry Smith (*B)->factortype = ftype; 4573423f386SBarry Smith (*B)->ops->ilufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal; 4583423f386SBarry Smith (*B)->ops->lufactorsymbolic = MatFactorSymbolic_LU_ConstantDiagonal; 4593423f386SBarry Smith (*B)->ops->iccfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal; 4603423f386SBarry Smith (*B)->ops->choleskyfactorsymbolic = MatFactorSymbolic_Cholesky_ConstantDiagonal; 4613423f386SBarry Smith 4623423f386SBarry Smith (*B)->ops->shift = NULL; 4633423f386SBarry Smith (*B)->ops->scale = NULL; 4643423f386SBarry Smith (*B)->ops->mult = NULL; 4653423f386SBarry Smith (*B)->ops->sor = NULL; 4663423f386SBarry Smith (*B)->ops->zeroentries = NULL; 4673423f386SBarry Smith 4689566063dSJacob Faibussowitsch PetscCall(PetscFree((*B)->solvertype)); 4699566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERPETSC, &(*B)->solvertype)); 4703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4713423f386SBarry Smith } 472