xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 0971522ceecba37df062696def42d3ac8bbed20b)
1*0971522cSBarry Smith /*
2*0971522cSBarry Smith 
3*0971522cSBarry Smith */
4*0971522cSBarry Smith #include "src/ksp/pc/pcimpl.h"     /*I "petscpc.h" I*/
5*0971522cSBarry Smith 
6*0971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
7*0971522cSBarry Smith struct _PC_FieldSplitLink {
8*0971522cSBarry Smith   PC                pc;
9*0971522cSBarry Smith   Vec               x,y;
10*0971522cSBarry Smith   PetscInt          nfields;
11*0971522cSBarry Smith   PetscInt          *fields;
12*0971522cSBarry Smith   PC_FieldSplitLink next;
13*0971522cSBarry Smith };
14*0971522cSBarry Smith 
15*0971522cSBarry Smith typedef struct {
16*0971522cSBarry Smith   PetscInt          bs;
17*0971522cSBarry Smith   PetscInt          nsplits;
18*0971522cSBarry Smith   PC_FieldSplitLink head;
19*0971522cSBarry Smith   Vec               *x,*y;
20*0971522cSBarry Smith } PC_FieldSplit;
21*0971522cSBarry Smith 
22*0971522cSBarry Smith #undef __FUNCT__
23*0971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
24*0971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
25*0971522cSBarry Smith {
26*0971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
27*0971522cSBarry Smith   PetscErrorCode    ierr;
28*0971522cSBarry Smith   PetscTruth        iascii;
29*0971522cSBarry Smith   PetscInt          i,j;
30*0971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
31*0971522cSBarry Smith 
32*0971522cSBarry Smith   PetscFunctionBegin;
33*0971522cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
34*0971522cSBarry Smith   if (iascii) {
35*0971522cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit: total splits = %D",jac->nsplits);CHKERRQ(ierr);
36*0971522cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following PC objects:\n");CHKERRQ(ierr);
37*0971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
38*0971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
39*0971522cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
40*0971522cSBarry Smith       for (j=0; j<link->nfields; j++) {
41*0971522cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"%D \n",link->fields[j]);CHKERRQ(ierr);
42*0971522cSBarry Smith       }
43*0971522cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
44*0971522cSBarry Smith       ierr = PCView(link->pc,viewer);CHKERRQ(ierr);
45*0971522cSBarry Smith       link = link->next;
46*0971522cSBarry Smith     }
47*0971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
48*0971522cSBarry Smith   } else {
49*0971522cSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
50*0971522cSBarry Smith   }
51*0971522cSBarry Smith   PetscFunctionReturn(0);
52*0971522cSBarry Smith }
53*0971522cSBarry Smith 
54*0971522cSBarry Smith #undef __FUNCT__
55*0971522cSBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
56*0971522cSBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
57*0971522cSBarry Smith {
58*0971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
59*0971522cSBarry Smith   PetscErrorCode    ierr;
60*0971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
61*0971522cSBarry Smith 
62*0971522cSBarry Smith   PetscFunctionBegin;
63*0971522cSBarry Smith   while (link) {
64*0971522cSBarry Smith     ierr = PCSetUp(link->pc);CHKERRQ(ierr);
65*0971522cSBarry Smith     link = link->next;
66*0971522cSBarry Smith   }
67*0971522cSBarry Smith   PetscFunctionReturn(0);
68*0971522cSBarry Smith }
69*0971522cSBarry Smith 
70*0971522cSBarry Smith #undef __FUNCT__
71*0971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
72*0971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
73*0971522cSBarry Smith {
74*0971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
75*0971522cSBarry Smith   PetscErrorCode    ierr;
76*0971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
77*0971522cSBarry Smith 
78*0971522cSBarry Smith   PetscFunctionBegin;
79*0971522cSBarry Smith   ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
80*0971522cSBarry Smith   while (link) {
81*0971522cSBarry Smith     ierr = PCApply(link->pc,link->x,link->y);CHKERRQ(ierr);
82*0971522cSBarry Smith     link = link->next;
83*0971522cSBarry Smith   }
84*0971522cSBarry Smith   ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
85*0971522cSBarry Smith   PetscFunctionReturn(0);
86*0971522cSBarry Smith }
87*0971522cSBarry Smith 
88*0971522cSBarry Smith #undef __FUNCT__
89*0971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
90*0971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
91*0971522cSBarry Smith {
92*0971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
93*0971522cSBarry Smith   PetscErrorCode    ierr;
94*0971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
95*0971522cSBarry Smith 
96*0971522cSBarry Smith   PetscFunctionBegin;
97*0971522cSBarry Smith   while (link) {
98*0971522cSBarry Smith     ierr = PCDestroy(link->pc);CHKERRQ(ierr);
99*0971522cSBarry Smith     ierr = PetscFree2(link,link->fields);CHKERRQ(ierr);
100*0971522cSBarry Smith   }
101*0971522cSBarry Smith   if (jac->x) {ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);}
102*0971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
103*0971522cSBarry Smith   PetscFunctionReturn(0);
104*0971522cSBarry Smith }
105*0971522cSBarry Smith 
106*0971522cSBarry Smith #undef __FUNCT__
107*0971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
108*0971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
109*0971522cSBarry Smith {
110*0971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
111*0971522cSBarry Smith   PetscErrorCode    ierr;
112*0971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
113*0971522cSBarry Smith 
114*0971522cSBarry Smith   PetscFunctionBegin;
115*0971522cSBarry Smith   if (!link) { /* user never set fields, so set them ourselves */
116*0971522cSBarry Smith     link = jac->head;
117*0971522cSBarry Smith   }
118*0971522cSBarry Smith   while (link) {
119*0971522cSBarry Smith     ierr = PCApply(link->pc,link->x,link->y);CHKERRQ(ierr);
120*0971522cSBarry Smith     link = link->next;
121*0971522cSBarry Smith   }
122*0971522cSBarry Smith   PetscFunctionReturn(0);
123*0971522cSBarry Smith }
124*0971522cSBarry Smith 
125*0971522cSBarry Smith /*------------------------------------------------------------------------------------*/
126*0971522cSBarry Smith 
127*0971522cSBarry Smith EXTERN_C_BEGIN
128*0971522cSBarry Smith #undef __FUNCT__
129*0971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
130*0971522cSBarry Smith PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
131*0971522cSBarry Smith {
132*0971522cSBarry Smith   PC_FieldSplit     *jac;
133*0971522cSBarry Smith   PetscErrorCode    ierr;
134*0971522cSBarry Smith   PetscInt          *myfields;
135*0971522cSBarry Smith   PC_FieldSplitLink link,next = jac->head;
136*0971522cSBarry Smith 
137*0971522cSBarry Smith   PetscFunctionBegin;
138*0971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
139*0971522cSBarry Smith   ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&link,n,PetscInt,&myfields);CHKERRQ(ierr);
140*0971522cSBarry Smith   ierr = PetscMemcpy(myfields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
141*0971522cSBarry Smith   link->fields  = myfields;
142*0971522cSBarry Smith   link->nfields = n;
143*0971522cSBarry Smith 
144*0971522cSBarry Smith   if (!next) {
145*0971522cSBarry Smith     jac->head = link;
146*0971522cSBarry Smith   } else {
147*0971522cSBarry Smith     while (next->next) {
148*0971522cSBarry Smith       next = next->next;
149*0971522cSBarry Smith     }
150*0971522cSBarry Smith     next->next = link;
151*0971522cSBarry Smith   }
152*0971522cSBarry Smith   jac->nsplits++;
153*0971522cSBarry Smith   PetscFunctionReturn(0);
154*0971522cSBarry Smith }
155*0971522cSBarry Smith EXTERN_C_END
156*0971522cSBarry Smith 
157*0971522cSBarry Smith 
158*0971522cSBarry Smith EXTERN_C_BEGIN
159*0971522cSBarry Smith #undef __FUNCT__
160*0971522cSBarry Smith #define __FUNCT__ "PCFieldSplitGetSubPC_FieldSplit"
161*0971522cSBarry Smith PetscErrorCode PCFieldSplitGetSubPC_FieldSplit(PC pc,PetscInt *n,PC **subpc)
162*0971522cSBarry Smith {
163*0971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
164*0971522cSBarry Smith   PetscErrorCode    ierr;
165*0971522cSBarry Smith   PetscInt          cnt = 0;
166*0971522cSBarry Smith   PC_FieldSplitLink link;
167*0971522cSBarry Smith 
168*0971522cSBarry Smith   PetscFunctionBegin;
169*0971522cSBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(PC*),subpc);CHKERRQ(ierr);
170*0971522cSBarry Smith   while (link) {
171*0971522cSBarry Smith     (*subpc)[cnt++] = link->pc;
172*0971522cSBarry Smith     link = link->next;
173*0971522cSBarry Smith   }
174*0971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
175*0971522cSBarry Smith   *n = jac->nsplits;
176*0971522cSBarry Smith   PetscFunctionReturn(0);
177*0971522cSBarry Smith }
178*0971522cSBarry Smith EXTERN_C_END
179*0971522cSBarry Smith 
180*0971522cSBarry Smith #undef __FUNCT__
181*0971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
182*0971522cSBarry Smith /*@
183*0971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
184*0971522cSBarry Smith 
185*0971522cSBarry Smith     Collective on PC
186*0971522cSBarry Smith 
187*0971522cSBarry Smith     Input Parameters:
188*0971522cSBarry Smith +   pc  - the preconditioner context
189*0971522cSBarry Smith .   n - the number of fields in this split
190*0971522cSBarry Smith .   fields - the fields in this split
191*0971522cSBarry Smith 
192*0971522cSBarry Smith     Level: intermediate
193*0971522cSBarry Smith 
194*0971522cSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT
195*0971522cSBarry Smith 
196*0971522cSBarry Smith @*/
197*0971522cSBarry Smith PetscErrorCode PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
198*0971522cSBarry Smith {
199*0971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
200*0971522cSBarry Smith 
201*0971522cSBarry Smith   PetscFunctionBegin;
202*0971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
203*0971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
204*0971522cSBarry Smith   if (f) {
205*0971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
206*0971522cSBarry Smith   }
207*0971522cSBarry Smith   PetscFunctionReturn(0);
208*0971522cSBarry Smith }
209*0971522cSBarry Smith 
210*0971522cSBarry Smith #undef __FUNCT__
211*0971522cSBarry Smith #define __FUNCT__ "PCFieldSplitGetSubPC"
212*0971522cSBarry Smith /*@C
213*0971522cSBarry Smith    PCFieldSplitGetSubPC - Gets the PC contexts for all splits
214*0971522cSBarry Smith 
215*0971522cSBarry Smith    Collective on PC
216*0971522cSBarry Smith 
217*0971522cSBarry Smith    Input Parameter:
218*0971522cSBarry Smith .  pc - the preconditioner context
219*0971522cSBarry Smith 
220*0971522cSBarry Smith    Output Parameters:
221*0971522cSBarry Smith +  n - the number of split
222*0971522cSBarry Smith -  pc - the array of PC contexts
223*0971522cSBarry Smith 
224*0971522cSBarry Smith    Note:
225*0971522cSBarry Smith    After PCFieldSplitGetSubPC() the array of PCs IS to be freed
226*0971522cSBarry Smith 
227*0971522cSBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubPC().
228*0971522cSBarry Smith 
229*0971522cSBarry Smith    Level: advanced
230*0971522cSBarry Smith 
231*0971522cSBarry Smith .keywords: PC,
232*0971522cSBarry Smith 
233*0971522cSBarry Smith .seealso: PCFIELDSPLIT
234*0971522cSBarry Smith @*/
235*0971522cSBarry Smith PetscErrorCode PCFieldSplitGetSubPC(PC pc,PetscInt *n,PC *subpc[])
236*0971522cSBarry Smith {
237*0971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,PC **);
238*0971522cSBarry Smith 
239*0971522cSBarry Smith   PetscFunctionBegin;
240*0971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
241*0971522cSBarry Smith   PetscValidIntPointer(n,2);
242*0971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubPC_C",(void (**)(void))&f);CHKERRQ(ierr);
243*0971522cSBarry Smith   if (f) {
244*0971522cSBarry Smith     ierr = (*f)(pc,n,subpc);CHKERRQ(ierr);
245*0971522cSBarry Smith   } else {
246*0971522cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subpc for this type of PC");
247*0971522cSBarry Smith   }
248*0971522cSBarry Smith 
249*0971522cSBarry Smith  PetscFunctionReturn(0);
250*0971522cSBarry Smith }
251*0971522cSBarry Smith 
252*0971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
253*0971522cSBarry Smith /*MC
254*0971522cSBarry Smith    PCFieldSplit - Preconditioner created by combining seperate preconditioners for individual
255*0971522cSBarry Smith                   fields or groups of fields
256*0971522cSBarry Smith 
257*0971522cSBarry Smith 
258*0971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
259*0971522cSBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_ilu_levels 1
260*0971522cSBarry Smith 
261*0971522cSBarry Smith      To set the options on the solvers seperate for each block call PCFieldSplitGetSubPC()
262*0971522cSBarry Smith          and set the options directly on the resulting PC object
263*0971522cSBarry Smith 
264*0971522cSBarry Smith    Level: intermediate
265*0971522cSBarry Smith 
266*0971522cSBarry Smith    Concepts: physics based preconditioners
267*0971522cSBarry Smith 
268*0971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
269*0971522cSBarry Smith            PCFieldSplitGetSubPC(), PCFieldSplitSetFields()
270*0971522cSBarry Smith M*/
271*0971522cSBarry Smith 
272*0971522cSBarry Smith EXTERN_C_BEGIN
273*0971522cSBarry Smith #undef __FUNCT__
274*0971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
275*0971522cSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc)
276*0971522cSBarry Smith {
277*0971522cSBarry Smith   PetscErrorCode ierr;
278*0971522cSBarry Smith   PC_FieldSplit  *jac;
279*0971522cSBarry Smith 
280*0971522cSBarry Smith   PetscFunctionBegin;
281*0971522cSBarry Smith   ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr);
282*0971522cSBarry Smith   PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));
283*0971522cSBarry Smith   ierr = PetscMemzero(jac,sizeof(PC_FieldSplit));CHKERRQ(ierr);
284*0971522cSBarry Smith   jac->bs      = -1;
285*0971522cSBarry Smith   jac->nsplits = 0;
286*0971522cSBarry Smith   pc->data     = (void*)jac;
287*0971522cSBarry Smith 
288*0971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
289*0971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
290*0971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
291*0971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
292*0971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
293*0971522cSBarry Smith   pc->ops->applyrichardson   = 0;
294*0971522cSBarry Smith 
295*0971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubPC_C","PCFieldSplitGetSubPC_FieldSplit",
296*0971522cSBarry Smith                     PCFieldSplitGetSubPC_FieldSplit);CHKERRQ(ierr);
297*0971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
298*0971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
299*0971522cSBarry Smith   PetscFunctionReturn(0);
300*0971522cSBarry Smith }
301*0971522cSBarry Smith EXTERN_C_END
302*0971522cSBarry Smith 
303*0971522cSBarry Smith 
304