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