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