xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision d07a9264daf19612ed0b9951a5a70bfd0f76d6cd)
10971522cSBarry Smith /*
20971522cSBarry Smith 
30971522cSBarry Smith */
40971522cSBarry Smith #include "src/ksp/pc/pcimpl.h"     /*I "petscpc.h" I*/
50971522cSBarry Smith 
60971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
70971522cSBarry Smith struct _PC_FieldSplitLink {
869a612a9SBarry Smith   KSP               ksp;
90971522cSBarry Smith   Vec               x,y;
100971522cSBarry Smith   PetscInt          nfields;
110971522cSBarry Smith   PetscInt          *fields;
121b9fc7fcSBarry Smith   VecScatter        sctx;
130971522cSBarry Smith   PC_FieldSplitLink next;
140971522cSBarry Smith };
150971522cSBarry Smith 
160971522cSBarry Smith typedef struct {
1779416396SBarry Smith   PCCompositeType   type;              /* additive or multiplicative */
1897bbdb24SBarry Smith   PetscTruth        defaultsplit;
190971522cSBarry Smith   PetscInt          bs;
200971522cSBarry Smith   PetscInt          nsplits;
2179416396SBarry Smith   Vec               *x,*y,w1,w2;
2297bbdb24SBarry Smith   Mat               *pmat;
2397bbdb24SBarry Smith   IS                *is;
2497bbdb24SBarry Smith   PC_FieldSplitLink head;
250971522cSBarry Smith } PC_FieldSplit;
260971522cSBarry Smith 
270971522cSBarry Smith #undef __FUNCT__
280971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
290971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
300971522cSBarry Smith {
310971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
320971522cSBarry Smith   PetscErrorCode    ierr;
330971522cSBarry Smith   PetscTruth        iascii;
340971522cSBarry Smith   PetscInt          i,j;
350971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
3679416396SBarry Smith   const char        *types[] = {"additive","multiplicative"};
370971522cSBarry Smith 
380971522cSBarry Smith   PetscFunctionBegin;
390971522cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
400971522cSBarry Smith   if (iascii) {
4179416396SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D",types[jac->type],jac->nsplits);CHKERRQ(ierr);
4269a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
430971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
440971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
450971522cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
4679416396SBarry Smith       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
470971522cSBarry Smith       for (j=0; j<link->nfields; j++) {
4879416396SBarry Smith         if (j > 0) {
4979416396SBarry Smith           ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
5079416396SBarry Smith         }
5179416396SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer," %D",link->fields[j]);CHKERRQ(ierr);
520971522cSBarry Smith       }
530971522cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
5479416396SBarry Smith       ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
5569a612a9SBarry Smith       ierr = KSPView(link->ksp,viewer);CHKERRQ(ierr);
560971522cSBarry Smith       link = link->next;
570971522cSBarry Smith     }
580971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
590971522cSBarry Smith   } else {
600971522cSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
610971522cSBarry Smith   }
620971522cSBarry Smith   PetscFunctionReturn(0);
630971522cSBarry Smith }
640971522cSBarry Smith 
650971522cSBarry Smith #undef __FUNCT__
6669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
6769a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
680971522cSBarry Smith {
690971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
700971522cSBarry Smith   PetscErrorCode    ierr;
710971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
7269a612a9SBarry Smith   PetscInt          i;
7379416396SBarry Smith   PetscTruth        flg = PETSC_FALSE,*fields;
740971522cSBarry Smith 
750971522cSBarry Smith   PetscFunctionBegin;
7679416396SBarry Smith   ierr = PetscOptionsGetLogical(pc->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
7779416396SBarry Smith   if (!link || flg) {
781b9fc7fcSBarry Smith     ierr = PetscLogInfo(pc,"PCFieldSplitSetDefaults: Using default splitting of fields");CHKERRQ(ierr);
79521d7252SBarry Smith     if (jac->bs <= 0) {
80521d7252SBarry Smith       ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
81521d7252SBarry Smith     }
8279416396SBarry Smith     ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
8379416396SBarry Smith     ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
8479416396SBarry Smith     while (link) {
8579416396SBarry Smith       for (i=0; i<link->nfields; i++) {
8679416396SBarry Smith         fields[link->fields[i]] = PETSC_TRUE;
87521d7252SBarry Smith       }
8879416396SBarry Smith       link = link->next;
8979416396SBarry Smith     }
9097bbdb24SBarry Smith     jac->defaultsplit = PETSC_TRUE;
9179416396SBarry Smith     for (i=0; i<jac->bs; i++) {
9279416396SBarry Smith       if (!fields[i]) {
9379416396SBarry Smith 	ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
9479416396SBarry Smith       } else {
9579416396SBarry Smith         jac->defaultsplit = PETSC_FALSE;
9679416396SBarry Smith       }
9779416396SBarry Smith     }
98521d7252SBarry Smith   }
9969a612a9SBarry Smith   PetscFunctionReturn(0);
10069a612a9SBarry Smith }
10169a612a9SBarry Smith 
10269a612a9SBarry Smith 
10369a612a9SBarry Smith #undef __FUNCT__
10469a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
10569a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
10669a612a9SBarry Smith {
10769a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
10869a612a9SBarry Smith   PetscErrorCode    ierr;
10969a612a9SBarry Smith   PC_FieldSplitLink link;
11069a612a9SBarry Smith   PetscInt          i,nsplit;
11169a612a9SBarry Smith   MatStructure      flag = pc->flag;
11269a612a9SBarry Smith 
11369a612a9SBarry Smith   PetscFunctionBegin;
11469a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
11597bbdb24SBarry Smith   nsplit = jac->nsplits;
11669a612a9SBarry Smith   link   = jac->head;
11797bbdb24SBarry Smith 
11897bbdb24SBarry Smith   /* get the matrices for each split */
11997bbdb24SBarry Smith   if (!jac->is) {
1201b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
12197bbdb24SBarry Smith 
1221b9fc7fcSBarry Smith     ierr   = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);
12397bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
1241b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
1251b9fc7fcSBarry Smith     ierr   = PetscMalloc(nsplit*sizeof(IS),&jac->is);CHKERRQ(ierr);
1261b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
1271b9fc7fcSBarry Smith       if (jac->defaultsplit) {
1281b9fc7fcSBarry Smith 	ierr = ISCreateStride(pc->comm,nslots,rstart+i,nsplit,&jac->is[i]);CHKERRQ(ierr);
12997bbdb24SBarry Smith       } else {
1301b9fc7fcSBarry Smith         PetscInt   *ii,j,k,nfields = link->nfields,*fields = link->fields;
1311b9fc7fcSBarry Smith         PetscTruth sorted;
1321b9fc7fcSBarry Smith         ierr = PetscMalloc(link->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
1331b9fc7fcSBarry Smith         for (j=0; j<nslots; j++) {
1341b9fc7fcSBarry Smith           for (k=0; k<nfields; k++) {
1351b9fc7fcSBarry Smith             ii[nfields*j + k] = rstart + bs*j + fields[k];
13697bbdb24SBarry Smith           }
13797bbdb24SBarry Smith         }
1381b9fc7fcSBarry Smith 	ierr = ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);CHKERRQ(ierr);
1391b9fc7fcSBarry Smith         ierr = ISSorted(jac->is[i],&sorted);CHKERRQ(ierr);
1401b9fc7fcSBarry Smith         if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
1411b9fc7fcSBarry Smith         ierr = PetscFree(ii);CHKERRQ(ierr);
1421b9fc7fcSBarry Smith         link = link->next;
1431b9fc7fcSBarry Smith       }
1441b9fc7fcSBarry Smith     }
1451b9fc7fcSBarry Smith   }
1461b9fc7fcSBarry Smith 
14797bbdb24SBarry Smith   if (!jac->pmat) {
14897bbdb24SBarry Smith     ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_INITIAL_MATRIX,&jac->pmat);CHKERRQ(ierr);
14997bbdb24SBarry Smith   } else {
15097bbdb24SBarry Smith     ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_REUSE_MATRIX,&jac->pmat);CHKERRQ(ierr);
15197bbdb24SBarry Smith   }
15297bbdb24SBarry Smith 
15397bbdb24SBarry Smith   /* set up the individual PCs */
15497bbdb24SBarry Smith   i    = 0;
1551b9fc7fcSBarry Smith   link = jac->head;
1560971522cSBarry Smith   while (link) {
15769a612a9SBarry Smith     ierr = KSPSetOperators(link->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
15869a612a9SBarry Smith     ierr = KSPSetFromOptions(link->ksp);CHKERRQ(ierr);
15969a612a9SBarry Smith     ierr = KSPSetUp(link->ksp);CHKERRQ(ierr);
16097bbdb24SBarry Smith     i++;
1610971522cSBarry Smith     link = link->next;
1620971522cSBarry Smith   }
16397bbdb24SBarry Smith 
16497bbdb24SBarry Smith   /* create work vectors for each split */
1651b9fc7fcSBarry Smith   if (!jac->x) {
16679416396SBarry Smith     Vec xtmp;
16797bbdb24SBarry Smith     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
16897bbdb24SBarry Smith     link = jac->head;
16997bbdb24SBarry Smith     for (i=0; i<nsplit; i++) {
17097bbdb24SBarry Smith       Mat A;
17169a612a9SBarry Smith       ierr      = KSPGetOperators(link->ksp,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr);
17297bbdb24SBarry Smith       ierr      = MatGetVecs(A,&link->x,&link->y);CHKERRQ(ierr);
17397bbdb24SBarry Smith       jac->x[i] = link->x;
17497bbdb24SBarry Smith       jac->y[i] = link->y;
17597bbdb24SBarry Smith       link      = link->next;
17697bbdb24SBarry Smith     }
17779416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
1781b9fc7fcSBarry Smith 
1791b9fc7fcSBarry Smith     link = jac->head;
1801b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
1811b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
1821b9fc7fcSBarry Smith       ierr = VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&link->sctx);CHKERRQ(ierr);
1831b9fc7fcSBarry Smith       link = link->next;
1841b9fc7fcSBarry Smith     }
1851b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
1861b9fc7fcSBarry Smith   }
1870971522cSBarry Smith   PetscFunctionReturn(0);
1880971522cSBarry Smith }
1890971522cSBarry Smith 
19079416396SBarry Smith #define FieldSplitSplitSolveAdd(link,xx,yy) \
19179416396SBarry Smith     (VecScatterBegin(xx,link->x,INSERT_VALUES,SCATTER_FORWARD,link->sctx) || \
19279416396SBarry Smith      VecScatterEnd(xx,link->x,INSERT_VALUES,SCATTER_FORWARD,link->sctx) || \
19379416396SBarry Smith      KSPSolve(link->ksp,link->x,link->y) || \
19479416396SBarry Smith      VecScatterBegin(link->y,yy,ADD_VALUES,SCATTER_REVERSE,link->sctx) || \
19579416396SBarry Smith      VecScatterEnd(link->y,yy,ADD_VALUES,SCATTER_REVERSE,link->sctx))
19679416396SBarry Smith 
1970971522cSBarry Smith #undef __FUNCT__
1980971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
1990971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
2000971522cSBarry Smith {
2010971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2020971522cSBarry Smith   PetscErrorCode    ierr;
2030971522cSBarry Smith   PC_FieldSplitLink link = jac->head;
20479416396SBarry Smith   PetscScalar       zero = 0.0,mone = -1.0;
2050971522cSBarry Smith 
2060971522cSBarry Smith   PetscFunctionBegin;
20779416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
2081b9fc7fcSBarry Smith     if (jac->defaultsplit) {
2090971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
2100971522cSBarry Smith       while (link) {
21169a612a9SBarry Smith 	ierr = KSPSolve(link->ksp,link->x,link->y);CHKERRQ(ierr);
2120971522cSBarry Smith 	link = link->next;
2130971522cSBarry Smith       }
2140971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
2151b9fc7fcSBarry Smith     } else {
2161b9fc7fcSBarry Smith       PetscInt    i = 0;
2171b9fc7fcSBarry Smith 
2181b9fc7fcSBarry Smith       ierr = VecSet(&zero,y);CHKERRQ(ierr);
2191b9fc7fcSBarry Smith       while (link) {
22079416396SBarry Smith         ierr = FieldSplitSplitSolveAdd(link,x,y);CHKERRQ(ierr);
2211b9fc7fcSBarry Smith 	link = link->next;
2221b9fc7fcSBarry Smith 	i++;
2231b9fc7fcSBarry Smith       }
2241b9fc7fcSBarry Smith     }
22579416396SBarry Smith   } else {
22679416396SBarry Smith     if (!jac->w1) {
22779416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
22879416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
22979416396SBarry Smith     }
23079416396SBarry Smith     ierr = VecSet(&zero,y);CHKERRQ(ierr);
23179416396SBarry Smith     ierr = FieldSplitSplitSolveAdd(link,x,y);CHKERRQ(ierr);
23279416396SBarry Smith     while (link->next) {
23379416396SBarry Smith       link = link->next;
23479416396SBarry Smith       ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr);
23579416396SBarry Smith       ierr = VecWAXPY(&mone,jac->w1,x,jac->w2);CHKERRQ(ierr);
23679416396SBarry Smith       ierr = FieldSplitSplitSolveAdd(link,jac->w2,y);CHKERRQ(ierr);
23779416396SBarry Smith     }
23879416396SBarry Smith   }
2390971522cSBarry Smith   PetscFunctionReturn(0);
2400971522cSBarry Smith }
2410971522cSBarry Smith 
2420971522cSBarry Smith #undef __FUNCT__
2430971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
2440971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
2450971522cSBarry Smith {
2460971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2470971522cSBarry Smith   PetscErrorCode    ierr;
24897bbdb24SBarry Smith   PC_FieldSplitLink link = jac->head,next;
2490971522cSBarry Smith 
2500971522cSBarry Smith   PetscFunctionBegin;
2510971522cSBarry Smith   while (link) {
25269a612a9SBarry Smith     ierr = KSPDestroy(link->ksp);CHKERRQ(ierr);
25369a612a9SBarry Smith     if (link->x) {ierr = VecDestroy(link->x);CHKERRQ(ierr);}
25469a612a9SBarry Smith     if (link->y) {ierr = VecDestroy(link->y);CHKERRQ(ierr);}
2551b9fc7fcSBarry Smith     if (link->sctx) {ierr = VecScatterDestroy(link->sctx);CHKERRQ(ierr);}
25697bbdb24SBarry Smith     next = link->next;
2570971522cSBarry Smith     ierr = PetscFree2(link,link->fields);CHKERRQ(ierr);
25897bbdb24SBarry Smith     link = next;
2590971522cSBarry Smith   }
2600971522cSBarry Smith   if (jac->x) {ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);}
26197bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
26297bbdb24SBarry Smith   if (jac->is) {
26397bbdb24SBarry Smith     PetscInt i;
26497bbdb24SBarry Smith     for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);}
26597bbdb24SBarry Smith     ierr = PetscFree(jac->is);CHKERRQ(ierr);
26697bbdb24SBarry Smith   }
26779416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
26879416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
2690971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
2700971522cSBarry Smith   PetscFunctionReturn(0);
2710971522cSBarry Smith }
2720971522cSBarry Smith 
2730971522cSBarry Smith #undef __FUNCT__
2740971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
2750971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
2761b9fc7fcSBarry Smith /*   This does not call KSPSetFromOptions() on the subksp's, see PCSetFromOptionsBJacobi/ASM() */
2770971522cSBarry Smith {
2781b9fc7fcSBarry Smith   PetscErrorCode ierr;
27979416396SBarry Smith   PetscInt       i = 0,nfields,fields[12],indx;
2801b9fc7fcSBarry Smith   PetscTruth     flg;
2811b9fc7fcSBarry Smith   char           optionname[128];
28279416396SBarry Smith   const char     *types[] = {"additive","multiplicative"};
2831b9fc7fcSBarry Smith 
2840971522cSBarry Smith   PetscFunctionBegin;
2851b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
28679416396SBarry Smith   ierr = PetscOptionsEList("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",types,2,types[0],&indx,&flg);CHKERRQ(ierr);
28779416396SBarry Smith   if (flg) {
28879416396SBarry Smith     ierr = PCFieldSplitSetType(pc,(PCCompositeType)indx);CHKERRQ(ierr);
28979416396SBarry Smith   }
2901b9fc7fcSBarry Smith   while (PETSC_TRUE) {
2911b9fc7fcSBarry Smith     sprintf(optionname,"-pc_fieldsplit_%d_fields",i++);
2921b9fc7fcSBarry Smith     nfields = 12;
2931b9fc7fcSBarry Smith     ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
2941b9fc7fcSBarry Smith     if (!flg) break;
2951b9fc7fcSBarry Smith     if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
2961b9fc7fcSBarry Smith     ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
2971b9fc7fcSBarry Smith   }
2981b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
2990971522cSBarry Smith   PetscFunctionReturn(0);
3000971522cSBarry Smith }
3010971522cSBarry Smith 
3020971522cSBarry Smith /*------------------------------------------------------------------------------------*/
3030971522cSBarry Smith 
3040971522cSBarry Smith EXTERN_C_BEGIN
3050971522cSBarry Smith #undef __FUNCT__
3060971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
3070971522cSBarry Smith PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
3080971522cSBarry Smith {
30997bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3100971522cSBarry Smith   PetscErrorCode    ierr;
3110971522cSBarry Smith   PC_FieldSplitLink link,next = jac->head;
31269a612a9SBarry Smith   char              prefix[128];
3130971522cSBarry Smith 
3140971522cSBarry Smith   PetscFunctionBegin;
3150971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
31697bbdb24SBarry Smith   ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&link,n,PetscInt,&link->fields);CHKERRQ(ierr);
31797bbdb24SBarry Smith   ierr = PetscMemcpy(link->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
3180971522cSBarry Smith   link->nfields = n;
31997bbdb24SBarry Smith   link->next    = PETSC_NULL;
32069a612a9SBarry Smith   ierr          = KSPCreate(pc->comm,&link->ksp);CHKERRQ(ierr);
32169a612a9SBarry Smith   ierr          = KSPSetType(link->ksp,KSPPREONLY);CHKERRQ(ierr);
32269a612a9SBarry Smith 
32369a612a9SBarry Smith   if (pc->prefix) {
3242e70eadcSBarry Smith     sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,jac->nsplits);
32569a612a9SBarry Smith   } else {
32669a612a9SBarry Smith     sprintf(prefix,"fieldsplit_%d_",jac->nsplits);
32769a612a9SBarry Smith   }
32869a612a9SBarry Smith   ierr = KSPSetOptionsPrefix(link->ksp,prefix);CHKERRQ(ierr);
3290971522cSBarry Smith 
3300971522cSBarry Smith   if (!next) {
3310971522cSBarry Smith     jac->head = link;
3320971522cSBarry Smith   } else {
3330971522cSBarry Smith     while (next->next) {
3340971522cSBarry Smith       next = next->next;
3350971522cSBarry Smith     }
3360971522cSBarry Smith     next->next = link;
3370971522cSBarry Smith   }
3380971522cSBarry Smith   jac->nsplits++;
3390971522cSBarry Smith   PetscFunctionReturn(0);
3400971522cSBarry Smith }
3410971522cSBarry Smith EXTERN_C_END
3420971522cSBarry Smith 
3430971522cSBarry Smith 
3440971522cSBarry Smith EXTERN_C_BEGIN
3450971522cSBarry Smith #undef __FUNCT__
34669a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
34769a612a9SBarry Smith PetscErrorCode PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
3480971522cSBarry Smith {
3490971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
3500971522cSBarry Smith   PetscErrorCode    ierr;
3510971522cSBarry Smith   PetscInt          cnt = 0;
352*d07a9264SSatish Balay   PC_FieldSplitLink link = jac->head;
3530971522cSBarry Smith 
3540971522cSBarry Smith   PetscFunctionBegin;
35569a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
3560971522cSBarry Smith   while (link) {
35769a612a9SBarry Smith     (*subksp)[cnt++] = link->ksp;
3580971522cSBarry Smith     link = link->next;
3590971522cSBarry Smith   }
3600971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
3610971522cSBarry Smith   *n = jac->nsplits;
3620971522cSBarry Smith   PetscFunctionReturn(0);
3630971522cSBarry Smith }
3640971522cSBarry Smith EXTERN_C_END
3650971522cSBarry Smith 
3660971522cSBarry Smith #undef __FUNCT__
3670971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
3680971522cSBarry Smith /*@
3690971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
3700971522cSBarry Smith 
3710971522cSBarry Smith     Collective on PC
3720971522cSBarry Smith 
3730971522cSBarry Smith     Input Parameters:
3740971522cSBarry Smith +   pc  - the preconditioner context
3750971522cSBarry Smith .   n - the number of fields in this split
3760971522cSBarry Smith .   fields - the fields in this split
3770971522cSBarry Smith 
3780971522cSBarry Smith     Level: intermediate
3790971522cSBarry Smith 
38069a612a9SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT
3810971522cSBarry Smith 
3820971522cSBarry Smith @*/
3830971522cSBarry Smith PetscErrorCode PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
3840971522cSBarry Smith {
3850971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
3860971522cSBarry Smith 
3870971522cSBarry Smith   PetscFunctionBegin;
3880971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3890971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
3900971522cSBarry Smith   if (f) {
3910971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
3920971522cSBarry Smith   }
3930971522cSBarry Smith   PetscFunctionReturn(0);
3940971522cSBarry Smith }
3950971522cSBarry Smith 
3960971522cSBarry Smith #undef __FUNCT__
39769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
3980971522cSBarry Smith /*@C
39969a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
4000971522cSBarry Smith 
40169a612a9SBarry Smith    Collective on KSP
4020971522cSBarry Smith 
4030971522cSBarry Smith    Input Parameter:
4040971522cSBarry Smith .  pc - the preconditioner context
4050971522cSBarry Smith 
4060971522cSBarry Smith    Output Parameters:
4070971522cSBarry Smith +  n - the number of split
40869a612a9SBarry Smith -  pc - the array of KSP contexts
4090971522cSBarry Smith 
4100971522cSBarry Smith    Note:
41169a612a9SBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed
4120971522cSBarry Smith 
41369a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
4140971522cSBarry Smith 
4150971522cSBarry Smith    Level: advanced
4160971522cSBarry Smith 
4170971522cSBarry Smith .seealso: PCFIELDSPLIT
4180971522cSBarry Smith @*/
41969a612a9SBarry Smith PetscErrorCode PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
4200971522cSBarry Smith {
42169a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
4220971522cSBarry Smith 
4230971522cSBarry Smith   PetscFunctionBegin;
4240971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4250971522cSBarry Smith   PetscValidIntPointer(n,2);
42669a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
4270971522cSBarry Smith   if (f) {
42869a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
4290971522cSBarry Smith   } else {
43069a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
4310971522cSBarry Smith   }
4320971522cSBarry Smith   PetscFunctionReturn(0);
4330971522cSBarry Smith }
4340971522cSBarry Smith 
43579416396SBarry Smith EXTERN_C_BEGIN
43679416396SBarry Smith #undef __FUNCT__
43779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
43879416396SBarry Smith PetscErrorCode PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
43979416396SBarry Smith {
44079416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
44179416396SBarry Smith 
44279416396SBarry Smith   PetscFunctionBegin;
44379416396SBarry Smith   jac->type = type;
44479416396SBarry Smith   PetscFunctionReturn(0);
44579416396SBarry Smith }
44679416396SBarry Smith EXTERN_C_END
44779416396SBarry Smith 
44879416396SBarry Smith #undef __FUNCT__
44979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
45079416396SBarry Smith /*@C
45179416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
45279416396SBarry Smith 
45379416396SBarry Smith    Collective on PC
45479416396SBarry Smith 
45579416396SBarry Smith    Input Parameter:
45679416396SBarry Smith .  pc - the preconditioner context
45779416396SBarry Smith .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE
45879416396SBarry Smith 
45979416396SBarry Smith    Options Database Key:
46079416396SBarry Smith .  -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type
46179416396SBarry Smith 
46279416396SBarry Smith    Level: Developer
46379416396SBarry Smith 
46479416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
46579416396SBarry Smith 
46679416396SBarry Smith .seealso: PCCompositeSetType()
46779416396SBarry Smith 
46879416396SBarry Smith @*/
46979416396SBarry Smith PetscErrorCode PCFieldSplitSetType(PC pc,PCCompositeType type)
47079416396SBarry Smith {
47179416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
47279416396SBarry Smith 
47379416396SBarry Smith   PetscFunctionBegin;
47479416396SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
47579416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
47679416396SBarry Smith   if (f) {
47779416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
47879416396SBarry Smith   }
47979416396SBarry Smith   PetscFunctionReturn(0);
48079416396SBarry Smith }
48179416396SBarry Smith 
4820971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
4830971522cSBarry Smith /*MC
4840971522cSBarry Smith    PCFieldSplit - Preconditioner created by combining seperate preconditioners for individual
4850971522cSBarry Smith                   fields or groups of fields
4860971522cSBarry Smith 
4870971522cSBarry Smith 
4880971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
4890971522cSBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_ilu_levels 1
4900971522cSBarry Smith 
49169a612a9SBarry Smith      To set the options on the solvers seperate for each block call PCFieldSplitGetSubKSP()
49269a612a9SBarry Smith          and set the options directly on the resulting KSP object
4930971522cSBarry Smith 
4940971522cSBarry Smith    Level: intermediate
4950971522cSBarry Smith 
49679416396SBarry Smith    Options Database Keys:
4972e70eadcSBarry Smith +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
4982e70eadcSBarry Smith .   -pc_splitfield_default - automatically add any fields to additional splits that have not
4992e70eadcSBarry Smith                               been supplied explicitly by -pc_splitfield_%d_fields
5002e70eadcSBarry Smith -   -pc_splitfield_type <additive,multiplicative>
50179416396SBarry Smith 
5020971522cSBarry Smith    Concepts: physics based preconditioners
5030971522cSBarry Smith 
5040971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
50569a612a9SBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields()
5060971522cSBarry Smith M*/
5070971522cSBarry Smith 
5080971522cSBarry Smith EXTERN_C_BEGIN
5090971522cSBarry Smith #undef __FUNCT__
5100971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
5110971522cSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc)
5120971522cSBarry Smith {
5130971522cSBarry Smith   PetscErrorCode ierr;
5140971522cSBarry Smith   PC_FieldSplit  *jac;
5150971522cSBarry Smith 
5160971522cSBarry Smith   PetscFunctionBegin;
5170971522cSBarry Smith   ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr);
5180971522cSBarry Smith   PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));
5190971522cSBarry Smith   ierr = PetscMemzero(jac,sizeof(PC_FieldSplit));CHKERRQ(ierr);
5200971522cSBarry Smith   jac->bs      = -1;
5210971522cSBarry Smith   jac->nsplits = 0;
52279416396SBarry Smith   jac->type    = PC_COMPOSITE_ADDITIVE;
5230971522cSBarry Smith   pc->data     = (void*)jac;
5240971522cSBarry Smith 
5250971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
5260971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
5270971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
5280971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
5290971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
5300971522cSBarry Smith   pc->ops->applyrichardson   = 0;
5310971522cSBarry Smith 
53269a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
53369a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
5340971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
5350971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
53679416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
53779416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
5380971522cSBarry Smith   PetscFunctionReturn(0);
5390971522cSBarry Smith }
5400971522cSBarry Smith EXTERN_C_END
5410971522cSBarry Smith 
5420971522cSBarry Smith 
543