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