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