Actual source code: nn.c

petsc-3.4.5 2014-06-29
  2: #include <../src/ksp/pc/impls/is/nn/nn.h>

  4: /* -------------------------------------------------------------------------- */
  5: /*
  6:    PCSetUp_NN - Prepares for the use of the NN preconditioner
  7:                     by setting data structures and options.

  9:    Input Parameter:
 10: .  pc - the preconditioner context

 12:    Application Interface Routine: PCSetUp()

 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: */
 20: static PetscErrorCode PCSetUp_NN(PC pc)
 21: {

 25:   if (!pc->setupcalled) {
 26:     /* Set up all the "iterative substructuring" common block */
 27:     PCISSetUp(pc);
 28:     /* Create the coarse matrix. */
 29:     PCNNCreateCoarseMatrix(pc);
 30:   }
 31:   return(0);
 32: }

 34: /* -------------------------------------------------------------------------- */
 35: /*
 36:    PCApply_NN - Applies the NN preconditioner to a vector.

 38:    Input Parameters:
 39: .  pc - the preconditioner context
 40: .  r - input vector (global)

 42:    Output Parameter:
 43: .  z - output vector (global)

 45:    Application Interface Routine: PCApply()
 46:  */
 49: static PetscErrorCode PCApply_NN(PC pc,Vec r,Vec z)
 50: {
 51:   PC_IS          *pcis = (PC_IS*)(pc->data);
 53:   PetscScalar    m_one = -1.0;
 54:   Vec            w     = pcis->vec1_global;

 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:   VecScatterBegin(pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
 63:   VecScatterEnd  (pcis->global_to_D,r,pcis->vec1_D,INSERT_VALUES,SCATTER_FORWARD);
 64:   KSPSolve(pcis->ksp_D,pcis->vec1_D,pcis->vec2_D);

 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:   MatMult(pcis->A_BI,pcis->vec2_D,pcis->vec1_B);
 71:   VecScale(pcis->vec1_B,m_one);
 72:   VecCopy(r,w);
 73:   VecScatterBegin(pcis->global_to_B,pcis->vec1_B,w,ADD_VALUES,SCATTER_REVERSE);
 74:   VecScatterEnd  (pcis->global_to_B,pcis->vec1_B,w,ADD_VALUES,SCATTER_REVERSE);

 76:   /*
 77:     Apply the interface preconditioner
 78:   */
 79:   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);

 82:   /*
 83:     Computing $ t_I^{(i)} = A_{IB}^{(i)} \tilde R_i z_B $
 84:     The result is stored in vec1_D.
 85:   */
 86:   VecScatterBegin(pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
 87:   VecScatterEnd  (pcis->global_to_B,z,pcis->vec1_B,INSERT_VALUES,SCATTER_FORWARD);
 88:   MatMult(pcis->A_IB,pcis->vec1_B,pcis->vec1_D);

 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:   VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
 96:   VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,INSERT_VALUES,SCATTER_REVERSE);
 97:   KSPSolve(pcis->ksp_D,pcis->vec1_D,pcis->vec2_D);
 98:   VecScale(pcis->vec2_D,m_one);
 99:   VecScatterBegin(pcis->global_to_D,pcis->vec2_D,z,ADD_VALUES,SCATTER_REVERSE);
100:   VecScatterEnd  (pcis->global_to_D,pcis->vec2_D,z,ADD_VALUES,SCATTER_REVERSE);
101:   return(0);
102: }

104: /* -------------------------------------------------------------------------- */
105: /*
106:    PCDestroy_NN - Destroys the private context for the NN preconditioner
107:    that was created with PCCreate_NN().

109:    Input Parameter:
110: .  pc - the preconditioner context

112:    Application Interface Routine: PCDestroy()
113: */
116: static PetscErrorCode PCDestroy_NN(PC pc)
117: {
118:   PC_NN          *pcnn = (PC_NN*)pc->data;

122:   PCISDestroy(pc);

124:   MatDestroy(&pcnn->coarse_mat);
125:   VecDestroy(&pcnn->coarse_x);
126:   VecDestroy(&pcnn->coarse_b);
127:   KSPDestroy(&pcnn->ksp_coarse);
128:   if (pcnn->DZ_IN) {
129:     PetscFree(pcnn->DZ_IN[0]);
130:     PetscFree(pcnn->DZ_IN);
131:   }

133:   /*
134:       Free the private data structure that was hanging off the PC
135:   */
136:   PetscFree(pc->data);
137:   return(0);
138: }

140: /* -------------------------------------------------------------------------- */
141: /*MC
142:    PCNN - Balancing Neumann-Neumann for scalar elliptic PDEs.

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 -

155:    Level: intermediate

157:    Notes: The matrix used with this preconditioner must be of type MATIS

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.).

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

167:    Contributed by Paulo Goldfeld

169: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,  MATIS
170: M*/

174: PETSC_EXTERN PetscErrorCode PCCreate_NN(PC pc)
175: {
177:   PC_NN          *pcnn;

180:   /*
181:      Creates the private data structure for this preconditioner and
182:      attach it to the PC object.
183:   */
184:   PetscNewLog(pc,PC_NN,&pcnn);
185:   pc->data = (void*)pcnn;

187:   PCISCreate(pc);
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;

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:   return(0);
210: }

212: /* -------------------------------------------------------------------------- */
213: /*
214:    PCNNCreateCoarseMatrix -
215: */
218: PetscErrorCode PCNNCreateCoarseMatrix(PC pc)
219: {
220:   MPI_Request    *send_request, *recv_request;
222:   PetscInt       i, j, k;
223:   PetscScalar    *mat;     /* Sub-matrix with this subdomain's contribution to the coarse matrix             */
224:   PetscScalar    **DZ_OUT; /* proc[k].DZ_OUT[i][] = bit of vector to be sent from processor k to processor i */

226:   /* aliasing some names */
227:   PC_IS       *pcis     = (PC_IS*)(pc->data);
228:   PC_NN       *pcnn     = (PC_NN*)pc->data;
229:   PetscInt    n_neigh   = pcis->n_neigh;
230:   PetscInt    *neigh    = pcis->neigh;
231:   PetscInt    *n_shared = pcis->n_shared;
232:   PetscInt    **shared  = pcis->shared;
233:   PetscScalar **DZ_IN;   /* Must be initialized after memory allocation. */

236:   /* Allocate memory for mat (the +1 is to handle the case n_neigh equal to zero) */
237:   PetscMalloc((n_neigh*n_neigh+1)*sizeof(PetscScalar),&mat);

239:   /* Allocate memory for DZ */
240:   /* Notice that DZ_OUT[0] is allocated some space that is never used. */
241:   /* This is just in order to DZ_OUT and DZ_IN to have exactly the same form. */
242:   {
243:     PetscInt size_of_Z = 0;
244:     PetscMalloc ((n_neigh+1)*sizeof(PetscScalar*),&pcnn->DZ_IN);
245:     DZ_IN = pcnn->DZ_IN;
246:     PetscMalloc ((n_neigh+1)*sizeof(PetscScalar*),&DZ_OUT);
247:     for (i=0; i<n_neigh; i++) size_of_Z += n_shared[i];
248:     PetscMalloc ((size_of_Z+1)*sizeof(PetscScalar),&DZ_IN[0]);
249:     PetscMalloc ((size_of_Z+1)*sizeof(PetscScalar),&DZ_OUT[0]);
250:   }
251:   for (i=1; i<n_neigh; i++) {
252:     DZ_IN[i]  = DZ_IN [i-1] + n_shared[i-1];
253:     DZ_OUT[i] = DZ_OUT[i-1] + n_shared[i-1];
254:   }

256:   /* Set the values of DZ_OUT, in order to send this info to the neighbours */
257:   /* First, set the auxiliary array pcis->work_N. */
258:   PCISScatterArrayNToVecB(pcis->work_N,pcis->D,INSERT_VALUES,SCATTER_REVERSE,pc);
259:   for (i=1; i<n_neigh; i++) {
260:     for (j=0; j<n_shared[i]; j++) {
261:       DZ_OUT[i][j] = pcis->work_N[shared[i][j]];
262:     }
263:   }

265:   /* Non-blocking send/receive the common-interface chunks of scaled nullspaces */
266:   /* Notice that send_request[] and recv_request[] could have one less element. */
267:   /* We make them longer to have request[i] corresponding to neigh[i].          */
268:   {
269:     PetscMPIInt tag;
270:     PetscObjectGetNewTag((PetscObject)pc,&tag);
271:     PetscMalloc((2*(n_neigh)+1)*sizeof(MPI_Request),&send_request);
272:     recv_request = send_request + (n_neigh);
273:     for (i=1; i<n_neigh; i++) {
274:       MPI_Isend((void*)(DZ_OUT[i]),n_shared[i],MPIU_SCALAR,neigh[i],tag,PetscObjectComm((PetscObject)pc),&(send_request[i]));
275:       MPI_Irecv((void*)(DZ_IN [i]),n_shared[i],MPIU_SCALAR,neigh[i],tag,PetscObjectComm((PetscObject)pc),&(recv_request[i]));
276:     }
277:   }

279:   /* Set DZ_IN[0][] (recall that neigh[0]==rank, always) */
280:   for (j=0; j<n_shared[0]; j++) DZ_IN[0][j] = pcis->work_N[shared[0][j]];

282:   /* Start computing with local D*Z while communication goes on.    */
283:   /* Apply Schur complement. The result is "stored" in vec (more    */
284:   /* precisely, vec points to the result, stored in pc_nn->vec1_B)  */
285:   /* and also scattered to pcnn->work_N.                            */
286:   PCNNApplySchurToChunk(pc,n_shared[0],shared[0],DZ_IN[0],pcis->work_N,pcis->vec1_B,
287:                                pcis->vec2_B,pcis->vec1_D,pcis->vec2_D);

289:   /* Compute the first column, while completing the receiving. */
290:   for (i=0; i<n_neigh; i++) {
291:     MPI_Status  stat;
292:     PetscMPIInt ind=0;
293:     if (i>0) { MPI_Waitany(n_neigh-1,recv_request+1,&ind,&stat); ind++;}
294:     mat[ind*n_neigh+0] = 0.0;
295:     for (k=0; k<n_shared[ind]; k++) mat[ind*n_neigh+0] += DZ_IN[ind][k] * pcis->work_N[shared[ind][k]];
296:   }

298:   /* Compute the remaining of the columns */
299:   for (j=1; j<n_neigh; j++) {
300:     PCNNApplySchurToChunk(pc,n_shared[j],shared[j],DZ_IN[j],pcis->work_N,pcis->vec1_B,
301:                                  pcis->vec2_B,pcis->vec1_D,pcis->vec2_D);
302:     for (i=0; i<n_neigh; i++) {
303:       mat[i*n_neigh+j] = 0.0;
304:       for (k=0; k<n_shared[i]; k++) mat[i*n_neigh+j] += DZ_IN[i][k] * pcis->work_N[shared[i][k]];
305:     }
306:   }

308:   /* Complete the sending. */
309:   if (n_neigh>1) {
310:     MPI_Status *stat;
311:     PetscMalloc((n_neigh-1)*sizeof(MPI_Status),&stat);
312:     if (n_neigh-1) {MPI_Waitall(n_neigh-1,&(send_request[1]),stat);}
313:     PetscFree(stat);
314:   }

316:   /* Free the memory for the MPI requests */
317:   PetscFree(send_request);

319:   /* Free the memory for DZ_OUT */
320:   if (DZ_OUT) {
321:     PetscFree(DZ_OUT[0]);
322:     PetscFree(DZ_OUT);
323:   }

325:   {
326:     PetscMPIInt size;
327:     MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
328:     /* Create the global coarse vectors (rhs and solution). */
329:     VecCreateMPI(PetscObjectComm((PetscObject)pc),1,size,&(pcnn->coarse_b));
330:     VecDuplicate(pcnn->coarse_b,&(pcnn->coarse_x));
331:     /* Create and set the global coarse AIJ matrix. */
332:     MatCreate(PetscObjectComm((PetscObject)pc),&(pcnn->coarse_mat));
333:     MatSetSizes(pcnn->coarse_mat,1,1,size,size);
334:     MatSetType(pcnn->coarse_mat,MATAIJ);
335:     MatSeqAIJSetPreallocation(pcnn->coarse_mat,1,NULL);
336:     MatMPIAIJSetPreallocation(pcnn->coarse_mat,1,NULL,1,NULL);
337:     MatSetValues(pcnn->coarse_mat,n_neigh,neigh,n_neigh,neigh,mat,ADD_VALUES);
338:     MatAssemblyBegin(pcnn->coarse_mat,MAT_FINAL_ASSEMBLY);
339:     MatAssemblyEnd  (pcnn->coarse_mat,MAT_FINAL_ASSEMBLY);
340:   }

342:   {
343:     PetscMPIInt rank;
344:     PetscScalar one = 1.0;
345:     MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
346:     /* "Zero out" rows of not-purely-Neumann subdomains */
347:     if (pcis->pure_neumann) {  /* does NOT zero the row; create an empty index set. The reason is that MatZeroRows() is collective. */
348:       MatZeroRows(pcnn->coarse_mat,0,NULL,one,0,0);
349:     } else { /* here it DOES zero the row, since it's not a floating subdomain. */
350:       PetscInt row = (PetscInt) rank;
351:       MatZeroRows(pcnn->coarse_mat,1,&row,one,0,0);
352:     }
353:   }

355:   /* Create the coarse linear solver context */
356:   {
357:     PC  pc_ctx, inner_pc;
358:     KSP inner_ksp;

360:     KSPCreate(PetscObjectComm((PetscObject)pc),&pcnn->ksp_coarse);
361:     PetscObjectIncrementTabLevel((PetscObject)pcnn->ksp_coarse,(PetscObject)pc,2);
362:     KSPSetOperators(pcnn->ksp_coarse,pcnn->coarse_mat,pcnn->coarse_mat,SAME_PRECONDITIONER);
363:     KSPGetPC(pcnn->ksp_coarse,&pc_ctx);
364:     PCSetType(pc_ctx,PCREDUNDANT);
365:     KSPSetType(pcnn->ksp_coarse,KSPPREONLY);
366:     PCRedundantGetKSP(pc_ctx,&inner_ksp);
367:     KSPGetPC(inner_ksp,&inner_pc);
368:     PCSetType(inner_pc,PCLU);
369:     KSPSetOptionsPrefix(pcnn->ksp_coarse,"nn_coarse_");
370:     KSPSetFromOptions(pcnn->ksp_coarse);
371:     /* the vectors in the following line are dummy arguments, just telling the KSP the vector size. Values are not used */
372:     KSPSetUp(pcnn->ksp_coarse);
373:   }

375:   /* Free the memory for mat */
376:   PetscFree(mat);

378:   /* for DEBUGGING, save the coarse matrix to a file. */
379:   {
380:     PetscBool flg = PETSC_FALSE;
381:     PetscOptionsGetBool(NULL,"-pc_nn_save_coarse_matrix",&flg,NULL);
382:     if (flg) {
383:       PetscViewer viewer;
384:       PetscViewerASCIIOpen(PETSC_COMM_WORLD,"coarse.m",&viewer);
385:       PetscViewerSetFormat(viewer,PETSC_VIEWER_ASCII_MATLAB);
386:       MatView(pcnn->coarse_mat,viewer);
387:       PetscViewerDestroy(&viewer);
388:     }
389:   }

391:   /*  Set the variable pcnn->factor_coarse_rhs. */
392:   pcnn->factor_coarse_rhs = (pcis->pure_neumann) ? 1.0 : 0.0;

394:   /* See historical note 02, at the bottom of this file. */
395:   return(0);
396: }

398: /* -------------------------------------------------------------------------- */
399: /*
400:    PCNNApplySchurToChunk -

402:    Input parameters:
403: .  pcnn
404: .  n - size of chunk
405: .  idx - indices of chunk
406: .  chunk - values

408:    Output parameters:
409: .  array_N - result of Schur complement applied to chunk, scattered to big array
410: .  vec1_B  - result of Schur complement applied to chunk
411: .  vec2_B  - garbage (used as work space)
412: .  vec1_D  - garbage (used as work space)
413: .  vec2_D  - garbage (used as work space)

415: */
418: 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)
419: {
421:   PetscInt       i;
422:   PC_IS          *pcis = (PC_IS*)(pc->data);

425:   PetscMemzero((void*)array_N, pcis->n*sizeof(PetscScalar));
426:   for (i=0; i<n; i++) array_N[idx[i]] = chunk[i];
427:   PCISScatterArrayNToVecB(array_N,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pc);
428:   PCISApplySchur(pc,vec2_B,vec1_B,(Vec)0,vec1_D,vec2_D);
429:   PCISScatterArrayNToVecB(array_N,vec1_B,INSERT_VALUES,SCATTER_REVERSE,pc);
430:   return(0);
431: }

433: /* -------------------------------------------------------------------------- */
434: /*
435:    PCNNApplyInterfacePreconditioner - Apply the interface preconditioner, i.e.,
436:                                       the preconditioner for the Schur complement.

438:    Input parameter:
439: .  r - global vector of interior and interface nodes. The values on the interior nodes are NOT used.

441:    Output parameters:
442: .  z - global vector of interior and interface nodes. The values on the interface are the result of
443:        the application of the interface preconditioner to the interface part of r. The values on the
444:        interior nodes are garbage.
445: .  work_N - array of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
446: .  vec1_B - vector of local interface nodes (including ghosts); returns garbage (used as work space)
447: .  vec2_B - vector of local interface nodes (including ghosts); returns garbage (used as work space)
448: .  vec3_B - vector of local interface nodes (including ghosts); returns garbage (used as work space)
449: .  vec1_D - vector of local interior nodes; returns garbage (used as work space)
450: .  vec2_D - vector of local interior nodes; returns garbage (used as work space)
451: .  vec1_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)
452: .  vec2_N - vector of local nodes (interior and interface, including ghosts); returns garbage (used as work space)

454: */
457: 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)
458: {
460:   PC_IS          *pcis = (PC_IS*)(pc->data);

463:   /*
464:     First balancing step.
465:   */
466:   {
467:     PetscBool flg = PETSC_FALSE;
468:     PetscOptionsGetBool(NULL,"-pc_nn_turn_off_first_balancing",&flg,NULL);
469:     if (!flg) {
470:       PCNNBalancing(pc,r,(Vec)0,z,vec1_B,vec2_B,(Vec)0,vec1_D,vec2_D,work_N);
471:     } else {
472:       VecCopy(r,z);
473:     }
474:   }

476:   /*
477:     Extract the local interface part of z and scale it by D
478:   */
479:   VecScatterBegin(pcis->global_to_B,z,vec1_B,INSERT_VALUES,SCATTER_FORWARD);
480:   VecScatterEnd  (pcis->global_to_B,z,vec1_B,INSERT_VALUES,SCATTER_FORWARD);
481:   VecPointwiseMult(vec2_B,pcis->D,vec1_B);

483:   /* Neumann Solver */
484:   PCISApplyInvSchur(pc,vec2_B,vec1_B,vec1_N,vec2_N);

486:   /*
487:     Second balancing step.
488:   */
489:   {
490:     PetscBool flg = PETSC_FALSE;
491:     PetscOptionsGetBool(NULL,"-pc_turn_off_second_balancing",&flg,NULL);
492:     if (!flg) {
493:       PCNNBalancing(pc,r,vec1_B,z,vec2_B,vec3_B,(Vec)0,vec1_D,vec2_D,work_N);
494:     } else {
495:       VecPointwiseMult(vec2_B,pcis->D,vec1_B);
496:       VecSet(z,0.0);
497:       VecScatterBegin(pcis->global_to_B,vec2_B,z,ADD_VALUES,SCATTER_REVERSE);
498:       VecScatterEnd  (pcis->global_to_B,vec2_B,z,ADD_VALUES,SCATTER_REVERSE);
499:     }
500:   }
501:   return(0);
502: }

504: /* -------------------------------------------------------------------------- */
505: /*
506:    PCNNBalancing - Computes z, as given in equations (15) and (16) (if the
507:                    input argument u is provided), or s, as given in equations
508:                    (12) and (13), if the input argument u is a null vector.
509:                    Notice that the input argument u plays the role of u_i in
510:                    equation (14). The equation numbers refer to [Man93].

512:    Input Parameters:
513: .  pcnn - NN preconditioner context.
514: .  r - MPI vector of all nodes (interior and interface). It's preserved.
515: .  u - (Optional) sequential vector of local interface nodes. It's preserved UNLESS vec3_B is null.

517:    Output Parameters:
518: .  z - MPI vector of interior and interface nodes. Returns s or z (see description above).
519: .  vec1_B - Sequential vector of local interface nodes. Workspace.
520: .  vec2_B - Sequential vector of local interface nodes. Workspace.
521: .  vec3_B - (Optional) sequential vector of local interface nodes. Workspace.
522: .  vec1_D - Sequential vector of local interior nodes. Workspace.
523: .  vec2_D - Sequential vector of local interior nodes. Workspace.
524: .  work_N - Array of all local nodes (interior and interface). Workspace.

526: */
529: 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)
530: {
532:   PetscInt       k;
533:   PetscScalar    value;
534:   PetscScalar    *lambda;
535:   PC_NN          *pcnn = (PC_NN*)(pc->data);
536:   PC_IS          *pcis = (PC_IS*)(pc->data);

539:   PetscLogEventBegin(PC_ApplyCoarse,0,0,0,0);
540:   if (u) {
541:     if (!vec3_B) vec3_B = u;
542:     VecPointwiseMult(vec1_B,pcis->D,u);
543:     VecSet(z,0.0);
544:     VecScatterBegin(pcis->global_to_B,vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
545:     VecScatterEnd  (pcis->global_to_B,vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
546:     VecScatterBegin(pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD);
547:     VecScatterEnd  (pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD);
548:     PCISApplySchur(pc,vec2_B,vec3_B,(Vec)0,vec1_D,vec2_D);
549:     VecScale(vec3_B,-1.0);
550:     VecCopy(r,z);
551:     VecScatterBegin(pcis->global_to_B,vec3_B,z,ADD_VALUES,SCATTER_REVERSE);
552:     VecScatterEnd  (pcis->global_to_B,vec3_B,z,ADD_VALUES,SCATTER_REVERSE);
553:   } else {
554:     VecCopy(r,z);
555:   }
556:   VecScatterBegin(pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD);
557:   VecScatterEnd  (pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD);
558:   PCISScatterArrayNToVecB(work_N,vec2_B,INSERT_VALUES,SCATTER_REVERSE,pc);
559:   for (k=0, value=0.0; k<pcis->n_shared[0]; k++) value += pcnn->DZ_IN[0][k] * work_N[pcis->shared[0][k]];
560:   value *= pcnn->factor_coarse_rhs;  /* This factor is set in CreateCoarseMatrix(). */
561:   {
562:     PetscMPIInt rank;
563:     MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
564:     VecSetValue(pcnn->coarse_b,rank,value,INSERT_VALUES);
565:     /*
566:        Since we are only inserting local values (one value actually) we don't need to do the
567:        reduction that tells us there is no data that needs to be moved. Hence we comment out these
568:        VecAssemblyBegin(pcnn->coarse_b);
569:        VecAssemblyEnd  (pcnn->coarse_b);
570:     */
571:   }
572:   KSPSolve(pcnn->ksp_coarse,pcnn->coarse_b,pcnn->coarse_x);
573:   if (!u) { VecScale(pcnn->coarse_x,-1.0); }
574:   VecGetArray(pcnn->coarse_x,&lambda);
575:   for (k=0; k<pcis->n_shared[0]; k++) work_N[pcis->shared[0][k]] = *lambda * pcnn->DZ_IN[0][k];
576:   VecRestoreArray(pcnn->coarse_x,&lambda);
577:   PCISScatterArrayNToVecB(work_N,vec2_B,INSERT_VALUES,SCATTER_FORWARD,pc);
578:   VecSet(z,0.0);
579:   VecScatterBegin(pcis->global_to_B,vec2_B,z,ADD_VALUES,SCATTER_REVERSE);
580:   VecScatterEnd  (pcis->global_to_B,vec2_B,z,ADD_VALUES,SCATTER_REVERSE);
581:   if (!u) {
582:     VecScatterBegin(pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD);
583:     VecScatterEnd  (pcis->global_to_B,z,vec2_B,INSERT_VALUES,SCATTER_FORWARD);
584:     PCISApplySchur(pc,vec2_B,vec1_B,(Vec)0,vec1_D,vec2_D);
585:     VecCopy(r,z);
586:   }
587:   VecScatterBegin(pcis->global_to_B,vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
588:   VecScatterEnd  (pcis->global_to_B,vec1_B,z,ADD_VALUES,SCATTER_REVERSE);
589:   PetscLogEventEnd(PC_ApplyCoarse,0,0,0,0);
590:   return(0);
591: }




597: /*  -------   E N D   O F   T H E   C O D E   -------  */
598: /*                                                     */
599: /*  From now on, "footnotes" (or "historical notes").  */
600: /*                                                     */
601: /*  -------------------------------------------------  */



605: /* --------------------------------------------------------------------------
606:    Historical note 01
607:    -------------------------------------------------------------------------- */
608: /*
609:    We considered the possibility of an alternative D_i that would still
610:    provide a partition of unity (i.e., $ \sum_i  N_i D_i N_i^T = I $).
611:    The basic principle was still the pseudo-inverse of the counting
612:    function; the difference was that we would not count subdomains
613:    that do not contribute to the coarse space (i.e., not pure-Neumann
614:    subdomains).

616:    This turned out to be a bad idea:  we would solve trivial Neumann
617:    problems in the not pure-Neumann subdomains, since we would be scaling
618:    the balanced residual by zero.
619: */




624: /* --------------------------------------------------------------------------
625:    Historical note 02
626:    -------------------------------------------------------------------------- */
627: /*
628:    We tried an alternative coarse problem, that would eliminate exactly a
629:    constant error. Turned out not to improve the overall convergence.
630: */