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 { 2179416396SBarry Smith PCCompositeType type; /* additive or multiplicative */ 2297bbdb24SBarry Smith PetscTruth defaultsplit; 230971522cSBarry Smith PetscInt bs; 24704ba839SBarry Smith PetscInt nsplits; 2579416396SBarry Smith Vec *x,*y,w1,w2; 2697bbdb24SBarry Smith Mat *pmat; 27704ba839SBarry Smith PetscTruth issetup; 2897bbdb24SBarry Smith PC_FieldSplitLink head; 290971522cSBarry Smith } PC_FieldSplit; 300971522cSBarry Smith 31*16913363SBarry Smith /* 32*16913363SBarry Smith Notes: there is no particular reason that pmat, x, and y are stored as arrays in PC_FieldSplit instead of 33*16913363SBarry Smith inside PC_FieldSplitLink, just historical. If you want to be able to add new fields after already using the 34*16913363SBarry Smith PC you could change this. 35*16913363SBarry Smith */ 360971522cSBarry Smith #undef __FUNCT__ 370971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 380971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 390971522cSBarry Smith { 400971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 410971522cSBarry Smith PetscErrorCode ierr; 420971522cSBarry Smith PetscTruth iascii; 430971522cSBarry Smith PetscInt i,j; 445a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 450971522cSBarry Smith 460971522cSBarry Smith PetscFunctionBegin; 470971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 480971522cSBarry Smith if (iascii) { 4951f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 5069a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 510971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 520971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 530971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 5479416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 555a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 5679416396SBarry Smith if (j > 0) { 5779416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 5879416396SBarry Smith } 595a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 600971522cSBarry Smith } 610971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 6279416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 635a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 645a9f2f41SSatish Balay ilink = ilink->next; 650971522cSBarry Smith } 660971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 670971522cSBarry Smith } else { 680971522cSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 690971522cSBarry Smith } 700971522cSBarry Smith PetscFunctionReturn(0); 710971522cSBarry Smith } 720971522cSBarry Smith 730971522cSBarry Smith #undef __FUNCT__ 7469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 7569a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 760971522cSBarry Smith { 770971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 780971522cSBarry Smith PetscErrorCode ierr; 795a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 8069a612a9SBarry Smith PetscInt i; 8179416396SBarry Smith PetscTruth flg = PETSC_FALSE,*fields; 820971522cSBarry Smith 830971522cSBarry Smith PetscFunctionBegin; 847adad957SLisandro Dalcin ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 855a9f2f41SSatish Balay if (!ilink || flg) { 86ae15b995SBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 87521d7252SBarry Smith if (jac->bs <= 0) { 88704ba839SBarry Smith if (pc->pmat) { 89521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 90704ba839SBarry Smith } else { 91704ba839SBarry Smith jac->bs = 1; 92704ba839SBarry Smith } 93521d7252SBarry Smith } 9479416396SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 9579416396SBarry Smith ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 965a9f2f41SSatish Balay while (ilink) { 975a9f2f41SSatish Balay for (i=0; i<ilink->nfields; i++) { 985a9f2f41SSatish Balay fields[ilink->fields[i]] = PETSC_TRUE; 99521d7252SBarry Smith } 1005a9f2f41SSatish Balay ilink = ilink->next; 10179416396SBarry Smith } 10297bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 10379416396SBarry Smith for (i=0; i<jac->bs; i++) { 10479416396SBarry Smith if (!fields[i]) { 10579416396SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 10679416396SBarry Smith } else { 10779416396SBarry Smith jac->defaultsplit = PETSC_FALSE; 10879416396SBarry Smith } 10979416396SBarry Smith } 110521d7252SBarry Smith } 11169a612a9SBarry Smith PetscFunctionReturn(0); 11269a612a9SBarry Smith } 11369a612a9SBarry Smith 11469a612a9SBarry Smith 11569a612a9SBarry Smith #undef __FUNCT__ 11669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 11769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 11869a612a9SBarry Smith { 11969a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 12069a612a9SBarry Smith PetscErrorCode ierr; 1215a9f2f41SSatish Balay PC_FieldSplitLink ilink; 122cf502942SBarry Smith PetscInt i,nsplit,ccsize; 12369a612a9SBarry Smith MatStructure flag = pc->flag; 124ccb205f8SBarry Smith PetscTruth sorted; 12569a612a9SBarry Smith 12669a612a9SBarry Smith PetscFunctionBegin; 12769a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 12897bbdb24SBarry Smith nsplit = jac->nsplits; 1295a9f2f41SSatish Balay ilink = jac->head; 13097bbdb24SBarry Smith 13197bbdb24SBarry Smith /* get the matrices for each split */ 132704ba839SBarry Smith if (!jac->issetup) { 1331b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 13497bbdb24SBarry Smith 135704ba839SBarry Smith jac->issetup = PETSC_TRUE; 136704ba839SBarry Smith 137704ba839SBarry Smith /* This is done here instead of in PCFieldSplitSetFields() because may not have matrix at that point */ 13851f519a2SBarry Smith bs = jac->bs; 13997bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 140cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 1411b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 1421b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 1431b9fc7fcSBarry Smith if (jac->defaultsplit) { 144704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+i,nsplit,&ilink->is);CHKERRQ(ierr); 145704ba839SBarry Smith ilink->csize = ccsize/nsplit; 146704ba839SBarry Smith } else if (!ilink->is) { 147ccb205f8SBarry Smith if (ilink->nfields > 1) { 1485a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 1495a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 1501b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 1511b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 1521b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 15397bbdb24SBarry Smith } 15497bbdb24SBarry Smith } 155704ba839SBarry Smith ierr = ISCreateGeneral(((PetscObject)pc)->comm,nslots*nfields,ii,&ilink->is);CHKERRQ(ierr); 156ccb205f8SBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 157ccb205f8SBarry Smith } else { 158704ba839SBarry Smith ierr = ISCreateStride(((PetscObject)pc)->comm,nslots,rstart+ilink->fields[0],bs,&ilink->is);CHKERRQ(ierr); 159ccb205f8SBarry Smith } 160704ba839SBarry Smith ilink->csize = (ccsize/bs)*ilink->nfields; 161704ba839SBarry Smith ierr = ISSorted(ilink->is,&sorted);CHKERRQ(ierr); 1621b9fc7fcSBarry Smith if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 1631b9fc7fcSBarry Smith } 164704ba839SBarry Smith ierr = ISAllGather(ilink->is,&ilink->cis);CHKERRQ(ierr); 165704ba839SBarry Smith ilink = ilink->next; 1661b9fc7fcSBarry Smith } 1671b9fc7fcSBarry Smith } 1681b9fc7fcSBarry Smith 169704ba839SBarry Smith ilink = jac->head; 17097bbdb24SBarry Smith if (!jac->pmat) { 171cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 172cf502942SBarry Smith for (i=0; i<nsplit; i++) { 173704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 174704ba839SBarry Smith ilink = ilink->next; 175cf502942SBarry Smith } 17697bbdb24SBarry Smith } else { 177cf502942SBarry Smith for (i=0; i<nsplit; i++) { 178704ba839SBarry Smith ierr = MatGetSubMatrix(pc->pmat,ilink->is,ilink->cis,ilink->csize,MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 179704ba839SBarry Smith ilink = ilink->next; 180cf502942SBarry Smith } 18197bbdb24SBarry Smith } 18297bbdb24SBarry Smith 18397bbdb24SBarry Smith /* set up the individual PCs */ 18497bbdb24SBarry Smith i = 0; 1855a9f2f41SSatish Balay ilink = jac->head; 1865a9f2f41SSatish Balay while (ilink) { 1875a9f2f41SSatish Balay ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 1885a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 1895a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 19097bbdb24SBarry Smith i++; 1915a9f2f41SSatish Balay ilink = ilink->next; 1920971522cSBarry Smith } 19397bbdb24SBarry Smith 19497bbdb24SBarry Smith /* create work vectors for each split */ 1951b9fc7fcSBarry Smith if (!jac->x) { 19679416396SBarry Smith Vec xtmp; 19797bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 1985a9f2f41SSatish Balay ilink = jac->head; 19997bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 200906ed7ccSBarry Smith Vec *vl,*vr; 201906ed7ccSBarry Smith 202906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 203906ed7ccSBarry Smith ilink->x = *vr; 204906ed7ccSBarry Smith ilink->y = *vl; 205906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 206906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 2075a9f2f41SSatish Balay jac->x[i] = ilink->x; 2085a9f2f41SSatish Balay jac->y[i] = ilink->y; 2095a9f2f41SSatish Balay ilink = ilink->next; 21097bbdb24SBarry Smith } 21179416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 2121b9fc7fcSBarry Smith 2135a9f2f41SSatish Balay ilink = jac->head; 2141b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 2151b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 216704ba839SBarry Smith ierr = VecScatterCreate(xtmp,ilink->is,jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 2175a9f2f41SSatish Balay ilink = ilink->next; 2181b9fc7fcSBarry Smith } 2191b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 2201b9fc7fcSBarry Smith } 2210971522cSBarry Smith PetscFunctionReturn(0); 2220971522cSBarry Smith } 2230971522cSBarry Smith 2245a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 225ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 226ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD) || \ 2275a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 228ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE) || \ 229ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->y,yy,ADD_VALUES,SCATTER_REVERSE)) 23079416396SBarry Smith 2310971522cSBarry Smith #undef __FUNCT__ 2320971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 2330971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 2340971522cSBarry Smith { 2350971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2360971522cSBarry Smith PetscErrorCode ierr; 2375a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 238e25d487eSBarry Smith PetscInt bs; 2390971522cSBarry Smith 2400971522cSBarry Smith PetscFunctionBegin; 24151f519a2SBarry Smith CHKMEMQ; 242e25d487eSBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 24351f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 24451f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 24551f519a2SBarry Smith 24679416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 2471b9fc7fcSBarry Smith if (jac->defaultsplit) { 2480971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 2495a9f2f41SSatish Balay while (ilink) { 2505a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 2515a9f2f41SSatish Balay ilink = ilink->next; 2520971522cSBarry Smith } 2530971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 2541b9fc7fcSBarry Smith } else { 255efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 2565a9f2f41SSatish Balay while (ilink) { 2575a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 2585a9f2f41SSatish Balay ilink = ilink->next; 2591b9fc7fcSBarry Smith } 2601b9fc7fcSBarry Smith } 261*16913363SBarry Smith } else if (jac->type == PC_COMPOSITE_MULTIPLICATIVE || jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 26279416396SBarry Smith if (!jac->w1) { 26379416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 26479416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 26579416396SBarry Smith } 266efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 2675a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 2685a9f2f41SSatish Balay while (ilink->next) { 2695a9f2f41SSatish Balay ilink = ilink->next; 27079416396SBarry Smith ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr); 271efb30889SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 2725a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 27379416396SBarry Smith } 27451f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 27551f519a2SBarry Smith while (ilink->previous) { 27651f519a2SBarry Smith ilink = ilink->previous; 27751f519a2SBarry Smith ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr); 27851f519a2SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 27951f519a2SBarry Smith ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 28079416396SBarry Smith } 28151f519a2SBarry Smith } 282*16913363SBarry Smith } else SETERRQ1(PETSC_ERR_SUP,"Unsupported or unknown composition",(int) jac->type); 28351f519a2SBarry Smith CHKMEMQ; 2840971522cSBarry Smith PetscFunctionReturn(0); 2850971522cSBarry Smith } 2860971522cSBarry Smith 287421e10b8SBarry Smith #define FieldSplitSplitSolveAddTranspose(ilink,xx,yy) \ 288ca9f406cSSatish Balay (VecScatterBegin(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 289ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,xx,ilink->y,INSERT_VALUES,SCATTER_FORWARD) || \ 290421e10b8SBarry Smith KSPSolveTranspose(ilink->ksp,ilink->y,ilink->x) || \ 291ca9f406cSSatish Balay VecScatterBegin(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE) || \ 292ca9f406cSSatish Balay VecScatterEnd(ilink->sctx,ilink->x,yy,ADD_VALUES,SCATTER_REVERSE)) 293421e10b8SBarry Smith 294421e10b8SBarry Smith #undef __FUNCT__ 295421e10b8SBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 296421e10b8SBarry Smith static PetscErrorCode PCApplyTranspose_FieldSplit(PC pc,Vec x,Vec y) 297421e10b8SBarry Smith { 298421e10b8SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 299421e10b8SBarry Smith PetscErrorCode ierr; 300421e10b8SBarry Smith PC_FieldSplitLink ilink = jac->head; 301421e10b8SBarry Smith PetscInt bs; 302421e10b8SBarry Smith 303421e10b8SBarry Smith PetscFunctionBegin; 304421e10b8SBarry Smith CHKMEMQ; 305421e10b8SBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 306421e10b8SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 307421e10b8SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 308421e10b8SBarry Smith 309421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 310421e10b8SBarry Smith if (jac->defaultsplit) { 311421e10b8SBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 312421e10b8SBarry Smith while (ilink) { 313421e10b8SBarry Smith ierr = KSPSolveTranspose(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 314421e10b8SBarry Smith ilink = ilink->next; 315421e10b8SBarry Smith } 316421e10b8SBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 317421e10b8SBarry Smith } else { 318421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 319421e10b8SBarry Smith while (ilink) { 320421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 321421e10b8SBarry Smith ilink = ilink->next; 322421e10b8SBarry Smith } 323421e10b8SBarry Smith } 324421e10b8SBarry Smith } else { 325421e10b8SBarry Smith if (!jac->w1) { 326421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 327421e10b8SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 328421e10b8SBarry Smith } 329421e10b8SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 330421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 331421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 332421e10b8SBarry Smith while (ilink->next) { 333421e10b8SBarry Smith ilink = ilink->next; 334421e10b8SBarry Smith ierr = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr); 335421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 336421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 337421e10b8SBarry Smith } 338421e10b8SBarry Smith while (ilink->previous) { 339421e10b8SBarry Smith ilink = ilink->previous; 340421e10b8SBarry Smith ierr = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr); 341421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 342421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 343421e10b8SBarry Smith } 344421e10b8SBarry Smith } else { 345421e10b8SBarry Smith while (ilink->next) { /* get to last entry in linked list */ 346421e10b8SBarry Smith ilink = ilink->next; 347421e10b8SBarry Smith } 348421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,x,y);CHKERRQ(ierr); 349421e10b8SBarry Smith while (ilink->previous) { 350421e10b8SBarry Smith ilink = ilink->previous; 351421e10b8SBarry Smith ierr = MatMultTranspose(pc->pmat,y,jac->w1);CHKERRQ(ierr); 352421e10b8SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 353421e10b8SBarry Smith ierr = FieldSplitSplitSolveAddTranspose(ilink,jac->w2,y);CHKERRQ(ierr); 354421e10b8SBarry Smith } 355421e10b8SBarry Smith } 356421e10b8SBarry Smith } 357421e10b8SBarry Smith CHKMEMQ; 358421e10b8SBarry Smith PetscFunctionReturn(0); 359421e10b8SBarry Smith } 360421e10b8SBarry Smith 3610971522cSBarry Smith #undef __FUNCT__ 3620971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 3630971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 3640971522cSBarry Smith { 3650971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3660971522cSBarry Smith PetscErrorCode ierr; 3675a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 3680971522cSBarry Smith 3690971522cSBarry Smith PetscFunctionBegin; 3705a9f2f41SSatish Balay while (ilink) { 3715a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 3725a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 3735a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 3745a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 375704ba839SBarry Smith if (ilink->is) {ierr = ISDestroy(ilink->is);CHKERRQ(ierr);} 376704ba839SBarry Smith if (ilink->cis) {ierr = ISDestroy(ilink->cis);CHKERRQ(ierr);} 3775a9f2f41SSatish Balay next = ilink->next; 378704ba839SBarry Smith ierr = PetscFree(ilink->fields);CHKERRQ(ierr); 379704ba839SBarry Smith ierr = PetscFree(ilink);CHKERRQ(ierr); 3805a9f2f41SSatish Balay ilink = next; 3810971522cSBarry Smith } 38205b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 38397bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 38479416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 38579416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 3860971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 3870971522cSBarry Smith PetscFunctionReturn(0); 3880971522cSBarry Smith } 3890971522cSBarry Smith 3900971522cSBarry Smith #undef __FUNCT__ 3910971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 3920971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 3930971522cSBarry Smith { 3941b9fc7fcSBarry Smith PetscErrorCode ierr; 39551f519a2SBarry Smith PetscInt i = 0,nfields,*fields,bs; 3961b9fc7fcSBarry Smith PetscTruth flg; 3971b9fc7fcSBarry Smith char optionname[128]; 3989dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3991b9fc7fcSBarry Smith 4000971522cSBarry Smith PetscFunctionBegin; 4011b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 40251f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 40351f519a2SBarry Smith if (flg) { 40451f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 40551f519a2SBarry Smith } 40651f519a2SBarry Smith if (jac->bs <= 0) { 407704ba839SBarry Smith 40851f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,1);CHKERRQ(ierr); 40951f519a2SBarry Smith } 4109dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 41151f519a2SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 4121b9fc7fcSBarry Smith while (PETSC_TRUE) { 41313f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 41451f519a2SBarry Smith nfields = jac->bs; 4151b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 4161b9fc7fcSBarry Smith if (!flg) break; 4171b9fc7fcSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 4181b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 4191b9fc7fcSBarry Smith } 42051f519a2SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 4211b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 4220971522cSBarry Smith PetscFunctionReturn(0); 4230971522cSBarry Smith } 4240971522cSBarry Smith 4250971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 4260971522cSBarry Smith 4270971522cSBarry Smith EXTERN_C_BEGIN 4280971522cSBarry Smith #undef __FUNCT__ 4290971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 430dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 4310971522cSBarry Smith { 43297bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4330971522cSBarry Smith PetscErrorCode ierr; 4345a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 43569a612a9SBarry Smith char prefix[128]; 43651f519a2SBarry Smith PetscInt i; 4370971522cSBarry Smith 4380971522cSBarry Smith PetscFunctionBegin; 4390971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 44051f519a2SBarry Smith for (i=0; i<n; i++) { 44151f519a2SBarry Smith if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 44251f519a2SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 44351f519a2SBarry Smith } 444704ba839SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 445704ba839SBarry Smith ierr = PetscMalloc(n*sizeof(PetscInt),&ilink->fields);CHKERRQ(ierr); 4465a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 4475a9f2f41SSatish Balay ilink->nfields = n; 4485a9f2f41SSatish Balay ilink->next = PETSC_NULL; 4497adad957SLisandro Dalcin ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 4505a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 45169a612a9SBarry Smith 4527adad957SLisandro Dalcin if (((PetscObject)pc)->prefix) { 4537adad957SLisandro Dalcin sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 45469a612a9SBarry Smith } else { 45513f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 45669a612a9SBarry Smith } 4575a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 4580971522cSBarry Smith 4590971522cSBarry Smith if (!next) { 4605a9f2f41SSatish Balay jac->head = ilink; 46151f519a2SBarry Smith ilink->previous = PETSC_NULL; 4620971522cSBarry Smith } else { 4630971522cSBarry Smith while (next->next) { 4640971522cSBarry Smith next = next->next; 4650971522cSBarry Smith } 4665a9f2f41SSatish Balay next->next = ilink; 46751f519a2SBarry Smith ilink->previous = next; 4680971522cSBarry Smith } 4690971522cSBarry Smith jac->nsplits++; 4700971522cSBarry Smith PetscFunctionReturn(0); 4710971522cSBarry Smith } 4720971522cSBarry Smith EXTERN_C_END 4730971522cSBarry Smith 4740971522cSBarry Smith 4750971522cSBarry Smith EXTERN_C_BEGIN 4760971522cSBarry Smith #undef __FUNCT__ 47769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 478dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 4790971522cSBarry Smith { 4800971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 4810971522cSBarry Smith PetscErrorCode ierr; 4820971522cSBarry Smith PetscInt cnt = 0; 4835a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 4840971522cSBarry Smith 4850971522cSBarry Smith PetscFunctionBegin; 48669a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 4875a9f2f41SSatish Balay while (ilink) { 4885a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 4895a9f2f41SSatish Balay ilink = ilink->next; 4900971522cSBarry Smith } 4910971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 4920971522cSBarry Smith *n = jac->nsplits; 4930971522cSBarry Smith PetscFunctionReturn(0); 4940971522cSBarry Smith } 4950971522cSBarry Smith EXTERN_C_END 4960971522cSBarry Smith 497704ba839SBarry Smith EXTERN_C_BEGIN 498704ba839SBarry Smith #undef __FUNCT__ 499704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS_FieldSplit" 500704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS_FieldSplit(PC pc,IS is) 501704ba839SBarry Smith { 502704ba839SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 503704ba839SBarry Smith PetscErrorCode ierr; 504704ba839SBarry Smith PC_FieldSplitLink ilink, next = jac->head; 505704ba839SBarry Smith char prefix[128]; 506704ba839SBarry Smith 507704ba839SBarry Smith PetscFunctionBegin; 508*16913363SBarry Smith ierr = PetscNew(struct _PC_FieldSplitLink,&ilink);CHKERRQ(ierr); 509704ba839SBarry Smith ilink->next = PETSC_NULL; 510704ba839SBarry Smith ierr = KSPCreate(((PetscObject)pc)->comm,&ilink->ksp);CHKERRQ(ierr); 511704ba839SBarry Smith ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 512704ba839SBarry Smith 513704ba839SBarry Smith if (((PetscObject)pc)->prefix) { 514704ba839SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",((PetscObject)pc)->prefix,(int)jac->nsplits); 515704ba839SBarry Smith } else { 516704ba839SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 517704ba839SBarry Smith } 518704ba839SBarry Smith ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 519704ba839SBarry Smith 520704ba839SBarry Smith if (!next) { 521704ba839SBarry Smith jac->head = ilink; 522704ba839SBarry Smith ilink->previous = PETSC_NULL; 523704ba839SBarry Smith } else { 524704ba839SBarry Smith while (next->next) { 525704ba839SBarry Smith next = next->next; 526704ba839SBarry Smith } 527704ba839SBarry Smith next->next = ilink; 528704ba839SBarry Smith ilink->previous = next; 529704ba839SBarry Smith } 530704ba839SBarry Smith jac->nsplits++; 531704ba839SBarry Smith 532704ba839SBarry Smith PetscFunctionReturn(0); 533704ba839SBarry Smith } 534704ba839SBarry Smith EXTERN_C_END 535704ba839SBarry Smith 5360971522cSBarry Smith #undef __FUNCT__ 5370971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 5380971522cSBarry Smith /*@ 5390971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 5400971522cSBarry Smith 5410971522cSBarry Smith Collective on PC 5420971522cSBarry Smith 5430971522cSBarry Smith Input Parameters: 5440971522cSBarry Smith + pc - the preconditioner context 5450971522cSBarry Smith . n - the number of fields in this split 5460971522cSBarry Smith . fields - the fields in this split 5470971522cSBarry Smith 5480971522cSBarry Smith Level: intermediate 5490971522cSBarry Smith 55051f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 5510971522cSBarry Smith 5520971522cSBarry Smith @*/ 553dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 5540971522cSBarry Smith { 5550971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 5560971522cSBarry Smith 5570971522cSBarry Smith PetscFunctionBegin; 5580971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5590971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 5600971522cSBarry Smith if (f) { 5610971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 5620971522cSBarry Smith } 5630971522cSBarry Smith PetscFunctionReturn(0); 5640971522cSBarry Smith } 5650971522cSBarry Smith 5660971522cSBarry Smith #undef __FUNCT__ 567704ba839SBarry Smith #define __FUNCT__ "PCFieldSplitSetIS" 568704ba839SBarry Smith /*@ 569704ba839SBarry Smith PCFieldSplitSetIS - Sets the exact elements for field 570704ba839SBarry Smith 571704ba839SBarry Smith Collective on PC 572704ba839SBarry Smith 573704ba839SBarry Smith Input Parameters: 574704ba839SBarry Smith + pc - the preconditioner context 575704ba839SBarry Smith . is - the index set that defines the vector elements in this field 576704ba839SBarry Smith 577704ba839SBarry Smith Level: intermediate 578704ba839SBarry Smith 579704ba839SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 580704ba839SBarry Smith 581704ba839SBarry Smith @*/ 582704ba839SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetIS(PC pc,IS is) 583704ba839SBarry Smith { 584704ba839SBarry Smith PetscErrorCode ierr,(*f)(PC,IS); 585704ba839SBarry Smith 586704ba839SBarry Smith PetscFunctionBegin; 587704ba839SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 588704ba839SBarry Smith PetscValidHeaderSpecific(is,IS_COOKIE,1); 589704ba839SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetIS_C",(void (**)(void))&f);CHKERRQ(ierr); 590704ba839SBarry Smith if (f) { 591704ba839SBarry Smith ierr = (*f)(pc,is);CHKERRQ(ierr); 592704ba839SBarry Smith } 593704ba839SBarry Smith PetscFunctionReturn(0); 594704ba839SBarry Smith } 595704ba839SBarry Smith 596704ba839SBarry Smith #undef __FUNCT__ 59751f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 59851f519a2SBarry Smith /*@ 59951f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 60051f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 60151f519a2SBarry Smith 60251f519a2SBarry Smith Collective on PC 60351f519a2SBarry Smith 60451f519a2SBarry Smith Input Parameters: 60551f519a2SBarry Smith + pc - the preconditioner context 60651f519a2SBarry Smith - bs - the block size 60751f519a2SBarry Smith 60851f519a2SBarry Smith Level: intermediate 60951f519a2SBarry Smith 61051f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 61151f519a2SBarry Smith 61251f519a2SBarry Smith @*/ 61351f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 61451f519a2SBarry Smith { 61551f519a2SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 61651f519a2SBarry Smith 61751f519a2SBarry Smith PetscFunctionBegin; 61851f519a2SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 61951f519a2SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 62051f519a2SBarry Smith if (f) { 62151f519a2SBarry Smith ierr = (*f)(pc,bs);CHKERRQ(ierr); 62251f519a2SBarry Smith } 62351f519a2SBarry Smith PetscFunctionReturn(0); 62451f519a2SBarry Smith } 62551f519a2SBarry Smith 62651f519a2SBarry Smith #undef __FUNCT__ 62769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 6280971522cSBarry Smith /*@C 62969a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 6300971522cSBarry Smith 63169a612a9SBarry Smith Collective on KSP 6320971522cSBarry Smith 6330971522cSBarry Smith Input Parameter: 6340971522cSBarry Smith . pc - the preconditioner context 6350971522cSBarry Smith 6360971522cSBarry Smith Output Parameters: 6370971522cSBarry Smith + n - the number of split 63869a612a9SBarry Smith - pc - the array of KSP contexts 6390971522cSBarry Smith 6400971522cSBarry Smith Note: 64169a612a9SBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed 6420971522cSBarry Smith 64369a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 6440971522cSBarry Smith 6450971522cSBarry Smith Level: advanced 6460971522cSBarry Smith 6470971522cSBarry Smith .seealso: PCFIELDSPLIT 6480971522cSBarry Smith @*/ 649dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 6500971522cSBarry Smith { 65169a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 6520971522cSBarry Smith 6530971522cSBarry Smith PetscFunctionBegin; 6540971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6550971522cSBarry Smith PetscValidIntPointer(n,2); 65669a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 6570971522cSBarry Smith if (f) { 65869a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 6590971522cSBarry Smith } else { 66069a612a9SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 6610971522cSBarry Smith } 6620971522cSBarry Smith PetscFunctionReturn(0); 6630971522cSBarry Smith } 6640971522cSBarry Smith 66579416396SBarry Smith EXTERN_C_BEGIN 66679416396SBarry Smith #undef __FUNCT__ 66779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 668dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 66979416396SBarry Smith { 67079416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 67179416396SBarry Smith 67279416396SBarry Smith PetscFunctionBegin; 67379416396SBarry Smith jac->type = type; 67479416396SBarry Smith PetscFunctionReturn(0); 67579416396SBarry Smith } 67679416396SBarry Smith EXTERN_C_END 67779416396SBarry Smith 67851f519a2SBarry Smith EXTERN_C_BEGIN 67951f519a2SBarry Smith #undef __FUNCT__ 68051f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 68151f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 68251f519a2SBarry Smith { 68351f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 68451f519a2SBarry Smith 68551f519a2SBarry Smith PetscFunctionBegin; 68651f519a2SBarry Smith if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 68751f519a2SBarry 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); 68851f519a2SBarry Smith jac->bs = bs; 68951f519a2SBarry Smith PetscFunctionReturn(0); 69051f519a2SBarry Smith } 69151f519a2SBarry Smith EXTERN_C_END 69251f519a2SBarry Smith 69379416396SBarry Smith #undef __FUNCT__ 69479416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 695bc08b0f1SBarry Smith /*@ 69679416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 69779416396SBarry Smith 69879416396SBarry Smith Collective on PC 69979416396SBarry Smith 70079416396SBarry Smith Input Parameter: 70179416396SBarry Smith . pc - the preconditioner context 702ce9499c7SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE 70379416396SBarry Smith 70479416396SBarry Smith Options Database Key: 705ce9499c7SBarry Smith . -pc_fieldsplit_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets fieldsplit preconditioner type 70679416396SBarry Smith 70779416396SBarry Smith Level: Developer 70879416396SBarry Smith 70979416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 71079416396SBarry Smith 71179416396SBarry Smith .seealso: PCCompositeSetType() 71279416396SBarry Smith 71379416396SBarry Smith @*/ 714dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 71579416396SBarry Smith { 71679416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 71779416396SBarry Smith 71879416396SBarry Smith PetscFunctionBegin; 71979416396SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 72079416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 72179416396SBarry Smith if (f) { 72279416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 72379416396SBarry Smith } 72479416396SBarry Smith PetscFunctionReturn(0); 72579416396SBarry Smith } 72679416396SBarry Smith 7270971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 7280971522cSBarry Smith /*MC 729a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 7300971522cSBarry Smith fields or groups of fields 7310971522cSBarry Smith 7320971522cSBarry Smith 7330971522cSBarry Smith To set options on the solvers for each block append -sub_ to all the PC 734d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 7350971522cSBarry Smith 736a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 73769a612a9SBarry Smith and set the options directly on the resulting KSP object 7380971522cSBarry Smith 7390971522cSBarry Smith Level: intermediate 7400971522cSBarry Smith 74179416396SBarry Smith Options Database Keys: 7422e70eadcSBarry Smith + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 7432e70eadcSBarry Smith . -pc_splitfield_default - automatically add any fields to additional splits that have not 7442e70eadcSBarry Smith been supplied explicitly by -pc_splitfield_%d_fields 74551f519a2SBarry Smith . -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 7462e70eadcSBarry Smith - -pc_splitfield_type <additive,multiplicative> 74779416396SBarry Smith 7480971522cSBarry Smith Concepts: physics based preconditioners 7490971522cSBarry Smith 7500971522cSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 751c8d321e0SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType() 7520971522cSBarry Smith M*/ 7530971522cSBarry Smith 7540971522cSBarry Smith EXTERN_C_BEGIN 7550971522cSBarry Smith #undef __FUNCT__ 7560971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 757dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 7580971522cSBarry Smith { 7590971522cSBarry Smith PetscErrorCode ierr; 7600971522cSBarry Smith PC_FieldSplit *jac; 7610971522cSBarry Smith 7620971522cSBarry Smith PetscFunctionBegin; 76338f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_FieldSplit,&jac);CHKERRQ(ierr); 7640971522cSBarry Smith jac->bs = -1; 7650971522cSBarry Smith jac->nsplits = 0; 76679416396SBarry Smith jac->type = PC_COMPOSITE_ADDITIVE; 76751f519a2SBarry Smith 7680971522cSBarry Smith pc->data = (void*)jac; 7690971522cSBarry Smith 7700971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 771421e10b8SBarry Smith pc->ops->applytranspose = PCApplyTranspose_FieldSplit; 7720971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 7730971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 7740971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 7750971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 7760971522cSBarry Smith pc->ops->applyrichardson = 0; 7770971522cSBarry Smith 77869a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 77969a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 7800971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 7810971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 782704ba839SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetIS_C","PCFieldSplitSetIS_FieldSplit", 783704ba839SBarry Smith PCFieldSplitSetIS_FieldSplit);CHKERRQ(ierr); 78479416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 78579416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 78651f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 78751f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 7880971522cSBarry Smith PetscFunctionReturn(0); 7890971522cSBarry Smith } 7900971522cSBarry Smith EXTERN_C_END 7910971522cSBarry Smith 7920971522cSBarry Smith 793