1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris 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; 380971522cSBarry Smith 390971522cSBarry Smith PetscFunctionBegin; 400971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 410971522cSBarry Smith if (iascii) { 429dcbbd2bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit with %s composition: total splits = %D",PCCompositeTypes[jac->type],jac->nsplits);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; 7779416396SBarry Smith ierr = PetscOptionsGetLogical(pc->prefix,"-pc_fieldsplit_default",&flg,PETSC_NULL);CHKERRQ(ierr); 785a9f2f41SSatish Balay if (!ilink || flg) { 7963ba0a88SBarry Smith ierr = PetscLogInfo((pc,"PCFieldSplitSetDefaults: 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; 11169a612a9SBarry Smith PetscInt i,nsplit; 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 1231b9fc7fcSBarry Smith ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr); 12497bbdb24SBarry Smith ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 1251b9fc7fcSBarry Smith nslots = (rend - rstart)/bs; 1261b9fc7fcSBarry Smith ierr = PetscMalloc(nsplit*sizeof(IS),&jac->is);CHKERRQ(ierr); 1271b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 1281b9fc7fcSBarry Smith if (jac->defaultsplit) { 1291b9fc7fcSBarry Smith ierr = ISCreateStride(pc->comm,nslots,rstart+i,nsplit,&jac->is[i]);CHKERRQ(ierr); 13097bbdb24SBarry Smith } else { 1315a9f2f41SSatish Balay PetscInt *ii,j,k,nfields = ilink->nfields,*fields = ilink->fields; 1321b9fc7fcSBarry Smith PetscTruth sorted; 1335a9f2f41SSatish Balay ierr = PetscMalloc(ilink->nfields*nslots*sizeof(PetscInt),&ii);CHKERRQ(ierr); 1341b9fc7fcSBarry Smith for (j=0; j<nslots; j++) { 1351b9fc7fcSBarry Smith for (k=0; k<nfields; k++) { 1361b9fc7fcSBarry Smith ii[nfields*j + k] = rstart + bs*j + fields[k]; 13797bbdb24SBarry Smith } 13897bbdb24SBarry Smith } 1391b9fc7fcSBarry Smith ierr = ISCreateGeneral(pc->comm,nslots*nfields,ii,&jac->is[i]);CHKERRQ(ierr); 1401b9fc7fcSBarry Smith ierr = ISSorted(jac->is[i],&sorted);CHKERRQ(ierr); 1411b9fc7fcSBarry Smith if (!sorted) SETERRQ(PETSC_ERR_USER,"Fields must be sorted when creating split"); 1421b9fc7fcSBarry Smith ierr = PetscFree(ii);CHKERRQ(ierr); 1435a9f2f41SSatish Balay ilink = ilink->next; 1441b9fc7fcSBarry Smith } 1451b9fc7fcSBarry Smith } 1461b9fc7fcSBarry Smith } 1471b9fc7fcSBarry Smith 14897bbdb24SBarry Smith if (!jac->pmat) { 14997bbdb24SBarry Smith ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_INITIAL_MATRIX,&jac->pmat);CHKERRQ(ierr); 15097bbdb24SBarry Smith } else { 15197bbdb24SBarry Smith ierr = MatGetSubMatrices(pc->pmat,nsplit,jac->is,jac->is,MAT_REUSE_MATRIX,&jac->pmat);CHKERRQ(ierr); 15297bbdb24SBarry Smith } 15397bbdb24SBarry Smith 15497bbdb24SBarry Smith /* set up the individual PCs */ 15597bbdb24SBarry Smith i = 0; 1565a9f2f41SSatish Balay ilink = jac->head; 1575a9f2f41SSatish Balay while (ilink) { 1585a9f2f41SSatish Balay ierr = KSPSetOperators(ilink->ksp,jac->pmat[i],jac->pmat[i],flag);CHKERRQ(ierr); 1595a9f2f41SSatish Balay ierr = KSPSetFromOptions(ilink->ksp);CHKERRQ(ierr); 1605a9f2f41SSatish Balay ierr = KSPSetUp(ilink->ksp);CHKERRQ(ierr); 16197bbdb24SBarry Smith i++; 1625a9f2f41SSatish Balay ilink = ilink->next; 1630971522cSBarry Smith } 16497bbdb24SBarry Smith 16597bbdb24SBarry Smith /* create work vectors for each split */ 1661b9fc7fcSBarry Smith if (!jac->x) { 16779416396SBarry Smith Vec xtmp; 16897bbdb24SBarry Smith ierr = PetscMalloc2(nsplit,Vec,&jac->x,nsplit,Vec,&jac->y);CHKERRQ(ierr); 1695a9f2f41SSatish Balay ilink = jac->head; 17097bbdb24SBarry Smith for (i=0; i<nsplit; i++) { 17197bbdb24SBarry Smith Mat A; 1725a9f2f41SSatish Balay ierr = KSPGetOperators(ilink->ksp,PETSC_NULL,&A,PETSC_NULL);CHKERRQ(ierr); 1735a9f2f41SSatish Balay ierr = MatGetVecs(A,&ilink->x,&ilink->y);CHKERRQ(ierr); 1745a9f2f41SSatish Balay jac->x[i] = ilink->x; 1755a9f2f41SSatish Balay jac->y[i] = ilink->y; 1765a9f2f41SSatish Balay ilink = ilink->next; 17797bbdb24SBarry Smith } 17879416396SBarry Smith /* compute scatter contexts needed by multiplicative versions and non-default splits */ 1791b9fc7fcSBarry Smith 1805a9f2f41SSatish Balay ilink = jac->head; 1811b9fc7fcSBarry Smith ierr = MatGetVecs(pc->pmat,&xtmp,PETSC_NULL);CHKERRQ(ierr); 1821b9fc7fcSBarry Smith for (i=0; i<nsplit; i++) { 1835a9f2f41SSatish Balay ierr = VecScatterCreate(xtmp,jac->is[i],jac->x[i],PETSC_NULL,&ilink->sctx);CHKERRQ(ierr); 1845a9f2f41SSatish Balay ilink = ilink->next; 1851b9fc7fcSBarry Smith } 1861b9fc7fcSBarry Smith ierr = VecDestroy(xtmp);CHKERRQ(ierr); 1871b9fc7fcSBarry Smith } 1880971522cSBarry Smith PetscFunctionReturn(0); 1890971522cSBarry Smith } 1900971522cSBarry Smith 1915a9f2f41SSatish Balay #define FieldSplitSplitSolveAdd(ilink,xx,yy) \ 1925a9f2f41SSatish Balay (VecScatterBegin(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \ 1935a9f2f41SSatish Balay VecScatterEnd(xx,ilink->x,INSERT_VALUES,SCATTER_FORWARD,ilink->sctx) || \ 1945a9f2f41SSatish Balay KSPSolve(ilink->ksp,ilink->x,ilink->y) || \ 1955a9f2f41SSatish Balay VecScatterBegin(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx) || \ 1965a9f2f41SSatish Balay VecScatterEnd(ilink->y,yy,ADD_VALUES,SCATTER_REVERSE,ilink->sctx)) 19779416396SBarry Smith 1980971522cSBarry Smith #undef __FUNCT__ 1990971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 2000971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 2010971522cSBarry Smith { 2020971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2030971522cSBarry Smith PetscErrorCode ierr; 2045a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 20579416396SBarry Smith PetscScalar zero = 0.0,mone = -1.0; 2060971522cSBarry Smith 2070971522cSBarry Smith PetscFunctionBegin; 20879416396SBarry Smith if (jac->type == PC_COMPOSITE_ADDITIVE) { 2091b9fc7fcSBarry Smith if (jac->defaultsplit) { 2100971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 2115a9f2f41SSatish Balay while (ilink) { 2125a9f2f41SSatish Balay ierr = KSPSolve(ilink->ksp,ilink->x,ilink->y);CHKERRQ(ierr); 2135a9f2f41SSatish Balay ilink = ilink->next; 2140971522cSBarry Smith } 2150971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 2161b9fc7fcSBarry Smith } else { 2171b9fc7fcSBarry Smith PetscInt i = 0; 2181b9fc7fcSBarry Smith 219*2dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 2205a9f2f41SSatish Balay while (ilink) { 2215a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 2225a9f2f41SSatish Balay ilink = ilink->next; 2231b9fc7fcSBarry Smith i++; 2241b9fc7fcSBarry Smith } 2251b9fc7fcSBarry Smith } 22679416396SBarry Smith } else { 22779416396SBarry Smith if (!jac->w1) { 22879416396SBarry Smith ierr = VecDuplicate(x,&jac->w1);CHKERRQ(ierr); 22979416396SBarry Smith ierr = VecDuplicate(x,&jac->w2);CHKERRQ(ierr); 23079416396SBarry Smith } 231*2dcb1b2aSMatthew Knepley ierr = VecSet(y,zero);CHKERRQ(ierr); 2325a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,x,y);CHKERRQ(ierr); 2335a9f2f41SSatish Balay while (ilink->next) { 2345a9f2f41SSatish Balay ilink = ilink->next; 23579416396SBarry Smith ierr = MatMult(pc->pmat,y,jac->w1);CHKERRQ(ierr); 236*2dcb1b2aSMatthew Knepley ierr = VecWAXPY(jac->w2,mone,jac->w1,x);CHKERRQ(ierr); 2375a9f2f41SSatish Balay ierr = FieldSplitSplitSolveAdd(ilink,jac->w2,y);CHKERRQ(ierr); 23879416396SBarry Smith } 23979416396SBarry Smith } 2400971522cSBarry Smith PetscFunctionReturn(0); 2410971522cSBarry Smith } 2420971522cSBarry Smith 2430971522cSBarry Smith #undef __FUNCT__ 2440971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 2450971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 2460971522cSBarry Smith { 2470971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2480971522cSBarry Smith PetscErrorCode ierr; 2495a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head,next; 2500971522cSBarry Smith 2510971522cSBarry Smith PetscFunctionBegin; 2525a9f2f41SSatish Balay while (ilink) { 2535a9f2f41SSatish Balay ierr = KSPDestroy(ilink->ksp);CHKERRQ(ierr); 2545a9f2f41SSatish Balay if (ilink->x) {ierr = VecDestroy(ilink->x);CHKERRQ(ierr);} 2555a9f2f41SSatish Balay if (ilink->y) {ierr = VecDestroy(ilink->y);CHKERRQ(ierr);} 2565a9f2f41SSatish Balay if (ilink->sctx) {ierr = VecScatterDestroy(ilink->sctx);CHKERRQ(ierr);} 2575a9f2f41SSatish Balay next = ilink->next; 2585a9f2f41SSatish Balay ierr = PetscFree2(ilink,ilink->fields);CHKERRQ(ierr); 2595a9f2f41SSatish Balay ilink = next; 2600971522cSBarry Smith } 2610971522cSBarry Smith if (jac->x) {ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);} 26297bbdb24SBarry Smith if (jac->pmat) {ierr = MatDestroyMatrices(jac->nsplits,&jac->pmat);CHKERRQ(ierr);} 26397bbdb24SBarry Smith if (jac->is) { 26497bbdb24SBarry Smith PetscInt i; 26597bbdb24SBarry Smith for (i=0; i<jac->nsplits; i++) {ierr = ISDestroy(jac->is[i]);CHKERRQ(ierr);} 26697bbdb24SBarry Smith ierr = PetscFree(jac->is);CHKERRQ(ierr); 26797bbdb24SBarry Smith } 26879416396SBarry Smith if (jac->w1) {ierr = VecDestroy(jac->w1);CHKERRQ(ierr);} 26979416396SBarry Smith if (jac->w2) {ierr = VecDestroy(jac->w2);CHKERRQ(ierr);} 2700971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 2710971522cSBarry Smith PetscFunctionReturn(0); 2720971522cSBarry Smith } 2730971522cSBarry Smith 2740971522cSBarry Smith #undef __FUNCT__ 2750971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 2760971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 2771b9fc7fcSBarry Smith /* This does not call KSPSetFromOptions() on the subksp's, see PCSetFromOptionsBJacobi/ASM() */ 2780971522cSBarry Smith { 2791b9fc7fcSBarry Smith PetscErrorCode ierr; 2809dcbbd2bSBarry Smith PetscInt i = 0,nfields,fields[12]; 2811b9fc7fcSBarry Smith PetscTruth flg; 2821b9fc7fcSBarry Smith char optionname[128]; 2839dcbbd2bSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 2841b9fc7fcSBarry Smith 2850971522cSBarry Smith PetscFunctionBegin; 2861b9fc7fcSBarry Smith ierr = PetscOptionsHead("FieldSplit options");CHKERRQ(ierr); 2879dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_fieldsplit_type","Type of composition","PCFieldSplitSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 2881b9fc7fcSBarry Smith while (PETSC_TRUE) { 28913f74950SBarry Smith sprintf(optionname,"-pc_fieldsplit_%d_fields",(int)i++); 2901b9fc7fcSBarry Smith nfields = 12; 2911b9fc7fcSBarry Smith ierr = PetscOptionsIntArray(optionname,"Fields in this split","PCFieldSplitSetFields",fields,&nfields,&flg);CHKERRQ(ierr); 2921b9fc7fcSBarry Smith if (!flg) break; 2931b9fc7fcSBarry Smith if (!nfields) SETERRQ(PETSC_ERR_USER,"Cannot list zero fields"); 2941b9fc7fcSBarry Smith ierr = PCFieldSplitSetFields(pc,nfields,fields);CHKERRQ(ierr); 2951b9fc7fcSBarry Smith } 2961b9fc7fcSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 2970971522cSBarry Smith PetscFunctionReturn(0); 2980971522cSBarry Smith } 2990971522cSBarry Smith 3000971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 3010971522cSBarry Smith 3020971522cSBarry Smith EXTERN_C_BEGIN 3030971522cSBarry Smith #undef __FUNCT__ 3040971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 305dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 3060971522cSBarry Smith { 30797bbdb24SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3080971522cSBarry Smith PetscErrorCode ierr; 3095a9f2f41SSatish Balay PC_FieldSplitLink ilink,next = jac->head; 31069a612a9SBarry Smith char prefix[128]; 3110971522cSBarry Smith 3120971522cSBarry Smith PetscFunctionBegin; 3130971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 3145a9f2f41SSatish Balay ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&ilink,n,PetscInt,&ilink->fields);CHKERRQ(ierr); 3155a9f2f41SSatish Balay ierr = PetscMemcpy(ilink->fields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 3165a9f2f41SSatish Balay ilink->nfields = n; 3175a9f2f41SSatish Balay ilink->next = PETSC_NULL; 3185a9f2f41SSatish Balay ierr = KSPCreate(pc->comm,&ilink->ksp);CHKERRQ(ierr); 3195a9f2f41SSatish Balay ierr = KSPSetType(ilink->ksp,KSPPREONLY);CHKERRQ(ierr); 32069a612a9SBarry Smith 32169a612a9SBarry Smith if (pc->prefix) { 32213f74950SBarry Smith sprintf(prefix,"%sfieldsplit_%d_",pc->prefix,(int)jac->nsplits); 32369a612a9SBarry Smith } else { 32413f74950SBarry Smith sprintf(prefix,"fieldsplit_%d_",(int)jac->nsplits); 32569a612a9SBarry Smith } 3265a9f2f41SSatish Balay ierr = KSPSetOptionsPrefix(ilink->ksp,prefix);CHKERRQ(ierr); 3270971522cSBarry Smith 3280971522cSBarry Smith if (!next) { 3295a9f2f41SSatish Balay jac->head = ilink; 3300971522cSBarry Smith } else { 3310971522cSBarry Smith while (next->next) { 3320971522cSBarry Smith next = next->next; 3330971522cSBarry Smith } 3345a9f2f41SSatish Balay next->next = ilink; 3350971522cSBarry Smith } 3360971522cSBarry Smith jac->nsplits++; 3370971522cSBarry Smith PetscFunctionReturn(0); 3380971522cSBarry Smith } 3390971522cSBarry Smith EXTERN_C_END 3400971522cSBarry Smith 3410971522cSBarry Smith 3420971522cSBarry Smith EXTERN_C_BEGIN 3430971522cSBarry Smith #undef __FUNCT__ 34469a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP_FieldSplit" 345dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP_FieldSplit(PC pc,PetscInt *n,KSP **subksp) 3460971522cSBarry Smith { 3470971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 3480971522cSBarry Smith PetscErrorCode ierr; 3490971522cSBarry Smith PetscInt cnt = 0; 3505a9f2f41SSatish Balay PC_FieldSplitLink ilink = jac->head; 3510971522cSBarry Smith 3520971522cSBarry Smith PetscFunctionBegin; 35369a612a9SBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(KSP*),subksp);CHKERRQ(ierr); 3545a9f2f41SSatish Balay while (ilink) { 3555a9f2f41SSatish Balay (*subksp)[cnt++] = ilink->ksp; 3565a9f2f41SSatish Balay ilink = ilink->next; 3570971522cSBarry Smith } 3580971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 3590971522cSBarry Smith *n = jac->nsplits; 3600971522cSBarry Smith PetscFunctionReturn(0); 3610971522cSBarry Smith } 3620971522cSBarry Smith EXTERN_C_END 3630971522cSBarry Smith 3640971522cSBarry Smith #undef __FUNCT__ 3650971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 3660971522cSBarry Smith /*@ 3670971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 3680971522cSBarry Smith 3690971522cSBarry Smith Collective on PC 3700971522cSBarry Smith 3710971522cSBarry Smith Input Parameters: 3720971522cSBarry Smith + pc - the preconditioner context 3730971522cSBarry Smith . n - the number of fields in this split 3740971522cSBarry Smith . fields - the fields in this split 3750971522cSBarry Smith 3760971522cSBarry Smith Level: intermediate 3770971522cSBarry Smith 37869a612a9SBarry Smith .seealso: PCFieldSplitGetSubKSP(), PCFIELDSPLIT 3790971522cSBarry Smith 3800971522cSBarry Smith @*/ 381dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 3820971522cSBarry Smith { 3830971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 3840971522cSBarry Smith 3850971522cSBarry Smith PetscFunctionBegin; 3860971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3870971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 3880971522cSBarry Smith if (f) { 3890971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 3900971522cSBarry Smith } 3910971522cSBarry Smith PetscFunctionReturn(0); 3920971522cSBarry Smith } 3930971522cSBarry Smith 3940971522cSBarry Smith #undef __FUNCT__ 39569a612a9SBarry Smith #define __FUNCT__ "PCFieldSplitGetSubKSP" 3960971522cSBarry Smith /*@C 39769a612a9SBarry Smith PCFieldSplitGetSubKSP - Gets the KSP contexts for all splits 3980971522cSBarry Smith 39969a612a9SBarry Smith Collective on KSP 4000971522cSBarry Smith 4010971522cSBarry Smith Input Parameter: 4020971522cSBarry Smith . pc - the preconditioner context 4030971522cSBarry Smith 4040971522cSBarry Smith Output Parameters: 4050971522cSBarry Smith + n - the number of split 40669a612a9SBarry Smith - pc - the array of KSP contexts 4070971522cSBarry Smith 4080971522cSBarry Smith Note: 40969a612a9SBarry Smith After PCFieldSplitGetSubKSP() the array of KSPs IS to be freed 4100971522cSBarry Smith 41169a612a9SBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubKSP(). 4120971522cSBarry Smith 4130971522cSBarry Smith Level: advanced 4140971522cSBarry Smith 4150971522cSBarry Smith .seealso: PCFIELDSPLIT 4160971522cSBarry Smith @*/ 417dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitGetSubKSP(PC pc,PetscInt *n,KSP *subksp[]) 4180971522cSBarry Smith { 41969a612a9SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,KSP **); 4200971522cSBarry Smith 4210971522cSBarry Smith PetscFunctionBegin; 4220971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4230971522cSBarry Smith PetscValidIntPointer(n,2); 42469a612a9SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 4250971522cSBarry Smith if (f) { 42669a612a9SBarry Smith ierr = (*f)(pc,n,subksp);CHKERRQ(ierr); 4270971522cSBarry Smith } else { 42869a612a9SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subksp for this type of PC"); 4290971522cSBarry Smith } 4300971522cSBarry Smith PetscFunctionReturn(0); 4310971522cSBarry Smith } 4320971522cSBarry Smith 43379416396SBarry Smith EXTERN_C_BEGIN 43479416396SBarry Smith #undef __FUNCT__ 43579416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType_FieldSplit" 436dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType_FieldSplit(PC pc,PCCompositeType type) 43779416396SBarry Smith { 43879416396SBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 43979416396SBarry Smith 44079416396SBarry Smith PetscFunctionBegin; 44179416396SBarry Smith jac->type = type; 44279416396SBarry Smith PetscFunctionReturn(0); 44379416396SBarry Smith } 44479416396SBarry Smith EXTERN_C_END 44579416396SBarry Smith 44679416396SBarry Smith #undef __FUNCT__ 44779416396SBarry Smith #define __FUNCT__ "PCFieldSplitSetType" 44879416396SBarry Smith /*@C 44979416396SBarry Smith PCFieldSplitSetType - Sets the type of fieldsplit preconditioner. 45079416396SBarry Smith 45179416396SBarry Smith Collective on PC 45279416396SBarry Smith 45379416396SBarry Smith Input Parameter: 45479416396SBarry Smith . pc - the preconditioner context 45579416396SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE 45679416396SBarry Smith 45779416396SBarry Smith Options Database Key: 45879416396SBarry Smith . -pc_fieldsplit_type <type: one of multiplicative, additive, special> - Sets fieldsplit preconditioner type 45979416396SBarry Smith 46079416396SBarry Smith Level: Developer 46179416396SBarry Smith 46279416396SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 46379416396SBarry Smith 46479416396SBarry Smith .seealso: PCCompositeSetType() 46579416396SBarry Smith 46679416396SBarry Smith @*/ 467dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFieldSplitSetType(PC pc,PCCompositeType type) 46879416396SBarry Smith { 46979416396SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 47079416396SBarry Smith 47179416396SBarry Smith PetscFunctionBegin; 47279416396SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 47379416396SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 47479416396SBarry Smith if (f) { 47579416396SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 47679416396SBarry Smith } 47779416396SBarry Smith PetscFunctionReturn(0); 47879416396SBarry Smith } 47979416396SBarry Smith 4800971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 4810971522cSBarry Smith /*MC 4827cdd61b2SBarry Smith PCFIELDSPLIT - Preconditioner created by combining seperate preconditioners for individual 4830971522cSBarry Smith fields or groups of fields 4840971522cSBarry Smith 4850971522cSBarry Smith 4860971522cSBarry Smith To set options on the solvers for each block append -sub_ to all the PC 4870971522cSBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_ilu_levels 1 4880971522cSBarry Smith 48969a612a9SBarry Smith To set the options on the solvers seperate for each block call PCFieldSplitGetSubKSP() 49069a612a9SBarry Smith and set the options directly on the resulting KSP object 4910971522cSBarry Smith 4920971522cSBarry Smith Level: intermediate 4930971522cSBarry Smith 49479416396SBarry Smith Options Database Keys: 4952e70eadcSBarry Smith + -pc_splitfield_%d_fields <a,b,..> - indicates the fields to be used in the %d'th split 4962e70eadcSBarry Smith . -pc_splitfield_default - automatically add any fields to additional splits that have not 4972e70eadcSBarry Smith been supplied explicitly by -pc_splitfield_%d_fields 4982e70eadcSBarry Smith - -pc_splitfield_type <additive,multiplicative> 49979416396SBarry Smith 5000971522cSBarry Smith Concepts: physics based preconditioners 5010971522cSBarry Smith 5020971522cSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 503c8d321e0SBarry Smith PCFieldSplitGetSubKSP(), PCFieldSplitSetFields(),PCFieldSplitSetType() 5040971522cSBarry Smith M*/ 5050971522cSBarry Smith 5060971522cSBarry Smith EXTERN_C_BEGIN 5070971522cSBarry Smith #undef __FUNCT__ 5080971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 509dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_FieldSplit(PC pc) 5100971522cSBarry Smith { 5110971522cSBarry Smith PetscErrorCode ierr; 5120971522cSBarry Smith PC_FieldSplit *jac; 5130971522cSBarry Smith 5140971522cSBarry Smith PetscFunctionBegin; 5150971522cSBarry Smith ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr); 51652e6d16bSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PC_FieldSplit));CHKERRQ(ierr); 5170971522cSBarry Smith jac->bs = -1; 5180971522cSBarry Smith jac->nsplits = 0; 51979416396SBarry Smith jac->type = PC_COMPOSITE_ADDITIVE; 5200971522cSBarry Smith pc->data = (void*)jac; 5210971522cSBarry Smith 5220971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 5230971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 5240971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 5250971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 5260971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 5270971522cSBarry Smith pc->ops->applyrichardson = 0; 5280971522cSBarry Smith 52969a612a9SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubKSP_C","PCFieldSplitGetSubKSP_FieldSplit", 53069a612a9SBarry Smith PCFieldSplitGetSubKSP_FieldSplit);CHKERRQ(ierr); 5310971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 5320971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 53379416396SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetType_C","PCFieldSplitSetType_FieldSplit", 53479416396SBarry Smith PCFieldSplitSetType_FieldSplit);CHKERRQ(ierr); 5350971522cSBarry Smith PetscFunctionReturn(0); 5360971522cSBarry Smith } 5370971522cSBarry Smith EXTERN_C_END 5380971522cSBarry Smith 5390971522cSBarry Smith 540