xref: /petsc/src/ksp/pc/impls/redundant/redundant.c (revision 41417bb58ccc2dfdcf1f98b94d1e30a9a42a4648)
1 
2 /*
3   This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner.
4 */
5 #include <petsc/private/pcimpl.h>
6 #include <petscksp.h>           /*I "petscksp.h" I*/
7 
8 typedef struct {
9   KSP                ksp;
10   PC                 pc;                   /* actual preconditioner used on each processor */
11   Vec                xsub,ysub;            /* vectors of a subcommunicator to hold parallel vectors of PetscObjectComm((PetscObject)pc) */
12   Vec                xdup,ydup;            /* parallel vector that congregates xsub or ysub facilitating vector scattering */
13   Mat                pmats;                /* matrix and optional preconditioner matrix belong to a subcommunicator */
14   VecScatter         scatterin,scatterout; /* scatter used to move all values to each processor group (subcommunicator) */
15   PetscBool          useparallelmat;
16   PetscSubcomm       psubcomm;
17   PetscInt           nsubcomm;             /* num of data structure PetscSubcomm */
18   PetscBool          shifttypeset;
19   MatFactorShiftType shifttype;
20 } PC_Redundant;
21 
22 PetscErrorCode  PCFactorSetShiftType_Redundant(PC pc,MatFactorShiftType shifttype)
23 {
24   PC_Redundant   *red = (PC_Redundant*)pc->data;
25 
26   PetscFunctionBegin;
27   if (red->ksp) {
28     PC pc;
29     PetscCall(KSPGetPC(red->ksp,&pc));
30     PetscCall(PCFactorSetShiftType(pc,shifttype));
31   } else {
32     red->shifttypeset = PETSC_TRUE;
33     red->shifttype    = shifttype;
34   }
35   PetscFunctionReturn(0);
36 }
37 
38 static PetscErrorCode PCView_Redundant(PC pc,PetscViewer viewer)
39 {
40   PC_Redundant   *red = (PC_Redundant*)pc->data;
41   PetscBool      iascii,isstring;
42   PetscViewer    subviewer;
43 
44   PetscFunctionBegin;
45   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
46   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring));
47   if (iascii) {
48     if (!red->psubcomm) {
49       PetscCall(PetscViewerASCIIPrintf(viewer,"  Not yet setup\n"));
50     } else {
51       PetscCall(PetscViewerASCIIPrintf(viewer,"  First (color=0) of %" PetscInt_FMT " PCs follows\n",red->nsubcomm));
52       PetscCall(PetscViewerGetSubViewer(viewer,((PetscObject)red->pc)->comm,&subviewer));
53       if (!red->psubcomm->color) { /* only view first redundant pc */
54         PetscCall(PetscViewerASCIIPushTab(subviewer));
55         PetscCall(KSPView(red->ksp,subviewer));
56         PetscCall(PetscViewerASCIIPopTab(subviewer));
57       }
58       PetscCall(PetscViewerRestoreSubViewer(viewer,((PetscObject)red->pc)->comm,&subviewer));
59     }
60   } else if (isstring) {
61     PetscCall(PetscViewerStringSPrintf(viewer," Redundant solver preconditioner"));
62   }
63   PetscFunctionReturn(0);
64 }
65 
66 #include <../src/mat/impls/aij/mpi/mpiaij.h>
67 static PetscErrorCode PCSetUp_Redundant(PC pc)
68 {
69   PC_Redundant   *red = (PC_Redundant*)pc->data;
70   PetscInt       mstart,mend,mlocal,M;
71   PetscMPIInt    size;
72   MPI_Comm       comm,subcomm;
73   Vec            x;
74 
75   PetscFunctionBegin;
76   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
77 
78   /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */
79   PetscCallMPI(MPI_Comm_size(comm,&size));
80   if (size == 1) red->useparallelmat = PETSC_FALSE;
81 
82   if (!pc->setupcalled) {
83     PetscInt mloc_sub;
84     if (!red->psubcomm) { /* create red->psubcomm, new ksp and pc over subcomm */
85       KSP ksp;
86       PetscCall(PCRedundantGetKSP(pc,&ksp));
87     }
88     subcomm = PetscSubcommChild(red->psubcomm);
89 
90     if (red->useparallelmat) {
91       /* grab the parallel matrix and put it into processors of a subcomminicator */
92       PetscCall(MatCreateRedundantMatrix(pc->pmat,red->psubcomm->n,subcomm,MAT_INITIAL_MATRIX,&red->pmats));
93 
94       PetscCallMPI(MPI_Comm_size(subcomm,&size));
95       if (size > 1) {
96         PetscBool foundpack,issbaij;
97         PetscCall(PetscObjectTypeCompare((PetscObject)red->pmats,MATMPISBAIJ,&issbaij));
98         if (!issbaij) {
99           PetscCall(MatGetFactorAvailable(red->pmats,NULL,MAT_FACTOR_LU,&foundpack));
100         } else {
101           PetscCall(MatGetFactorAvailable(red->pmats,NULL,MAT_FACTOR_CHOLESKY,&foundpack));
102         }
103         if (!foundpack) { /* reset default ksp and pc */
104           PetscCall(KSPSetType(red->ksp,KSPGMRES));
105           PetscCall(PCSetType(red->pc,PCBJACOBI));
106         } else {
107           PetscCall(PCFactorSetMatSolverType(red->pc,NULL));
108         }
109       }
110 
111       PetscCall(KSPSetOperators(red->ksp,red->pmats,red->pmats));
112 
113       /* get working vectors xsub and ysub */
114       PetscCall(MatCreateVecs(red->pmats,&red->xsub,&red->ysub));
115 
116       /* create working vectors xdup and ydup.
117        xdup concatenates all xsub's contigously to form a mpi vector over dupcomm  (see PetscSubcommCreate_interlaced())
118        ydup concatenates all ysub and has empty local arrays because ysub's arrays will be place into it.
119        Note: we use communicator dupcomm, not PetscObjectComm((PetscObject)pc)! */
120       PetscCall(MatGetLocalSize(red->pmats,&mloc_sub,NULL));
121       PetscCall(VecCreateMPI(PetscSubcommContiguousParent(red->psubcomm),mloc_sub,PETSC_DECIDE,&red->xdup));
122       PetscCall(VecCreateMPIWithArray(PetscSubcommContiguousParent(red->psubcomm),1,mloc_sub,PETSC_DECIDE,NULL,&red->ydup));
123 
124       /* create vecscatters */
125       if (!red->scatterin) { /* efficiency of scatterin is independent from psubcomm_type! */
126         IS       is1,is2;
127         PetscInt *idx1,*idx2,i,j,k;
128 
129         PetscCall(MatCreateVecs(pc->pmat,&x,NULL));
130         PetscCall(VecGetSize(x,&M));
131         PetscCall(VecGetOwnershipRange(x,&mstart,&mend));
132         mlocal = mend - mstart;
133         PetscCall(PetscMalloc2(red->psubcomm->n*mlocal,&idx1,red->psubcomm->n*mlocal,&idx2));
134         j    = 0;
135         for (k=0; k<red->psubcomm->n; k++) {
136           for (i=mstart; i<mend; i++) {
137             idx1[j]   = i;
138             idx2[j++] = i + M*k;
139           }
140         }
141         PetscCall(ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1));
142         PetscCall(ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2));
143         PetscCall(VecScatterCreate(x,is1,red->xdup,is2,&red->scatterin));
144         PetscCall(ISDestroy(&is1));
145         PetscCall(ISDestroy(&is2));
146 
147         /* Impl below is good for PETSC_SUBCOMM_INTERLACED (no inter-process communication) and PETSC_SUBCOMM_CONTIGUOUS (communication within subcomm) */
148         PetscCall(ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*M,1,&is1));
149         PetscCall(ISCreateStride(comm,mlocal,mstart,1,&is2));
150         PetscCall(VecScatterCreate(red->xdup,is1,x,is2,&red->scatterout));
151         PetscCall(ISDestroy(&is1));
152         PetscCall(ISDestroy(&is2));
153         PetscCall(PetscFree2(idx1,idx2));
154         PetscCall(VecDestroy(&x));
155       }
156     } else { /* !red->useparallelmat */
157       PetscCall(KSPSetOperators(red->ksp,pc->mat,pc->pmat));
158     }
159   } else { /* pc->setupcalled */
160     if (red->useparallelmat) {
161       MatReuse       reuse;
162       /* grab the parallel matrix and put it into processors of a subcomminicator */
163       /*--------------------------------------------------------------------------*/
164       if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
165         /* destroy old matrices */
166         PetscCall(MatDestroy(&red->pmats));
167         reuse = MAT_INITIAL_MATRIX;
168       } else {
169         reuse = MAT_REUSE_MATRIX;
170       }
171       PetscCall(MatCreateRedundantMatrix(pc->pmat,red->psubcomm->n,PetscSubcommChild(red->psubcomm),reuse,&red->pmats));
172       PetscCall(KSPSetOperators(red->ksp,red->pmats,red->pmats));
173     } else { /* !red->useparallelmat */
174       PetscCall(KSPSetOperators(red->ksp,pc->mat,pc->pmat));
175     }
176   }
177 
178   if (pc->setfromoptionscalled) {
179     PetscCall(KSPSetFromOptions(red->ksp));
180   }
181   PetscCall(KSPSetUp(red->ksp));
182   PetscFunctionReturn(0);
183 }
184 
185 static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y)
186 {
187   PC_Redundant   *red = (PC_Redundant*)pc->data;
188   PetscScalar    *array;
189 
190   PetscFunctionBegin;
191   if (!red->useparallelmat) {
192     PetscCall(KSPSolve(red->ksp,x,y));
193     PetscCall(KSPCheckSolve(red->ksp,pc,y));
194     PetscFunctionReturn(0);
195   }
196 
197   /* scatter x to xdup */
198   PetscCall(VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD));
199   PetscCall(VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD));
200 
201   /* place xdup's local array into xsub */
202   PetscCall(VecGetArray(red->xdup,&array));
203   PetscCall(VecPlaceArray(red->xsub,(const PetscScalar*)array));
204 
205   /* apply preconditioner on each processor */
206   PetscCall(KSPSolve(red->ksp,red->xsub,red->ysub));
207   PetscCall(KSPCheckSolve(red->ksp,pc,red->ysub));
208   PetscCall(VecResetArray(red->xsub));
209   PetscCall(VecRestoreArray(red->xdup,&array));
210 
211   /* place ysub's local array into ydup */
212   PetscCall(VecGetArray(red->ysub,&array));
213   PetscCall(VecPlaceArray(red->ydup,(const PetscScalar*)array));
214 
215   /* scatter ydup to y */
216   PetscCall(VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD));
217   PetscCall(VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD));
218   PetscCall(VecResetArray(red->ydup));
219   PetscCall(VecRestoreArray(red->ysub,&array));
220   PetscFunctionReturn(0);
221 }
222 
223 static PetscErrorCode PCApplyTranspose_Redundant(PC pc,Vec x,Vec y)
224 {
225   PC_Redundant   *red = (PC_Redundant*)pc->data;
226   PetscScalar    *array;
227 
228   PetscFunctionBegin;
229   if (!red->useparallelmat) {
230     PetscCall(KSPSolveTranspose(red->ksp,x,y));
231     PetscCall(KSPCheckSolve(red->ksp,pc,y));
232     PetscFunctionReturn(0);
233   }
234 
235   /* scatter x to xdup */
236   PetscCall(VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD));
237   PetscCall(VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD));
238 
239   /* place xdup's local array into xsub */
240   PetscCall(VecGetArray(red->xdup,&array));
241   PetscCall(VecPlaceArray(red->xsub,(const PetscScalar*)array));
242 
243   /* apply preconditioner on each processor */
244   PetscCall(KSPSolveTranspose(red->ksp,red->xsub,red->ysub));
245   PetscCall(KSPCheckSolve(red->ksp,pc,red->ysub));
246   PetscCall(VecResetArray(red->xsub));
247   PetscCall(VecRestoreArray(red->xdup,&array));
248 
249   /* place ysub's local array into ydup */
250   PetscCall(VecGetArray(red->ysub,&array));
251   PetscCall(VecPlaceArray(red->ydup,(const PetscScalar*)array));
252 
253   /* scatter ydup to y */
254   PetscCall(VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD));
255   PetscCall(VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD));
256   PetscCall(VecResetArray(red->ydup));
257   PetscCall(VecRestoreArray(red->ysub,&array));
258   PetscFunctionReturn(0);
259 }
260 
261 static PetscErrorCode PCReset_Redundant(PC pc)
262 {
263   PC_Redundant   *red = (PC_Redundant*)pc->data;
264 
265   PetscFunctionBegin;
266   if (red->useparallelmat) {
267     PetscCall(VecScatterDestroy(&red->scatterin));
268     PetscCall(VecScatterDestroy(&red->scatterout));
269     PetscCall(VecDestroy(&red->ysub));
270     PetscCall(VecDestroy(&red->xsub));
271     PetscCall(VecDestroy(&red->xdup));
272     PetscCall(VecDestroy(&red->ydup));
273   }
274   PetscCall(MatDestroy(&red->pmats));
275   PetscCall(KSPReset(red->ksp));
276   PetscFunctionReturn(0);
277 }
278 
279 static PetscErrorCode PCDestroy_Redundant(PC pc)
280 {
281   PC_Redundant   *red = (PC_Redundant*)pc->data;
282 
283   PetscFunctionBegin;
284   PetscCall(PCReset_Redundant(pc));
285   PetscCall(KSPDestroy(&red->ksp));
286   PetscCall(PetscSubcommDestroy(&red->psubcomm));
287   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetScatter_C",NULL));
288   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetNumber_C",NULL));
289   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetKSP_C",NULL));
290   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetOperators_C",NULL));
291   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",NULL));
292   PetscCall(PetscFree(pc->data));
293   PetscFunctionReturn(0);
294 }
295 
296 static PetscErrorCode PCSetFromOptions_Redundant(PetscOptionItems *PetscOptionsObject,PC pc)
297 {
298   PC_Redundant   *red = (PC_Redundant*)pc->data;
299 
300   PetscFunctionBegin;
301   PetscOptionsHeadBegin(PetscOptionsObject,"Redundant options");
302   PetscCall(PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,NULL));
303   PetscOptionsHeadEnd();
304   PetscFunctionReturn(0);
305 }
306 
307 static PetscErrorCode PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds)
308 {
309   PC_Redundant *red = (PC_Redundant*)pc->data;
310 
311   PetscFunctionBegin;
312   red->nsubcomm = nreds;
313   PetscFunctionReturn(0);
314 }
315 
316 /*@
317    PCRedundantSetNumber - Sets the number of redundant preconditioner contexts.
318 
319    Logically Collective on PC
320 
321    Input Parameters:
322 +  pc - the preconditioner context
323 -  nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and
324                               use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
325 
326    Level: advanced
327 
328 @*/
329 PetscErrorCode PCRedundantSetNumber(PC pc,PetscInt nredundant)
330 {
331   PetscFunctionBegin;
332   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
333   PetscCheck(nredundant > 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "num of redundant pc %" PetscInt_FMT " must be positive",nredundant);
334   PetscTryMethod(pc,"PCRedundantSetNumber_C",(PC,PetscInt),(pc,nredundant));
335   PetscFunctionReturn(0);
336 }
337 
338 static PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out)
339 {
340   PC_Redundant   *red = (PC_Redundant*)pc->data;
341 
342   PetscFunctionBegin;
343   PetscCall(PetscObjectReference((PetscObject)in));
344   PetscCall(VecScatterDestroy(&red->scatterin));
345 
346   red->scatterin  = in;
347 
348   PetscCall(PetscObjectReference((PetscObject)out));
349   PetscCall(VecScatterDestroy(&red->scatterout));
350   red->scatterout = out;
351   PetscFunctionReturn(0);
352 }
353 
354 /*@
355    PCRedundantSetScatter - Sets the scatter used to copy values into the
356      redundant local solve and the scatter to move them back into the global
357      vector.
358 
359    Logically Collective on PC
360 
361    Input Parameters:
362 +  pc - the preconditioner context
363 .  in - the scatter to move the values in
364 -  out - the scatter to move them out
365 
366    Level: advanced
367 
368 @*/
369 PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out)
370 {
371   PetscFunctionBegin;
372   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
373   PetscValidHeaderSpecific(in,PETSCSF_CLASSID,2);
374   PetscValidHeaderSpecific(out,PETSCSF_CLASSID,3);
375   PetscTryMethod(pc,"PCRedundantSetScatter_C",(PC,VecScatter,VecScatter),(pc,in,out));
376   PetscFunctionReturn(0);
377 }
378 
379 static PetscErrorCode PCRedundantGetKSP_Redundant(PC pc,KSP *innerksp)
380 {
381   PC_Redundant   *red = (PC_Redundant*)pc->data;
382   MPI_Comm       comm,subcomm;
383   const char     *prefix;
384   PetscBool      issbaij;
385 
386   PetscFunctionBegin;
387   if (!red->psubcomm) {
388     PetscCall(PCGetOptionsPrefix(pc,&prefix));
389 
390     PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
391     PetscCall(PetscSubcommCreate(comm,&red->psubcomm));
392     PetscCall(PetscSubcommSetNumber(red->psubcomm,red->nsubcomm));
393     PetscCall(PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_CONTIGUOUS));
394 
395     PetscCall(PetscSubcommSetOptionsPrefix(red->psubcomm,prefix));
396     PetscCall(PetscSubcommSetFromOptions(red->psubcomm));
397     PetscCall(PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm)));
398 
399     /* create a new PC that processors in each subcomm have copy of */
400     subcomm = PetscSubcommChild(red->psubcomm);
401 
402     PetscCall(KSPCreate(subcomm,&red->ksp));
403     PetscCall(KSPSetErrorIfNotConverged(red->ksp,pc->erroriffailure));
404     PetscCall(PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1));
405     PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)red->ksp));
406     PetscCall(KSPSetType(red->ksp,KSPPREONLY));
407     PetscCall(KSPGetPC(red->ksp,&red->pc));
408     PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQSBAIJ,&issbaij));
409     if (!issbaij) {
410       PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPISBAIJ,&issbaij));
411     }
412     if (!issbaij) {
413       PetscCall(PCSetType(red->pc,PCLU));
414     } else {
415       PetscCall(PCSetType(red->pc,PCCHOLESKY));
416     }
417     if (red->shifttypeset) {
418       PetscCall(PCFactorSetShiftType(red->pc,red->shifttype));
419       red->shifttypeset = PETSC_FALSE;
420     }
421     PetscCall(KSPSetOptionsPrefix(red->ksp,prefix));
422     PetscCall(KSPAppendOptionsPrefix(red->ksp,"redundant_"));
423   }
424   *innerksp = red->ksp;
425   PetscFunctionReturn(0);
426 }
427 
428 /*@
429    PCRedundantGetKSP - Gets the less parallel KSP created by the redundant PC.
430 
431    Not Collective
432 
433    Input Parameter:
434 .  pc - the preconditioner context
435 
436    Output Parameter:
437 .  innerksp - the KSP on the smaller set of processes
438 
439    Level: advanced
440 
441 @*/
442 PetscErrorCode PCRedundantGetKSP(PC pc,KSP *innerksp)
443 {
444   PetscFunctionBegin;
445   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
446   PetscValidPointer(innerksp,2);
447   PetscUseMethod(pc,"PCRedundantGetKSP_C",(PC,KSP*),(pc,innerksp));
448   PetscFunctionReturn(0);
449 }
450 
451 static PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat)
452 {
453   PC_Redundant *red = (PC_Redundant*)pc->data;
454 
455   PetscFunctionBegin;
456   if (mat)  *mat  = red->pmats;
457   if (pmat) *pmat = red->pmats;
458   PetscFunctionReturn(0);
459 }
460 
461 /*@
462    PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix
463 
464    Not Collective
465 
466    Input Parameter:
467 .  pc - the preconditioner context
468 
469    Output Parameters:
470 +  mat - the matrix
471 -  pmat - the (possibly different) preconditioner matrix
472 
473    Level: advanced
474 
475 @*/
476 PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat)
477 {
478   PetscFunctionBegin;
479   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
480   if (mat)  PetscValidPointer(mat,2);
481   if (pmat) PetscValidPointer(pmat,3);
482   PetscUseMethod(pc,"PCRedundantGetOperators_C",(PC,Mat*,Mat*),(pc,mat,pmat));
483   PetscFunctionReturn(0);
484 }
485 
486 /* -------------------------------------------------------------------------------------*/
487 /*MC
488      PCREDUNDANT - Runs a KSP solver with preconditioner for the entire problem on subgroups of processors
489 
490      Options for the redundant preconditioners can be set with -redundant_pc_xxx for the redundant KSP with -redundant_ksp_xxx
491 
492   Options Database:
493 .  -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and
494                               use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
495 
496    Level: intermediate
497 
498    Notes:
499     The default KSP is preonly and the default PC is LU or CHOLESKY if Pmat is of type MATSBAIJ.
500 
501    PCFactorSetShiftType() applied to this PC will convey they shift type into the inner PC if it is factorization based.
502 
503    Developer Notes:
504     Note that PCSetInitialGuessNonzero()  is not used by this class but likely should be.
505 
506 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PCRedundantSetScatter()`,
507           `PCRedundantGetKSP()`, `PCRedundantGetOperators()`, `PCRedundantSetNumber()`
508 M*/
509 
510 PETSC_EXTERN PetscErrorCode PCCreate_Redundant(PC pc)
511 {
512   PC_Redundant   *red;
513   PetscMPIInt    size;
514 
515   PetscFunctionBegin;
516   PetscCall(PetscNewLog(pc,&red));
517   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size));
518 
519   red->nsubcomm       = size;
520   red->useparallelmat = PETSC_TRUE;
521   pc->data            = (void*)red;
522 
523   pc->ops->apply          = PCApply_Redundant;
524   pc->ops->applytranspose = PCApplyTranspose_Redundant;
525   pc->ops->setup          = PCSetUp_Redundant;
526   pc->ops->destroy        = PCDestroy_Redundant;
527   pc->ops->reset          = PCReset_Redundant;
528   pc->ops->setfromoptions = PCSetFromOptions_Redundant;
529   pc->ops->view           = PCView_Redundant;
530 
531   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetScatter_C",PCRedundantSetScatter_Redundant));
532   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetNumber_C",PCRedundantSetNumber_Redundant));
533   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetKSP_C",PCRedundantGetKSP_Redundant));
534   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetOperators_C",PCRedundantGetOperators_Redundant));
535   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Redundant));
536   PetscFunctionReturn(0);
537 }
538