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