xref: /petsc/src/ksp/pc/tutorials/ex4.c (revision 78a3110e6e5b4f375af22805359ffc6186cef0f0)
1*78a3110eSAlexander static char help[] = "Applies the 2023 preconditioner of Benzi and Faccio\n\n";
2*78a3110eSAlexander 
3*78a3110eSAlexander #include <petscmat.h>
4*78a3110eSAlexander #include <petscviewer.h>
5*78a3110eSAlexander #include <petscvec.h>
6*78a3110eSAlexander #include <petscis.h>
7*78a3110eSAlexander #include <petscksp.h>
8*78a3110eSAlexander 
9*78a3110eSAlexander /*
10*78a3110eSAlexander  * This example reproduces the preconditioner outlined in Benzi's paper
11*78a3110eSAlexander  * https://doi.org/10.1137/22M1505529. The problem considered is:
12*78a3110eSAlexander  *
13*78a3110eSAlexander  * (A + gamma UU^T)x = b
14*78a3110eSAlexander  *
15*78a3110eSAlexander  * whose structure arises from, for example, grad-div stabilization in the
16*78a3110eSAlexander  * Navier-Stokes momentum equation. In the code we will also refer to
17*78a3110eSAlexander  * gamma UU^T as J. The preconditioner developed by Benzi is:
18*78a3110eSAlexander  *
19*78a3110eSAlexander  * P_alpha = (A + alpha I)(alpha I + gamma UU^T)
20*78a3110eSAlexander  *
21*78a3110eSAlexander  * Another variant which may yield better convergence depending on the specific
22*78a3110eSAlexander  * problem is
23*78a3110eSAlexander  *
24*78a3110eSAlexander  * P_alpha = (A + alpha D) D^-1 (alpha D + gamma UU^T)
25*78a3110eSAlexander  *
26*78a3110eSAlexander  * where D = diag(A + gamma UU^T). This is the variant implemented
27*78a3110eSAlexander  * here. Application of the preconditioner involves (approximate) solution of
28*78a3110eSAlexander  * two systems, one with (A + alpha D), and another with (alpha D + gamma
29*78a3110eSAlexander  * UU^T). For small alpha (which generally yields the best overall
30*78a3110eSAlexander  * preconditioner), (alpha D + gamma UU^T) is ill-conditioned. To combat this we
31*78a3110eSAlexander  * solve (alpha D + gamma UU^T) using the Sherman-Morrison-Woodbury (SMW) matrix
32*78a3110eSAlexander  * identity, which effectively converts the grad-div structure to a much nicer
33*78a3110eSAlexander  * div-grad (laplacian) structure.
34*78a3110eSAlexander  *
35*78a3110eSAlexander  * The matrices used as input can be generated by running the matlab/octave
36*78a3110eSAlexander  * program IFISS. The particular matrices checked into the datafiles repository
37*78a3110eSAlexander  * and used in testing of this example correspond to a leaky lid-driven cavity
38*78a3110eSAlexander  * with a stretched grid and Q2-Q1 finite elements. The matrices are taken from
39*78a3110eSAlexander  * the last iteration of a Picard solve with tolerance 1e-8 with a viscosity of
40*78a3110eSAlexander  * 0.1 and a 32x32 grid. We summarize below iteration counts from running this
41*78a3110eSAlexander  * preconditioner for different grids and viscosity with a KSP tolerance of 1e-6.
42*78a3110eSAlexander  *
43*78a3110eSAlexander  *       32x32 64x64 128x128
44*78a3110eSAlexander  * 0.1   28    36    43
45*78a3110eSAlexander  * 0.01  59    75    73
46*78a3110eSAlexander  * 0.002 136   161   167
47*78a3110eSAlexander  *
48*78a3110eSAlexander  * A reader of Benzi's paper will note that the performance shown above with
49*78a3110eSAlexander  * respect to decreasing viscosity is significantly worse than in the
50*78a3110eSAlexander  * paper. This is actually because of the choice of RHS. In Benzi's work, the
51*78a3110eSAlexander  * RHS was generated by multiplying the operator with a vector of 1s whereas
52*78a3110eSAlexander  * here we generate the RHS using a random vector. The iteration counts from the
53*78a3110eSAlexander  * Benzi paper can be reproduced by changing the RHS generation in this example,
54*78a3110eSAlexander  * but we choose to use the more difficult RHS as the resulting performance may
55*78a3110eSAlexander  * more closely match what users experience in "physical" contexts.
56*78a3110eSAlexander  */
57*78a3110eSAlexander 
58*78a3110eSAlexander PetscErrorCode CreateAndLoadMat(const char *mat_name, Mat *mat)
59*78a3110eSAlexander {
60*78a3110eSAlexander   PetscViewer viewer;
61*78a3110eSAlexander   char        file[PETSC_MAX_PATH_LEN];
62*78a3110eSAlexander   char        flag_name[10] = "-f";
63*78a3110eSAlexander   PetscBool   flg;
64*78a3110eSAlexander 
65*78a3110eSAlexander   PetscFunctionBeginUser;
66*78a3110eSAlexander   PetscCall(PetscOptionsGetString(NULL, NULL, strcat(flag_name, mat_name), file, sizeof(file), &flg));
67*78a3110eSAlexander   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate file with the -f<mat_name> option");
68*78a3110eSAlexander   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &viewer));
69*78a3110eSAlexander   PetscCall(MatCreate(PETSC_COMM_WORLD, mat));
70*78a3110eSAlexander   PetscCall(MatSetType(*mat, MATAIJ));
71*78a3110eSAlexander   PetscCall(PetscObjectSetName((PetscObject)*mat, mat_name));
72*78a3110eSAlexander   PetscCall(MatSetFromOptions(*mat));
73*78a3110eSAlexander   PetscCall(MatLoad(*mat, viewer));
74*78a3110eSAlexander   PetscCall(PetscViewerDestroy(&viewer));
75*78a3110eSAlexander   PetscFunctionReturn(PETSC_SUCCESS);
76*78a3110eSAlexander }
77*78a3110eSAlexander 
78*78a3110eSAlexander typedef struct {
79*78a3110eSAlexander   Mat       U, UT, D, aD, aDinv, I_plus_gammaUTaDinvU;
80*78a3110eSAlexander   PC        smw_cholesky;
81*78a3110eSAlexander   PetscReal gamma, alpha;
82*78a3110eSAlexander   PetscBool setup_called;
83*78a3110eSAlexander } SmwPCCtx;
84*78a3110eSAlexander 
85*78a3110eSAlexander PetscErrorCode SmwSetup(PC pc)
86*78a3110eSAlexander {
87*78a3110eSAlexander   SmwPCCtx *ctx;
88*78a3110eSAlexander   Vec       aDVec;
89*78a3110eSAlexander 
90*78a3110eSAlexander   PetscFunctionBeginUser;
91*78a3110eSAlexander   PetscCall(PCShellGetContext(pc, &ctx));
92*78a3110eSAlexander 
93*78a3110eSAlexander   if (ctx->setup_called) PetscFunctionReturn(PETSC_SUCCESS);
94*78a3110eSAlexander 
95*78a3110eSAlexander   // Create aD
96*78a3110eSAlexander   PetscCall(MatDuplicate(ctx->D, MAT_COPY_VALUES, &ctx->aD));
97*78a3110eSAlexander   PetscCall(MatScale(ctx->aD, ctx->alpha));
98*78a3110eSAlexander 
99*78a3110eSAlexander   // Create aDinv
100*78a3110eSAlexander   PetscCall(MatDuplicate(ctx->aD, MAT_DO_NOT_COPY_VALUES, &ctx->aDinv));
101*78a3110eSAlexander   PetscCall(MatCreateVecs(ctx->aD, &aDVec, NULL));
102*78a3110eSAlexander   PetscCall(MatGetDiagonal(ctx->aD, aDVec));
103*78a3110eSAlexander   PetscCall(VecReciprocal(aDVec));
104*78a3110eSAlexander   PetscCall(MatDiagonalSet(ctx->aDinv, aDVec, INSERT_VALUES));
105*78a3110eSAlexander 
106*78a3110eSAlexander   // Create UT
107*78a3110eSAlexander   PetscCall(MatTranspose(ctx->U, MAT_INITIAL_MATRIX, &ctx->UT));
108*78a3110eSAlexander 
109*78a3110eSAlexander   // Create sum Mat
110*78a3110eSAlexander   PetscCall(MatMatMatMult(ctx->UT, ctx->aDinv, ctx->U, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &ctx->I_plus_gammaUTaDinvU));
111*78a3110eSAlexander   PetscCall(MatScale(ctx->I_plus_gammaUTaDinvU, ctx->gamma));
112*78a3110eSAlexander   PetscCall(MatShift(ctx->I_plus_gammaUTaDinvU, 1.));
113*78a3110eSAlexander 
114*78a3110eSAlexander   PetscCall(PCCreate(PETSC_COMM_WORLD, &ctx->smw_cholesky));
115*78a3110eSAlexander   PetscCall(PCSetType(ctx->smw_cholesky, PCCHOLESKY));
116*78a3110eSAlexander   PetscCall(PCSetOperators(ctx->smw_cholesky, ctx->I_plus_gammaUTaDinvU, ctx->I_plus_gammaUTaDinvU));
117*78a3110eSAlexander   PetscCall(PCSetOptionsPrefix(ctx->smw_cholesky, "smw_"));
118*78a3110eSAlexander   PetscCall(PCSetFromOptions(ctx->smw_cholesky));
119*78a3110eSAlexander   PetscCall(PCSetUp(ctx->smw_cholesky));
120*78a3110eSAlexander 
121*78a3110eSAlexander   PetscCall(VecDestroy(&aDVec));
122*78a3110eSAlexander 
123*78a3110eSAlexander   ctx->setup_called = PETSC_TRUE;
124*78a3110eSAlexander   PetscFunctionReturn(PETSC_SUCCESS);
125*78a3110eSAlexander }
126*78a3110eSAlexander 
127*78a3110eSAlexander PetscErrorCode SmwApply(PC pc, Vec x, Vec y)
128*78a3110eSAlexander {
129*78a3110eSAlexander   SmwPCCtx *ctx;
130*78a3110eSAlexander   Vec       vel0, pressure0, pressure1;
131*78a3110eSAlexander 
132*78a3110eSAlexander   PetscFunctionBeginUser;
133*78a3110eSAlexander   PetscCall(PCShellGetContext(pc, &ctx));
134*78a3110eSAlexander 
135*78a3110eSAlexander   PetscCall(MatCreateVecs(ctx->UT, &vel0, &pressure0));
136*78a3110eSAlexander   PetscCall(VecDuplicate(pressure0, &pressure1));
137*78a3110eSAlexander 
138*78a3110eSAlexander   // First term
139*78a3110eSAlexander   PetscCall(MatMult(ctx->aDinv, x, vel0));
140*78a3110eSAlexander   PetscCall(MatMult(ctx->UT, vel0, pressure0));
141*78a3110eSAlexander   PetscCall(PCApply(ctx->smw_cholesky, pressure0, pressure1));
142*78a3110eSAlexander   PetscCall(MatMult(ctx->U, pressure1, vel0));
143*78a3110eSAlexander   PetscCall(MatMult(ctx->aDinv, vel0, y));
144*78a3110eSAlexander   PetscCall(VecScale(y, -ctx->gamma));
145*78a3110eSAlexander 
146*78a3110eSAlexander   // Second term
147*78a3110eSAlexander   PetscCall(MatMult(ctx->aDinv, x, vel0));
148*78a3110eSAlexander 
149*78a3110eSAlexander   PetscCall(VecAXPY(y, 1, vel0));
150*78a3110eSAlexander 
151*78a3110eSAlexander   PetscCall(VecDestroy(&vel0));
152*78a3110eSAlexander   PetscCall(VecDestroy(&pressure0));
153*78a3110eSAlexander   PetscCall(VecDestroy(&pressure1));
154*78a3110eSAlexander   PetscFunctionReturn(PETSC_SUCCESS);
155*78a3110eSAlexander }
156*78a3110eSAlexander 
157*78a3110eSAlexander int main(int argc, char **args)
158*78a3110eSAlexander {
159*78a3110eSAlexander   Mat               A, B, Q, Acondensed, Bcondensed, BT, J, AplusJ, QInv, D, AplusD, JplusD, U;
160*78a3110eSAlexander   Mat               AplusJarray[2];
161*78a3110eSAlexander   Vec               bound, x, b, Qdiag, DVec;
162*78a3110eSAlexander   PetscBool         flg;
163*78a3110eSAlexander   PetscViewer       viewer;
164*78a3110eSAlexander   char              file[PETSC_MAX_PATH_LEN];
165*78a3110eSAlexander   PetscInt         *boundary_indices;
166*78a3110eSAlexander   PetscInt          boundary_indices_size, am, an, bm, bn, condensed_am, astart, aend, Dstart, Dend, num_local_bnd_dofs = 0;
167*78a3110eSAlexander   const PetscScalar zero = 0;
168*78a3110eSAlexander   IS                boundary_is, bulk_is;
169*78a3110eSAlexander   KSP               ksp;
170*78a3110eSAlexander   PC                pc, pcA, pcJ;
171*78a3110eSAlexander   PetscRandom       rctx;
172*78a3110eSAlexander   PetscReal        *boundary_indices_values;
173*78a3110eSAlexander   PetscReal         gamma = 100, alpha = .01;
174*78a3110eSAlexander   PetscMPIInt       rank;
175*78a3110eSAlexander   SmwPCCtx          ctx;
176*78a3110eSAlexander 
177*78a3110eSAlexander   PetscFunctionBeginUser;
178*78a3110eSAlexander   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
179*78a3110eSAlexander   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
180*78a3110eSAlexander 
181*78a3110eSAlexander   PetscCall(CreateAndLoadMat("A", &A));
182*78a3110eSAlexander   PetscCall(CreateAndLoadMat("B", &B));
183*78a3110eSAlexander   PetscCall(CreateAndLoadMat("Q", &Q));
184*78a3110eSAlexander 
185*78a3110eSAlexander   PetscCall(PetscOptionsGetString(NULL, NULL, "-fbound", file, sizeof(file), &flg));
186*78a3110eSAlexander   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate file with the -fbound option");
187*78a3110eSAlexander 
188*78a3110eSAlexander   if (rank == 0) {
189*78a3110eSAlexander     PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, file, FILE_MODE_READ, &viewer));
190*78a3110eSAlexander     PetscCall(VecCreate(PETSC_COMM_SELF, &bound));
191*78a3110eSAlexander     PetscCall(PetscObjectSetName((PetscObject)bound, "bound"));
192*78a3110eSAlexander     PetscCall(VecSetType(bound, VECSEQ));
193*78a3110eSAlexander     PetscCall(VecLoad(bound, viewer));
194*78a3110eSAlexander     PetscCall(PetscViewerDestroy(&viewer));
195*78a3110eSAlexander     PetscCall(VecGetLocalSize(bound, &boundary_indices_size));
196*78a3110eSAlexander     PetscCall(VecGetArray(bound, &boundary_indices_values));
197*78a3110eSAlexander   }
198*78a3110eSAlexander   PetscCallMPI(MPI_Bcast(&boundary_indices_size, 1, MPIU_INT, 0, PETSC_COMM_WORLD));
199*78a3110eSAlexander   if (rank != 0) PetscCall(PetscMalloc1(boundary_indices_size, &boundary_indices_values));
200*78a3110eSAlexander   PetscCallMPI(MPI_Bcast(boundary_indices_values, boundary_indices_size, MPIU_SCALAR, 0, PETSC_COMM_WORLD));
201*78a3110eSAlexander 
202*78a3110eSAlexander   PetscCall(MatGetSize(A, &am, NULL));
203*78a3110eSAlexander   // The total number of dofs for a given velocity component
204*78a3110eSAlexander   const PetscInt nc = am / 2;
205*78a3110eSAlexander   PetscCall(MatGetOwnershipRange(A, &astart, &aend));
206*78a3110eSAlexander 
207*78a3110eSAlexander   PetscCall(PetscMalloc1(2 * boundary_indices_size, &boundary_indices));
208*78a3110eSAlexander 
209*78a3110eSAlexander   //
210*78a3110eSAlexander   // The dof index ordering appears to be all vx dofs and then all vy dofs.
211*78a3110eSAlexander   //
212*78a3110eSAlexander 
213*78a3110eSAlexander   // First do vx
214*78a3110eSAlexander   for (PetscInt i = 0; i < boundary_indices_size; ++i) {
215*78a3110eSAlexander     // MATLAB uses 1-based indexing
216*78a3110eSAlexander     const PetscInt bnd_dof = (PetscInt)boundary_indices_values[i] - 1;
217*78a3110eSAlexander     if ((bnd_dof >= astart) && (bnd_dof < aend)) boundary_indices[num_local_bnd_dofs++] = bnd_dof;
218*78a3110eSAlexander   }
219*78a3110eSAlexander 
220*78a3110eSAlexander   // Now vy
221*78a3110eSAlexander   for (PetscInt i = 0; i < boundary_indices_size; ++i) {
222*78a3110eSAlexander     // MATLAB uses 1-based indexing
223*78a3110eSAlexander     const PetscInt bnd_dof = ((PetscInt)boundary_indices_values[i] - 1) + nc;
224*78a3110eSAlexander     if ((bnd_dof >= astart) && (bnd_dof < aend)) boundary_indices[num_local_bnd_dofs++] = bnd_dof;
225*78a3110eSAlexander   }
226*78a3110eSAlexander   if (rank == 0) PetscCall(VecRestoreArray(bound, &boundary_indices_values));
227*78a3110eSAlexander   else PetscCall(PetscFree(boundary_indices_values));
228*78a3110eSAlexander 
229*78a3110eSAlexander   PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, num_local_bnd_dofs, boundary_indices, PETSC_USE_POINTER, &boundary_is));
230*78a3110eSAlexander   PetscCall(ISComplement(boundary_is, astart, aend, &bulk_is));
231*78a3110eSAlexander 
232*78a3110eSAlexander   PetscCall(MatCreateSubMatrix(A, bulk_is, bulk_is, MAT_INITIAL_MATRIX, &Acondensed));
233*78a3110eSAlexander   // Can't pass null for row index set :-(
234*78a3110eSAlexander   PetscCall(MatTranspose(B, MAT_INPLACE_MATRIX, &B));
235*78a3110eSAlexander   PetscCall(MatCreateSubMatrix(B, bulk_is, NULL, MAT_INITIAL_MATRIX, &Bcondensed));
236*78a3110eSAlexander   PetscCall(MatGetLocalSize(Acondensed, &am, &an));
237*78a3110eSAlexander   PetscCall(MatGetLocalSize(Bcondensed, &bm, &bn));
238*78a3110eSAlexander 
239*78a3110eSAlexander   // Create QInv
240*78a3110eSAlexander   PetscCall(MatCreateVecs(Q, &Qdiag, NULL));
241*78a3110eSAlexander   PetscCall(MatGetDiagonal(Q, Qdiag));
242*78a3110eSAlexander   PetscCall(VecReciprocal(Qdiag));
243*78a3110eSAlexander   PetscCall(MatDuplicate(Q, MAT_DO_NOT_COPY_VALUES, &QInv));
244*78a3110eSAlexander   PetscCall(MatDiagonalSet(QInv, Qdiag, INSERT_VALUES));
245*78a3110eSAlexander   PetscCall(MatAssemblyBegin(QInv, MAT_FINAL_ASSEMBLY));
246*78a3110eSAlexander   PetscCall(MatAssemblyEnd(QInv, MAT_FINAL_ASSEMBLY));
247*78a3110eSAlexander 
248*78a3110eSAlexander   // Create J
249*78a3110eSAlexander   PetscCall(MatTranspose(Bcondensed, MAT_INITIAL_MATRIX, &BT));
250*78a3110eSAlexander   PetscCall(MatMatMatMult(Bcondensed, QInv, BT, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &J));
251*78a3110eSAlexander   PetscCall(MatScale(J, gamma));
252*78a3110eSAlexander 
253*78a3110eSAlexander   // Create sum of A + J
254*78a3110eSAlexander   AplusJarray[0] = Acondensed;
255*78a3110eSAlexander   AplusJarray[1] = J;
256*78a3110eSAlexander   PetscCall(MatCreateComposite(PETSC_COMM_WORLD, 2, AplusJarray, &AplusJ));
257*78a3110eSAlexander 
258*78a3110eSAlexander   // Create decomposition matrices
259*78a3110eSAlexander   // We've already used Qdiag, which currently represents Q^-1,  for its necessary purposes. Let's
260*78a3110eSAlexander   // convert it to represent Q^(-1/2)
261*78a3110eSAlexander   PetscCall(VecSqrtAbs(Qdiag));
262*78a3110eSAlexander   // We can similarly reuse Qinv
263*78a3110eSAlexander   PetscCall(MatDiagonalSet(QInv, Qdiag, INSERT_VALUES));
264*78a3110eSAlexander   PetscCall(MatAssemblyBegin(QInv, MAT_FINAL_ASSEMBLY));
265*78a3110eSAlexander   PetscCall(MatAssemblyEnd(QInv, MAT_FINAL_ASSEMBLY));
266*78a3110eSAlexander   // Create U
267*78a3110eSAlexander   PetscCall(MatMatMult(Bcondensed, QInv, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &U));
268*78a3110eSAlexander 
269*78a3110eSAlexander   // Create x and b
270*78a3110eSAlexander   PetscCall(MatCreateVecs(AplusJ, &x, &b));
271*78a3110eSAlexander   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &rctx));
272*78a3110eSAlexander   PetscCall(VecSetRandom(x, rctx));
273*78a3110eSAlexander   PetscCall(PetscRandomDestroy(&rctx));
274*78a3110eSAlexander   PetscCall(MatMult(AplusJ, x, b));
275*78a3110eSAlexander 
276*78a3110eSAlexander   // Compute preconditioner operators
277*78a3110eSAlexander   PetscCall(MatGetLocalSize(Acondensed, &condensed_am, NULL));
278*78a3110eSAlexander   PetscCall(MatCreate(PETSC_COMM_WORLD, &D));
279*78a3110eSAlexander   PetscCall(MatSetType(D, MATAIJ));
280*78a3110eSAlexander   PetscCall(MatSetSizes(D, condensed_am, condensed_am, PETSC_DETERMINE, PETSC_DETERMINE));
281*78a3110eSAlexander   PetscCall(MatGetOwnershipRange(D, &Dstart, &Dend));
282*78a3110eSAlexander   for (PetscInt i = Dstart; i < Dend; ++i) PetscCall(MatSetValues(D, 1, &i, 1, &i, &zero, INSERT_VALUES));
283*78a3110eSAlexander   PetscCall(MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY));
284*78a3110eSAlexander   PetscCall(MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY));
285*78a3110eSAlexander   PetscCall(MatCreateVecs(D, &DVec, NULL));
286*78a3110eSAlexander   PetscCall(MatGetDiagonal(AplusJ, DVec));
287*78a3110eSAlexander   PetscCall(MatDiagonalSet(D, DVec, INSERT_VALUES));
288*78a3110eSAlexander   PetscCall(MatDuplicate(Acondensed, MAT_COPY_VALUES, &AplusD));
289*78a3110eSAlexander   PetscCall(MatAXPY(AplusD, alpha, D, SUBSET_NONZERO_PATTERN));
290*78a3110eSAlexander   PetscCall(MatDuplicate(J, MAT_COPY_VALUES, &JplusD));
291*78a3110eSAlexander   PetscCall(MatAXPY(JplusD, alpha, D, SUBSET_NONZERO_PATTERN));
292*78a3110eSAlexander 
293*78a3110eSAlexander   // Initialize our SMW context
294*78a3110eSAlexander   ctx.U            = U;
295*78a3110eSAlexander   ctx.D            = D;
296*78a3110eSAlexander   ctx.gamma        = gamma;
297*78a3110eSAlexander   ctx.alpha        = alpha;
298*78a3110eSAlexander   ctx.setup_called = PETSC_FALSE;
299*78a3110eSAlexander 
300*78a3110eSAlexander   // Set preconditioner operators
301*78a3110eSAlexander   PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
302*78a3110eSAlexander   PetscCall(KSPSetType(ksp, KSPFGMRES));
303*78a3110eSAlexander   PetscCall(KSPSetOperators(ksp, AplusJ, AplusJ));
304*78a3110eSAlexander   PetscCall(KSPSetNormType(ksp, KSP_NORM_UNPRECONDITIONED));
305*78a3110eSAlexander   PetscCall(KSPGMRESSetRestart(ksp, 300));
306*78a3110eSAlexander   PetscCall(KSPGetPC(ksp, &pc));
307*78a3110eSAlexander   PetscCall(PCSetType(pc, PCCOMPOSITE));
308*78a3110eSAlexander   PetscCall(PCCompositeSetType(pc, PC_COMPOSITE_SPECIAL));
309*78a3110eSAlexander   PetscCall(PCCompositeAddPCType(pc, PCILU));
310*78a3110eSAlexander   PetscCall(PCCompositeAddPCType(pc, PCSHELL));
311*78a3110eSAlexander   PetscCall(PCCompositeGetPC(pc, 0, &pcA));
312*78a3110eSAlexander   PetscCall(PCCompositeGetPC(pc, 1, &pcJ));
313*78a3110eSAlexander   PetscCall(PCSetOperators(pcA, AplusD, AplusD));
314*78a3110eSAlexander   PetscCall(PCSetOperators(pcJ, JplusD, JplusD));
315*78a3110eSAlexander   PetscCall(PCShellSetContext(pcJ, &ctx));
316*78a3110eSAlexander   PetscCall(PCShellSetApply(pcJ, SmwApply));
317*78a3110eSAlexander   PetscCall(PCShellSetSetUp(pcJ, SmwSetup));
318*78a3110eSAlexander   PetscCall(PCCompositeSpecialSetAlphaMat(pc, D));
319*78a3110eSAlexander 
320*78a3110eSAlexander   // Solve
321*78a3110eSAlexander   PetscCall(KSPSetFromOptions(ksp));
322*78a3110eSAlexander   PetscCall(KSPSolve(ksp, b, x));
323*78a3110eSAlexander 
324*78a3110eSAlexander   PetscCall(MatDestroy(&A));
325*78a3110eSAlexander   PetscCall(MatDestroy(&B));
326*78a3110eSAlexander   PetscCall(MatDestroy(&Q));
327*78a3110eSAlexander   PetscCall(MatDestroy(&Acondensed));
328*78a3110eSAlexander   PetscCall(MatDestroy(&Bcondensed));
329*78a3110eSAlexander   PetscCall(MatDestroy(&BT));
330*78a3110eSAlexander   PetscCall(MatDestroy(&J));
331*78a3110eSAlexander   PetscCall(MatDestroy(&AplusJ));
332*78a3110eSAlexander   PetscCall(MatDestroy(&QInv));
333*78a3110eSAlexander   PetscCall(MatDestroy(&D));
334*78a3110eSAlexander   PetscCall(MatDestroy(&AplusD));
335*78a3110eSAlexander   PetscCall(MatDestroy(&JplusD));
336*78a3110eSAlexander   PetscCall(MatDestroy(&U));
337*78a3110eSAlexander   if (rank == 0) PetscCall(VecDestroy(&bound));
338*78a3110eSAlexander   PetscCall(VecDestroy(&x));
339*78a3110eSAlexander   PetscCall(VecDestroy(&b));
340*78a3110eSAlexander   PetscCall(VecDestroy(&Qdiag));
341*78a3110eSAlexander   PetscCall(VecDestroy(&DVec));
342*78a3110eSAlexander   PetscCall(ISDestroy(&boundary_is));
343*78a3110eSAlexander   PetscCall(ISDestroy(&bulk_is));
344*78a3110eSAlexander   PetscCall(KSPDestroy(&ksp));
345*78a3110eSAlexander   PetscCall(PetscFree(boundary_indices));
346*78a3110eSAlexander   PetscCall(MatDestroy(&ctx.UT));
347*78a3110eSAlexander   PetscCall(MatDestroy(&ctx.I_plus_gammaUTaDinvU));
348*78a3110eSAlexander   PetscCall(MatDestroy(&ctx.aD));
349*78a3110eSAlexander   PetscCall(MatDestroy(&ctx.aDinv));
350*78a3110eSAlexander   PetscCall(PCDestroy(&ctx.smw_cholesky));
351*78a3110eSAlexander 
352*78a3110eSAlexander   PetscCall(PetscFinalize());
353*78a3110eSAlexander   return 0;
354*78a3110eSAlexander }
355*78a3110eSAlexander 
356*78a3110eSAlexander /*TEST
357*78a3110eSAlexander 
358*78a3110eSAlexander    build:
359*78a3110eSAlexander       requires: !complex
360*78a3110eSAlexander 
361*78a3110eSAlexander    test:
362*78a3110eSAlexander       args: -fA ${DATAFILESPATH}/matrices/ifiss/A -fB ${DATAFILESPATH}/matrices/ifiss/B -fQ ${DATAFILESPATH}/matrices/ifiss/Q -fbound ${DATAFILESPATH}/is/ifiss/bound -ksp_monitor
363*78a3110eSAlexander       requires: datafilespath defined(PETSC_USE_64BIT_INDICES) !complex double
364*78a3110eSAlexander 
365*78a3110eSAlexander    test:
366*78a3110eSAlexander       suffix: 2
367*78a3110eSAlexander       nsize: 2
368*78a3110eSAlexander       args: -fA ${DATAFILESPATH}/matrices/ifiss/A -fB ${DATAFILESPATH}/matrices/ifiss/B -fQ ${DATAFILESPATH}/matrices/ifiss/Q -fbound ${DATAFILESPATH}/is/ifiss/bound -ksp_monitor
369*78a3110eSAlexander       requires: datafilespath defined(PETSC_USE_64BIT_INDICES) !complex double strumpack
370*78a3110eSAlexander 
371*78a3110eSAlexander TEST*/
372