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; 15*51f519a2SBarry Smith PC_FieldSplitLink next,previous; 160971522cSBarry Smith }; 170971522cSBarry Smith 180971522cSBarry Smith typedef struct { 1979416396SBarry Smith PCCompositeType type; /* additive or multiplicative */ 2097bbdb24SBarry Smith PetscTruth defaultsplit; 210971522cSBarry Smith PetscInt bs; 22cf502942SBarry Smith PetscInt nsplits,*csize; 2379416396SBarry Smith Vec *x,*y,w1,w2; 2497bbdb24SBarry Smith Mat *pmat; 25cf502942SBarry Smith IS *is,*cis; 2697bbdb24SBarry Smith PC_FieldSplitLink head; 270971522cSBarry Smith } PC_FieldSplit; 280971522cSBarry Smith 290971522cSBarry Smith #undef __FUNCT__ 300971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 310971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 320971522cSBarry Smith { 330971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 340971522cSBarry Smith PetscErrorCode ierr; 350971522cSBarry Smith PetscTruth iascii; 360971522cSBarry Smith PetscInt i,j; 375a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 380971522cSBarry Smith 390971522cSBarry Smith PetscFunctionBegin; 400971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 410971522cSBarry Smith if (iascii) { 42*51f519a2SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D, blocksize = %D\n",PCCompositeTypes[jac->type],jac->nsplits,jac->bs);CHKERRQ(ierr); 4369a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 440971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 450971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 460971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 4779416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 485a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 4979416396SBarry Smith if (j > 0) { 5079416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 5179416396SBarry Smith } 525a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 530971522cSBarry Smith } 540971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 5579416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 565a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 575a9f2f41SSatish Balay ilink = ilink->next; 580971522cSBarry Smith } 590971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 600971522cSBarry Smith } else { 610971522cSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 620971522cSBarry Smith } 630971522cSBarry Smith PetscFunctionReturn(0); 640971522cSBarry Smith } 650971522cSBarry Smith 660971522cSBarry Smith #undef __FUNCT__ 6769a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 6869a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 690971522cSBarry Smith { 700971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 710971522cSBarry Smith PetscErrorCode ierr; 725a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 7369a612a9SBarry Smith PetscInt i; 7479416396SBarry Smith PetscTruth flg = PETSC_FALSE,*fields; 750971522cSBarry Smith 760971522cSBarry Smith PetscFunctionBegin; 774dc4c822SBarry Smith ierr = PetscOptionsGetTruth(pc->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 785a9f2f41SSatish Balay if (!ilink || flg) { 79ae15b995SBarry Smith ierr = PetscInfo(pc,"Using default splitting of fields\n");CHKERRQ(ierr); 80521d7252SBarry Smith if (jac->bs <= 0) { 81521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 82521d7252SBarry Smith } 8379416396SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 8479416396SBarry Smith ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 855a9f2f41SSatish Balay while (ilink) { 865a9f2f41SSatish Balay for (i=0; i<ilink->nfields; i++) { 875a9f2f41SSatish Balay fields[ilink->fields[i]] = PETSC_TRUE; 88521d7252SBarry Smith } 895a9f2f41SSatish Balay ilink = ilink->next; 9079416396SBarry Smith } 9197bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 9279416396SBarry Smith for (i=0; i<jac->bs; i++) { 9379416396SBarry Smith if (!fields[i]) { 9479416396SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 9579416396SBarry Smith } else { 9679416396SBarry Smith jac->defaultsplit = PETSC_FALSE; 9779416396SBarry Smith } 9879416396SBarry Smith } 99521d7252SBarry Smith } 10069a612a9SBarry Smith PetscFunctionReturn(0); 10169a612a9SBarry Smith } 10269a612a9SBarry Smith 10369a612a9SBarry Smith 10469a612a9SBarry Smith #undef __FUNCT__ 10569a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 10669a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 10769a612a9SBarry Smith { 10869a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 10969a612a9SBarry Smith PetscErrorCode ierr; 1105a9f2f41SSatish Balay PC_FieldSplitLink ilink; 111cf502942SBarry Smith PetscInt i,nsplit,ccsize; 11269a612a9SBarry Smith MatStructure flag = pc->flag; 11369a612a9SBarry Smith 11469a612a9SBarry Smith PetscFunctionBegin; 11569a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 11697bbdb24SBarry Smith nsplit = jac->nsplits; 1175a9f2f41SSatish Balay ilink = jac->head; 11897bbdb24SBarry Smith 11997bbdb24SBarry Smith /* get the matrices for each split */ 12097bbdb24SBarry Smith if (!jac->is) { 1211b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 12297bbdb24SBarry Smith 123*51f519a2SBarry Smith bs = jac->bs; 12497bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 125cf502942SBarry Smith ierr = MatGetLocalSize(pc->pmat,PETSC_NULL,&ccsize);CHKERRQ(ierr); 1261b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 1271b9fc7fcSBarry Smith ierr = PetscMalloc(nsplit*sizeof(IS),&jac->is);CHKERRQ(ierr); 128cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(IS),&jac->cis);CHKERRQ(ierr); 129cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(PetscInt),&jac->csize);CHKERRQ(ierr); 1301b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 1311b9fc7fcSBarry Smith if (jac->defaultsplit) { 1321b9fc7fcSBarry Smith ierr = ISCreateStride(pc->comm,nslots,rstart+i,nsplit,&jac->is[i]);CHKERRQ(ierr); 133cf502942SBarry Smith jac->csize[i] = ccsize/nsplit; 13497bbdb24SBarry Smith } else { 1355a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 1361b9fc7fcSBarry Smith PetscTruth sorted; 1375a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 1381b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 1391b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 1401b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 14197bbdb24SBarry Smith } 14297bbdb24SBarry Smith } 1431b9fc7fcSBarry Smith ierr = ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);CHKERRQ(ierr); 144cf502942SBarry Smith jac->csize[i] = (ccsize/bs)*ilink->nfields; 1451b9fc7fcSBarry Smith ierr = ISSorted(jac->is[i],&sorted);CHKERRQ(ierr); 1461b9fc7fcSBarry Smith if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 1471b9fc7fcSBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 1485a9f2f41SSatish Balay ilink = ilink->next; 1491b9fc7fcSBarry Smith } 150cf502942SBarry Smith ierr = ISAllGather(jac->is[i],&jac->cis[i]);CHKERRQ(ierr); 1511b9fc7fcSBarry Smith } 1521b9fc7fcSBarry Smith } 1531b9fc7fcSBarry Smith 15497bbdb24SBarry Smith if (!jac->pmat) { 155cf502942SBarry Smith ierr = PetscMalloc(nsplit*sizeof(Mat),&jac->pmat);CHKERRQ(ierr); 156cf502942SBarry Smith for (i=0; i<nsplit; i++) { 157cf502942SBarry Smith ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_INITIAL_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 158cf502942SBarry Smith } 15997bbdb24SBarry Smith } else { 160cf502942SBarry Smith for (i=0; i<nsplit; i++) { 161cf502942SBarry Smith ierr = MatGetSubMatrix(pc->pmat,jac->is[i],jac->cis[i],jac->csize[i],MAT_REUSE_MATRIX,&jac->pmat[i]);CHKERRQ(ierr); 162cf502942SBarry Smith } 16397bbdb24SBarry Smith } 16497bbdb24SBarry Smith 16597bbdb24SBarry Smith /* set up the individual PCs */ 16697bbdb24SBarry Smith i = 0; 1675a9f2f41SSatish Balay ilink = jac->head; 1685a9f2f41SSatish Balay while (ilink) { 1695a9f2f41SSatish Balay ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 1705a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 1715a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 17297bbdb24SBarry Smith i++; 1735a9f2f41SSatish Balay ilink = ilink->next; 1740971522cSBarry Smith } 17597bbdb24SBarry Smith 17697bbdb24SBarry Smith /* create work vectors for each split */ 1771b9fc7fcSBarry Smith if (!jac->x) { 17879416396SBarry Smith Vec xtmp; 17997bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 1805a9f2f41SSatish Balay ilink = jac->head; 18197bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 182906ed7ccSBarry Smith Vec *vl,*vr; 183906ed7ccSBarry Smith 184906ed7ccSBarry Smith ierr = KSPGetVecs(ilink->ksp,1,&vr,1,&vl);CHKERRQ(ierr); 185906ed7ccSBarry Smith ilink->x = *vr; 186906ed7ccSBarry Smith ilink->y = *vl; 187906ed7ccSBarry Smith ierr = PetscFree(vr);CHKERRQ(ierr); 188906ed7ccSBarry Smith ierr = PetscFree(vl);CHKERRQ(ierr); 1895a9f2f41SSatish Balay jac->x[i] = ilink->x; 1905a9f2f41SSatish Balay jac->y[i] = ilink->y; 1915a9f2f41SSatish Balay ilink = ilink->next; 19297bbdb24SBarry Smith } 19379416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 1941b9fc7fcSBarry Smith 1955a9f2f41SSatish Balay ilink = jac->head; 1961b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 1971b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 1985a9f2f41SSatish Balay ierr = VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 1995a9f2f41SSatish Balay ilink = ilink->next; 2001b9fc7fcSBarry Smith } 2011b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 2021b9fc7fcSBarry Smith } 2030971522cSBarry Smith PetscFunctionReturn(0); 2040971522cSBarry Smith } 2050971522cSBarry Smith 2065a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 2075a9f2f41SSatish Balay (VecScatterBegin(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \ 2085a9f2f41SSatish Balay VecScatterEnd(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \ 2095a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 2105a9f2f41SSatish Balay VecScatterBegin(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \ 2115a9f2f41SSatish Balay VecScatterEnd(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx)) 21279416396SBarry Smith 2130971522cSBarry Smith #undef __FUNCT__ 2140971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 2150971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 2160971522cSBarry Smith { 2170971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2180971522cSBarry Smith PetscErrorCode ierr; 2195a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 220e25d487eSBarry Smith PetscInt bs; 2210971522cSBarry Smith 2220971522cSBarry Smith PetscFunctionBegin; 223*51f519a2SBarry Smith CHKMEMQ; 224e25d487eSBarry Smith ierr = VecGetBlockSize(x,&bs);CHKERRQ(ierr); 225*51f519a2SBarry Smith ierr = VecSetBlockSize(x,jac->bs);CHKERRQ(ierr); 226*51f519a2SBarry Smith ierr = VecSetBlockSize(y,jac->bs);CHKERRQ(ierr); 227*51f519a2SBarry Smith 22879416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 2291b9fc7fcSBarry Smith if (jac->defaultsplit) { 2300971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 2315a9f2f41SSatish Balay while (ilink) { 2325a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 2335a9f2f41SSatish Balay ilink = ilink->next; 2340971522cSBarry Smith } 2350971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 2361b9fc7fcSBarry Smith } else { 2371b9fc7fcSBarry Smith PetscInt i = 0; 2381b9fc7fcSBarry Smith 239efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 2405a9f2f41SSatish Balay while (ilink) { 2415a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 2425a9f2f41SSatish Balay ilink = ilink->next; 2431b9fc7fcSBarry Smith i++; 2441b9fc7fcSBarry Smith } 2451b9fc7fcSBarry Smith } 24679416396SBarry Smith } else { 24779416396SBarry Smith if (!jac->w1) { 24879416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 24979416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 25079416396SBarry Smith } 251efb30889SBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 2525a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 2535a9f2f41SSatish Balay while (ilink->next) { 2545a9f2f41SSatish Balay ilink = ilink->next; 25579416396SBarry Smith ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr); 256efb30889SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 2575a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 25879416396SBarry Smith } 259*51f519a2SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 260*51f519a2SBarry Smith while (ilink->previous) { 261*51f519a2SBarry Smith ilink = ilink->previous; 262*51f519a2SBarry Smith ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr); 263*51f519a2SBarry Smith ierr = VecWAXPY(jac->w2,-1.0,jac->w1,x);CHKERRQ(ierr); 264*51f519a2SBarry Smith ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 26579416396SBarry Smith } 266*51f519a2SBarry Smith } 267*51f519a2SBarry Smith } 268*51f519a2SBarry Smith CHKMEMQ; 2690971522cSBarry Smith PetscFunctionReturn(0); 2700971522cSBarry Smith } 2710971522cSBarry Smith 2720971522cSBarry Smith #undef __FUNCT__ 2730971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 2740971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 2750971522cSBarry Smith { 2760971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2770971522cSBarry Smith PetscErrorCode ierr; 2785a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 2790971522cSBarry Smith 2800971522cSBarry Smith PetscFunctionBegin; 2815a9f2f41SSatish Balay while (ilink) { 2825a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 2835a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 2845a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 2855a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 2865a9f2f41SSatish Balay next = ilink->next; 2875a9f2f41SSatish Balay ierr = PetscFree2(ilink,ilink->fields);CHKERRQ(ierr); 2885a9f2f41SSatish Balay ilink = next; 2890971522cSBarry Smith } 29005b42c5fSBarry Smith ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr); 29197bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 29297bbdb24SBarry Smith if (jac->is) { 29397bbdb24SBarry Smith PetscInt i; 29497bbdb24SBarry Smith for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);} 29597bbdb24SBarry Smith ierr = PetscFree(jac->is);CHKERRQ(ierr); 29697bbdb24SBarry Smith } 297cf502942SBarry Smith if (jac->cis) { 298cf502942SBarry Smith PetscInt i; 299cf502942SBarry Smith for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->cis[i]);CHKERRQ(ierr);} 300cf502942SBarry Smith ierr = PetscFree(jac->cis);CHKERRQ(ierr); 301cf502942SBarry Smith } 30279416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 30379416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 304cf502942SBarry Smith ierr = PetscFree(jac->csize);CHKERRQ(ierr); 3050971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 3060971522cSBarry Smith PetscFunctionReturn(0); 3070971522cSBarry Smith } 3080971522cSBarry Smith 3090971522cSBarry Smith #undef __FUNCT__ 3100971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 3110971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 3120971522cSBarry Smith { 3131b9fc7fcSBarry Smith PetscErrorCode ierr; 314*51f519a2SBarry Smith PetscInt i = 0,nfields,*fields,bs; 3151b9fc7fcSBarry Smith PetscTruth flg; 3161b9fc7fcSBarry Smith char optionname[128]; 3179dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3181b9fc7fcSBarry Smith 3190971522cSBarry Smith PetscFunctionBegin; 3201b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 321*51f519a2SBarry Smith ierr = PetscOptionsInt("-pc_fieldsplit_block_size","Blocksize that defines number of fields","PCFieldSplitSetBlockSize",jac->bs,&bs,&flg);CHKERRQ(ierr); 322*51f519a2SBarry Smith if (flg) { 323*51f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,bs);CHKERRQ(ierr); 324*51f519a2SBarry Smith } 325*51f519a2SBarry Smith if (jac->bs <= 0) { 326*51f519a2SBarry Smith ierr = PCFieldSplitSetBlockSize(pc,1);CHKERRQ(ierr); 327*51f519a2SBarry Smith } 3289dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 329*51f519a2SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscInt),&fields);CHKERRQ(ierr); 3301b9fc7fcSBarry Smith while (PETSC_TRUE) { 33113f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 332*51f519a2SBarry Smith nfields = jac->bs; 3331b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 3341b9fc7fcSBarry Smith if (!flg) break; 3351b9fc7fcSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 3361b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 3371b9fc7fcSBarry Smith } 338*51f519a2SBarry Smith ierr = PetscFree(fields);CHKERRQ(ierr); 3391b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3400971522cSBarry Smith PetscFunctionReturn(0); 3410971522cSBarry Smith } 3420971522cSBarry Smith 3430971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 3440971522cSBarry Smith 3450971522cSBarry Smith EXTERN_C_BEGIN 3460971522cSBarry Smith #undef __FUNCT__ 3470971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 348dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 3490971522cSBarry Smith { 35097bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3510971522cSBarry Smith PetscErrorCode ierr; 3525a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 35369a612a9SBarry Smith char prefix[128]; 354*51f519a2SBarry Smith PetscInt i; 3550971522cSBarry Smith 3560971522cSBarry Smith PetscFunctionBegin; 3570971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 358*51f519a2SBarry Smith for (i=0; i<n; i++) { 359*51f519a2SBarry Smith if (fields[i] >= jac->bs) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Field %D requested but only %D exist",fields[i],jac->bs); 360*51f519a2SBarry Smith if (fields[i] < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Negative field %D requested",fields[i]); 361*51f519a2SBarry Smith } 3625a9f2f41SSatish Balay ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);CHKERRQ(ierr); 3635a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 3645a9f2f41SSatish Balay ilink->nfields = n; 3655a9f2f41SSatish Balay ilink->next = PETSC_NULL; 3665a9f2f41SSatish Balay ierr = KSPCreate(pc->comm,&ilink->ksp);CHKERRQ(ierr); 3675a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 36869a612a9SBarry Smith 36969a612a9SBarry Smith if (pc->prefix) { 37013f74950SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits); 37169a612a9SBarry Smith } else { 37213f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 37369a612a9SBarry Smith } 3745a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 3750971522cSBarry Smith 3760971522cSBarry Smith if (!next) { 3775a9f2f41SSatish Balay jac->head = ilink; 378*51f519a2SBarry Smith ilink->previous = PETSC_NULL; 3790971522cSBarry Smith } else { 3800971522cSBarry Smith while (next->next) { 3810971522cSBarry Smith next = next->next; 3820971522cSBarry Smith } 3835a9f2f41SSatish Balay next->next = ilink; 384*51f519a2SBarry Smith ilink->previous = next; 3850971522cSBarry Smith } 3860971522cSBarry Smith jac->nsplits++; 3870971522cSBarry Smith PetscFunctionReturn(0); 3880971522cSBarry Smith } 3890971522cSBarry Smith EXTERN_C_END 3900971522cSBarry Smith 3910971522cSBarry Smith 3920971522cSBarry Smith EXTERN_C_BEGIN 3930971522cSBarry Smith #undef __FUNCT__ 39469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 395dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 3960971522cSBarry Smith { 3970971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3980971522cSBarry Smith PetscErrorCode ierr; 3990971522cSBarry Smith PetscInt cnt = 0; 4005a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 4010971522cSBarry Smith 4020971522cSBarry Smith PetscFunctionBegin; 40369a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 4045a9f2f41SSatish Balay while (ilink) { 4055a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 4065a9f2f41SSatish Balay ilink = ilink->next; 4070971522cSBarry Smith } 4080971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 4090971522cSBarry Smith *n = jac->nsplits; 4100971522cSBarry Smith PetscFunctionReturn(0); 4110971522cSBarry Smith } 4120971522cSBarry Smith EXTERN_C_END 4130971522cSBarry Smith 4140971522cSBarry Smith #undef __FUNCT__ 4150971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 4160971522cSBarry Smith /*@ 4170971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 4180971522cSBarry Smith 4190971522cSBarry Smith Collective on PC 4200971522cSBarry Smith 4210971522cSBarry Smith Input Parameters: 4220971522cSBarry Smith + pc - the preconditioner context 4230971522cSBarry Smith . n - the number of fields in this split 4240971522cSBarry Smith . fields - the fields in this split 4250971522cSBarry Smith 4260971522cSBarry Smith Level: intermediate 4270971522cSBarry Smith 428*51f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetBlockSize() 4290971522cSBarry Smith 4300971522cSBarry Smith @*/ 431dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 4320971522cSBarry Smith { 4330971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 4340971522cSBarry Smith 4350971522cSBarry Smith PetscFunctionBegin; 4360971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4370971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 4380971522cSBarry Smith if (f) { 4390971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 4400971522cSBarry Smith } 4410971522cSBarry Smith PetscFunctionReturn(0); 4420971522cSBarry Smith } 4430971522cSBarry Smith 4440971522cSBarry Smith #undef __FUNCT__ 445*51f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize" 446*51f519a2SBarry Smith /*@ 447*51f519a2SBarry Smith PCFieldSplitSetBlockSize - Sets the block size for defining where fields start in the 448*51f519a2SBarry Smith fieldsplit preconditioner. If not set the matrix block size is used. 449*51f519a2SBarry Smith 450*51f519a2SBarry Smith Collective on PC 451*51f519a2SBarry Smith 452*51f519a2SBarry Smith Input Parameters: 453*51f519a2SBarry Smith + pc - the preconditioner context 454*51f519a2SBarry Smith - bs - the block size 455*51f519a2SBarry Smith 456*51f519a2SBarry Smith Level: intermediate 457*51f519a2SBarry Smith 458*51f519a2SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT, PCFieldSplitSetFields() 459*51f519a2SBarry Smith 460*51f519a2SBarry Smith @*/ 461*51f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize(PC pc,PetscInt bs) 462*51f519a2SBarry Smith { 463*51f519a2SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt); 464*51f519a2SBarry Smith 465*51f519a2SBarry Smith PetscFunctionBegin; 466*51f519a2SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 467*51f519a2SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 468*51f519a2SBarry Smith if (f) { 469*51f519a2SBarry Smith ierr = (*f)(pc,bs);CHKERRQ(ierr); 470*51f519a2SBarry Smith } 471*51f519a2SBarry Smith PetscFunctionReturn(0); 472*51f519a2SBarry Smith } 473*51f519a2SBarry Smith 474*51f519a2SBarry Smith #undef __FUNCT__ 47569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 4760971522cSBarry Smith /*@C 47769a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 4780971522cSBarry Smith 47969a612a9SBarry Smith Collective on KSP 4800971522cSBarry Smith 4810971522cSBarry Smith Input Parameter: 4820971522cSBarry Smith . pc - the preconditioner context 4830971522cSBarry Smith 4840971522cSBarry Smith Output Parameters: 4850971522cSBarry Smith + n - the number of split 48669a612a9SBarry Smith - pc - the array of KSP contexts 4870971522cSBarry Smith 4880971522cSBarry Smith Note: 48969a612a9SBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed 4900971522cSBarry Smith 49169a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 4920971522cSBarry Smith 4930971522cSBarry Smith Level: advanced 4940971522cSBarry Smith 4950971522cSBarry Smith .seealso: PCFIELDSPLIT 4960971522cSBarry Smith @*/ 497dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 4980971522cSBarry Smith { 49969a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 5000971522cSBarry Smith 5010971522cSBarry Smith PetscFunctionBegin; 5020971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5030971522cSBarry Smith PetscValidIntPointer(n,2); 50469a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 5050971522cSBarry Smith if (f) { 50669a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 5070971522cSBarry Smith } else { 50869a612a9SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 5090971522cSBarry Smith } 5100971522cSBarry Smith PetscFunctionReturn(0); 5110971522cSBarry Smith } 5120971522cSBarry Smith 51379416396SBarry Smith EXTERN_C_BEGIN 51479416396SBarry Smith #undef __FUNCT__ 51579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 516dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 51779416396SBarry Smith { 51879416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 51979416396SBarry Smith 52079416396SBarry Smith PetscFunctionBegin; 52179416396SBarry Smith jac->type = type; 52279416396SBarry Smith PetscFunctionReturn(0); 52379416396SBarry Smith } 52479416396SBarry Smith EXTERN_C_END 52579416396SBarry Smith 526*51f519a2SBarry Smith EXTERN_C_BEGIN 527*51f519a2SBarry Smith #undef __FUNCT__ 528*51f519a2SBarry Smith #define __FUNCT__ "PCFieldSplitSetBlockSize_FieldSplit" 529*51f519a2SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetBlockSize_FieldSplit(PC pc,PetscInt bs) 530*51f519a2SBarry Smith { 531*51f519a2SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 532*51f519a2SBarry Smith 533*51f519a2SBarry Smith PetscFunctionBegin; 534*51f519a2SBarry Smith if (bs < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Blocksize must be positive, you gave %D",bs); 535*51f519a2SBarry 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); 536*51f519a2SBarry Smith jac->bs = bs; 537*51f519a2SBarry Smith PetscFunctionReturn(0); 538*51f519a2SBarry Smith } 539*51f519a2SBarry Smith EXTERN_C_END 540*51f519a2SBarry Smith 54179416396SBarry Smith #undef __FUNCT__ 54279416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 54379416396SBarry Smith /*@C 54479416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 54579416396SBarry Smith 54679416396SBarry Smith Collective on PC 54779416396SBarry Smith 54879416396SBarry Smith Input Parameter: 54979416396SBarry Smith . pc - the preconditioner context 55079416396SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE 55179416396SBarry Smith 55279416396SBarry Smith Options Database Key: 55379416396SBarry Smith . -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type 55479416396SBarry Smith 55579416396SBarry Smith Level: Developer 55679416396SBarry Smith 55779416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 55879416396SBarry Smith 55979416396SBarry Smith .seealso: PCCompositeSetType() 56079416396SBarry Smith 56179416396SBarry Smith @*/ 562dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 56379416396SBarry Smith { 56479416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 56579416396SBarry Smith 56679416396SBarry Smith PetscFunctionBegin; 56779416396SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 56879416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 56979416396SBarry Smith if (f) { 57079416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 57179416396SBarry Smith } 57279416396SBarry Smith PetscFunctionReturn(0); 57379416396SBarry Smith } 57479416396SBarry Smith 5750971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 5760971522cSBarry Smith /*MC 577a8c7a070SBarry Smith PCFIELDSPLIT - Preconditioner created by combining separate preconditioners for individual 5780971522cSBarry Smith fields or groups of fields 5790971522cSBarry Smith 5800971522cSBarry Smith 5810971522cSBarry Smith To set options on the solvers for each block append -sub_ to all the PC 582d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 5830971522cSBarry Smith 584a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCFieldSplitGetSubKSP() 58569a612a9SBarry Smith and set the options directly on the resulting KSP object 5860971522cSBarry Smith 5870971522cSBarry Smith Level: intermediate 5880971522cSBarry Smith 58979416396SBarry Smith Options Database Keys: 5902e70eadcSBarry Smith + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 5912e70eadcSBarry Smith . -pc_splitfield_default - automatically add any fields to additional splits that have not 5922e70eadcSBarry Smith been supplied explicitly by -pc_splitfield_%d_fields 593*51f519a2SBarry Smith . -pc_splitfield_block_size <bs> - size of block that defines fields (i.e. there are bs fields) 5942e70eadcSBarry Smith - -pc_splitfield_type <additive,multiplicative> 59579416396SBarry Smith 5960971522cSBarry Smith Concepts: physics based preconditioners 5970971522cSBarry Smith 5980971522cSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 599c8d321e0SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType() 6000971522cSBarry Smith M*/ 6010971522cSBarry Smith 6020971522cSBarry Smith EXTERN_C_BEGIN 6030971522cSBarry Smith #undef __FUNCT__ 6040971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 605dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 6060971522cSBarry Smith { 6070971522cSBarry Smith PetscErrorCode ierr; 6080971522cSBarry Smith PC_FieldSplit *jac; 6090971522cSBarry Smith 6100971522cSBarry Smith PetscFunctionBegin; 6110971522cSBarry Smith ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr); 61252e6d16bSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));CHKERRQ(ierr); 6130971522cSBarry Smith jac->bs = -1; 6140971522cSBarry Smith jac->nsplits = 0; 61579416396SBarry Smith jac->type = PC_COMPOSITE_ADDITIVE; 616*51f519a2SBarry Smith 6170971522cSBarry Smith pc->data = (void*)jac; 6180971522cSBarry Smith 6190971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 6200971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 6210971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 6220971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 6230971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 6240971522cSBarry Smith pc->ops->applyrichardson = 0; 6250971522cSBarry Smith 62669a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 62769a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 6280971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 6290971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 63079416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 63179416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 632*51f519a2SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetBlockSize_C","PCFieldSplitSetBlockSize_FieldSplit", 633*51f519a2SBarry Smith PCFieldSplitSetBlockSize_FieldSplit);CHKERRQ(ierr); 6340971522cSBarry Smith PetscFunctionReturn(0); 6350971522cSBarry Smith } 6360971522cSBarry Smith EXTERN_C_END 6370971522cSBarry Smith 6380971522cSBarry Smith 639