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