xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 668f157ea6d169bb50bcdb9ebcdd418abd089fa7)
1 #define PETSCKSP_DLL
2 
3 #include "private/pcimpl.h"     /*I "petscpc.h" I*/
4 
5 const char *PCFieldSplitSchurPreTypes[] = {"self","diag","user","PCFieldSplitSchurPreType","PC_FIELDSPLIT_SCHUR_PRE_",0};
6 
7 typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
8 struct _PC_FieldSplitLink {
9   KSP               ksp;
10   Vec               x,y;
11   PetscInt          nfields;
12   PetscInt          *fields;
13   VecScatter        sctx;
14   IS                is;
15   PC_FieldSplitLink next,previous;
16 };
17 
18 typedef struct {
19   PCCompositeType   type;
20   PetscTruth        defaultsplit; /* Flag for a system with a set of 'k' scalar fields with the same layout (and bs = k) */
21   PetscTruth        realdiagonal; /* Flag to use the diagonal blocks of mat preconditioned by pmat, instead of just pmat */
22   PetscInt          bs;           /* Block size for IS and Mat structures */
23   PetscInt          nsplits;      /* Number of field divisions defined */
24   Vec               *x,*y,w1,w2;
25   Mat               *mat;         /* The diagonal block for each split */
26   Mat               *pmat;        /* The preconditioning diagonal block for each split */
27   Mat               *Afield;      /* The rows of the matrix associated with each split */
28   PetscTruth        issetup;
29   /* Only used when Schur complement preconditioning is used */
30   Mat               B;            /* The (0,1) block */
31   Mat               C;            /* The (1,0) block */
32   Mat               schur;        /* The Schur complement S = D - C A^{-1} B */
33   Mat               schur_user;   /* User-provided preconditioning matrix for the Schur complement */
34   PCFieldSplitSchurPreType schurpre; /* Determines which preconditioning matrix is used for the Schur complement */
35   KSP               kspschur;     /* The solver for S */
36   PC_FieldSplitLink head;
37 } PC_FieldSplit;
38 
39 /*
40     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
41    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
42    PC you could change this.
43 */
44 
45 /* This helper is so that setting a user-provided preconditioning matrix is orthogonal to choosing to use it.  This way the
46 * application-provided FormJacobian can provide this matrix without interfering with the user's (command-line) choices. */
47 static Mat FieldSplitSchurPre(PC_FieldSplit *jac)
48 {
49   switch (jac->schurpre) {
50     case PC_FIELDSPLIT_SCHUR_PRE_SELF: return jac->schur;
51     case PC_FIELDSPLIT_SCHUR_PRE_DIAG: return jac->pmat[1];
52     case PC_FIELDSPLIT_SCHUR_PRE_USER: /* Use a user-provided matrix if it is given, otherwise diagonal block */
53     default:
54       return jac->schur_user ? jac->schur_user : jac->pmat[1];
55   }
56 }
57 
58 
59 #undef __FUNCT__
60 #define __FUNCT__ "PCView_FieldSplit"
61 static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
62 {
63   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
64   PetscErrorCode    ierr;
65   PetscTruth        iascii;
66   PetscInt          i,j;
67   PC_FieldSplitLink ilink = jac->head;
68 
69   PetscFunctionBegin;
70   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
71   if (iascii) {
72     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
73     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
74     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
75     for (i=0; i<jac->nsplits; i++) {
76       if (ilink->fields) {
77 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
78 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
79 	for (j=0; j<ilink->nfields; j++) {
80 	  if (j > 0) {
81 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
82 	  }
83 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
84 	}
85 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
86         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
87       } else {
88 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
89       }
90       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
91       ilink = ilink->next;
92     }
93     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
94   } else {
95     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
96   }
97   PetscFunctionReturn(0);
98 }
99 
100 #undef __FUNCT__
101 #define __FUNCT__ "PCView_FieldSplit_Schur"
102 static PetscErrorCode PCView_FieldSplit_Schur(PC pc,PetscViewer viewer)
103 {
104   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
105   PetscErrorCode    ierr;
106   PetscTruth        iascii;
107   PetscInt          i,j;
108   PC_FieldSplitLink ilink = jac->head;
109   KSP               ksp;
110 
111   PetscFunctionBegin;
112   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
113   if (iascii) {
114     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with Schur preconditioner, blocksize = %D\n",jac->bs);CHKERRQ(ierr);
115     ierr = PetscViewerASCIIPrintf(viewer,"  Split info:\n");CHKERRQ(ierr);
116     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
117     for (i=0; i<jac->nsplits; i++) {
118       if (ilink->fields) {
119 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
120 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
121 	for (j=0; j<ilink->nfields; j++) {
122 	  if (j > 0) {
123 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
124 	  }
125 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
126 	}
127 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
128         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
129       } else {
130 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
131       }
132       ilink = ilink->next;
133     }
134     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for A block \n");CHKERRQ(ierr);
135     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
136     if (jac->schur) {
137       ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
138       ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
139     } else {
140       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
141     }
142     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
143     ierr = PetscViewerASCIIPrintf(viewer,"KSP solver for S = D - C inv(A) B \n");CHKERRQ(ierr);
144     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
145     if (jac->kspschur) {
146       ierr = KSPView(jac->kspschur,viewer);CHKERRQ(ierr);
147     } else {
148       ierr = PetscViewerASCIIPrintf(viewer,"  not yet available\n");CHKERRQ(ierr);
149     }
150     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
151     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
152   } else {
153     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
154   }
155   PetscFunctionReturn(0);
156 }
157 
158 #undef __FUNCT__
159 #define __FUNCT__ "PCFieldSplitSetDefaults"
160 static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
161 {
162   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
163   PetscErrorCode    ierr;
164   PC_FieldSplitLink ilink = jac->head;
165   PetscInt          i = 0,*ifields,nfields;
166   PetscTruth        flg = PETSC_FALSE,*fields,flg2;
167   char              optionname[128];
168 
169   PetscFunctionBegin;
170   if (!ilink) {
171 
172     if (jac->bs <= 0) {
173       if (pc->pmat) {
174         ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
175       } else {
176         jac->bs = 1;
177       }
178     }
179 
180     ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
181     if (!flg) {
182       /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
183          then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
184       flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */
185       ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
186       while (PETSC_TRUE) {
187         sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
188         nfields = jac->bs;
189         ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr);
190         if (!flg2) break;
191         if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
192         flg = PETSC_FALSE;
193         ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr);
194       }
195       ierr = PetscFree(ifields);CHKERRQ(ierr);
196     }
197 
198     if (flg) {
199       ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
200       ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
201       ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
202       while (ilink) {
203 	for (i=0; i<ilink->nfields; i++) {
204 	  fields[ilink->fields[i]] = PETSC_TRUE;
205 	}
206 	ilink = ilink->next;
207       }
208       jac->defaultsplit = PETSC_TRUE;
209       for (i=0; i<jac->bs; i++) {
210 	if (!fields[i]) {
211 	  ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
212 	} else {
213 	  jac->defaultsplit = PETSC_FALSE;
214 	}
215       }
216       ierr = PetscFree(fields);CHKERRQ(ierr);
217     }
218   } else if (jac->nsplits == 1) {
219     if (ilink->is) {
220       IS       is2;
221       PetscInt nmin,nmax;
222 
223       ierr = MatGetOwnershipRange(pc->mat,&nmin,&nmax);CHKERRQ(ierr);
224       ierr = ISComplement(ilink->is,nmin,nmax,&is2);CHKERRQ(ierr);
225       ierr = PCFieldSplitSetIS(pc,is2);CHKERRQ(ierr);
226       ierr = ISDestroy(is2);CHKERRQ(ierr);
227     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must provide at least two sets of fields to PCFieldSplit()");
228   }
229   if (jac->nsplits < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Unhandled case, must have at least two fields");
230   PetscFunctionReturn(0);
231 }
232 
233 
234 #undef __FUNCT__
235 #define __FUNCT__ "PCSetUp_FieldSplit"
236 static PetscErrorCode PCSetUp_FieldSplit(PC pc)
237 {
238   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
239   PetscErrorCode    ierr;
240   PC_FieldSplitLink ilink;
241   PetscInt          i,nsplit,ccsize;
242   MatStructure      flag = pc->flag;
243   PetscTruth        sorted;
244 
245   PetscFunctionBegin;
246   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
247   nsplit = jac->nsplits;
248   ilink  = jac->head;
249 
250   /* get the matrices for each split */
251   if (!jac->issetup) {
252     PetscInt rstart,rend,nslots,bs;
253 
254     jac->issetup = PETSC_TRUE;
255 
256     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
257     bs     = jac->bs;
258     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
259     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
260     nslots = (rend - rstart)/bs;
261     for (i=0; i<nsplit; i++) {
262       if (jac->defaultsplit) {
263         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
264       } else if (!ilink->is) {
265         if (ilink->nfields > 1) {
266           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
267           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
268           for (j=0; j<nslots; j++) {
269             for (k=0; k<nfields; k++) {
270               ii[nfields*j + k] = rstart + bs*j + fields[k];
271             }
272           }
273           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr);
274           ierr = PetscFree(ii);CHKERRQ(ierr);
275         } else {
276           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
277         }
278       }
279       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
280       if (!sorted) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Fields must be sorted when creating split");
281       ilink = ilink->next;
282     }
283   }
284 
285   ilink  = jac->head;
286   if (!jac->pmat) {
287     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
288     for (i=0; i<nsplit; i++) {
289       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
290       ilink = ilink->next;
291     }
292   } else {
293     for (i=0; i<nsplit; i++) {
294       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
295       ilink = ilink->next;
296     }
297   }
298   if (jac->realdiagonal) {
299     ilink = jac->head;
300     if (!jac->mat) {
301       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->mat);CHKERRQ(ierr);
302       for (i=0; i<nsplit; i++) {
303         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_INITIAL_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
304         ilink = ilink->next;
305       }
306     } else {
307       for (i=0; i<nsplit; i++) {
308         ierr = MatGetSubMatrix(pc->mat,ilink->is,ilink->is,MAT_REUSE_MATRIX,&jac->mat[i]);CHKERRQ(ierr);
309         ilink = ilink->next;
310       }
311     }
312   } else {
313     jac->mat = jac->pmat;
314   }
315 
316   if (jac->type != PC_COMPOSITE_ADDITIVE  && jac->type != PC_COMPOSITE_SCHUR) {
317     /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
318     ilink  = jac->head;
319     if (!jac->Afield) {
320       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
321       for (i=0; i<nsplit; i++) {
322         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
323         ilink = ilink->next;
324       }
325     } else {
326       for (i=0; i<nsplit; i++) {
327         ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
328         ilink = ilink->next;
329       }
330     }
331   }
332 
333   if (jac->type == PC_COMPOSITE_SCHUR) {
334     IS       ccis;
335     PetscInt rstart,rend;
336     if (nsplit != 2) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"To use Schur complement preconditioner you must have exactly 2 fields");
337 
338     /* When extracting off-diagonal submatrices, we take complements from this range */
339     ierr  = MatGetOwnershipRangeColumn(pc->mat,&rstart,&rend);CHKERRQ(ierr);
340 
341     /* need to handle case when one is resetting up the preconditioner */
342     if (jac->schur) {
343       ilink = jac->head;
344       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
345       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->B);CHKERRQ(ierr);
346       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
347       ilink = ilink->next;
348       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
349       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_REUSE_MATRIX,&jac->C);CHKERRQ(ierr);
350       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
351       ierr  = MatSchurComplementUpdate(jac->schur,jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->pmat[1],pc->flag);CHKERRQ(ierr);
352       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),pc->flag);CHKERRQ(ierr);
353 
354      } else {
355       KSP ksp;
356 
357       /* extract the B and C matrices */
358       ilink = jac->head;
359       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
360       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->B);CHKERRQ(ierr);
361       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
362       ilink = ilink->next;
363       ierr  = ISComplement(ilink->is,rstart,rend,&ccis);CHKERRQ(ierr);
364       ierr  = MatGetSubMatrix(pc->mat,ilink->is,ccis,MAT_INITIAL_MATRIX,&jac->C);CHKERRQ(ierr);
365       ierr  = ISDestroy(ccis);CHKERRQ(ierr);
366       /* Better would be to use 'mat[0]' (diagonal block of the real matrix) preconditioned by pmat[0] */
367       ierr  = MatCreateSchurComplement(jac->mat[0],jac->pmat[0],jac->B,jac->C,jac->mat[1],&jac->schur);CHKERRQ(ierr);
368       ierr  = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
369       ierr  = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,2);CHKERRQ(ierr);
370       ierr  = MatSetFromOptions(jac->schur);CHKERRQ(ierr);
371 
372       ierr  = KSPCreate(((PetscObject)pc)->comm,&jac->kspschur);CHKERRQ(ierr);
373       ierr  = PetscObjectIncrementTabLevel((PetscObject)jac->kspschur,(PetscObject)pc,1);CHKERRQ(ierr);
374       ierr  = KSPSetOperators(jac->kspschur,jac->schur,FieldSplitSchurPre(jac),DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
375       if (jac->schurpre == PC_FIELDSPLIT_SCHUR_PRE_SELF) {
376         PC pc;
377         ierr = KSPGetPC(jac->kspschur,&pc);CHKERRQ(ierr);
378         ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
379         /* Note: This is bad if there exist preconditioners for MATSCHURCOMPLEMENT */
380       }
381       ierr  = KSPSetOptionsPrefix(jac->kspschur,((PetscObject)pc)->prefix);CHKERRQ(ierr);
382       ierr  = KSPAppendOptionsPrefix(jac->kspschur,"fieldsplit_1_");CHKERRQ(ierr);
383       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
384       ierr = KSPSetFromOptions(jac->kspschur);CHKERRQ(ierr);
385 
386       ierr = PetscMalloc2(2,Vec,&jac->x,2,Vec,&jac->y);CHKERRQ(ierr);
387       ierr = MatGetVecs(jac->pmat[0],&jac->x[0],&jac->y[0]);CHKERRQ(ierr);
388       ierr = MatGetVecs(jac->pmat[1],&jac->x[1],&jac->y[1]);CHKERRQ(ierr);
389       ilink = jac->head;
390       ilink->x = jac->x[0]; ilink->y = jac->y[0];
391       ilink = ilink->next;
392       ilink->x = jac->x[1]; ilink->y = jac->y[1];
393     }
394   } else {
395     /* set up the individual PCs */
396     i    = 0;
397     ilink = jac->head;
398     while (ilink) {
399       ierr = KSPSetOperators(ilink->ksp,jac->mat[i],jac->pmat[i],flag);CHKERRQ(ierr);
400       /* really want setfromoptions called in PCSetFromOptions_FieldSplit(), but it is not ready yet */
401       ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
402       ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
403       i++;
404       ilink = ilink->next;
405     }
406 
407     /* create work vectors for each split */
408     if (!jac->x) {
409       ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
410       ilink = jac->head;
411       for (i=0; i<nsplit; i++) {
412         Vec *vl,*vr;
413 
414         ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
415         ilink->x  = *vr;
416         ilink->y  = *vl;
417         ierr      = PetscFree(vr);CHKERRQ(ierr);
418         ierr      = PetscFree(vl);CHKERRQ(ierr);
419         jac->x[i] = ilink->x;
420         jac->y[i] = ilink->y;
421         ilink     = ilink->next;
422       }
423     }
424   }
425 
426 
427   if (!jac->head->sctx) {
428     Vec xtmp;
429 
430     /* compute scatter contexts needed by multiplicative versions and non-default splits */
431 
432     ilink = jac->head;
433     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
434     for (i=0; i<nsplit; i++) {
435       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
436       ilink = ilink->next;
437     }
438     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
439   }
440   PetscFunctionReturn(0);
441 }
442 
443 #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
444     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
445      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
446      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
447      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
448      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
449 
450 #undef __FUNCT__
451 #define __FUNCT__ "PCApply_FieldSplit_Schur"
452 static PetscErrorCode PCApply_FieldSplit_Schur(PC pc,Vec x,Vec y)
453 {
454   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
455   PetscErrorCode    ierr;
456   KSP               ksp;
457   PC_FieldSplitLink ilinkA = jac->head, ilinkD = ilinkA->next;
458 
459   PetscFunctionBegin;
460   ierr = MatSchurComplementGetKSP(jac->schur,&ksp);CHKERRQ(ierr);
461 
462   ierr = VecScatterBegin(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
463   ierr = VecScatterEnd(ilinkA->sctx,x,ilinkA->x,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
464   ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
465   ierr = MatMult(jac->C,ilinkA->y,ilinkD->x);CHKERRQ(ierr);
466   ierr = VecScale(ilinkD->x,-1.0);CHKERRQ(ierr);
467   ierr = VecScatterBegin(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
468   ierr = VecScatterEnd(ilinkD->sctx,x,ilinkD->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
469 
470   ierr = KSPSolve(jac->kspschur,ilinkD->x,ilinkD->y);CHKERRQ(ierr);
471   ierr = VecScatterBegin(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
472   ierr = VecScatterEnd(ilinkD->sctx,ilinkD->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
473 
474   ierr = MatMult(jac->B,ilinkD->y,ilinkA->y);CHKERRQ(ierr);
475   ierr = VecAXPY(ilinkA->x,-1.0,ilinkA->y);CHKERRQ(ierr);
476   ierr = KSPSolve(ksp,ilinkA->x,ilinkA->y);CHKERRQ(ierr);
477   ierr = VecScatterBegin(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
478   ierr = VecScatterEnd(ilinkA->sctx,ilinkA->y,y,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
479 
480   PetscFunctionReturn(0);
481 }
482 
483 #undef __FUNCT__
484 #define __FUNCT__ "PCApply_FieldSplit"
485 static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
486 {
487   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
488   PetscErrorCode    ierr;
489   PC_FieldSplitLink ilink = jac->head;
490   PetscInt          cnt;
491 
492   PetscFunctionBegin;
493   CHKMEMQ;
494   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
495   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
496 
497   if (jac->type == PC_COMPOSITE_ADDITIVE) {
498     if (jac->defaultsplit) {
499       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
500       while (ilink) {
501         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
502         ilink = ilink->next;
503       }
504       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
505     } else {
506       ierr = VecSet(y,0.0);CHKERRQ(ierr);
507       while (ilink) {
508         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
509         ilink = ilink->next;
510       }
511     }
512   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
513     if (!jac->w1) {
514       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
515       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
516     }
517     ierr = VecSet(y,0.0);CHKERRQ(ierr);
518     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
519     cnt = 1;
520     while (ilink->next) {
521       ilink = ilink->next;
522       /* compute the residual only over the part of the vector needed */
523       ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
524       ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
525       ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
526       ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
527       ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
528       ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
529       ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
530     }
531     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
532       cnt -= 2;
533       while (ilink->previous) {
534         ilink = ilink->previous;
535         /* compute the residual only over the part of the vector needed */
536         ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
537         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
538         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
539         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
540         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
541         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
542         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
543       }
544     }
545   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
546   CHKMEMQ;
547   PetscFunctionReturn(0);
548 }
549 
550 #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
551     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
552      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
553      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
554      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
555      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
556 
557 #undef __FUNCT__
558 #define __FUNCT__ "PCApply_FieldSplit"
559 static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
560 {
561   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
562   PetscErrorCode    ierr;
563   PC_FieldSplitLink ilink = jac->head;
564 
565   PetscFunctionBegin;
566   CHKMEMQ;
567   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
568   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
569 
570   if (jac->type == PC_COMPOSITE_ADDITIVE) {
571     if (jac->defaultsplit) {
572       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
573       while (ilink) {
574 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
575 	ilink = ilink->next;
576       }
577       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
578     } else {
579       ierr = VecSet(y,0.0);CHKERRQ(ierr);
580       while (ilink) {
581         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
582 	ilink = ilink->next;
583       }
584     }
585   } else {
586     if (!jac->w1) {
587       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
588       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
589     }
590     ierr = VecSet(y,0.0);CHKERRQ(ierr);
591     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
592       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
593       while (ilink->next) {
594         ilink = ilink->next;
595         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
596         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
597         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
598       }
599       while (ilink->previous) {
600         ilink = ilink->previous;
601         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
602         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
603         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
604       }
605     } else {
606       while (ilink->next) {   /* get to last entry in linked list */
607 	ilink = ilink->next;
608       }
609       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
610       while (ilink->previous) {
611 	ilink = ilink->previous;
612 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
613 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
614 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
615       }
616     }
617   }
618   CHKMEMQ;
619   PetscFunctionReturn(0);
620 }
621 
622 #undef __FUNCT__
623 #define __FUNCT__ "PCDestroy_FieldSplit"
624 static PetscErrorCode PCDestroy_FieldSplit(PC pc)
625 {
626   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
627   PetscErrorCode    ierr;
628   PC_FieldSplitLink ilink = jac->head,next;
629 
630   PetscFunctionBegin;
631   while (ilink) {
632     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
633     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
634     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
635     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
636     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
637     next = ilink->next;
638     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
639     ierr = PetscFree(ilink);CHKERRQ(ierr);
640     ilink = next;
641   }
642   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
643   if (jac->mat && jac->mat != jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->mat);CHKERRQ(ierr);}
644   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
645   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
646   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
647   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
648   if (jac->schur) {ierr = MatDestroy(jac->schur);CHKERRQ(ierr);}
649   if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
650   if (jac->kspschur) {ierr = KSPDestroy(jac->kspschur);CHKERRQ(ierr);}
651   if (jac->B) {ierr = MatDestroy(jac->B);CHKERRQ(ierr);}
652   if (jac->C) {ierr = MatDestroy(jac->C);CHKERRQ(ierr);}
653   ierr = PetscFree(jac);CHKERRQ(ierr);
654   PetscFunctionReturn(0);
655 }
656 
657 #undef __FUNCT__
658 #define __FUNCT__ "PCSetFromOptions_FieldSplit"
659 static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
660 {
661   PetscErrorCode  ierr;
662   PetscInt        i = 0,nfields,*fields,bs;
663   PetscTruth      flg;
664   char            optionname[128];
665   PC_FieldSplit   *jac = (PC_FieldSplit*)pc->data;
666   PCCompositeType ctype;
667 
668   PetscFunctionBegin;
669   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
670   ierr = PetscOptionsTruth("-pc_fieldsplit_real_diagonal","Use diagonal blocks of the operator","PCFieldSplitSetRealDiagonal",jac->realdiagonal,&jac->realdiagonal,PETSC_NULL);CHKERRQ(ierr);
671   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
672   if (flg) {
673     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
674   }
675 
676   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&ctype,&flg);CHKERRQ(ierr);
677   if (flg) {
678     ierr = PCFieldSplitSetType(pc,ctype);CHKERRQ(ierr);
679   }
680 
681   /* Only setup fields once */
682   if ((jac->bs > 0) && (jac->nsplits == 0)) {
683     /* only allow user to set fields from command line if bs is already known.
684        otherwise user can set them in PCFieldSplitSetDefaults() */
685     ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
686     while (PETSC_TRUE) {
687       sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
688       nfields = jac->bs;
689       ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
690       if (!flg) break;
691       if (!nfields) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot list zero fields");
692       ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
693     }
694     ierr = PetscFree(fields);CHKERRQ(ierr);
695   }
696   ierr = PetscOptionsEnum("-pc_fieldsplit_schur_precondition","How to build preconditioner for Schur complement","PCFieldSplitSchurPrecondition",PCFieldSplitSchurPreTypes,(PetscEnum)jac->schurpre,(PetscEnum*)&jac->schurpre,PETSC_NULL);CHKERRQ(ierr);
697   ierr = PetscOptionsTail();CHKERRQ(ierr);
698   PetscFunctionReturn(0);
699 }
700 
701 /*------------------------------------------------------------------------------------*/
702 
703 EXTERN_C_BEGIN
704 #undef __FUNCT__
705 #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
706 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
707 {
708   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
709   PetscErrorCode    ierr;
710   PC_FieldSplitLink ilink,next = jac->head;
711   char              prefix[128];
712   PetscInt          i;
713 
714   PetscFunctionBegin;
715   if (n <= 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
716   for (i=0; i<n; i++) {
717     if (fields[i] >= jac->bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
718     if (fields[i] < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
719   }
720   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
721   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
722   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
723   ilink->nfields = n;
724   ilink->next    = PETSC_NULL;
725   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
726   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
727   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
728 
729   if (((PetscObject)pc)->prefix) {
730     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
731   } else {
732     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
733   }
734   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
735 
736   if (!next) {
737     jac->head       = ilink;
738     ilink->previous = PETSC_NULL;
739   } else {
740     while (next->next) {
741       next = next->next;
742     }
743     next->next      = ilink;
744     ilink->previous = next;
745   }
746   jac->nsplits++;
747   PetscFunctionReturn(0);
748 }
749 EXTERN_C_END
750 
751 EXTERN_C_BEGIN
752 #undef __FUNCT__
753 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit_Schur"
754 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit_Schur(PC pc,PetscInt *n,KSP **subksp)
755 {
756   PC_FieldSplit *jac = (PC_FieldSplit*)pc->data;
757   PetscErrorCode ierr;
758 
759   PetscFunctionBegin;
760   ierr = PetscMalloc(jac->nsplits*sizeof(KSP),subksp);CHKERRQ(ierr);
761   ierr = MatSchurComplementGetKSP(jac->schur,*subksp);CHKERRQ(ierr);
762   (*subksp)[1] = jac->kspschur;
763   *n = jac->nsplits;
764   PetscFunctionReturn(0);
765 }
766 EXTERN_C_END
767 
768 EXTERN_C_BEGIN
769 #undef __FUNCT__
770 #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
771 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
772 {
773   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
774   PetscErrorCode    ierr;
775   PetscInt          cnt = 0;
776   PC_FieldSplitLink ilink = jac->head;
777 
778   PetscFunctionBegin;
779   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
780   while (ilink) {
781     (*subksp)[cnt++] = ilink->ksp;
782     ilink = ilink->next;
783   }
784   if (cnt != jac->nsplits) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
785   *n = jac->nsplits;
786   PetscFunctionReturn(0);
787 }
788 EXTERN_C_END
789 
790 EXTERN_C_BEGIN
791 #undef __FUNCT__
792 #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
793 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is)
794 {
795   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
796   PetscErrorCode    ierr;
797   PC_FieldSplitLink ilink, next = jac->head;
798   char              prefix[128];
799 
800   PetscFunctionBegin;
801   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
802   ilink->is      = is;
803   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
804   ilink->next    = PETSC_NULL;
805   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
806   ierr           = PetscObjectIncrementTabLevel((PetscObject)ilink->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
807   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
808 
809   if (((PetscObject)pc)->prefix) {
810     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
811   } else {
812     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
813   }
814   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
815 
816   if (!next) {
817     jac->head       = ilink;
818     ilink->previous = PETSC_NULL;
819   } else {
820     while (next->next) {
821       next = next->next;
822     }
823     next->next      = ilink;
824     ilink->previous = next;
825   }
826   jac->nsplits++;
827 
828   PetscFunctionReturn(0);
829 }
830 EXTERN_C_END
831 
832 #undef __FUNCT__
833 #define __FUNCT__ "PCFieldSplitSetFields"
834 /*@
835     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
836 
837     Collective on PC
838 
839     Input Parameters:
840 +   pc  - the preconditioner context
841 .   n - the number of fields in this split
842 .   fields - the fields in this split
843 
844     Level: intermediate
845 
846     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
847 
848      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
849      size is three then one can define a field as 0, or 1 or 2 or 0,1 or 0,2 or 1,2 which mean
850      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
851      where the numbered entries indicate what is in the field.
852 
853 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
854 
855 @*/
856 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
857 {
858   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
859 
860   PetscFunctionBegin;
861   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
862   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
863   if (f) {
864     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
865   }
866   PetscFunctionReturn(0);
867 }
868 
869 #undef __FUNCT__
870 #define __FUNCT__ "PCFieldSplitSetIS"
871 /*@
872     PCFieldSplitSetIS - Sets the exact elements for field
873 
874     Collective on PC
875 
876     Input Parameters:
877 +   pc  - the preconditioner context
878 .   is - the index set that defines the vector elements in this field
879 
880 
881     Notes:
882     Use PCFieldSplitSetFields(), for fields defined by strided types.
883 
884     This function is called once per split (it creates a new split each time).
885 
886     Level: intermediate
887 
888 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
889 
890 @*/
891 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is)
892 {
893   PetscErrorCode ierr,(*f)(PC,IS);
894 
895   PetscFunctionBegin;
896   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
897   PetscValidHeaderSpecific(is,IS_CLASSID,1);
898   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr);
899   if (f) {
900     ierr = (*f)(pc,is);CHKERRQ(ierr);
901   }
902   PetscFunctionReturn(0);
903 }
904 
905 #undef __FUNCT__
906 #define __FUNCT__ "PCFieldSplitSetBlockSize"
907 /*@
908     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
909       fieldsplit preconditioner. If not set the matrix block size is used.
910 
911     Collective on PC
912 
913     Input Parameters:
914 +   pc  - the preconditioner context
915 -   bs - the block size
916 
917     Level: intermediate
918 
919 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
920 
921 @*/
922 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
923 {
924   PetscErrorCode ierr,(*f)(PC,PetscInt);
925 
926   PetscFunctionBegin;
927   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
928   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
929   if (f) {
930     ierr = (*f)(pc,bs);CHKERRQ(ierr);
931   }
932   PetscFunctionReturn(0);
933 }
934 
935 #undef __FUNCT__
936 #define __FUNCT__ "PCFieldSplitGetSubKSP"
937 /*@C
938    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
939 
940    Collective on KSP
941 
942    Input Parameter:
943 .  pc - the preconditioner context
944 
945    Output Parameters:
946 +  n - the number of split
947 -  pc - the array of KSP contexts
948 
949    Note:
950    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
951    (not the KSP just the array that contains them).
952 
953    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
954 
955    Level: advanced
956 
957 .seealso: PCFIELDSPLIT
958 @*/
959 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
960 {
961   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
962 
963   PetscFunctionBegin;
964   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
965   PetscValidIntPointer(n,2);
966   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
967   if (f) {
968     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
969   } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
970   PetscFunctionReturn(0);
971 }
972 
973 #undef __FUNCT__
974 #define __FUNCT__ "PCFieldSplitSchurPrecondition"
975 /*@
976     PCFieldSplitSchurPrecondition -  Indicates if the Schur complement is preconditioned by a preconditioner constructed by the
977       D matrix. Otherwise no preconditioner is used.
978 
979     Collective on PC
980 
981     Input Parameters:
982 +   pc  - the preconditioner context
983 .   ptype - which matrix to use for preconditioning the Schur complement
984 -   userpre - matrix to use for preconditioning, or PETSC_NULL
985 
986     Notes:
987     The default is to use the block on the diagonal of the preconditioning matrix.  This is D, in the (1,1) position.
988     There are currently no preconditioners that work directly with the Schur complement so setting
989     PC_FIELDSPLIT_SCHUR_PRE_SELF is observationally equivalent to -fieldsplit_1_pc_type none.
990 
991     Options Database:
992 .     -pc_fieldsplit_schur_precondition <self,user,diag> default is diag
993 
994     Level: intermediate
995 
996 .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields(), PCFieldSplitSchurPreType
997 
998 @*/
999 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1000 {
1001   PetscErrorCode ierr,(*f)(PC,PCFieldSplitSchurPreType,Mat);
1002 
1003   PetscFunctionBegin;
1004   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1005   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSchurPrecondition_C",(void (**)(void))&f);CHKERRQ(ierr);
1006   if (f) {
1007     ierr = (*f)(pc,ptype,pre);CHKERRQ(ierr);
1008   }
1009   PetscFunctionReturn(0);
1010 }
1011 
1012 EXTERN_C_BEGIN
1013 #undef __FUNCT__
1014 #define __FUNCT__ "PCFieldSplitSchurPrecondition_FieldSplit"
1015 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSchurPrecondition_FieldSplit(PC pc,PCFieldSplitSchurPreType ptype,Mat pre)
1016 {
1017   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1018   PetscErrorCode  ierr;
1019 
1020   PetscFunctionBegin;
1021   jac->schurpre = ptype;
1022   if (pre) {
1023     if (jac->schur_user) {ierr = MatDestroy(jac->schur_user);CHKERRQ(ierr);}
1024     jac->schur_user = pre;
1025     ierr = PetscObjectReference((PetscObject)jac->schur_user);CHKERRQ(ierr);
1026   }
1027   PetscFunctionReturn(0);
1028 }
1029 EXTERN_C_END
1030 
1031 #undef __FUNCT__
1032 #define __FUNCT__ "PCFieldSplitGetSchurBlocks"
1033 /*@C
1034    PCFieldSplitGetSchurBlocks - Gets the all matrix blocks for the Schur complement
1035 
1036    Collective on KSP
1037 
1038    Input Parameter:
1039 .  pc - the preconditioner context
1040 
1041    Output Parameters:
1042 +  A - the (0,0) block
1043 .  B - the (0,1) block
1044 .  C - the (1,0) block
1045 -  D - the (1,1) block
1046 
1047    Level: advanced
1048 
1049 .seealso: PCFIELDSPLIT
1050 @*/
1051 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSchurBlocks(PC pc,Mat *A,Mat *B,Mat *C, Mat *D)
1052 {
1053   PC_FieldSplit *jac = (PC_FieldSplit *) pc->data;
1054 
1055   PetscFunctionBegin;
1056   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1057   if (jac->type != PC_COMPOSITE_SCHUR) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "FieldSplit is not using a Schur complement approach.");
1058   if (A) *A = jac->pmat[0];
1059   if (B) *B = jac->B;
1060   if (C) *C = jac->C;
1061   if (D) *D = jac->pmat[1];
1062   PetscFunctionReturn(0);
1063 }
1064 
1065 EXTERN_C_BEGIN
1066 #undef __FUNCT__
1067 #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
1068 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
1069 {
1070   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1071   PetscErrorCode ierr;
1072 
1073   PetscFunctionBegin;
1074   jac->type = type;
1075   if (type == PC_COMPOSITE_SCHUR) {
1076     pc->ops->apply = PCApply_FieldSplit_Schur;
1077     pc->ops->view  = PCView_FieldSplit_Schur;
1078     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit_Schur",PCFieldSplitGetSubKSP_FieldSplit_Schur);CHKERRQ(ierr);
1079     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","PCFieldSplitSchurPrecondition_FieldSplit",PCFieldSplitSchurPrecondition_FieldSplit);CHKERRQ(ierr);
1080 
1081   } else {
1082     pc->ops->apply = PCApply_FieldSplit;
1083     pc->ops->view  = PCView_FieldSplit;
1084     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
1085     ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSchurPrecondition_C","",0);CHKERRQ(ierr);
1086   }
1087   PetscFunctionReturn(0);
1088 }
1089 EXTERN_C_END
1090 
1091 EXTERN_C_BEGIN
1092 #undef __FUNCT__
1093 #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
1094 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
1095 {
1096   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
1097 
1098   PetscFunctionBegin;
1099   if (bs < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
1100   if (jac->bs > 0 && jac->bs != bs) SETERRQ2(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs);
1101   jac->bs = bs;
1102   PetscFunctionReturn(0);
1103 }
1104 EXTERN_C_END
1105 
1106 #undef __FUNCT__
1107 #define __FUNCT__ "PCFieldSplitSetType"
1108 /*@
1109    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
1110 
1111    Collective on PC
1112 
1113    Input Parameter:
1114 .  pc - the preconditioner context
1115 .  type - PC_COMPOSITE_ADDITIVE, PC_COMPOSITE_MULTIPLICATIVE (default), PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL, PC_COMPOSITE_SCHUR
1116 
1117    Options Database Key:
1118 .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special, schur> - Sets fieldsplit preconditioner type
1119 
1120    Level: Developer
1121 
1122 .keywords: PC, set, type, composite preconditioner, additive, multiplicative
1123 
1124 .seealso: PCCompositeSetType()
1125 
1126 @*/
1127 PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
1128 {
1129   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
1130 
1131   PetscFunctionBegin;
1132   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1133   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
1134   if (f) {
1135     ierr = (*f)(pc,type);CHKERRQ(ierr);
1136   }
1137   PetscFunctionReturn(0);
1138 }
1139 
1140 /* -------------------------------------------------------------------------------------*/
1141 /*MC
1142    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
1143                   fields or groups of fields
1144 
1145 
1146      To set options on the solvers for each block append -fieldsplit_ to all the PC
1147         options database keys. For example, -fieldsplit_pc_type ilu -fieldsplit_pc_factor_levels 1
1148 
1149      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
1150          and set the options directly on the resulting KSP object
1151 
1152    Level: intermediate
1153 
1154    Options Database Keys:
1155 +   -pc_fieldsplit_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
1156 .   -pc_fieldsplit_default - automatically add any fields to additional splits that have not
1157                               been supplied explicitly by -pc_fieldsplit_%d_fields
1158 .   -pc_fieldsplit_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
1159 .    -pc_fieldsplit_type <additive,multiplicative,schur,symmetric_multiplicative>
1160 .    -pc_fieldsplit_schur_precondition <true,false> default is true
1161 
1162 -    Options prefix for inner solvers when using Schur complement preconditioner are -fieldsplit_0_ and -fieldsplit_1_
1163      for all other solvers they are -fieldsplit_%d_ for the dth field, use -fieldsplit_ for all fields
1164 
1165 
1166    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
1167      to define a field by an arbitrary collection of entries.
1168 
1169       If no fields are set the default is used. The fields are defined by entries strided by bs,
1170       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
1171       if this is not called the block size defaults to the blocksize of the second matrix passed
1172       to KSPSetOperators()/PCSetOperators().
1173 
1174       Currently for the multiplicative version, the updated residual needed for the next field
1175      solve is computed via a matrix vector product over the entire array. An optimization would be
1176      to update the residual only for the part of the right hand side associated with the next field
1177      solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the
1178      part of the matrix needed to just update part of the residual).
1179 
1180      For the Schur complement preconditioner if J = ( A B )
1181                                                     ( C D )
1182      the preconditioner is
1183               (I   -B inv(A)) ( inv(A)   0    ) (I         0  )
1184               (0    I       ) (   0    inv(S) ) (-C inv(A) I  )
1185      where the action of inv(A) is applied using the KSP solver with prefix -fieldsplit_0_. The action of
1186      inv(S) is computed using the KSP solver with prefix -schur_. For PCFieldSplitGetKSP() when field number is
1187      0 it returns the KSP associated with -fieldsplit_0_ while field number 1 gives -fieldsplit_1_ KSP. By default
1188      D is used to construct a preconditioner for S, use PCFieldSplitSchurPrecondition() to turn on or off this
1189      option.
1190 
1191      If only one set of indices (one IS) is provided with PCFieldSplitSetIS() then the complement of that IS
1192      is used automatically for a second block.
1193 
1194    Concepts: physics based preconditioners, block preconditioners
1195 
1196 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, Block_Preconditioners
1197            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
1198 M*/
1199 
1200 EXTERN_C_BEGIN
1201 #undef __FUNCT__
1202 #define __FUNCT__ "PCCreate_FieldSplit"
1203 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
1204 {
1205   PetscErrorCode ierr;
1206   PC_FieldSplit  *jac;
1207 
1208   PetscFunctionBegin;
1209   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
1210   jac->bs        = -1;
1211   jac->nsplits   = 0;
1212   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
1213   jac->schurpre  = PC_FIELDSPLIT_SCHUR_PRE_USER; /* Try user preconditioner first, fall back on diagonal */
1214 
1215   pc->data     = (void*)jac;
1216 
1217   pc->ops->apply             = PCApply_FieldSplit;
1218   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
1219   pc->ops->setup             = PCSetUp_FieldSplit;
1220   pc->ops->destroy           = PCDestroy_FieldSplit;
1221   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
1222   pc->ops->view              = PCView_FieldSplit;
1223   pc->ops->applyrichardson   = 0;
1224 
1225   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
1226                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
1227   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
1228                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
1229   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
1230                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
1231   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
1232                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
1233   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
1234                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
1235   PetscFunctionReturn(0);
1236 }
1237 EXTERN_C_END
1238 
1239 
1240 /*MC
1241    Block_Preconditioners - PETSc provides support for a variety of "block preconditioners", this provides an
1242           overview of these methods.
1243 
1244       Consider the solution to ( A B ) (x_1)  =  (b_1)
1245                                ( C D ) (x_2)     (b_2)
1246 
1247       Important special cases, the Stokes equation: C = B' and D = 0  (A   B) (x_1) = (b_1)
1248                                                                        B'  0) (x_2)   (b_2)
1249 
1250       One of the goals of the PCFieldSplit preconditioner in PETSc is to provide a variety of preconditioners
1251       for this block system.
1252 
1253       Consider an additional matrix (Ap  Bp)
1254                                     (Cp  Dp) where some or all of the entries may be the same as
1255       in the original matrix (for example Ap == A).
1256 
1257       In the following, A^ denotes the approximate application of the inverse of A, possibly using Ap in the
1258       approximation. In PETSc this simply means one has called KSPSetOperators(ksp,A,Ap,...) or KSPSetOperators(ksp,Ap,Ap,...)
1259 
1260       Block Jacobi:   x_1 = A^ b_1
1261                       x_2 = D^ b_2
1262 
1263       Lower block Gauss-Seidel:   x_1 = A^ b_1
1264                             x_2 = D^ (b_2 - C x_1)       variant x_2 = D^ (b_2 - Cp x_1)
1265 
1266       Symmetric Gauss-Seidel:  x_1 = x_1 + A^(b_1 - A x_1 - B x_2)    variant  x_1 = x_1 + A^(b_1 - Ap x_1 - Bp x_2)
1267           Interestingly this form is not actually a symmetric matrix, the symmetric version is
1268                               x_1 = A^(b_1 - B x_2)      variant x_1 = A^(b_1 - Bp x_2)
1269 
1270    Level: intermediate
1271 
1272 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCFIELDSPLIT
1273            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS(), PCFieldSplitSchurPrecondition()
1274 M*/
1275