xref: /petsc/src/ksp/pc/impls/is/nn/nn.c (revision d2522c19e8fa9bca20aaca277941d9a63e71db6a)
1 
2 #include <../src/ksp/pc/impls/is/nn/nn.h>
3 
4 /* -------------------------------------------------------------------------- */
5 /*
6    PCSetUp_NN - Prepares for the use of the NN preconditioner
7                     by setting data structures and options.
8 
9    Input Parameter:
10 .  pc - the preconditioner context
11 
12    Application Interface Routine: PCSetUp()
13 
14    Notes:
15    The interface routine PCSetUp() is not usually called directly by
16    the user, but instead is called by PCApply() if necessary.
17 */
18 static PetscErrorCode PCSetUp_NN(PC pc) {
19   PetscFunctionBegin;
20   if (!pc->setupcalled) {
21     /* Set up all the "iterative substructuring" common block */
22     PetscCall(PCISSetUp(pc, PETSC_TRUE, PETSC_TRUE));
23     /* Create the coarse matrix. */
24     PetscCall(PCNNCreateCoarseMatrix(pc));
25   }
26   PetscFunctionReturn(0);
27 }
28 
29 /* -------------------------------------------------------------------------- */
30 /*
31    PCApply_NN - Applies the NN preconditioner to a vector.
32 
33    Input Parameters:
34 .  pc - the preconditioner context
35 .  r - input vector (global)
36 
37    Output Parameter:
38 .  z - output vector (global)
39 
40    Application Interface Routine: PCApply()
41  */
42 static PetscErrorCode PCApply_NN(PC pc, Vec r, Vec z) {
43   PC_IS      *pcis  = (PC_IS *)(pc->data);
44   PetscScalar m_one = -1.0;
45   Vec         w     = pcis->vec1_global;
46 
47   PetscFunctionBegin;
48   /*
49     Dirichlet solvers.
50     Solving $ B_I^{(i)}r_I^{(i)} $ at each processor.
51     Storing the local results at vec2_D
52   */
53   PetscCall(VecScatterBegin(pcis->global_to_D, r, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD));
54   PetscCall(VecScatterEnd(pcis->global_to_D, r, pcis->vec1_D, INSERT_VALUES, SCATTER_FORWARD));
55   PetscCall(KSPSolve(pcis->ksp_D, pcis->vec1_D, pcis->vec2_D));
56 
57   /*
58     Computing $ r_B - \sum_j \tilde R_j^T A_{BI}^{(j)} (B_I^{(j)}r_I^{(j)}) $ .
59     Storing the result in the interface portion of the global vector w.
60   */
61   PetscCall(MatMult(pcis->A_BI, pcis->vec2_D, pcis->vec1_B));
62   PetscCall(VecScale(pcis->vec1_B, m_one));
63   PetscCall(VecCopy(r, w));
64   PetscCall(VecScatterBegin(pcis->global_to_B, pcis->vec1_B, w, ADD_VALUES, SCATTER_REVERSE));
65   PetscCall(VecScatterEnd(pcis->global_to_B, pcis->vec1_B, w, ADD_VALUES, SCATTER_REVERSE));
66 
67   /*
68     Apply the interface preconditioner
69   */
70   PetscCall(PCNNApplyInterfacePreconditioner(pc, w, z, pcis->work_N, pcis->vec1_B, pcis->vec2_B, pcis->vec3_B, pcis->vec1_D, pcis->vec3_D, pcis->vec1_N, pcis->vec2_N));
71 
72   /*
73     Computing $ t_I^{(i)} = A_{IB}^{(i)} \tilde R_i z_B $
74     The result is stored in vec1_D.
75   */
76   PetscCall(VecScatterBegin(pcis->global_to_B, z, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
77   PetscCall(VecScatterEnd(pcis->global_to_B, z, pcis->vec1_B, INSERT_VALUES, SCATTER_FORWARD));
78   PetscCall(MatMult(pcis->A_IB, pcis->vec1_B, pcis->vec1_D));
79 
80   /*
81     Dirichlet solvers.
82     Computing $ B_I^{(i)}t_I^{(i)} $ and sticking into the global vector the blocks
83     $ B_I^{(i)}r_I^{(i)} - B_I^{(i)}t_I^{(i)} $.
84   */
85   PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec2_D, z, INSERT_VALUES, SCATTER_REVERSE));
86   PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec2_D, z, INSERT_VALUES, SCATTER_REVERSE));
87   PetscCall(KSPSolve(pcis->ksp_D, pcis->vec1_D, pcis->vec2_D));
88   PetscCall(VecScale(pcis->vec2_D, m_one));
89   PetscCall(VecScatterBegin(pcis->global_to_D, pcis->vec2_D, z, ADD_VALUES, SCATTER_REVERSE));
90   PetscCall(VecScatterEnd(pcis->global_to_D, pcis->vec2_D, z, ADD_VALUES, SCATTER_REVERSE));
91   PetscFunctionReturn(0);
92 }
93 
94 /* -------------------------------------------------------------------------- */
95 /*
96    PCDestroy_NN - Destroys the private context for the NN preconditioner
97    that was created with PCCreate_NN().
98 
99    Input Parameter:
100 .  pc - the preconditioner context
101 
102    Application Interface Routine: PCDestroy()
103 */
104 static PetscErrorCode PCDestroy_NN(PC pc) {
105   PC_NN *pcnn = (PC_NN *)pc->data;
106 
107   PetscFunctionBegin;
108   PetscCall(PCISDestroy(pc));
109 
110   PetscCall(MatDestroy(&pcnn->coarse_mat));
111   PetscCall(VecDestroy(&pcnn->coarse_x));
112   PetscCall(VecDestroy(&pcnn->coarse_b));
113   PetscCall(KSPDestroy(&pcnn->ksp_coarse));
114   if (pcnn->DZ_IN) {
115     PetscCall(PetscFree(pcnn->DZ_IN[0]));
116     PetscCall(PetscFree(pcnn->DZ_IN));
117   }
118 
119   /*
120       Free the private data structure that was hanging off the PC
121   */
122   PetscCall(PetscFree(pc->data));
123   PetscFunctionReturn(0);
124 }
125 
126 /* -------------------------------------------------------------------------- */
127 /*MC
128    PCNN - Balancing Neumann-Neumann for scalar elliptic PDEs.
129 
130    Options Database Keys:
131 +    -pc_nn_turn_off_first_balancing - do not balance the residual before solving the local Neumann problems
132                                        (this skips the first coarse grid solve in the preconditioner)
133 .    -pc_nn_turn_off_second_balancing - do not balance the solution solving the local Neumann problems
134                                        (this skips the second coarse grid solve in the preconditioner)
135 .    -pc_is_damp_fixed <fact> -
136 .    -pc_is_remove_nullspace_fixed -
137 .    -pc_is_set_damping_factor_floating <fact> -
138 .    -pc_is_not_damp_floating -
139 -    -pc_is_not_remove_nullspace_floating -
140 
141    Level: intermediate
142 
143    Notes:
144     The matrix used with this preconditioner must be of type MATIS
145 
146           Unlike more 'conventional' Neumann-Neumann preconditioners this iterates over ALL the
147           degrees of freedom, NOT just those on the interface (this allows the use of approximate solvers
148           on the subdomains; though in our experience using approximate solvers is slower.).
149 
150           Options for the coarse grid preconditioner can be set with -nn_coarse_pc_xxx
151           Options for the Dirichlet subproblem preconditioner can be set with -is_localD_pc_xxx
152           Options for the Neumann subproblem preconditioner can be set with -is_localN_pc_xxx
153 
154    Contributed by Paulo Goldfeld
155 
156 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `MATIS`
157 M*/
158 
159 PETSC_EXTERN PetscErrorCode PCCreate_NN(PC pc) {
160   PC_NN *pcnn;
161 
162   PetscFunctionBegin;
163   /*
164      Creates the private data structure for this preconditioner and
165      attach it to the PC object.
166   */
167   PetscCall(PetscNewLog(pc, &pcnn));
168   pc->data = (void *)pcnn;
169 
170   PetscCall(PCISCreate(pc));
171   pcnn->coarse_mat = NULL;
172   pcnn->coarse_x   = NULL;
173   pcnn->coarse_b   = NULL;
174   pcnn->ksp_coarse = NULL;
175   pcnn->DZ_IN      = NULL;
176 
177   /*
178       Set the pointers for the functions that are provided above.
179       Now when the user-level routines (such as PCApply(), PCDestroy(), etc.)
180       are called, they will automatically call these functions.  Note we
181       choose not to provide a couple of these functions since they are
182       not needed.
183   */
184   pc->ops->apply               = PCApply_NN;
185   pc->ops->applytranspose      = NULL;
186   pc->ops->setup               = PCSetUp_NN;
187   pc->ops->destroy             = PCDestroy_NN;
188   pc->ops->view                = NULL;
189   pc->ops->applyrichardson     = NULL;
190   pc->ops->applysymmetricleft  = NULL;
191   pc->ops->applysymmetricright = NULL;
192   PetscFunctionReturn(0);
193 }
194 
195 /* -------------------------------------------------------------------------- */
196 /*
197    PCNNCreateCoarseMatrix -
198 */
199 PetscErrorCode PCNNCreateCoarseMatrix(PC pc) {
200   MPI_Request  *send_request, *recv_request;
201   PetscInt      i, j, k;
202   PetscScalar  *mat;    /* Sub-matrix with this subdomain's contribution to the coarse matrix             */
203   PetscScalar **DZ_OUT; /* proc[k].DZ_OUT[i][] = bit of vector to be sent from processor k to processor i */
204 
205   /* aliasing some names */
206   PC_IS        *pcis     = (PC_IS *)(pc->data);
207   PC_NN        *pcnn     = (PC_NN *)pc->data;
208   PetscInt      n_neigh  = pcis->n_neigh;
209   PetscInt     *neigh    = pcis->neigh;
210   PetscInt     *n_shared = pcis->n_shared;
211   PetscInt    **shared   = pcis->shared;
212   PetscScalar **DZ_IN; /* Must be initialized after memory allocation. */
213 
214   PetscFunctionBegin;
215   /* Allocate memory for mat (the +1 is to handle the case n_neigh equal to zero) */
216   PetscCall(PetscMalloc1(n_neigh * n_neigh + 1, &mat));
217 
218   /* Allocate memory for DZ */
219   /* Notice that DZ_OUT[0] is allocated some space that is never used. */
220   /* This is just in order to DZ_OUT and DZ_IN to have exactly the same form. */
221   {
222     PetscInt size_of_Z = 0;
223     PetscCall(PetscMalloc((n_neigh + 1) * sizeof(PetscScalar *), &pcnn->DZ_IN));
224     DZ_IN = pcnn->DZ_IN;
225     PetscCall(PetscMalloc((n_neigh + 1) * sizeof(PetscScalar *), &DZ_OUT));
226     for (i = 0; i < n_neigh; i++) size_of_Z += n_shared[i];
227     PetscCall(PetscMalloc((size_of_Z + 1) * sizeof(PetscScalar), &DZ_IN[0]));
228     PetscCall(PetscMalloc((size_of_Z + 1) * sizeof(PetscScalar), &DZ_OUT[0]));
229   }
230   for (i = 1; i < n_neigh; i++) {
231     DZ_IN[i]  = DZ_IN[i - 1] + n_shared[i - 1];
232     DZ_OUT[i] = DZ_OUT[i - 1] + n_shared[i - 1];
233   }
234 
235   /* Set the values of DZ_OUT, in order to send this info to the neighbours */
236   /* First, set the auxiliary array pcis->work_N. */
237   PetscCall(PCISScatterArrayNToVecB(pcis->work_N, pcis->D, INSERT_VALUES, SCATTER_REVERSE, pc));
238   for (i = 1; i < n_neigh; i++) {
239     for (j = 0; j < n_shared[i]; j++) { DZ_OUT[i][j] = pcis->work_N[shared[i][j]]; }
240   }
241 
242   /* Non-blocking send/receive the common-interface chunks of scaled nullspaces */
243   /* Notice that send_request[] and recv_request[] could have one less element. */
244   /* We make them longer to have request[i] corresponding to neigh[i].          */
245   {
246     PetscMPIInt tag;
247     PetscCall(PetscObjectGetNewTag((PetscObject)pc, &tag));
248     PetscCall(PetscMalloc2(n_neigh + 1, &send_request, n_neigh + 1, &recv_request));
249     for (i = 1; i < n_neigh; i++) {
250       PetscCallMPI(MPI_Isend((void *)(DZ_OUT[i]), n_shared[i], MPIU_SCALAR, neigh[i], tag, PetscObjectComm((PetscObject)pc), &(send_request[i])));
251       PetscCallMPI(MPI_Irecv((void *)(DZ_IN[i]), n_shared[i], MPIU_SCALAR, neigh[i], tag, PetscObjectComm((PetscObject)pc), &(recv_request[i])));
252     }
253   }
254 
255   /* Set DZ_IN[0][] (recall that neigh[0]==rank, always) */
256   for (j = 0; j < n_shared[0]; j++) DZ_IN[0][j] = pcis->work_N[shared[0][j]];
257 
258   /* Start computing with local D*Z while communication goes on.    */
259   /* Apply Schur complement. The result is "stored" in vec (more    */
260   /* precisely, vec points to the result, stored in pc_nn->vec1_B)  */
261   /* and also scattered to pcnn->work_N.                            */
262   PetscCall(PCNNApplySchurToChunk(pc, n_shared[0], shared[0], DZ_IN[0], pcis->work_N, pcis->vec1_B, pcis->vec2_B, pcis->vec1_D, pcis->vec2_D));
263 
264   /* Compute the first column, while completing the receiving. */
265   for (i = 0; i < n_neigh; i++) {
266     MPI_Status  stat;
267     PetscMPIInt ind = 0;
268     if (i > 0) {
269       PetscCallMPI(MPI_Waitany(n_neigh - 1, recv_request + 1, &ind, &stat));
270       ind++;
271     }
272     mat[ind * n_neigh + 0] = 0.0;
273     for (k = 0; k < n_shared[ind]; k++) mat[ind * n_neigh + 0] += DZ_IN[ind][k] * pcis->work_N[shared[ind][k]];
274   }
275 
276   /* Compute the remaining of the columns */
277   for (j = 1; j < n_neigh; j++) {
278     PetscCall(PCNNApplySchurToChunk(pc, n_shared[j], shared[j], DZ_IN[j], pcis->work_N, pcis->vec1_B, pcis->vec2_B, pcis->vec1_D, pcis->vec2_D));
279     for (i = 0; i < n_neigh; i++) {
280       mat[i * n_neigh + j] = 0.0;
281       for (k = 0; k < n_shared[i]; k++) mat[i * n_neigh + j] += DZ_IN[i][k] * pcis->work_N[shared[i][k]];
282     }
283   }
284 
285   /* Complete the sending. */
286   if (n_neigh > 1) {
287     MPI_Status *stat;
288     PetscCall(PetscMalloc1(n_neigh - 1, &stat));
289     if (n_neigh - 1) PetscCallMPI(MPI_Waitall(n_neigh - 1, &(send_request[1]), stat));
290     PetscCall(PetscFree(stat));
291   }
292 
293   /* Free the memory for the MPI requests */
294   PetscCall(PetscFree2(send_request, recv_request));
295 
296   /* Free the memory for DZ_OUT */
297   if (DZ_OUT) {
298     PetscCall(PetscFree(DZ_OUT[0]));
299     PetscCall(PetscFree(DZ_OUT));
300   }
301 
302   {
303     PetscMPIInt size;
304     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
305     /* Create the global coarse vectors (rhs and solution). */
306     PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)pc), 1, size, &(pcnn->coarse_b)));
307     PetscCall(VecDuplicate(pcnn->coarse_b, &(pcnn->coarse_x)));
308     /* Create and set the global coarse AIJ matrix. */
309     PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &(pcnn->coarse_mat)));
310     PetscCall(MatSetSizes(pcnn->coarse_mat, 1, 1, size, size));
311     PetscCall(MatSetType(pcnn->coarse_mat, MATAIJ));
312     PetscCall(MatSeqAIJSetPreallocation(pcnn->coarse_mat, 1, NULL));
313     PetscCall(MatMPIAIJSetPreallocation(pcnn->coarse_mat, 1, NULL, n_neigh, NULL));
314     PetscCall(MatSetOption(pcnn->coarse_mat, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
315     PetscCall(MatSetOption(pcnn->coarse_mat, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE));
316     PetscCall(MatSetValues(pcnn->coarse_mat, n_neigh, neigh, n_neigh, neigh, mat, ADD_VALUES));
317     PetscCall(MatAssemblyBegin(pcnn->coarse_mat, MAT_FINAL_ASSEMBLY));
318     PetscCall(MatAssemblyEnd(pcnn->coarse_mat, MAT_FINAL_ASSEMBLY));
319   }
320 
321   {
322     PetscMPIInt rank;
323     PetscScalar one = 1.0;
324     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
325     /* "Zero out" rows of not-purely-Neumann subdomains */
326     if (pcis->pure_neumann) { /* does NOT zero the row; create an empty index set. The reason is that MatZeroRows() is collective. */
327       PetscCall(MatZeroRows(pcnn->coarse_mat, 0, NULL, one, NULL, NULL));
328     } else { /* here it DOES zero the row, since it's not a floating subdomain. */
329       PetscInt row = (PetscInt)rank;
330       PetscCall(MatZeroRows(pcnn->coarse_mat, 1, &row, one, NULL, NULL));
331     }
332   }
333 
334   /* Create the coarse linear solver context */
335   {
336     PC  pc_ctx, inner_pc;
337     KSP inner_ksp;
338 
339     PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &pcnn->ksp_coarse));
340     PetscCall(PetscObjectIncrementTabLevel((PetscObject)pcnn->ksp_coarse, (PetscObject)pc, 2));
341     PetscCall(KSPSetOperators(pcnn->ksp_coarse, pcnn->coarse_mat, pcnn->coarse_mat));
342     PetscCall(KSPGetPC(pcnn->ksp_coarse, &pc_ctx));
343     PetscCall(PCSetType(pc_ctx, PCREDUNDANT));
344     PetscCall(KSPSetType(pcnn->ksp_coarse, KSPPREONLY));
345     PetscCall(PCRedundantGetKSP(pc_ctx, &inner_ksp));
346     PetscCall(KSPGetPC(inner_ksp, &inner_pc));
347     PetscCall(PCSetType(inner_pc, PCLU));
348     PetscCall(KSPSetOptionsPrefix(pcnn->ksp_coarse, "nn_coarse_"));
349     PetscCall(KSPSetFromOptions(pcnn->ksp_coarse));
350     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
351     PetscCall(KSPSetUp(pcnn->ksp_coarse));
352   }
353 
354   /* Free the memory for mat */
355   PetscCall(PetscFree(mat));
356 
357   /* for DEBUGGING, save the coarse matrix to a file. */
358   {
359     PetscBool flg = PETSC_FALSE;
360     PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_nn_save_coarse_matrix", &flg, NULL));
361     if (flg) {
362       PetscViewer viewer;
363       PetscCall(PetscViewerASCIIOpen(PETSC_COMM_WORLD, "coarse.m", &viewer));
364       PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB));
365       PetscCall(MatView(pcnn->coarse_mat, viewer));
366       PetscCall(PetscViewerPopFormat(viewer));
367       PetscCall(PetscViewerDestroy(&viewer));
368     }
369   }
370 
371   /*  Set the variable pcnn->factor_coarse_rhs. */
372   pcnn->factor_coarse_rhs = (pcis->pure_neumann) ? 1.0 : 0.0;
373 
374   /* See historical note 02, at the bottom of this file. */
375   PetscFunctionReturn(0);
376 }
377 
378 /* -------------------------------------------------------------------------- */
379 /*
380    PCNNApplySchurToChunk -
381 
382    Input parameters:
383 .  pcnn
384 .  n - size of chunk
385 .  idx - indices of chunk
386 .  chunk - values
387 
388    Output parameters:
389 .  array_N - result of Schur complement applied to chunk, scattered to big array
390 .  vec1_B  - result of Schur complement applied to chunk
391 .  vec2_B  - garbage (used as work space)
392 .  vec1_D  - garbage (used as work space)
393 .  vec2_D  - garbage (used as work space)
394 
395 */
396 PetscErrorCode PCNNApplySchurToChunk(PC pc, PetscInt n, PetscInt *idx, PetscScalar *chunk, PetscScalar *array_N, Vec vec1_B, Vec vec2_B, Vec vec1_D, Vec vec2_D) {
397   PetscInt i;
398   PC_IS   *pcis = (PC_IS *)(pc->data);
399 
400   PetscFunctionBegin;
401   PetscCall(PetscArrayzero(array_N, pcis->n));
402   for (i = 0; i < n; i++) array_N[idx[i]] = chunk[i];
403   PetscCall(PCISScatterArrayNToVecB(array_N, vec2_B, INSERT_VALUES, SCATTER_FORWARD, pc));
404   PetscCall(PCISApplySchur(pc, vec2_B, vec1_B, (Vec)0, vec1_D, vec2_D));
405   PetscCall(PCISScatterArrayNToVecB(array_N, vec1_B, INSERT_VALUES, SCATTER_REVERSE, pc));
406   PetscFunctionReturn(0);
407 }
408 
409 /* -------------------------------------------------------------------------- */
410 /*
411    PCNNApplyInterfacePreconditioner - Apply the interface preconditioner, i.e.,
412                                       the preconditioner for the Schur complement.
413 
414    Input parameter:
415 .  r - global vector of interior and interface nodes. The values on the interior nodes are NOT used.
416 
417    Output parameters:
418 .  z - global vector of interior and interface nodes. The values on the interface are the result of
419        the application of the interface preconditioner to the interface part of r. The values on the
420        interior nodes are garbage.
421 .  work_N - array of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
422 .  vec1_B - vector of local interface nodes (including ghosts); returns garbage (used as work space)
423 .  vec2_B - vector of local interface nodes (including ghosts); returns garbage (used as work space)
424 .  vec3_B - vector of local interface nodes (including ghosts); returns garbage (used as work space)
425 .  vec1_D - vector of local interior nodes; returns garbage (used as work space)
426 .  vec2_D - vector of local interior nodes; returns garbage (used as work space)
427 .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
428 .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
429 
430 */
431 PetscErrorCode PCNNApplyInterfacePreconditioner(PC pc, Vec r, Vec z, PetscScalar *work_N, Vec vec1_B, Vec vec2_B, Vec vec3_B, Vec vec1_D, Vec vec2_D, Vec vec1_N, Vec vec2_N) {
432   PC_IS *pcis = (PC_IS *)(pc->data);
433 
434   PetscFunctionBegin;
435   /*
436     First balancing step.
437   */
438   {
439     PetscBool flg = PETSC_FALSE;
440     PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_nn_turn_off_first_balancing", &flg, NULL));
441     if (!flg) {
442       PetscCall(PCNNBalancing(pc, r, (Vec)0, z, vec1_B, vec2_B, (Vec)0, vec1_D, vec2_D, work_N));
443     } else {
444       PetscCall(VecCopy(r, z));
445     }
446   }
447 
448   /*
449     Extract the local interface part of z and scale it by D
450   */
451   PetscCall(VecScatterBegin(pcis->global_to_B, z, vec1_B, INSERT_VALUES, SCATTER_FORWARD));
452   PetscCall(VecScatterEnd(pcis->global_to_B, z, vec1_B, INSERT_VALUES, SCATTER_FORWARD));
453   PetscCall(VecPointwiseMult(vec2_B, pcis->D, vec1_B));
454 
455   /* Neumann Solver */
456   PetscCall(PCISApplyInvSchur(pc, vec2_B, vec1_B, vec1_N, vec2_N));
457 
458   /*
459     Second balancing step.
460   */
461   {
462     PetscBool flg = PETSC_FALSE;
463     PetscCall(PetscOptionsGetBool(NULL, NULL, "-pc_turn_off_second_balancing", &flg, NULL));
464     if (!flg) {
465       PetscCall(PCNNBalancing(pc, r, vec1_B, z, vec2_B, vec3_B, (Vec)0, vec1_D, vec2_D, work_N));
466     } else {
467       PetscCall(VecPointwiseMult(vec2_B, pcis->D, vec1_B));
468       PetscCall(VecSet(z, 0.0));
469       PetscCall(VecScatterBegin(pcis->global_to_B, vec2_B, z, ADD_VALUES, SCATTER_REVERSE));
470       PetscCall(VecScatterEnd(pcis->global_to_B, vec2_B, z, ADD_VALUES, SCATTER_REVERSE));
471     }
472   }
473   PetscFunctionReturn(0);
474 }
475 
476 /* -------------------------------------------------------------------------- */
477 /*
478    PCNNBalancing - Computes z, as given in equations (15) and (16) (if the
479                    input argument u is provided), or s, as given in equations
480                    (12) and (13), if the input argument u is a null vector.
481                    Notice that the input argument u plays the role of u_i in
482                    equation (14). The equation numbers refer to [Man93].
483 
484    Input Parameters:
485 .  pcnn - NN preconditioner context.
486 .  r - MPI vector of all nodes (interior and interface). It's preserved.
487 .  u - (Optional) sequential vector of local interface nodes. It's preserved UNLESS vec3_B is null.
488 
489    Output Parameters:
490 .  z - MPI vector of interior and interface nodes. Returns s or z (see description above).
491 .  vec1_B - Sequential vector of local interface nodes. Workspace.
492 .  vec2_B - Sequential vector of local interface nodes. Workspace.
493 .  vec3_B - (Optional) sequential vector of local interface nodes. Workspace.
494 .  vec1_D - Sequential vector of local interior nodes. Workspace.
495 .  vec2_D - Sequential vector of local interior nodes. Workspace.
496 .  work_N - Array of all local nodes (interior and interface). Workspace.
497 
498 */
499 PetscErrorCode PCNNBalancing(PC pc, Vec r, Vec u, Vec z, Vec vec1_B, Vec vec2_B, Vec vec3_B, Vec vec1_D, Vec vec2_D, PetscScalar *work_N) {
500   PetscInt     k;
501   PetscScalar  value;
502   PetscScalar *lambda;
503   PC_NN       *pcnn = (PC_NN *)(pc->data);
504   PC_IS       *pcis = (PC_IS *)(pc->data);
505 
506   PetscFunctionBegin;
507   PetscCall(PetscLogEventBegin(PC_ApplyCoarse, pc, 0, 0, 0));
508   if (u) {
509     if (!vec3_B) vec3_B = u;
510     PetscCall(VecPointwiseMult(vec1_B, pcis->D, u));
511     PetscCall(VecSet(z, 0.0));
512     PetscCall(VecScatterBegin(pcis->global_to_B, vec1_B, z, ADD_VALUES, SCATTER_REVERSE));
513     PetscCall(VecScatterEnd(pcis->global_to_B, vec1_B, z, ADD_VALUES, SCATTER_REVERSE));
514     PetscCall(VecScatterBegin(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD));
515     PetscCall(VecScatterEnd(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD));
516     PetscCall(PCISApplySchur(pc, vec2_B, vec3_B, (Vec)0, vec1_D, vec2_D));
517     PetscCall(VecScale(vec3_B, -1.0));
518     PetscCall(VecCopy(r, z));
519     PetscCall(VecScatterBegin(pcis->global_to_B, vec3_B, z, ADD_VALUES, SCATTER_REVERSE));
520     PetscCall(VecScatterEnd(pcis->global_to_B, vec3_B, z, ADD_VALUES, SCATTER_REVERSE));
521   } else {
522     PetscCall(VecCopy(r, z));
523   }
524   PetscCall(VecScatterBegin(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD));
525   PetscCall(VecScatterEnd(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD));
526   PetscCall(PCISScatterArrayNToVecB(work_N, vec2_B, INSERT_VALUES, SCATTER_REVERSE, pc));
527   for (k = 0, value = 0.0; k < pcis->n_shared[0]; k++) value += pcnn->DZ_IN[0][k] * work_N[pcis->shared[0][k]];
528   value *= pcnn->factor_coarse_rhs; /* This factor is set in CreateCoarseMatrix(). */
529   {
530     PetscMPIInt rank;
531     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
532     PetscCall(VecSetValue(pcnn->coarse_b, rank, value, INSERT_VALUES));
533     /*
534        Since we are only inserting local values (one value actually) we don't need to do the
535        reduction that tells us there is no data that needs to be moved. Hence we comment out these
536        PetscCall(VecAssemblyBegin(pcnn->coarse_b));
537        PetscCall(VecAssemblyEnd  (pcnn->coarse_b));
538     */
539   }
540   PetscCall(KSPSolve(pcnn->ksp_coarse, pcnn->coarse_b, pcnn->coarse_x));
541   if (!u) PetscCall(VecScale(pcnn->coarse_x, -1.0));
542   PetscCall(VecGetArray(pcnn->coarse_x, &lambda));
543   for (k = 0; k < pcis->n_shared[0]; k++) work_N[pcis->shared[0][k]] = *lambda * pcnn->DZ_IN[0][k];
544   PetscCall(VecRestoreArray(pcnn->coarse_x, &lambda));
545   PetscCall(PCISScatterArrayNToVecB(work_N, vec2_B, INSERT_VALUES, SCATTER_FORWARD, pc));
546   PetscCall(VecSet(z, 0.0));
547   PetscCall(VecScatterBegin(pcis->global_to_B, vec2_B, z, ADD_VALUES, SCATTER_REVERSE));
548   PetscCall(VecScatterEnd(pcis->global_to_B, vec2_B, z, ADD_VALUES, SCATTER_REVERSE));
549   if (!u) {
550     PetscCall(VecScatterBegin(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD));
551     PetscCall(VecScatterEnd(pcis->global_to_B, z, vec2_B, INSERT_VALUES, SCATTER_FORWARD));
552     PetscCall(PCISApplySchur(pc, vec2_B, vec1_B, (Vec)0, vec1_D, vec2_D));
553     PetscCall(VecCopy(r, z));
554   }
555   PetscCall(VecScatterBegin(pcis->global_to_B, vec1_B, z, ADD_VALUES, SCATTER_REVERSE));
556   PetscCall(VecScatterEnd(pcis->global_to_B, vec1_B, z, ADD_VALUES, SCATTER_REVERSE));
557   PetscCall(PetscLogEventEnd(PC_ApplyCoarse, pc, 0, 0, 0));
558   PetscFunctionReturn(0);
559 }
560 
561 /*  -------   E N D   O F   T H E   C O D E   -------  */
562 /*                                                     */
563 /*  From now on, "footnotes" (or "historical notes").  */
564 /*                                                     */
565 /*  -------------------------------------------------  */
566 
567 /* --------------------------------------------------------------------------
568    Historical note 01
569    -------------------------------------------------------------------------- */
570 /*
571    We considered the possibility of an alternative D_i that would still
572    provide a partition of unity (i.e., $ \sum_i  N_i D_i N_i^T = I $).
573    The basic principle was still the pseudo-inverse of the counting
574    function; the difference was that we would not count subdomains
575    that do not contribute to the coarse space (i.e., not pure-Neumann
576    subdomains).
577 
578    This turned out to be a bad idea:  we would solve trivial Neumann
579    problems in the not pure-Neumann subdomains, since we would be scaling
580    the balanced residual by zero.
581 */
582 
583 /* --------------------------------------------------------------------------
584    Historical note 02
585    -------------------------------------------------------------------------- */
586 /*
587    We tried an alternative coarse problem, that would eliminate exactly a
588    constant error. Turned out not to improve the overall convergence.
589 */
590