xref: /petsc/src/ksp/pc/tutorials/ex4.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
178a3110eSAlexander static char help[] = "Applies the 2023 preconditioner of Benzi and Faccio\n\n";
278a3110eSAlexander 
378a3110eSAlexander #include <petscmat.h>
478a3110eSAlexander #include <petscviewer.h>
578a3110eSAlexander #include <petscvec.h>
678a3110eSAlexander #include <petscis.h>
778a3110eSAlexander #include <petscksp.h>
878a3110eSAlexander 
978a3110eSAlexander /*
1078a3110eSAlexander  * This example reproduces the preconditioner outlined in Benzi's paper
1178a3110eSAlexander  * https://doi.org/10.1137/22M1505529. The problem considered is:
1278a3110eSAlexander  *
1378a3110eSAlexander  * (A + gamma UU^T)x = b
1478a3110eSAlexander  *
1578a3110eSAlexander  * whose structure arises from, for example, grad-div stabilization in the
1678a3110eSAlexander  * Navier-Stokes momentum equation. In the code we will also refer to
1778a3110eSAlexander  * gamma UU^T as J. The preconditioner developed by Benzi is:
1878a3110eSAlexander  *
1978a3110eSAlexander  * P_alpha = (A + alpha I)(alpha I + gamma UU^T)
2078a3110eSAlexander  *
2178a3110eSAlexander  * Another variant which may yield better convergence depending on the specific
2278a3110eSAlexander  * problem is
2378a3110eSAlexander  *
2478a3110eSAlexander  * P_alpha = (A + alpha D) D^-1 (alpha D + gamma UU^T)
2578a3110eSAlexander  *
2678a3110eSAlexander  * where D = diag(A + gamma UU^T). This is the variant implemented
2778a3110eSAlexander  * here. Application of the preconditioner involves (approximate) solution of
2878a3110eSAlexander  * two systems, one with (A + alpha D), and another with (alpha D + gamma
2978a3110eSAlexander  * UU^T). For small alpha (which generally yields the best overall
3078a3110eSAlexander  * preconditioner), (alpha D + gamma UU^T) is ill-conditioned. To combat this we
3178a3110eSAlexander  * solve (alpha D + gamma UU^T) using the Sherman-Morrison-Woodbury (SMW) matrix
3278a3110eSAlexander  * identity, which effectively converts the grad-div structure to a much nicer
3378a3110eSAlexander  * div-grad (laplacian) structure.
3478a3110eSAlexander  *
3578a3110eSAlexander  * The matrices used as input can be generated by running the matlab/octave
3678a3110eSAlexander  * program IFISS. The particular matrices checked into the datafiles repository
3778a3110eSAlexander  * and used in testing of this example correspond to a leaky lid-driven cavity
3878a3110eSAlexander  * with a stretched grid and Q2-Q1 finite elements. The matrices are taken from
3978a3110eSAlexander  * the last iteration of a Picard solve with tolerance 1e-8 with a viscosity of
4078a3110eSAlexander  * 0.1 and a 32x32 grid. We summarize below iteration counts from running this
4178a3110eSAlexander  * preconditioner for different grids and viscosity with a KSP tolerance of 1e-6.
4278a3110eSAlexander  *
4378a3110eSAlexander  *       32x32 64x64 128x128
4478a3110eSAlexander  * 0.1   28    36    43
4578a3110eSAlexander  * 0.01  59    75    73
4678a3110eSAlexander  * 0.002 136   161   167
4778a3110eSAlexander  *
4878a3110eSAlexander  * A reader of Benzi's paper will note that the performance shown above with
4978a3110eSAlexander  * respect to decreasing viscosity is significantly worse than in the
5078a3110eSAlexander  * paper. This is actually because of the choice of RHS. In Benzi's work, the
5178a3110eSAlexander  * RHS was generated by multiplying the operator with a vector of 1s whereas
5278a3110eSAlexander  * here we generate the RHS using a random vector. The iteration counts from the
5378a3110eSAlexander  * Benzi paper can be reproduced by changing the RHS generation in this example,
5478a3110eSAlexander  * but we choose to use the more difficult RHS as the resulting performance may
5578a3110eSAlexander  * more closely match what users experience in "physical" contexts.
5678a3110eSAlexander  */
5778a3110eSAlexander 
CreateAndLoadMat(const char * mat_name,Mat * mat)5878a3110eSAlexander PetscErrorCode CreateAndLoadMat(const char *mat_name, Mat *mat)
5978a3110eSAlexander {
6078a3110eSAlexander   PetscViewer viewer;
6178a3110eSAlexander   char        file[PETSC_MAX_PATH_LEN];
6278a3110eSAlexander   char        flag_name[10] = "-f";
6378a3110eSAlexander   PetscBool   flg;
6478a3110eSAlexander 
6578a3110eSAlexander   PetscFunctionBeginUser;
6678a3110eSAlexander   PetscCall(PetscOptionsGetString(NULL, NULL, strcat(flag_name, mat_name), file, sizeof(file), &flg));
6778a3110eSAlexander   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate file with the -f<mat_name> option");
6878a3110eSAlexander   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &viewer));
6978a3110eSAlexander   PetscCall(MatCreate(PETSC_COMM_WORLD, mat));
7078a3110eSAlexander   PetscCall(MatSetType(*mat, MATAIJ));
7178a3110eSAlexander   PetscCall(PetscObjectSetName((PetscObject)*mat, mat_name));
7278a3110eSAlexander   PetscCall(MatSetFromOptions(*mat));
7378a3110eSAlexander   PetscCall(MatLoad(*mat, viewer));
7478a3110eSAlexander   PetscCall(PetscViewerDestroy(&viewer));
7578a3110eSAlexander   PetscFunctionReturn(PETSC_SUCCESS);
7678a3110eSAlexander }
7778a3110eSAlexander 
7878a3110eSAlexander typedef struct {
7978a3110eSAlexander   Mat       U, UT, D, aD, aDinv, I_plus_gammaUTaDinvU;
8078a3110eSAlexander   PC        smw_cholesky;
8178a3110eSAlexander   PetscReal gamma, alpha;
8278a3110eSAlexander   PetscBool setup_called;
8378a3110eSAlexander } SmwPCCtx;
8478a3110eSAlexander 
SmwSetup(PC pc)8578a3110eSAlexander PetscErrorCode SmwSetup(PC pc)
8678a3110eSAlexander {
8778a3110eSAlexander   SmwPCCtx *ctx;
8878a3110eSAlexander   Vec       aDVec;
8978a3110eSAlexander 
9078a3110eSAlexander   PetscFunctionBeginUser;
9178a3110eSAlexander   PetscCall(PCShellGetContext(pc, &ctx));
9278a3110eSAlexander 
9378a3110eSAlexander   if (ctx->setup_called) PetscFunctionReturn(PETSC_SUCCESS);
9478a3110eSAlexander 
9578a3110eSAlexander   // Create aD
9678a3110eSAlexander   PetscCall(MatDuplicate(ctx->D, MAT_COPY_VALUES, &ctx->aD));
9778a3110eSAlexander   PetscCall(MatScale(ctx->aD, ctx->alpha));
9878a3110eSAlexander 
9978a3110eSAlexander   // Create aDinv
10078a3110eSAlexander   PetscCall(MatDuplicate(ctx->aD, MAT_DO_NOT_COPY_VALUES, &ctx->aDinv));
10178a3110eSAlexander   PetscCall(MatCreateVecs(ctx->aD, &aDVec, NULL));
10278a3110eSAlexander   PetscCall(MatGetDiagonal(ctx->aD, aDVec));
10378a3110eSAlexander   PetscCall(VecReciprocal(aDVec));
10478a3110eSAlexander   PetscCall(MatDiagonalSet(ctx->aDinv, aDVec, INSERT_VALUES));
10578a3110eSAlexander 
10678a3110eSAlexander   // Create UT
10778a3110eSAlexander   PetscCall(MatTranspose(ctx->U, MAT_INITIAL_MATRIX, &ctx->UT));
10878a3110eSAlexander 
10978a3110eSAlexander   // Create sum Mat
110fb842aefSJose E. Roman   PetscCall(MatMatMatMult(ctx->UT, ctx->aDinv, ctx->U, MAT_INITIAL_MATRIX, PETSC_CURRENT, &ctx->I_plus_gammaUTaDinvU));
11178a3110eSAlexander   PetscCall(MatScale(ctx->I_plus_gammaUTaDinvU, ctx->gamma));
11278a3110eSAlexander   PetscCall(MatShift(ctx->I_plus_gammaUTaDinvU, 1.));
11378a3110eSAlexander 
11478a3110eSAlexander   PetscCall(PCCreate(PETSC_COMM_WORLD, &ctx->smw_cholesky));
11578a3110eSAlexander   PetscCall(PCSetType(ctx->smw_cholesky, PCCHOLESKY));
11678a3110eSAlexander   PetscCall(PCSetOperators(ctx->smw_cholesky, ctx->I_plus_gammaUTaDinvU, ctx->I_plus_gammaUTaDinvU));
11778a3110eSAlexander   PetscCall(PCSetOptionsPrefix(ctx->smw_cholesky, "smw_"));
11878a3110eSAlexander   PetscCall(PCSetFromOptions(ctx->smw_cholesky));
11978a3110eSAlexander   PetscCall(PCSetUp(ctx->smw_cholesky));
12078a3110eSAlexander 
12178a3110eSAlexander   PetscCall(VecDestroy(&aDVec));
12278a3110eSAlexander 
12378a3110eSAlexander   ctx->setup_called = PETSC_TRUE;
12478a3110eSAlexander   PetscFunctionReturn(PETSC_SUCCESS);
12578a3110eSAlexander }
12678a3110eSAlexander 
SmwApply(PC pc,Vec x,Vec y)12778a3110eSAlexander PetscErrorCode SmwApply(PC pc, Vec x, Vec y)
12878a3110eSAlexander {
12978a3110eSAlexander   SmwPCCtx *ctx;
13078a3110eSAlexander   Vec       vel0, pressure0, pressure1;
13178a3110eSAlexander 
13278a3110eSAlexander   PetscFunctionBeginUser;
13378a3110eSAlexander   PetscCall(PCShellGetContext(pc, &ctx));
13478a3110eSAlexander 
13578a3110eSAlexander   PetscCall(MatCreateVecs(ctx->UT, &vel0, &pressure0));
13678a3110eSAlexander   PetscCall(VecDuplicate(pressure0, &pressure1));
13778a3110eSAlexander 
13878a3110eSAlexander   // First term
13978a3110eSAlexander   PetscCall(MatMult(ctx->aDinv, x, vel0));
14078a3110eSAlexander   PetscCall(MatMult(ctx->UT, vel0, pressure0));
14178a3110eSAlexander   PetscCall(PCApply(ctx->smw_cholesky, pressure0, pressure1));
14278a3110eSAlexander   PetscCall(MatMult(ctx->U, pressure1, vel0));
14378a3110eSAlexander   PetscCall(MatMult(ctx->aDinv, vel0, y));
14478a3110eSAlexander   PetscCall(VecScale(y, -ctx->gamma));
14578a3110eSAlexander 
14678a3110eSAlexander   // Second term
14778a3110eSAlexander   PetscCall(MatMult(ctx->aDinv, x, vel0));
14878a3110eSAlexander 
14978a3110eSAlexander   PetscCall(VecAXPY(y, 1, vel0));
15078a3110eSAlexander 
15178a3110eSAlexander   PetscCall(VecDestroy(&vel0));
15278a3110eSAlexander   PetscCall(VecDestroy(&pressure0));
15378a3110eSAlexander   PetscCall(VecDestroy(&pressure1));
15478a3110eSAlexander   PetscFunctionReturn(PETSC_SUCCESS);
15578a3110eSAlexander }
15678a3110eSAlexander 
main(int argc,char ** args)15778a3110eSAlexander int main(int argc, char **args)
15878a3110eSAlexander {
15978a3110eSAlexander   Mat               A, B, Q, Acondensed, Bcondensed, BT, J, AplusJ, QInv, D, AplusD, JplusD, U;
16078a3110eSAlexander   Mat               AplusJarray[2];
16178a3110eSAlexander   Vec               bound, x, b, Qdiag, DVec;
16278a3110eSAlexander   PetscBool         flg;
16378a3110eSAlexander   PetscViewer       viewer;
16478a3110eSAlexander   char              file[PETSC_MAX_PATH_LEN];
16578a3110eSAlexander   PetscInt         *boundary_indices;
16678a3110eSAlexander   PetscInt          boundary_indices_size, am, an, bm, bn, condensed_am, astart, aend, Dstart, Dend, num_local_bnd_dofs = 0;
16778a3110eSAlexander   const PetscScalar zero = 0;
16878a3110eSAlexander   IS                boundary_is, bulk_is;
16978a3110eSAlexander   KSP               ksp;
17078a3110eSAlexander   PC                pc, pcA, pcJ;
17178a3110eSAlexander   PetscRandom       rctx;
17278a3110eSAlexander   PetscReal        *boundary_indices_values;
17378a3110eSAlexander   PetscReal         gamma = 100, alpha = .01;
17478a3110eSAlexander   PetscMPIInt       rank;
17578a3110eSAlexander   SmwPCCtx          ctx;
17678a3110eSAlexander 
17778a3110eSAlexander   PetscFunctionBeginUser;
178*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
17978a3110eSAlexander   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
18078a3110eSAlexander 
18178a3110eSAlexander   PetscCall(CreateAndLoadMat("A", &A));
18278a3110eSAlexander   PetscCall(CreateAndLoadMat("B", &B));
18378a3110eSAlexander   PetscCall(CreateAndLoadMat("Q", &Q));
18478a3110eSAlexander 
18578a3110eSAlexander   PetscCall(PetscOptionsGetString(NULL, NULL, "-fbound", file, sizeof(file), &flg));
18678a3110eSAlexander   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate file with the -fbound option");
18778a3110eSAlexander 
18878a3110eSAlexander   if (rank == 0) {
18978a3110eSAlexander     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, file, FILE_MODE_READ, &viewer));
19078a3110eSAlexander     PetscCall(VecCreate(PETSC_COMM_SELF, &bound));
19178a3110eSAlexander     PetscCall(PetscObjectSetName((PetscObject)bound, "bound"));
19278a3110eSAlexander     PetscCall(VecSetType(bound, VECSEQ));
19378a3110eSAlexander     PetscCall(VecLoad(bound, viewer));
19478a3110eSAlexander     PetscCall(PetscViewerDestroy(&viewer));
19578a3110eSAlexander     PetscCall(VecGetLocalSize(bound, &boundary_indices_size));
19678a3110eSAlexander     PetscCall(VecGetArray(bound, &boundary_indices_values));
19778a3110eSAlexander   }
19878a3110eSAlexander   PetscCallMPI(MPI_Bcast(&boundary_indices_size, 1, MPIU_INT, 0, PETSC_COMM_WORLD));
19978a3110eSAlexander   if (rank != 0) PetscCall(PetscMalloc1(boundary_indices_size, &boundary_indices_values));
2006497c311SBarry Smith   PetscCallMPI(MPI_Bcast(boundary_indices_values, (PetscMPIInt)boundary_indices_size, MPIU_SCALAR, 0, PETSC_COMM_WORLD));
20178a3110eSAlexander 
20278a3110eSAlexander   PetscCall(MatGetSize(A, &am, NULL));
20378a3110eSAlexander   // The total number of dofs for a given velocity component
20478a3110eSAlexander   const PetscInt nc = am / 2;
20578a3110eSAlexander   PetscCall(MatGetOwnershipRange(A, &astart, &aend));
20678a3110eSAlexander 
20778a3110eSAlexander   PetscCall(PetscMalloc1(2 * boundary_indices_size, &boundary_indices));
20878a3110eSAlexander 
20978a3110eSAlexander   //
21078a3110eSAlexander   // The dof index ordering appears to be all vx dofs and then all vy dofs.
21178a3110eSAlexander   //
21278a3110eSAlexander 
21378a3110eSAlexander   // First do vx
21478a3110eSAlexander   for (PetscInt i = 0; i < boundary_indices_size; ++i) {
21578a3110eSAlexander     // MATLAB uses 1-based indexing
21678a3110eSAlexander     const PetscInt bnd_dof = (PetscInt)boundary_indices_values[i] - 1;
21778a3110eSAlexander     if ((bnd_dof >= astart) && (bnd_dof < aend)) boundary_indices[num_local_bnd_dofs++] = bnd_dof;
21878a3110eSAlexander   }
21978a3110eSAlexander 
22078a3110eSAlexander   // Now vy
22178a3110eSAlexander   for (PetscInt i = 0; i < boundary_indices_size; ++i) {
22278a3110eSAlexander     // MATLAB uses 1-based indexing
22378a3110eSAlexander     const PetscInt bnd_dof = ((PetscInt)boundary_indices_values[i] - 1) + nc;
22478a3110eSAlexander     if ((bnd_dof >= astart) && (bnd_dof < aend)) boundary_indices[num_local_bnd_dofs++] = bnd_dof;
22578a3110eSAlexander   }
22678a3110eSAlexander   if (rank == 0) PetscCall(VecRestoreArray(bound, &boundary_indices_values));
22778a3110eSAlexander   else PetscCall(PetscFree(boundary_indices_values));
22878a3110eSAlexander 
22978a3110eSAlexander   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, num_local_bnd_dofs, boundary_indices, PETSC_USE_POINTER, &boundary_is));
23078a3110eSAlexander   PetscCall(ISComplement(boundary_is, astart, aend, &bulk_is));
23178a3110eSAlexander 
23278a3110eSAlexander   PetscCall(MatCreateSubMatrix(A, bulk_is, bulk_is, MAT_INITIAL_MATRIX, &Acondensed));
23378a3110eSAlexander   // Can't pass null for row index set :-(
23478a3110eSAlexander   PetscCall(MatTranspose(B, MAT_INPLACE_MATRIX, &B));
23578a3110eSAlexander   PetscCall(MatCreateSubMatrix(B, bulk_is, NULL, MAT_INITIAL_MATRIX, &Bcondensed));
23678a3110eSAlexander   PetscCall(MatGetLocalSize(Acondensed, &am, &an));
23778a3110eSAlexander   PetscCall(MatGetLocalSize(Bcondensed, &bm, &bn));
23878a3110eSAlexander 
23978a3110eSAlexander   // Create QInv
24078a3110eSAlexander   PetscCall(MatCreateVecs(Q, &Qdiag, NULL));
24178a3110eSAlexander   PetscCall(MatGetDiagonal(Q, Qdiag));
24278a3110eSAlexander   PetscCall(VecReciprocal(Qdiag));
24378a3110eSAlexander   PetscCall(MatDuplicate(Q, MAT_DO_NOT_COPY_VALUES, &QInv));
24478a3110eSAlexander   PetscCall(MatDiagonalSet(QInv, Qdiag, INSERT_VALUES));
24578a3110eSAlexander   PetscCall(MatAssemblyBegin(QInv, MAT_FINAL_ASSEMBLY));
24678a3110eSAlexander   PetscCall(MatAssemblyEnd(QInv, MAT_FINAL_ASSEMBLY));
24778a3110eSAlexander 
24878a3110eSAlexander   // Create J
24978a3110eSAlexander   PetscCall(MatTranspose(Bcondensed, MAT_INITIAL_MATRIX, &BT));
250fb842aefSJose E. Roman   PetscCall(MatMatMatMult(Bcondensed, QInv, BT, MAT_INITIAL_MATRIX, PETSC_CURRENT, &J));
25178a3110eSAlexander   PetscCall(MatScale(J, gamma));
25278a3110eSAlexander 
25378a3110eSAlexander   // Create sum of A + J
25478a3110eSAlexander   AplusJarray[0] = Acondensed;
25578a3110eSAlexander   AplusJarray[1] = J;
25678a3110eSAlexander   PetscCall(MatCreateComposite(PETSC_COMM_WORLD, 2, AplusJarray, &AplusJ));
25778a3110eSAlexander 
25878a3110eSAlexander   // Create decomposition matrices
25978a3110eSAlexander   // We've already used Qdiag, which currently represents Q^-1,  for its necessary purposes. Let's
26078a3110eSAlexander   // convert it to represent Q^(-1/2)
26178a3110eSAlexander   PetscCall(VecSqrtAbs(Qdiag));
26278a3110eSAlexander   // We can similarly reuse Qinv
26378a3110eSAlexander   PetscCall(MatDiagonalSet(QInv, Qdiag, INSERT_VALUES));
26478a3110eSAlexander   PetscCall(MatAssemblyBegin(QInv, MAT_FINAL_ASSEMBLY));
26578a3110eSAlexander   PetscCall(MatAssemblyEnd(QInv, MAT_FINAL_ASSEMBLY));
26678a3110eSAlexander   // Create U
267fb842aefSJose E. Roman   PetscCall(MatMatMult(Bcondensed, QInv, MAT_INITIAL_MATRIX, PETSC_CURRENT, &U));
26878a3110eSAlexander 
26978a3110eSAlexander   // Create x and b
27078a3110eSAlexander   PetscCall(MatCreateVecs(AplusJ, &x, &b));
27178a3110eSAlexander   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx));
27278a3110eSAlexander   PetscCall(VecSetRandom(x, rctx));
27378a3110eSAlexander   PetscCall(PetscRandomDestroy(&rctx));
27478a3110eSAlexander   PetscCall(MatMult(AplusJ, x, b));
27578a3110eSAlexander 
27678a3110eSAlexander   // Compute preconditioner operators
27778a3110eSAlexander   PetscCall(MatGetLocalSize(Acondensed, &condensed_am, NULL));
27878a3110eSAlexander   PetscCall(MatCreate(PETSC_COMM_WORLD, &D));
27978a3110eSAlexander   PetscCall(MatSetType(D, MATAIJ));
28078a3110eSAlexander   PetscCall(MatSetSizes(D, condensed_am, condensed_am, PETSC_DETERMINE, PETSC_DETERMINE));
28178a3110eSAlexander   PetscCall(MatGetOwnershipRange(D, &Dstart, &Dend));
28278a3110eSAlexander   for (PetscInt i = Dstart; i < Dend; ++i) PetscCall(MatSetValues(D, 1, &i, 1, &i, &zero, INSERT_VALUES));
28378a3110eSAlexander   PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
28478a3110eSAlexander   PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
28578a3110eSAlexander   PetscCall(MatCreateVecs(D, &DVec, NULL));
28678a3110eSAlexander   PetscCall(MatGetDiagonal(AplusJ, DVec));
28778a3110eSAlexander   PetscCall(MatDiagonalSet(D, DVec, INSERT_VALUES));
28878a3110eSAlexander   PetscCall(MatDuplicate(Acondensed, MAT_COPY_VALUES, &AplusD));
28978a3110eSAlexander   PetscCall(MatAXPY(AplusD, alpha, D, SUBSET_NONZERO_PATTERN));
29078a3110eSAlexander   PetscCall(MatDuplicate(J, MAT_COPY_VALUES, &JplusD));
29178a3110eSAlexander   PetscCall(MatAXPY(JplusD, alpha, D, SUBSET_NONZERO_PATTERN));
29278a3110eSAlexander 
29378a3110eSAlexander   // Initialize our SMW context
29478a3110eSAlexander   ctx.U            = U;
29578a3110eSAlexander   ctx.D            = D;
29678a3110eSAlexander   ctx.gamma        = gamma;
29778a3110eSAlexander   ctx.alpha        = alpha;
29878a3110eSAlexander   ctx.setup_called = PETSC_FALSE;
29978a3110eSAlexander 
30078a3110eSAlexander   // Set preconditioner operators
30178a3110eSAlexander   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
30278a3110eSAlexander   PetscCall(KSPSetType(ksp, KSPFGMRES));
30378a3110eSAlexander   PetscCall(KSPSetOperators(ksp, AplusJ, AplusJ));
30478a3110eSAlexander   PetscCall(KSPSetNormType(ksp, KSP_NORM_UNPRECONDITIONED));
30578a3110eSAlexander   PetscCall(KSPGMRESSetRestart(ksp, 300));
30678a3110eSAlexander   PetscCall(KSPGetPC(ksp, &pc));
30778a3110eSAlexander   PetscCall(PCSetType(pc, PCCOMPOSITE));
30878a3110eSAlexander   PetscCall(PCCompositeSetType(pc, PC_COMPOSITE_SPECIAL));
30978a3110eSAlexander   PetscCall(PCCompositeAddPCType(pc, PCILU));
31078a3110eSAlexander   PetscCall(PCCompositeAddPCType(pc, PCSHELL));
31178a3110eSAlexander   PetscCall(PCCompositeGetPC(pc, 0, &pcA));
31278a3110eSAlexander   PetscCall(PCCompositeGetPC(pc, 1, &pcJ));
31378a3110eSAlexander   PetscCall(PCSetOperators(pcA, AplusD, AplusD));
31478a3110eSAlexander   PetscCall(PCSetOperators(pcJ, JplusD, JplusD));
31578a3110eSAlexander   PetscCall(PCShellSetContext(pcJ, &ctx));
31678a3110eSAlexander   PetscCall(PCShellSetApply(pcJ, SmwApply));
31778a3110eSAlexander   PetscCall(PCShellSetSetUp(pcJ, SmwSetup));
31878a3110eSAlexander   PetscCall(PCCompositeSpecialSetAlphaMat(pc, D));
31978a3110eSAlexander 
32078a3110eSAlexander   // Solve
32178a3110eSAlexander   PetscCall(KSPSetFromOptions(ksp));
32278a3110eSAlexander   PetscCall(KSPSolve(ksp, b, x));
32378a3110eSAlexander 
32478a3110eSAlexander   PetscCall(MatDestroy(&A));
32578a3110eSAlexander   PetscCall(MatDestroy(&B));
32678a3110eSAlexander   PetscCall(MatDestroy(&Q));
32778a3110eSAlexander   PetscCall(MatDestroy(&Acondensed));
32878a3110eSAlexander   PetscCall(MatDestroy(&Bcondensed));
32978a3110eSAlexander   PetscCall(MatDestroy(&BT));
33078a3110eSAlexander   PetscCall(MatDestroy(&J));
33178a3110eSAlexander   PetscCall(MatDestroy(&AplusJ));
33278a3110eSAlexander   PetscCall(MatDestroy(&QInv));
33378a3110eSAlexander   PetscCall(MatDestroy(&D));
33478a3110eSAlexander   PetscCall(MatDestroy(&AplusD));
33578a3110eSAlexander   PetscCall(MatDestroy(&JplusD));
33678a3110eSAlexander   PetscCall(MatDestroy(&U));
33778a3110eSAlexander   if (rank == 0) PetscCall(VecDestroy(&bound));
33878a3110eSAlexander   PetscCall(VecDestroy(&x));
33978a3110eSAlexander   PetscCall(VecDestroy(&b));
34078a3110eSAlexander   PetscCall(VecDestroy(&Qdiag));
34178a3110eSAlexander   PetscCall(VecDestroy(&DVec));
34278a3110eSAlexander   PetscCall(ISDestroy(&boundary_is));
34378a3110eSAlexander   PetscCall(ISDestroy(&bulk_is));
34478a3110eSAlexander   PetscCall(KSPDestroy(&ksp));
34578a3110eSAlexander   PetscCall(PetscFree(boundary_indices));
34678a3110eSAlexander   PetscCall(MatDestroy(&ctx.UT));
34778a3110eSAlexander   PetscCall(MatDestroy(&ctx.I_plus_gammaUTaDinvU));
34878a3110eSAlexander   PetscCall(MatDestroy(&ctx.aD));
34978a3110eSAlexander   PetscCall(MatDestroy(&ctx.aDinv));
35078a3110eSAlexander   PetscCall(PCDestroy(&ctx.smw_cholesky));
35178a3110eSAlexander 
35278a3110eSAlexander   PetscCall(PetscFinalize());
35378a3110eSAlexander   return 0;
35478a3110eSAlexander }
35578a3110eSAlexander 
35678a3110eSAlexander /*TEST
35778a3110eSAlexander 
35878a3110eSAlexander    build:
35978a3110eSAlexander       requires: !complex
36078a3110eSAlexander 
36178a3110eSAlexander    test:
36278a3110eSAlexander       args: -fA ${DATAFILESPATH}/matrices/ifiss/A -fB ${DATAFILESPATH}/matrices/ifiss/B -fQ ${DATAFILESPATH}/matrices/ifiss/Q -fbound ${DATAFILESPATH}/is/ifiss/bound -ksp_monitor
36378a3110eSAlexander       requires: datafilespath defined(PETSC_USE_64BIT_INDICES) !complex double
36478a3110eSAlexander 
36578a3110eSAlexander    test:
36678a3110eSAlexander       suffix: 2
36778a3110eSAlexander       nsize: 2
36878a3110eSAlexander       args: -fA ${DATAFILESPATH}/matrices/ifiss/A -fB ${DATAFILESPATH}/matrices/ifiss/B -fQ ${DATAFILESPATH}/matrices/ifiss/Q -fbound ${DATAFILESPATH}/is/ifiss/bound -ksp_monitor
36978a3110eSAlexander       requires: datafilespath defined(PETSC_USE_64BIT_INDICES) !complex double strumpack
37078a3110eSAlexander 
37178a3110eSAlexander TEST*/
372