xref: /petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c (revision 11755939e33fc80afac7f43ea383c6942b9c52a8)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
30971522cSBarry Smith /*
40971522cSBarry Smith 
50971522cSBarry Smith */
66356e834SBarry Smith #include "private/pcimpl.h"     /*I "petscpc.h" I*/
70971522cSBarry Smith 
80971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink;
90971522cSBarry Smith struct _PC_FieldSplitLink {
1069a612a9SBarry Smith   KSP               ksp;
110971522cSBarry Smith   Vec               x,y;
120971522cSBarry Smith   PetscInt          nfields;
130971522cSBarry Smith   PetscInt          *fields;
141b9fc7fcSBarry Smith   VecScatter        sctx;
15704ba839SBarry Smith   IS                is,cis;
16704ba839SBarry Smith   PetscInt          csize;
1751f519a2SBarry Smith   PC_FieldSplitLink next,previous;
180971522cSBarry Smith };
190971522cSBarry Smith 
200971522cSBarry Smith typedef struct {
2168dd23aaSBarry Smith   PCCompositeType   type;
2297bbdb24SBarry Smith   PetscTruth        defaultsplit;
230971522cSBarry Smith   PetscInt          bs;
24704ba839SBarry Smith   PetscInt          nsplits;
2579416396SBarry Smith   Vec               *x,*y,w1,w2;
2697bbdb24SBarry Smith   Mat               *pmat;
2768dd23aaSBarry Smith   Mat               *Afield; /* the rows of the matrix associated with each field */
28704ba839SBarry Smith   PetscTruth        issetup;
2997bbdb24SBarry Smith   PC_FieldSplitLink head;
300971522cSBarry Smith } PC_FieldSplit;
310971522cSBarry Smith 
3216913363SBarry Smith /*
3316913363SBarry Smith     Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of
3416913363SBarry Smith    inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the
3516913363SBarry Smith    PC you could change this.
3616913363SBarry Smith */
370971522cSBarry Smith #undef __FUNCT__
380971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit"
390971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer)
400971522cSBarry Smith {
410971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
420971522cSBarry Smith   PetscErrorCode    ierr;
430971522cSBarry Smith   PetscTruth        iascii;
440971522cSBarry Smith   PetscInt          i,j;
455a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
460971522cSBarry Smith 
470971522cSBarry Smith   PetscFunctionBegin;
480971522cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
490971522cSBarry Smith   if (iascii) {
5051f519a2SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr);
5169a612a9SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr);
520971522cSBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
530971522cSBarry Smith     for (i=0; i<jac->nsplits; i++) {
541ab39975SBarry Smith       if (ilink->fields) {
550971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr);
5679416396SBarry Smith 	ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr);
575a9f2f41SSatish Balay 	for (j=0; j<ilink->nfields; j++) {
5879416396SBarry Smith 	  if (j > 0) {
5979416396SBarry Smith 	    ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr);
6079416396SBarry Smith 	  }
615a9f2f41SSatish Balay 	  ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr);
620971522cSBarry Smith 	}
630971522cSBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr);
6479416396SBarry Smith         ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr);
651ab39975SBarry Smith       } else {
661ab39975SBarry Smith 	ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Defined by IS\n",i);CHKERRQ(ierr);
671ab39975SBarry Smith       }
685a9f2f41SSatish Balay       ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr);
695a9f2f41SSatish Balay       ilink = ilink->next;
700971522cSBarry Smith     }
710971522cSBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
720971522cSBarry Smith   } else {
730971522cSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name);
740971522cSBarry Smith   }
750971522cSBarry Smith   PetscFunctionReturn(0);
760971522cSBarry Smith }
770971522cSBarry Smith 
780971522cSBarry Smith #undef __FUNCT__
7969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults"
8069a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc)
810971522cSBarry Smith {
820971522cSBarry Smith   PC_FieldSplit     *jac  = (PC_FieldSplit*)pc->data;
830971522cSBarry Smith   PetscErrorCode    ierr;
845a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
85d32f9abdSBarry Smith   PetscInt          i = 0,*ifields,nfields;
86d32f9abdSBarry Smith   PetscTruth        flg = PETSC_FALSE,*fields,flg2;
87d32f9abdSBarry Smith   char              optionname[128];
880971522cSBarry Smith 
890971522cSBarry Smith   PetscFunctionBegin;
90d32f9abdSBarry Smith   if (!ilink) {
91d32f9abdSBarry Smith 
92521d7252SBarry Smith     if (jac->bs <= 0) {
93704ba839SBarry Smith       if (pc->pmat) {
94521d7252SBarry Smith         ierr   = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr);
95704ba839SBarry Smith       } else {
96704ba839SBarry Smith         jac->bs = 1;
97704ba839SBarry Smith       }
98521d7252SBarry Smith     }
99d32f9abdSBarry Smith 
100d32f9abdSBarry Smith     ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr);
101d32f9abdSBarry Smith     if (!flg) {
102d32f9abdSBarry Smith       /* Allow user to set fields from command line,  if bs was known at the time of PCSetFromOptions_FieldSplit()
103d32f9abdSBarry Smith          then it is set there. This is not ideal because we should only have options set in XXSetFromOptions(). */
104d32f9abdSBarry Smith       flg = PETSC_TRUE; /* switched off automatically if user sets fields manually here */
105d32f9abdSBarry Smith       ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&ifields);CHKERRQ(ierr);
106d32f9abdSBarry Smith       while (PETSC_TRUE) {
107d32f9abdSBarry Smith         sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
108d32f9abdSBarry Smith         nfields = jac->bs;
109d32f9abdSBarry Smith         ierr    = PetscOptionsGetIntArray(((PetscObject)pc)->prefix,optionname,ifields,&nfields,&flg2);CHKERRQ(ierr);
110d32f9abdSBarry Smith         if (!flg2) break;
111d32f9abdSBarry Smith         if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
112d32f9abdSBarry Smith         flg = PETSC_FALSE;
113d32f9abdSBarry Smith         ierr = PCFieldSplitSetFields(pc,nfields,ifields);CHKERRQ(ierr);
114d32f9abdSBarry Smith       }
115d32f9abdSBarry Smith       ierr = PetscFree(ifields);CHKERRQ(ierr);
116d32f9abdSBarry Smith     }
117d32f9abdSBarry Smith 
118d32f9abdSBarry Smith     if (flg) {
119d32f9abdSBarry Smith       ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr);
12079416396SBarry Smith       ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr);
12179416396SBarry Smith       ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr);
1225a9f2f41SSatish Balay       while (ilink) {
1235a9f2f41SSatish Balay 	for (i=0; i<ilink->nfields; i++) {
1245a9f2f41SSatish Balay 	  fields[ilink->fields[i]] = PETSC_TRUE;
125521d7252SBarry Smith 	}
1265a9f2f41SSatish Balay 	ilink = ilink->next;
12779416396SBarry Smith       }
12897bbdb24SBarry Smith       jac->defaultsplit = PETSC_TRUE;
12979416396SBarry Smith       for (i=0; i<jac->bs; i++) {
13079416396SBarry Smith 	if (!fields[i]) {
13179416396SBarry Smith 	  ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr);
13279416396SBarry Smith 	} else {
13379416396SBarry Smith 	  jac->defaultsplit = PETSC_FALSE;
13479416396SBarry Smith 	}
13579416396SBarry Smith       }
13668dd23aaSBarry Smith       ierr = PetscFree(fields);CHKERRQ(ierr);
137521d7252SBarry Smith     }
138d32f9abdSBarry Smith   }
13969a612a9SBarry Smith   PetscFunctionReturn(0);
14069a612a9SBarry Smith }
14169a612a9SBarry Smith 
14269a612a9SBarry Smith 
14369a612a9SBarry Smith #undef __FUNCT__
14469a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit"
14569a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc)
14669a612a9SBarry Smith {
14769a612a9SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
14869a612a9SBarry Smith   PetscErrorCode    ierr;
1495a9f2f41SSatish Balay   PC_FieldSplitLink ilink;
150cf502942SBarry Smith   PetscInt          i,nsplit,ccsize;
15169a612a9SBarry Smith   MatStructure      flag = pc->flag;
15268dd23aaSBarry Smith   PetscTruth        sorted,getsub;
15369a612a9SBarry Smith 
15469a612a9SBarry Smith   PetscFunctionBegin;
15569a612a9SBarry Smith   ierr   = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr);
15697bbdb24SBarry Smith   nsplit = jac->nsplits;
1575a9f2f41SSatish Balay   ilink  = jac->head;
15897bbdb24SBarry Smith 
15997bbdb24SBarry Smith   /* get the matrices for each split */
160704ba839SBarry Smith   if (!jac->issetup) {
1611b9fc7fcSBarry Smith     PetscInt rstart,rend,nslots,bs;
16297bbdb24SBarry Smith 
163704ba839SBarry Smith     jac->issetup = PETSC_TRUE;
164704ba839SBarry Smith 
165704ba839SBarry Smith     /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */
16651f519a2SBarry Smith     bs     = jac->bs;
16797bbdb24SBarry Smith     ierr   = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr);
168cf502942SBarry Smith     ierr   = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr);
1691b9fc7fcSBarry Smith     nslots = (rend - rstart)/bs;
1701b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
1711b9fc7fcSBarry Smith       if (jac->defaultsplit) {
172704ba839SBarry Smith         ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr);
173704ba839SBarry Smith       } else if (!ilink->is) {
174ccb205f8SBarry Smith         if (ilink->nfields > 1) {
1755a9f2f41SSatish Balay           PetscInt   *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields;
1765a9f2f41SSatish Balay           ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr);
1771b9fc7fcSBarry Smith           for (j=0; j<nslots; j++) {
1781b9fc7fcSBarry Smith             for (k=0; k<nfields; k++) {
1791b9fc7fcSBarry Smith               ii[nfields*j + k] = rstart + bs*j + fields[k];
18097bbdb24SBarry Smith             }
18197bbdb24SBarry Smith           }
182704ba839SBarry Smith           ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr);
183ccb205f8SBarry Smith           ierr = PetscFree(ii);CHKERRQ(ierr);
184ccb205f8SBarry Smith         } else {
185704ba839SBarry Smith           ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr);
186ccb205f8SBarry Smith         }
1873e197d65SBarry Smith       }
1883e197d65SBarry Smith       ierr = ISGetLocalSize(ilink->is,&ilink->csize);CHKERRQ(ierr);
189704ba839SBarry Smith       ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr);
1901b9fc7fcSBarry Smith       if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split");
191704ba839SBarry Smith       ierr = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr);
192704ba839SBarry Smith       ilink = ilink->next;
1931b9fc7fcSBarry Smith     }
1941b9fc7fcSBarry Smith   }
1951b9fc7fcSBarry Smith 
196704ba839SBarry Smith   ilink  = jac->head;
19797bbdb24SBarry Smith   if (!jac->pmat) {
198cf502942SBarry Smith     ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr);
199cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
200704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
201704ba839SBarry Smith       ilink = ilink->next;
202cf502942SBarry Smith     }
20397bbdb24SBarry Smith   } else {
204cf502942SBarry Smith     for (i=0; i<nsplit; i++) {
205704ba839SBarry Smith       ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr);
206704ba839SBarry Smith       ilink = ilink->next;
207cf502942SBarry Smith     }
20897bbdb24SBarry Smith   }
20997bbdb24SBarry Smith 
21068dd23aaSBarry Smith   /* extract the rows of the matrix associated with each field: used for efficient computation of residual inside algorithm */
21168dd23aaSBarry Smith   ierr = MatHasOperation(pc->mat,MATOP_GET_SUBMATRIX,&getsub);CHKERRQ(ierr);
21268dd23aaSBarry Smith   if (getsub && jac->type != PC_COMPOSITE_ADDITIVE) {
21368dd23aaSBarry Smith     ilink  = jac->head;
21468dd23aaSBarry Smith     if (!jac->Afield) {
21568dd23aaSBarry Smith       ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->Afield);CHKERRQ(ierr);
21668dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
217*11755939SBarry Smith 	ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_INITIAL_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
21868dd23aaSBarry Smith 	ilink = ilink->next;
21968dd23aaSBarry Smith       }
22068dd23aaSBarry Smith     } else {
22168dd23aaSBarry Smith       for (i=0; i<nsplit; i++) {
222*11755939SBarry Smith 	ierr = MatGetSubMatrix(pc->mat,ilink->is,PETSC_NULL,PETSC_DECIDE,MAT_REUSE_MATRIX,&jac->Afield[i]);CHKERRQ(ierr);
22368dd23aaSBarry Smith 	ilink = ilink->next;
22468dd23aaSBarry Smith       }
22568dd23aaSBarry Smith     }
22668dd23aaSBarry Smith   }
22768dd23aaSBarry Smith 
22868dd23aaSBarry Smith 
22997bbdb24SBarry Smith   /* set up the individual PCs */
23097bbdb24SBarry Smith   i    = 0;
2315a9f2f41SSatish Balay   ilink = jac->head;
2325a9f2f41SSatish Balay   while (ilink) {
2335a9f2f41SSatish Balay     ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr);
2345a9f2f41SSatish Balay     ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr);
2355a9f2f41SSatish Balay     ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr);
23697bbdb24SBarry Smith     i++;
2375a9f2f41SSatish Balay     ilink = ilink->next;
2380971522cSBarry Smith   }
23997bbdb24SBarry Smith 
24097bbdb24SBarry Smith   /* create work vectors for each split */
2411b9fc7fcSBarry Smith   if (!jac->x) {
24279416396SBarry Smith     Vec xtmp;
24397bbdb24SBarry Smith     ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr);
2445a9f2f41SSatish Balay     ilink = jac->head;
24597bbdb24SBarry Smith     for (i=0; i<nsplit; i++) {
246906ed7ccSBarry Smith       Vec *vl,*vr;
247906ed7ccSBarry Smith 
248906ed7ccSBarry Smith       ierr      = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr);
249906ed7ccSBarry Smith       ilink->x  = *vr;
250906ed7ccSBarry Smith       ilink->y  = *vl;
251906ed7ccSBarry Smith       ierr      = PetscFree(vr);CHKERRQ(ierr);
252906ed7ccSBarry Smith       ierr      = PetscFree(vl);CHKERRQ(ierr);
2535a9f2f41SSatish Balay       jac->x[i] = ilink->x;
2545a9f2f41SSatish Balay       jac->y[i] = ilink->y;
2555a9f2f41SSatish Balay       ilink     = ilink->next;
25697bbdb24SBarry Smith     }
25779416396SBarry Smith     /* compute scatter contexts needed by multiplicative versions and non-default splits */
2581b9fc7fcSBarry Smith 
2595a9f2f41SSatish Balay     ilink = jac->head;
2601b9fc7fcSBarry Smith     ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr);
2611b9fc7fcSBarry Smith     for (i=0; i<nsplit; i++) {
262704ba839SBarry Smith       ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr);
2635a9f2f41SSatish Balay       ilink = ilink->next;
2641b9fc7fcSBarry Smith     }
2651b9fc7fcSBarry Smith     ierr = VecDestroy(xtmp);CHKERRQ(ierr);
2661b9fc7fcSBarry Smith   }
2670971522cSBarry Smith   PetscFunctionReturn(0);
2680971522cSBarry Smith }
2690971522cSBarry Smith 
2705a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \
271ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
272ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \
2735a9f2f41SSatish Balay      KSPSolve(ilink->ksp,ilink->x,ilink->y) || \
274ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \
275ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE))
27679416396SBarry Smith 
2770971522cSBarry Smith #undef __FUNCT__
2780971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
2790971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y)
2800971522cSBarry Smith {
2810971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
2820971522cSBarry Smith   PetscErrorCode    ierr;
2835a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
2843e197d65SBarry Smith   PetscInt          bs,cnt;
2850971522cSBarry Smith 
2860971522cSBarry Smith   PetscFunctionBegin;
28751f519a2SBarry Smith   CHKMEMQ;
288e25d487eSBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
28951f519a2SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
29051f519a2SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
29151f519a2SBarry Smith 
29279416396SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
2931b9fc7fcSBarry Smith     if (jac->defaultsplit) {
2940971522cSBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
2955a9f2f41SSatish Balay       while (ilink) {
2965a9f2f41SSatish Balay         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
2975a9f2f41SSatish Balay         ilink = ilink->next;
2980971522cSBarry Smith       }
2990971522cSBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
3001b9fc7fcSBarry Smith     } else {
301efb30889SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
3025a9f2f41SSatish Balay       while (ilink) {
3035a9f2f41SSatish Balay         ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
3045a9f2f41SSatish Balay         ilink = ilink->next;
3051b9fc7fcSBarry Smith       }
3061b9fc7fcSBarry Smith     }
30716913363SBarry Smith   } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
30879416396SBarry Smith     if (!jac->w1) {
30979416396SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
31079416396SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
31179416396SBarry Smith     }
312efb30889SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
3135a9f2f41SSatish Balay     ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr);
3143e197d65SBarry Smith     cnt = 1;
3155a9f2f41SSatish Balay     while (ilink->next) {
3165a9f2f41SSatish Balay       ilink = ilink->next;
3173e197d65SBarry Smith       if (jac->Afield) {
3183e197d65SBarry Smith         /* compute the residual only over the part of the vector needed */
3193e197d65SBarry Smith         ierr = MatMult(jac->Afield[cnt++],y,ilink->x);CHKERRQ(ierr);
3203e197d65SBarry Smith         ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
3213e197d65SBarry Smith         ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3223e197d65SBarry Smith         ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
3233e197d65SBarry Smith         ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
3243e197d65SBarry Smith         ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
3253e197d65SBarry Smith         ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
3263e197d65SBarry Smith       } else {
3273e197d65SBarry Smith         /* compute the residual over the entire vector */
3281ab39975SBarry Smith 	ierr = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
329efb30889SBarry Smith 	ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
3305a9f2f41SSatish Balay 	ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
33179416396SBarry Smith       }
3323e197d65SBarry Smith     }
33351f519a2SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
334*11755939SBarry Smith       cnt -= 2;
33551f519a2SBarry Smith       while (ilink->previous) {
33651f519a2SBarry Smith         ilink = ilink->previous;
337*11755939SBarry Smith         if (jac->Afield) {
338*11755939SBarry Smith 	  /* compute the residual only over the part of the vector needed */
339*11755939SBarry Smith 	  ierr = MatMult(jac->Afield[cnt--],y,ilink->x);CHKERRQ(ierr);
340*11755939SBarry Smith 	  ierr = VecScale(ilink->x,-1.0);CHKERRQ(ierr);
341*11755939SBarry Smith 	  ierr = VecScatterBegin(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
342*11755939SBarry Smith 	  ierr = VecScatterEnd(ilink->sctx,x,ilink->x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
343*11755939SBarry Smith 	  ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
344*11755939SBarry Smith 	  ierr = VecScatterBegin(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
345*11755939SBarry Smith 	  ierr = VecScatterEnd(ilink->sctx,ilink->y,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
346*11755939SBarry Smith         } else {
3471ab39975SBarry Smith 	  ierr  = MatMult(pc->mat,y,jac->w1);CHKERRQ(ierr);
34851f519a2SBarry Smith 	  ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
34951f519a2SBarry Smith 	  ierr  = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr);
35079416396SBarry Smith         }
35151f519a2SBarry Smith       }
352*11755939SBarry Smith     }
35316913363SBarry Smith   } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type);
35451f519a2SBarry Smith   CHKMEMQ;
3550971522cSBarry Smith   PetscFunctionReturn(0);
3560971522cSBarry Smith }
3570971522cSBarry Smith 
358421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \
359ca9f406cSSatish Balay     (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
360ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \
361421e10b8SBarry Smith      KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \
362ca9f406cSSatish Balay      VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \
363ca9f406cSSatish Balay      VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE))
364421e10b8SBarry Smith 
365421e10b8SBarry Smith #undef __FUNCT__
366421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit"
367421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y)
368421e10b8SBarry Smith {
369421e10b8SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
370421e10b8SBarry Smith   PetscErrorCode    ierr;
371421e10b8SBarry Smith   PC_FieldSplitLink ilink = jac->head;
372421e10b8SBarry Smith   PetscInt          bs;
373421e10b8SBarry Smith 
374421e10b8SBarry Smith   PetscFunctionBegin;
375421e10b8SBarry Smith   CHKMEMQ;
376421e10b8SBarry Smith   ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr);
377421e10b8SBarry Smith   ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr);
378421e10b8SBarry Smith   ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr);
379421e10b8SBarry Smith 
380421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_ADDITIVE) {
381421e10b8SBarry Smith     if (jac->defaultsplit) {
382421e10b8SBarry Smith       ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr);
383421e10b8SBarry Smith       while (ilink) {
384421e10b8SBarry Smith 	ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr);
385421e10b8SBarry Smith 	ilink = ilink->next;
386421e10b8SBarry Smith       }
387421e10b8SBarry Smith       ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr);
388421e10b8SBarry Smith     } else {
389421e10b8SBarry Smith       ierr = VecSet(y,0.0);CHKERRQ(ierr);
390421e10b8SBarry Smith       while (ilink) {
391421e10b8SBarry Smith         ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
392421e10b8SBarry Smith 	ilink = ilink->next;
393421e10b8SBarry Smith       }
394421e10b8SBarry Smith     }
395421e10b8SBarry Smith   } else {
396421e10b8SBarry Smith     if (!jac->w1) {
397421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr);
398421e10b8SBarry Smith       ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr);
399421e10b8SBarry Smith     }
400421e10b8SBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
401421e10b8SBarry Smith     if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
402421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
403421e10b8SBarry Smith       while (ilink->next) {
404421e10b8SBarry Smith         ilink = ilink->next;
4059989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
406421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
407421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
408421e10b8SBarry Smith       }
409421e10b8SBarry Smith       while (ilink->previous) {
410421e10b8SBarry Smith         ilink = ilink->previous;
4119989ab13SBarry Smith         ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
412421e10b8SBarry Smith         ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
413421e10b8SBarry Smith         ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
414421e10b8SBarry Smith       }
415421e10b8SBarry Smith     } else {
416421e10b8SBarry Smith       while (ilink->next) {   /* get to last entry in linked list */
417421e10b8SBarry Smith 	ilink = ilink->next;
418421e10b8SBarry Smith       }
419421e10b8SBarry Smith       ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr);
420421e10b8SBarry Smith       while (ilink->previous) {
421421e10b8SBarry Smith 	ilink = ilink->previous;
4229989ab13SBarry Smith 	ierr  = MatMultTranspose(pc->mat,y,jac->w1);CHKERRQ(ierr);
423421e10b8SBarry Smith 	ierr  = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr);
424421e10b8SBarry Smith 	ierr  = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr);
425421e10b8SBarry Smith       }
426421e10b8SBarry Smith     }
427421e10b8SBarry Smith   }
428421e10b8SBarry Smith   CHKMEMQ;
429421e10b8SBarry Smith   PetscFunctionReturn(0);
430421e10b8SBarry Smith }
431421e10b8SBarry Smith 
4320971522cSBarry Smith #undef __FUNCT__
4330971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit"
4340971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc)
4350971522cSBarry Smith {
4360971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
4370971522cSBarry Smith   PetscErrorCode    ierr;
4385a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head,next;
4390971522cSBarry Smith 
4400971522cSBarry Smith   PetscFunctionBegin;
4415a9f2f41SSatish Balay   while (ilink) {
4425a9f2f41SSatish Balay     ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr);
4435a9f2f41SSatish Balay     if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);}
4445a9f2f41SSatish Balay     if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);}
4455a9f2f41SSatish Balay     if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);}
446704ba839SBarry Smith     if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);}
447704ba839SBarry Smith     if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);}
4485a9f2f41SSatish Balay     next = ilink->next;
449704ba839SBarry Smith     ierr = PetscFree(ilink->fields);CHKERRQ(ierr);
450704ba839SBarry Smith     ierr = PetscFree(ilink);CHKERRQ(ierr);
4515a9f2f41SSatish Balay     ilink = next;
4520971522cSBarry Smith   }
45305b42c5fSBarry Smith   ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);
45497bbdb24SBarry Smith   if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);}
45568dd23aaSBarry Smith   if (jac->Afield) {ierr = MatDestroyMatrices(jac->nsplits,&jac->Afield);CHKERRQ(ierr);}
45679416396SBarry Smith   if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);}
45779416396SBarry Smith   if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);}
4580971522cSBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
4590971522cSBarry Smith   PetscFunctionReturn(0);
4600971522cSBarry Smith }
4610971522cSBarry Smith 
4620971522cSBarry Smith #undef __FUNCT__
4630971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit"
4640971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc)
4650971522cSBarry Smith {
4661b9fc7fcSBarry Smith   PetscErrorCode ierr;
46751f519a2SBarry Smith   PetscInt       i = 0,nfields,*fields,bs;
4681b9fc7fcSBarry Smith   PetscTruth     flg;
4691b9fc7fcSBarry Smith   char           optionname[128];
4709dcbbd2bSBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
4711b9fc7fcSBarry Smith 
4720971522cSBarry Smith   PetscFunctionBegin;
4731b9fc7fcSBarry Smith   ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr);
47451f519a2SBarry Smith   ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr);
47551f519a2SBarry Smith   if (flg) {
47651f519a2SBarry Smith     ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr);
47751f519a2SBarry Smith   }
478704ba839SBarry Smith 
4799dcbbd2bSBarry Smith   ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
480d32f9abdSBarry Smith 
481d32f9abdSBarry Smith   if (jac->bs > 0) {
482d32f9abdSBarry Smith     /* only allow user to set fields from command line if bs is already known.
483d32f9abdSBarry Smith        otherwise user can set them in PCFieldSplitSetDefaults() */
48451f519a2SBarry Smith     ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr);
4851b9fc7fcSBarry Smith     while (PETSC_TRUE) {
48613f74950SBarry Smith       sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++);
48751f519a2SBarry Smith       nfields = jac->bs;
4881b9fc7fcSBarry Smith       ierr    = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr);
4891b9fc7fcSBarry Smith       if (!flg) break;
4901b9fc7fcSBarry Smith       if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields");
4911b9fc7fcSBarry Smith       ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr);
4921b9fc7fcSBarry Smith     }
49351f519a2SBarry Smith     ierr = PetscFree(fields);CHKERRQ(ierr);
494d32f9abdSBarry Smith   }
4951b9fc7fcSBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
4960971522cSBarry Smith   PetscFunctionReturn(0);
4970971522cSBarry Smith }
4980971522cSBarry Smith 
4990971522cSBarry Smith /*------------------------------------------------------------------------------------*/
5000971522cSBarry Smith 
5010971522cSBarry Smith EXTERN_C_BEGIN
5020971522cSBarry Smith #undef __FUNCT__
5030971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit"
504dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields)
5050971522cSBarry Smith {
50697bbdb24SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5070971522cSBarry Smith   PetscErrorCode    ierr;
5085a9f2f41SSatish Balay   PC_FieldSplitLink ilink,next = jac->head;
50969a612a9SBarry Smith   char              prefix[128];
51051f519a2SBarry Smith   PetscInt          i;
5110971522cSBarry Smith 
5120971522cSBarry Smith   PetscFunctionBegin;
5130971522cSBarry Smith   if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested");
51451f519a2SBarry Smith   for (i=0; i<n; i++) {
51551f519a2SBarry Smith     if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs);
51651f519a2SBarry Smith     if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]);
51751f519a2SBarry Smith   }
518704ba839SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
519704ba839SBarry Smith   ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr);
5205a9f2f41SSatish Balay   ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr);
5215a9f2f41SSatish Balay   ilink->nfields = n;
5225a9f2f41SSatish Balay   ilink->next    = PETSC_NULL;
5237adad957SLisandro Dalcin   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
5245a9f2f41SSatish Balay   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
52569a612a9SBarry Smith 
5267adad957SLisandro Dalcin   if (((PetscObject)pc)->prefix) {
5277adad957SLisandro Dalcin     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
52869a612a9SBarry Smith   } else {
52913f74950SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
53069a612a9SBarry Smith   }
5315a9f2f41SSatish Balay   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
5320971522cSBarry Smith 
5330971522cSBarry Smith   if (!next) {
5345a9f2f41SSatish Balay     jac->head       = ilink;
53551f519a2SBarry Smith     ilink->previous = PETSC_NULL;
5360971522cSBarry Smith   } else {
5370971522cSBarry Smith     while (next->next) {
5380971522cSBarry Smith       next = next->next;
5390971522cSBarry Smith     }
5405a9f2f41SSatish Balay     next->next      = ilink;
54151f519a2SBarry Smith     ilink->previous = next;
5420971522cSBarry Smith   }
5430971522cSBarry Smith   jac->nsplits++;
5440971522cSBarry Smith   PetscFunctionReturn(0);
5450971522cSBarry Smith }
5460971522cSBarry Smith EXTERN_C_END
5470971522cSBarry Smith 
5480971522cSBarry Smith 
5490971522cSBarry Smith EXTERN_C_BEGIN
5500971522cSBarry Smith #undef __FUNCT__
55169a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit"
552dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp)
5530971522cSBarry Smith {
5540971522cSBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
5550971522cSBarry Smith   PetscErrorCode    ierr;
5560971522cSBarry Smith   PetscInt          cnt = 0;
5575a9f2f41SSatish Balay   PC_FieldSplitLink ilink = jac->head;
5580971522cSBarry Smith 
5590971522cSBarry Smith   PetscFunctionBegin;
56069a612a9SBarry Smith   ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr);
5615a9f2f41SSatish Balay   while (ilink) {
5625a9f2f41SSatish Balay     (*subksp)[cnt++] = ilink->ksp;
5635a9f2f41SSatish Balay     ilink = ilink->next;
5640971522cSBarry Smith   }
5650971522cSBarry Smith   if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits);
5660971522cSBarry Smith   *n = jac->nsplits;
5670971522cSBarry Smith   PetscFunctionReturn(0);
5680971522cSBarry Smith }
5690971522cSBarry Smith EXTERN_C_END
5700971522cSBarry Smith 
571704ba839SBarry Smith EXTERN_C_BEGIN
572704ba839SBarry Smith #undef __FUNCT__
573704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit"
574704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is)
575704ba839SBarry Smith {
576704ba839SBarry Smith   PC_FieldSplit     *jac = (PC_FieldSplit*)pc->data;
577704ba839SBarry Smith   PetscErrorCode    ierr;
578704ba839SBarry Smith   PC_FieldSplitLink ilink, next = jac->head;
579704ba839SBarry Smith   char              prefix[128];
580704ba839SBarry Smith 
581704ba839SBarry Smith   PetscFunctionBegin;
58216913363SBarry Smith   ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr);
5831ab39975SBarry Smith   ilink->is      = is;
5841ab39975SBarry Smith   ierr           = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
585704ba839SBarry Smith   ilink->next    = PETSC_NULL;
586704ba839SBarry Smith   ierr           = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr);
587704ba839SBarry Smith   ierr           = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr);
588704ba839SBarry Smith 
589704ba839SBarry Smith   if (((PetscObject)pc)->prefix) {
590704ba839SBarry Smith     sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits);
591704ba839SBarry Smith   } else {
592704ba839SBarry Smith     sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits);
593704ba839SBarry Smith   }
594704ba839SBarry Smith   ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr);
595704ba839SBarry Smith 
596704ba839SBarry Smith   if (!next) {
597704ba839SBarry Smith     jac->head       = ilink;
598704ba839SBarry Smith     ilink->previous = PETSC_NULL;
599704ba839SBarry Smith   } else {
600704ba839SBarry Smith     while (next->next) {
601704ba839SBarry Smith       next = next->next;
602704ba839SBarry Smith     }
603704ba839SBarry Smith     next->next      = ilink;
604704ba839SBarry Smith     ilink->previous = next;
605704ba839SBarry Smith   }
606704ba839SBarry Smith   jac->nsplits++;
607704ba839SBarry Smith 
608704ba839SBarry Smith   PetscFunctionReturn(0);
609704ba839SBarry Smith }
610704ba839SBarry Smith EXTERN_C_END
611704ba839SBarry Smith 
6120971522cSBarry Smith #undef __FUNCT__
6130971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields"
6140971522cSBarry Smith /*@
6150971522cSBarry Smith     PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner
6160971522cSBarry Smith 
6170971522cSBarry Smith     Collective on PC
6180971522cSBarry Smith 
6190971522cSBarry Smith     Input Parameters:
6200971522cSBarry Smith +   pc  - the preconditioner context
6210971522cSBarry Smith .   n - the number of fields in this split
6220971522cSBarry Smith .   fields - the fields in this split
6230971522cSBarry Smith 
6240971522cSBarry Smith     Level: intermediate
6250971522cSBarry Smith 
626d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetIS() to set a completely general set of indices as a field.
627d32f9abdSBarry Smith 
628d32f9abdSBarry Smith      The PCFieldSplitSetFields() is for defining fields as a strided blocks. For example, if the block
629d32f9abdSBarry Smith      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
630d32f9abdSBarry Smith      0xx3xx6xx9xx12 ... x1xx4xx7xx ... xx2xx5xx8xx.. 01x34x67x... 0x1x3x5x7.. x12x45x78x....
631d32f9abdSBarry Smith      where the numbered entries indicate what is in the field.
632d32f9abdSBarry Smith 
633d32f9abdSBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize(), PCFieldSplitSetIS()
6340971522cSBarry Smith 
6350971522cSBarry Smith @*/
636dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields)
6370971522cSBarry Smith {
6380971522cSBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *);
6390971522cSBarry Smith 
6400971522cSBarry Smith   PetscFunctionBegin;
6410971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6420971522cSBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr);
6430971522cSBarry Smith   if (f) {
6440971522cSBarry Smith     ierr = (*f)(pc,n,fields);CHKERRQ(ierr);
6450971522cSBarry Smith   }
6460971522cSBarry Smith   PetscFunctionReturn(0);
6470971522cSBarry Smith }
6480971522cSBarry Smith 
6490971522cSBarry Smith #undef __FUNCT__
650704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS"
651704ba839SBarry Smith /*@
652704ba839SBarry Smith     PCFieldSplitSetIS - Sets the exact elements for field
653704ba839SBarry Smith 
654704ba839SBarry Smith     Collective on PC
655704ba839SBarry Smith 
656704ba839SBarry Smith     Input Parameters:
657704ba839SBarry Smith +   pc  - the preconditioner context
658704ba839SBarry Smith .   is - the index set that defines the vector elements in this field
659704ba839SBarry Smith 
660d32f9abdSBarry Smith 
661d32f9abdSBarry Smith     Notes: Use PCFieldSplitSetFields(), for fields defined by strided types.
662d32f9abdSBarry Smith 
663704ba839SBarry Smith     Level: intermediate
664704ba839SBarry Smith 
665704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize()
666704ba839SBarry Smith 
667704ba839SBarry Smith @*/
668704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is)
669704ba839SBarry Smith {
670704ba839SBarry Smith   PetscErrorCode ierr,(*f)(PC,IS);
671704ba839SBarry Smith 
672704ba839SBarry Smith   PetscFunctionBegin;
673704ba839SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
674704ba839SBarry Smith   PetscValidHeaderSpecific(is,IS_COOKIE,1);
675704ba839SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr);
676704ba839SBarry Smith   if (f) {
677704ba839SBarry Smith     ierr = (*f)(pc,is);CHKERRQ(ierr);
678704ba839SBarry Smith   }
679704ba839SBarry Smith   PetscFunctionReturn(0);
680704ba839SBarry Smith }
681704ba839SBarry Smith 
682704ba839SBarry Smith #undef __FUNCT__
68351f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize"
68451f519a2SBarry Smith /*@
68551f519a2SBarry Smith     PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the
68651f519a2SBarry Smith       fieldsplit preconditioner. If not set the matrix block size is used.
68751f519a2SBarry Smith 
68851f519a2SBarry Smith     Collective on PC
68951f519a2SBarry Smith 
69051f519a2SBarry Smith     Input Parameters:
69151f519a2SBarry Smith +   pc  - the preconditioner context
69251f519a2SBarry Smith -   bs - the block size
69351f519a2SBarry Smith 
69451f519a2SBarry Smith     Level: intermediate
69551f519a2SBarry Smith 
69651f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields()
69751f519a2SBarry Smith 
69851f519a2SBarry Smith @*/
69951f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs)
70051f519a2SBarry Smith {
70151f519a2SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt);
70251f519a2SBarry Smith 
70351f519a2SBarry Smith   PetscFunctionBegin;
70451f519a2SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
70551f519a2SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
70651f519a2SBarry Smith   if (f) {
70751f519a2SBarry Smith     ierr = (*f)(pc,bs);CHKERRQ(ierr);
70851f519a2SBarry Smith   }
70951f519a2SBarry Smith   PetscFunctionReturn(0);
71051f519a2SBarry Smith }
71151f519a2SBarry Smith 
71251f519a2SBarry Smith #undef __FUNCT__
71369a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP"
7140971522cSBarry Smith /*@C
71569a612a9SBarry Smith    PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits
7160971522cSBarry Smith 
71769a612a9SBarry Smith    Collective on KSP
7180971522cSBarry Smith 
7190971522cSBarry Smith    Input Parameter:
7200971522cSBarry Smith .  pc - the preconditioner context
7210971522cSBarry Smith 
7220971522cSBarry Smith    Output Parameters:
7230971522cSBarry Smith +  n - the number of split
72469a612a9SBarry Smith -  pc - the array of KSP contexts
7250971522cSBarry Smith 
7260971522cSBarry Smith    Note:
727d32f9abdSBarry Smith    After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed by the user
728d32f9abdSBarry Smith    (not the KSP just the array that contains them).
7290971522cSBarry Smith 
73069a612a9SBarry Smith    You must call KSPSetUp() before calling PCFieldSplitGetSubKSP().
7310971522cSBarry Smith 
7320971522cSBarry Smith    Level: advanced
7330971522cSBarry Smith 
7340971522cSBarry Smith .seealso: PCFIELDSPLIT
7350971522cSBarry Smith @*/
736dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[])
7370971522cSBarry Smith {
73869a612a9SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **);
7390971522cSBarry Smith 
7400971522cSBarry Smith   PetscFunctionBegin;
7410971522cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7420971522cSBarry Smith   PetscValidIntPointer(n,2);
74369a612a9SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr);
7440971522cSBarry Smith   if (f) {
74569a612a9SBarry Smith     ierr = (*f)(pc,n,subksp);CHKERRQ(ierr);
7460971522cSBarry Smith   } else {
74769a612a9SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC");
7480971522cSBarry Smith   }
7490971522cSBarry Smith   PetscFunctionReturn(0);
7500971522cSBarry Smith }
7510971522cSBarry Smith 
75279416396SBarry Smith EXTERN_C_BEGIN
75379416396SBarry Smith #undef __FUNCT__
75479416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit"
755dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type)
75679416396SBarry Smith {
75779416396SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
75879416396SBarry Smith 
75979416396SBarry Smith   PetscFunctionBegin;
76079416396SBarry Smith   jac->type = type;
76179416396SBarry Smith   PetscFunctionReturn(0);
76279416396SBarry Smith }
76379416396SBarry Smith EXTERN_C_END
76479416396SBarry Smith 
76551f519a2SBarry Smith EXTERN_C_BEGIN
76651f519a2SBarry Smith #undef __FUNCT__
76751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit"
76851f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs)
76951f519a2SBarry Smith {
77051f519a2SBarry Smith   PC_FieldSplit  *jac = (PC_FieldSplit*)pc->data;
77151f519a2SBarry Smith 
77251f519a2SBarry Smith   PetscFunctionBegin;
77351f519a2SBarry Smith   if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs);
77451f519a2SBarry Smith   if (jac->bs > 0 && jac->bs != bs) SETERRQ2(PETSC_ERR_ARG_WRONGSTATE,"Cannot change fieldsplit blocksize from %D to %D after it has been set",jac->bs,bs);
77551f519a2SBarry Smith   jac->bs = bs;
77651f519a2SBarry Smith   PetscFunctionReturn(0);
77751f519a2SBarry Smith }
77851f519a2SBarry Smith EXTERN_C_END
77951f519a2SBarry Smith 
78079416396SBarry Smith #undef __FUNCT__
78179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType"
782bc08b0f1SBarry Smith /*@
78379416396SBarry Smith    PCFieldSplitSetType - Sets the type of fieldsplit preconditioner.
78479416396SBarry Smith 
78579416396SBarry Smith    Collective on PC
78679416396SBarry Smith 
78779416396SBarry Smith    Input Parameter:
78879416396SBarry Smith .  pc - the preconditioner context
789ce9499c7SBarry Smith .  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE
79079416396SBarry Smith 
79179416396SBarry Smith    Options Database Key:
792ce9499c7SBarry Smith .  -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type
79379416396SBarry Smith 
79479416396SBarry Smith    Level: Developer
79579416396SBarry Smith 
79679416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
79779416396SBarry Smith 
79879416396SBarry Smith .seealso: PCCompositeSetType()
79979416396SBarry Smith 
80079416396SBarry Smith @*/
801dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type)
80279416396SBarry Smith {
80379416396SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
80479416396SBarry Smith 
80579416396SBarry Smith   PetscFunctionBegin;
80679416396SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
80779416396SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
80879416396SBarry Smith   if (f) {
80979416396SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
81079416396SBarry Smith   }
81179416396SBarry Smith   PetscFunctionReturn(0);
81279416396SBarry Smith }
81379416396SBarry Smith 
8140971522cSBarry Smith /* -------------------------------------------------------------------------------------*/
8150971522cSBarry Smith /*MC
816a8c7a070SBarry Smith    PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual
8170971522cSBarry Smith                   fields or groups of fields
8180971522cSBarry Smith 
8190971522cSBarry Smith 
8200971522cSBarry Smith      To set options on the solvers for each block append -sub_ to all the PC
821d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1
8220971522cSBarry Smith 
823a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP()
82469a612a9SBarry Smith          and set the options directly on the resulting KSP object
8250971522cSBarry Smith 
8260971522cSBarry Smith    Level: intermediate
8270971522cSBarry Smith 
82879416396SBarry Smith    Options Database Keys:
8292e70eadcSBarry Smith +   -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split
8302e70eadcSBarry Smith .   -pc_splitfield_default - automatically add any fields to additional splits that have not
8312e70eadcSBarry Smith                               been supplied explicitly by -pc_splitfield_%d_fields
83251f519a2SBarry Smith .   -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields)
8332e70eadcSBarry Smith -   -pc_splitfield_type <additive,multiplicative>
83479416396SBarry Smith 
835d32f9abdSBarry Smith    Notes: use PCFieldSplitSetFields() to set fields defined by "strided" entries and PCFieldSplitSetIS()
836d32f9abdSBarry Smith      to define a field by an arbitrary collection of entries.
837d32f9abdSBarry Smith 
838d32f9abdSBarry Smith       If no fields are set the default is used. The fields are defined by entries strided by bs,
839d32f9abdSBarry Smith       beginning at 0 then 1, etc to bs-1. The block size can be set with PCFieldSplitSetBlockSize(),
840d32f9abdSBarry Smith       if this is not called the block size defaults to the blocksize of the second matrix passed
841d32f9abdSBarry Smith       to KSPSetOperators()/PCSetOperators().
842d32f9abdSBarry Smith 
843d32f9abdSBarry Smith       Currently for the multiplicative version, the updated residual needed for the next field
844d32f9abdSBarry Smith      solve is computed via a matrix vector product over the entire array. An optimization would be
845d32f9abdSBarry Smith      to update the residual only for the part of the right hand side associated with the next field
846d32f9abdSBarry Smith      solve. (This would involve more MatGetSubMatrix() calls or some other mechanism to compute the
847d32f9abdSBarry Smith      part of the matrix needed to just update part of the residual).
848d32f9abdSBarry Smith 
8490971522cSBarry Smith    Concepts: physics based preconditioners
8500971522cSBarry Smith 
8510971522cSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
852d32f9abdSBarry Smith            PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(), PCFieldSplitSetType(), PCFieldSplitSetIS()
8530971522cSBarry Smith M*/
8540971522cSBarry Smith 
8550971522cSBarry Smith EXTERN_C_BEGIN
8560971522cSBarry Smith #undef __FUNCT__
8570971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit"
858dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc)
8590971522cSBarry Smith {
8600971522cSBarry Smith   PetscErrorCode ierr;
8610971522cSBarry Smith   PC_FieldSplit  *jac;
8620971522cSBarry Smith 
8630971522cSBarry Smith   PetscFunctionBegin;
86438f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr);
8650971522cSBarry Smith   jac->bs        = -1;
8660971522cSBarry Smith   jac->nsplits   = 0;
8673e197d65SBarry Smith   jac->type      = PC_COMPOSITE_MULTIPLICATIVE;
86851f519a2SBarry Smith 
8690971522cSBarry Smith   pc->data     = (void*)jac;
8700971522cSBarry Smith 
8710971522cSBarry Smith   pc->ops->apply             = PCApply_FieldSplit;
872421e10b8SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_FieldSplit;
8730971522cSBarry Smith   pc->ops->setup             = PCSetUp_FieldSplit;
8740971522cSBarry Smith   pc->ops->destroy           = PCDestroy_FieldSplit;
8750971522cSBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_FieldSplit;
8760971522cSBarry Smith   pc->ops->view              = PCView_FieldSplit;
8770971522cSBarry Smith   pc->ops->applyrichardson   = 0;
8780971522cSBarry Smith 
87969a612a9SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit",
88069a612a9SBarry Smith                     PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr);
8810971522cSBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit",
8820971522cSBarry Smith                     PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr);
883704ba839SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit",
884704ba839SBarry Smith                     PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr);
88579416396SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit",
88679416396SBarry Smith                     PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr);
88751f519a2SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit",
88851f519a2SBarry Smith                     PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr);
8890971522cSBarry Smith   PetscFunctionReturn(0);
8900971522cSBarry Smith }
8910971522cSBarry Smith EXTERN_C_END
8920971522cSBarry Smith 
8930971522cSBarry Smith 
894