1*dba47a55SKris Buschelman #define PETSCKSP_DLL 2*dba47a55SKris Buschelman 30971522cSBarry Smith /* 40971522cSBarry Smith 50971522cSBarry Smith */ 60971522cSBarry Smith #include "src/ksp/pc/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; 150971522cSBarry Smith PC_FieldSplitLink next; 160971522cSBarry Smith }; 170971522cSBarry Smith 180971522cSBarry Smith typedef struct { 1979416396SBarry Smith PCCompositeType type; /* additive or multiplicative */ 2097bbdb24SBarry Smith PetscTruth defaultsplit; 210971522cSBarry Smith PetscInt bs; 220971522cSBarry Smith PetscInt nsplits; 2379416396SBarry Smith Vec *x,*y,w1,w2; 2497bbdb24SBarry Smith Mat *pmat; 2597bbdb24SBarry Smith IS *is; 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; 3879416396SBarry Smith const char *types[] = {"additive","multiplicative"}; 390971522cSBarry Smith 400971522cSBarry Smith PetscFunctionBegin; 410971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 420971522cSBarry Smith if (iascii) { 4379416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D",types[jac->type],jac->nsplits);CHKERRQ(ierr); 4469a612a9SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following KSP objects:\n");CHKERRQ(ierr); 450971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 460971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 470971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 4879416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 495a9f2f41SSatish Balay for (j=0; j<ilink->nfields; j++) { 5079416396SBarry Smith if (j > 0) { 5179416396SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,",");CHKERRQ(ierr); 5279416396SBarry Smith } 535a9f2f41SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," %D",ilink->fields[j]);CHKERRQ(ierr); 540971522cSBarry Smith } 550971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 5679416396SBarry Smith ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 575a9f2f41SSatish Balay ierr = KSPView(ilink->ksp,viewer);CHKERRQ(ierr); 585a9f2f41SSatish Balay ilink = ilink->next; 590971522cSBarry Smith } 600971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 610971522cSBarry Smith } else { 620971522cSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 630971522cSBarry Smith } 640971522cSBarry Smith PetscFunctionReturn(0); 650971522cSBarry Smith } 660971522cSBarry Smith 670971522cSBarry Smith #undef __FUNCT__ 6869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitSetDefaults" 6969a612a9SBarry Smith static PetscErrorCode PCFieldSplitSetDefaults(PC pc) 700971522cSBarry Smith { 710971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 720971522cSBarry Smith PetscErrorCode ierr; 735a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 7469a612a9SBarry Smith PetscInt i; 7579416396SBarry Smith PetscTruth flg = PETSC_FALSE,*fields; 760971522cSBarry Smith 770971522cSBarry Smith PetscFunctionBegin; 7879416396SBarry Smith ierr = PetscOptionsGetLogical(pc->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 795a9f2f41SSatish Balay if (!ilink || flg) { 8063ba0a88SBarry Smith ierr = PetscLogInfo((pc,"PCFieldSplitSetDefaults: Using default splitting of fields\n"));CHKERRQ(ierr); 81521d7252SBarry Smith if (jac->bs <= 0) { 82521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 83521d7252SBarry Smith } 8479416396SBarry Smith ierr = PetscMalloc(jac->bs*sizeof(PetscTruth),&fields);CHKERRQ(ierr); 8579416396SBarry Smith ierr = PetscMemzero(fields,jac->bs*sizeof(PetscTruth));CHKERRQ(ierr); 865a9f2f41SSatish Balay while (ilink) { 875a9f2f41SSatish Balay for (i=0; i<ilink->nfields; i++) { 885a9f2f41SSatish Balay fields[ilink->fields[i]] = PETSC_TRUE; 89521d7252SBarry Smith } 905a9f2f41SSatish Balay ilink = ilink->next; 9179416396SBarry Smith } 9297bbdb24SBarry Smith jac->defaultsplit = PETSC_TRUE; 9379416396SBarry Smith for (i=0; i<jac->bs; i++) { 9479416396SBarry Smith if (!fields[i]) { 9579416396SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 9679416396SBarry Smith } else { 9779416396SBarry Smith jac->defaultsplit = PETSC_FALSE; 9879416396SBarry Smith } 9979416396SBarry Smith } 100521d7252SBarry Smith } 10169a612a9SBarry Smith PetscFunctionReturn(0); 10269a612a9SBarry Smith } 10369a612a9SBarry Smith 10469a612a9SBarry Smith 10569a612a9SBarry Smith #undef __FUNCT__ 10669a612a9SBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 10769a612a9SBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 10869a612a9SBarry Smith { 10969a612a9SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 11069a612a9SBarry Smith PetscErrorCode ierr; 1115a9f2f41SSatish Balay PC_FieldSplitLink ilink; 11269a612a9SBarry Smith PetscInt i,nsplit; 11369a612a9SBarry Smith MatStructure flag = pc->flag; 11469a612a9SBarry Smith 11569a612a9SBarry Smith PetscFunctionBegin; 11669a612a9SBarry Smith ierr = PCFieldSplitSetDefaults(pc);CHKERRQ(ierr); 11797bbdb24SBarry Smith nsplit = jac->nsplits; 1185a9f2f41SSatish Balay ilink = jac->head; 11997bbdb24SBarry Smith 12097bbdb24SBarry Smith /* get the matrices for each split */ 12197bbdb24SBarry Smith if (!jac->is) { 1221b9fc7fcSBarry Smith PetscInt rstart,rend,nslots,bs; 12397bbdb24SBarry Smith 1241b9fc7fcSBarry Smith ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr); 12597bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 1261b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 1271b9fc7fcSBarry Smith ierr = PetscMalloc(nsplit*sizeof(IS),&jac->is);CHKERRQ(ierr); 1281b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 1291b9fc7fcSBarry Smith if (jac->defaultsplit) { 1301b9fc7fcSBarry Smith ierr = ISCreateStride(pc->comm,nslots,rstart+i,nsplit,&jac->is[i]);CHKERRQ(ierr); 13197bbdb24SBarry Smith } else { 1325a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 1331b9fc7fcSBarry Smith PetscTruth sorted; 1345a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 1351b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 1361b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 1371b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 13897bbdb24SBarry Smith } 13997bbdb24SBarry Smith } 1401b9fc7fcSBarry Smith ierr = ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);CHKERRQ(ierr); 1411b9fc7fcSBarry Smith ierr = ISSorted(jac->is[i],&sorted);CHKERRQ(ierr); 1421b9fc7fcSBarry Smith if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 1431b9fc7fcSBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 1445a9f2f41SSatish Balay ilink = ilink->next; 1451b9fc7fcSBarry Smith } 1461b9fc7fcSBarry Smith } 1471b9fc7fcSBarry Smith } 1481b9fc7fcSBarry Smith 14997bbdb24SBarry Smith if (!jac->pmat) { 15097bbdb24SBarry Smith ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_INITIAL_MATRIX,&jac->pmat);CHKERRQ(ierr); 15197bbdb24SBarry Smith } else { 15297bbdb24SBarry Smith ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_REUSE_MATRIX,&jac->pmat);CHKERRQ(ierr); 15397bbdb24SBarry Smith } 15497bbdb24SBarry Smith 15597bbdb24SBarry Smith /* set up the individual PCs */ 15697bbdb24SBarry Smith i = 0; 1575a9f2f41SSatish Balay ilink = jac->head; 1585a9f2f41SSatish Balay while (ilink) { 1595a9f2f41SSatish Balay ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 1605a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 1615a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 16297bbdb24SBarry Smith i++; 1635a9f2f41SSatish Balay ilink = ilink->next; 1640971522cSBarry Smith } 16597bbdb24SBarry Smith 16697bbdb24SBarry Smith /* create work vectors for each split */ 1671b9fc7fcSBarry Smith if (!jac->x) { 16879416396SBarry Smith Vec xtmp; 16997bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 1705a9f2f41SSatish Balay ilink = jac->head; 17197bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 17297bbdb24SBarry Smith Mat A; 1735a9f2f41SSatish Balay ierr = KSPGetOperators(ilink->ksp,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr); 1745a9f2f41SSatish Balay ierr = MatGetVecs(A,&ilink->x,&ilink->y);CHKERRQ(ierr); 1755a9f2f41SSatish Balay jac->x[i] = ilink->x; 1765a9f2f41SSatish Balay jac->y[i] = ilink->y; 1775a9f2f41SSatish Balay ilink = ilink->next; 17897bbdb24SBarry Smith } 17979416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 1801b9fc7fcSBarry Smith 1815a9f2f41SSatish Balay ilink = jac->head; 1821b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 1831b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 1845a9f2f41SSatish Balay ierr = VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 1855a9f2f41SSatish Balay ilink = ilink->next; 1861b9fc7fcSBarry Smith } 1871b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 1881b9fc7fcSBarry Smith } 1890971522cSBarry Smith PetscFunctionReturn(0); 1900971522cSBarry Smith } 1910971522cSBarry Smith 1925a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 1935a9f2f41SSatish Balay (VecScatterBegin(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \ 1945a9f2f41SSatish Balay VecScatterEnd(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \ 1955a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 1965a9f2f41SSatish Balay VecScatterBegin(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \ 1975a9f2f41SSatish Balay VecScatterEnd(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx)) 19879416396SBarry Smith 1990971522cSBarry Smith #undef __FUNCT__ 2000971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 2010971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 2020971522cSBarry Smith { 2030971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2040971522cSBarry Smith PetscErrorCode ierr; 2055a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 20679416396SBarry Smith PetscScalar zero = 0.0,mone = -1.0; 2070971522cSBarry Smith 2080971522cSBarry Smith PetscFunctionBegin; 20979416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 2101b9fc7fcSBarry Smith if (jac->defaultsplit) { 2110971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 2125a9f2f41SSatish Balay while (ilink) { 2135a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 2145a9f2f41SSatish Balay ilink = ilink->next; 2150971522cSBarry Smith } 2160971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 2171b9fc7fcSBarry Smith } else { 2181b9fc7fcSBarry Smith PetscInt i = 0; 2191b9fc7fcSBarry Smith 2201b9fc7fcSBarry Smith ierr = VecSet(&zero,y);CHKERRQ(ierr); 2215a9f2f41SSatish Balay while (ilink) { 2225a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 2235a9f2f41SSatish Balay ilink = ilink->next; 2241b9fc7fcSBarry Smith i++; 2251b9fc7fcSBarry Smith } 2261b9fc7fcSBarry Smith } 22779416396SBarry Smith } else { 22879416396SBarry Smith if (!jac->w1) { 22979416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 23079416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 23179416396SBarry Smith } 23279416396SBarry Smith ierr = VecSet(&zero,y);CHKERRQ(ierr); 2335a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 2345a9f2f41SSatish Balay while (ilink->next) { 2355a9f2f41SSatish Balay ilink = ilink->next; 23679416396SBarry Smith ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr); 23779416396SBarry Smith ierr = VecWAXPY(&mone,jac->w1,x,jac->w2);CHKERRQ(ierr); 2385a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 23979416396SBarry Smith } 24079416396SBarry Smith } 2410971522cSBarry Smith PetscFunctionReturn(0); 2420971522cSBarry Smith } 2430971522cSBarry Smith 2440971522cSBarry Smith #undef __FUNCT__ 2450971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 2460971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 2470971522cSBarry Smith { 2480971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2490971522cSBarry Smith PetscErrorCode ierr; 2505a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 2510971522cSBarry Smith 2520971522cSBarry Smith PetscFunctionBegin; 2535a9f2f41SSatish Balay while (ilink) { 2545a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 2555a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 2565a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 2575a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 2585a9f2f41SSatish Balay next = ilink->next; 2595a9f2f41SSatish Balay ierr = PetscFree2(ilink,ilink->fields);CHKERRQ(ierr); 2605a9f2f41SSatish Balay ilink = next; 2610971522cSBarry Smith } 2620971522cSBarry Smith if (jac->x) {ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);} 26397bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 26497bbdb24SBarry Smith if (jac->is) { 26597bbdb24SBarry Smith PetscInt i; 26697bbdb24SBarry Smith for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);} 26797bbdb24SBarry Smith ierr = PetscFree(jac->is);CHKERRQ(ierr); 26897bbdb24SBarry Smith } 26979416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 27079416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 2710971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 2720971522cSBarry Smith PetscFunctionReturn(0); 2730971522cSBarry Smith } 2740971522cSBarry Smith 2750971522cSBarry Smith #undef __FUNCT__ 2760971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 2770971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 2781b9fc7fcSBarry Smith /* This does not call KSPSetFromOptions() on the subksp's, see PCSetFromOptionsBJacobi/ASM() */ 2790971522cSBarry Smith { 2801b9fc7fcSBarry Smith PetscErrorCode ierr; 28179416396SBarry Smith PetscInt i = 0,nfields,fields[12],indx; 2821b9fc7fcSBarry Smith PetscTruth flg; 2831b9fc7fcSBarry Smith char optionname[128]; 28479416396SBarry Smith const char *types[] = {"additive","multiplicative"}; 2851b9fc7fcSBarry Smith 2860971522cSBarry Smith PetscFunctionBegin; 2871b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 28879416396SBarry Smith ierr = PetscOptionsEList("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",types,2,types[0],&indx,&flg);CHKERRQ(ierr); 28979416396SBarry Smith if (flg) { 29079416396SBarry Smith ierr = PCFieldSplitSetType(pc,(PCCompositeType)indx);CHKERRQ(ierr); 29179416396SBarry Smith } 2921b9fc7fcSBarry Smith while (PETSC_TRUE) { 29313f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 2941b9fc7fcSBarry Smith nfields = 12; 2951b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 2961b9fc7fcSBarry Smith if (!flg) break; 2971b9fc7fcSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 2981b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 2991b9fc7fcSBarry Smith } 3001b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3010971522cSBarry Smith PetscFunctionReturn(0); 3020971522cSBarry Smith } 3030971522cSBarry Smith 3040971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 3050971522cSBarry Smith 3060971522cSBarry Smith EXTERN_C_BEGIN 3070971522cSBarry Smith #undef __FUNCT__ 3080971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 309*dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 3100971522cSBarry Smith { 31197bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3120971522cSBarry Smith PetscErrorCode ierr; 3135a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 31469a612a9SBarry Smith char prefix[128]; 3150971522cSBarry Smith 3160971522cSBarry Smith PetscFunctionBegin; 3170971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 3185a9f2f41SSatish Balay ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);CHKERRQ(ierr); 3195a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 3205a9f2f41SSatish Balay ilink->nfields = n; 3215a9f2f41SSatish Balay ilink->next = PETSC_NULL; 3225a9f2f41SSatish Balay ierr = KSPCreate(pc->comm,&ilink->ksp);CHKERRQ(ierr); 3235a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 32469a612a9SBarry Smith 32569a612a9SBarry Smith if (pc->prefix) { 32613f74950SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits); 32769a612a9SBarry Smith } else { 32813f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 32969a612a9SBarry Smith } 3305a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 3310971522cSBarry Smith 3320971522cSBarry Smith if (!next) { 3335a9f2f41SSatish Balay jac->head = ilink; 3340971522cSBarry Smith } else { 3350971522cSBarry Smith while (next->next) { 3360971522cSBarry Smith next = next->next; 3370971522cSBarry Smith } 3385a9f2f41SSatish Balay next->next = ilink; 3390971522cSBarry Smith } 3400971522cSBarry Smith jac->nsplits++; 3410971522cSBarry Smith PetscFunctionReturn(0); 3420971522cSBarry Smith } 3430971522cSBarry Smith EXTERN_C_END 3440971522cSBarry Smith 3450971522cSBarry Smith 3460971522cSBarry Smith EXTERN_C_BEGIN 3470971522cSBarry Smith #undef __FUNCT__ 34869a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 349*dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 3500971522cSBarry Smith { 3510971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3520971522cSBarry Smith PetscErrorCode ierr; 3530971522cSBarry Smith PetscInt cnt = 0; 3545a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 3550971522cSBarry Smith 3560971522cSBarry Smith PetscFunctionBegin; 35769a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 3585a9f2f41SSatish Balay while (ilink) { 3595a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 3605a9f2f41SSatish Balay ilink = ilink->next; 3610971522cSBarry Smith } 3620971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 3630971522cSBarry Smith *n = jac->nsplits; 3640971522cSBarry Smith PetscFunctionReturn(0); 3650971522cSBarry Smith } 3660971522cSBarry Smith EXTERN_C_END 3670971522cSBarry Smith 3680971522cSBarry Smith #undef __FUNCT__ 3690971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 3700971522cSBarry Smith /*@ 3710971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 3720971522cSBarry Smith 3730971522cSBarry Smith Collective on PC 3740971522cSBarry Smith 3750971522cSBarry Smith Input Parameters: 3760971522cSBarry Smith + pc - the preconditioner context 3770971522cSBarry Smith . n - the number of fields in this split 3780971522cSBarry Smith . fields - the fields in this split 3790971522cSBarry Smith 3800971522cSBarry Smith Level: intermediate 3810971522cSBarry Smith 38269a612a9SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT 3830971522cSBarry Smith 3840971522cSBarry Smith @*/ 385*dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 3860971522cSBarry Smith { 3870971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 3880971522cSBarry Smith 3890971522cSBarry Smith PetscFunctionBegin; 3900971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3910971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 3920971522cSBarry Smith if (f) { 3930971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 3940971522cSBarry Smith } 3950971522cSBarry Smith PetscFunctionReturn(0); 3960971522cSBarry Smith } 3970971522cSBarry Smith 3980971522cSBarry Smith #undef __FUNCT__ 39969a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 4000971522cSBarry Smith /*@C 40169a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 4020971522cSBarry Smith 40369a612a9SBarry Smith Collective on KSP 4040971522cSBarry Smith 4050971522cSBarry Smith Input Parameter: 4060971522cSBarry Smith . pc - the preconditioner context 4070971522cSBarry Smith 4080971522cSBarry Smith Output Parameters: 4090971522cSBarry Smith + n - the number of split 41069a612a9SBarry Smith - pc - the array of KSP contexts 4110971522cSBarry Smith 4120971522cSBarry Smith Note: 41369a612a9SBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed 4140971522cSBarry Smith 41569a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 4160971522cSBarry Smith 4170971522cSBarry Smith Level: advanced 4180971522cSBarry Smith 4190971522cSBarry Smith .seealso: PCFIELDSPLIT 4200971522cSBarry Smith @*/ 421*dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 4220971522cSBarry Smith { 42369a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 4240971522cSBarry Smith 4250971522cSBarry Smith PetscFunctionBegin; 4260971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4270971522cSBarry Smith PetscValidIntPointer(n,2); 42869a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 4290971522cSBarry Smith if (f) { 43069a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 4310971522cSBarry Smith } else { 43269a612a9SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 4330971522cSBarry Smith } 4340971522cSBarry Smith PetscFunctionReturn(0); 4350971522cSBarry Smith } 4360971522cSBarry Smith 43779416396SBarry Smith EXTERN_C_BEGIN 43879416396SBarry Smith #undef __FUNCT__ 43979416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 440*dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 44179416396SBarry Smith { 44279416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 44379416396SBarry Smith 44479416396SBarry Smith PetscFunctionBegin; 44579416396SBarry Smith jac->type = type; 44679416396SBarry Smith PetscFunctionReturn(0); 44779416396SBarry Smith } 44879416396SBarry Smith EXTERN_C_END 44979416396SBarry Smith 45079416396SBarry Smith #undef __FUNCT__ 45179416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 45279416396SBarry Smith /*@C 45379416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 45479416396SBarry Smith 45579416396SBarry Smith Collective on PC 45679416396SBarry Smith 45779416396SBarry Smith Input Parameter: 45879416396SBarry Smith . pc - the preconditioner context 45979416396SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE 46079416396SBarry Smith 46179416396SBarry Smith Options Database Key: 46279416396SBarry Smith . -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type 46379416396SBarry Smith 46479416396SBarry Smith Level: Developer 46579416396SBarry Smith 46679416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 46779416396SBarry Smith 46879416396SBarry Smith .seealso: PCCompositeSetType() 46979416396SBarry Smith 47079416396SBarry Smith @*/ 471*dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 47279416396SBarry Smith { 47379416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 47479416396SBarry Smith 47579416396SBarry Smith PetscFunctionBegin; 47679416396SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 47779416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 47879416396SBarry Smith if (f) { 47979416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 48079416396SBarry Smith } 48179416396SBarry Smith PetscFunctionReturn(0); 48279416396SBarry Smith } 48379416396SBarry Smith 4840971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 4850971522cSBarry Smith /*MC 4860971522cSBarry Smith PCFieldSplit - Preconditioner created by combining seperate preconditioners for individual 4870971522cSBarry Smith fields or groups of fields 4880971522cSBarry Smith 4890971522cSBarry Smith 4900971522cSBarry Smith To set options on the solvers for each block append -sub_ to all the PC 4910971522cSBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_ilu_levels 1 4920971522cSBarry Smith 49369a612a9SBarry Smith To set the options on the solvers seperate for each block call PCFieldSplitGetSubKSP() 49469a612a9SBarry Smith and set the options directly on the resulting KSP object 4950971522cSBarry Smith 4960971522cSBarry Smith Level: intermediate 4970971522cSBarry Smith 49879416396SBarry Smith Options Database Keys: 4992e70eadcSBarry Smith + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 5002e70eadcSBarry Smith . -pc_splitfield_default - automatically add any fields to additional splits that have not 5012e70eadcSBarry Smith been supplied explicitly by -pc_splitfield_%d_fields 5022e70eadcSBarry Smith - -pc_splitfield_type <additive,multiplicative> 50379416396SBarry Smith 5040971522cSBarry Smith Concepts: physics based preconditioners 5050971522cSBarry Smith 5060971522cSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 507c8d321e0SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType() 5080971522cSBarry Smith M*/ 5090971522cSBarry Smith 5100971522cSBarry Smith EXTERN_C_BEGIN 5110971522cSBarry Smith #undef __FUNCT__ 5120971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 513*dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 5140971522cSBarry Smith { 5150971522cSBarry Smith PetscErrorCode ierr; 5160971522cSBarry Smith PC_FieldSplit *jac; 5170971522cSBarry Smith 5180971522cSBarry Smith PetscFunctionBegin; 5190971522cSBarry Smith ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr); 52052e6d16bSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));CHKERRQ(ierr); 5210971522cSBarry Smith jac->bs = -1; 5220971522cSBarry Smith jac->nsplits = 0; 52379416396SBarry Smith jac->type = PC_COMPOSITE_ADDITIVE; 5240971522cSBarry Smith pc->data = (void*)jac; 5250971522cSBarry Smith 5260971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 5270971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 5280971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 5290971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 5300971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 5310971522cSBarry Smith pc->ops->applyrichardson = 0; 5320971522cSBarry Smith 53369a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 53469a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 5350971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 5360971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 53779416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 53879416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 5390971522cSBarry Smith PetscFunctionReturn(0); 5400971522cSBarry Smith } 5410971522cSBarry Smith EXTERN_C_END 5420971522cSBarry Smith 5430971522cSBarry Smith 544