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