10971522cSBarry Smith /* 20971522cSBarry Smith 30971522cSBarry Smith */ 40971522cSBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 50971522cSBarry Smith 60971522cSBarry Smith typedef struct _PC_FieldSplitLink *PC_FieldSplitLink; 70971522cSBarry Smith struct _PC_FieldSplitLink { 80971522cSBarry Smith PC pc; 90971522cSBarry Smith Vec x,y; 100971522cSBarry Smith PetscInt nfields; 110971522cSBarry Smith PetscInt *fields; 120971522cSBarry Smith PC_FieldSplitLink next; 130971522cSBarry Smith }; 140971522cSBarry Smith 150971522cSBarry Smith typedef struct { 160971522cSBarry Smith PetscInt bs; 170971522cSBarry Smith PetscInt nsplits; 180971522cSBarry Smith PC_FieldSplitLink head; 190971522cSBarry Smith Vec *x,*y; 200971522cSBarry Smith } PC_FieldSplit; 210971522cSBarry Smith 220971522cSBarry Smith #undef __FUNCT__ 230971522cSBarry Smith #define __FUNCT__ "PCView_FieldSplit" 240971522cSBarry Smith static PetscErrorCode PCView_FieldSplit(PC pc,PetscViewer viewer) 250971522cSBarry Smith { 260971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 270971522cSBarry Smith PetscErrorCode ierr; 280971522cSBarry Smith PetscTruth iascii; 290971522cSBarry Smith PetscInt i,j; 300971522cSBarry Smith PC_FieldSplitLink link = jac->head; 310971522cSBarry Smith 320971522cSBarry Smith PetscFunctionBegin; 330971522cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 340971522cSBarry Smith if (iascii) { 350971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," FieldSplit: total splits = %D",jac->nsplits);CHKERRQ(ierr); 360971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Solver info for each split is in the following PC objects:\n");CHKERRQ(ierr); 370971522cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 380971522cSBarry Smith for (i=0; i<jac->nsplits; i++) { 390971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Split number %D Fields ",i);CHKERRQ(ierr); 400971522cSBarry Smith for (j=0; j<link->nfields; j++) { 410971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"%D \n",link->fields[j]);CHKERRQ(ierr); 420971522cSBarry Smith } 430971522cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 440971522cSBarry Smith ierr = PCView(link->pc,viewer);CHKERRQ(ierr); 450971522cSBarry Smith link = link->next; 460971522cSBarry Smith } 470971522cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 480971522cSBarry Smith } else { 490971522cSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCFieldSplit",((PetscObject)viewer)->type_name); 500971522cSBarry Smith } 510971522cSBarry Smith PetscFunctionReturn(0); 520971522cSBarry Smith } 530971522cSBarry Smith 540971522cSBarry Smith #undef __FUNCT__ 550971522cSBarry Smith #define __FUNCT__ "PCSetUp_FieldSplit" 560971522cSBarry Smith static PetscErrorCode PCSetUp_FieldSplit(PC pc) 570971522cSBarry Smith { 580971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 590971522cSBarry Smith PetscErrorCode ierr; 600971522cSBarry Smith PC_FieldSplitLink link = jac->head; 61521d7252SBarry Smith PetscInt i; 620971522cSBarry Smith 630971522cSBarry Smith PetscFunctionBegin; 64521d7252SBarry Smith if (jac->bs <= 0) { 65521d7252SBarry Smith ierr = MatGetBlockSize(pc->pmat,&jac->bs);CHKERRQ(ierr); 66521d7252SBarry Smith } 67521d7252SBarry Smith if (!link) { /* user has not split fields so use default */ 68521d7252SBarry Smith for (i=0; i<jac->bs; i++) { 69521d7252SBarry Smith ierr = PCFieldSplitSetFields(pc,1,&i);CHKERRQ(ierr); 70521d7252SBarry Smith } 71521d7252SBarry Smith link = jac->head; 72521d7252SBarry Smith } 730971522cSBarry Smith while (link) { 740971522cSBarry Smith ierr = PCSetUp(link->pc);CHKERRQ(ierr); 750971522cSBarry Smith link = link->next; 760971522cSBarry Smith } 770971522cSBarry Smith PetscFunctionReturn(0); 780971522cSBarry Smith } 790971522cSBarry Smith 800971522cSBarry Smith #undef __FUNCT__ 810971522cSBarry Smith #define __FUNCT__ "PCApply_FieldSplit" 820971522cSBarry Smith static PetscErrorCode PCApply_FieldSplit(PC pc,Vec x,Vec y) 830971522cSBarry Smith { 840971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 850971522cSBarry Smith PetscErrorCode ierr; 860971522cSBarry Smith PC_FieldSplitLink link = jac->head; 870971522cSBarry Smith 880971522cSBarry Smith PetscFunctionBegin; 890971522cSBarry Smith ierr = VecStrideGatherAll(x,jac->x,INSERT_VALUES);CHKERRQ(ierr); 900971522cSBarry Smith while (link) { 910971522cSBarry Smith ierr = PCApply(link->pc,link->x,link->y);CHKERRQ(ierr); 920971522cSBarry Smith link = link->next; 930971522cSBarry Smith } 940971522cSBarry Smith ierr = VecStrideScatterAll(jac->y,y,INSERT_VALUES);CHKERRQ(ierr); 950971522cSBarry Smith PetscFunctionReturn(0); 960971522cSBarry Smith } 970971522cSBarry Smith 980971522cSBarry Smith #undef __FUNCT__ 990971522cSBarry Smith #define __FUNCT__ "PCDestroy_FieldSplit" 1000971522cSBarry Smith static PetscErrorCode PCDestroy_FieldSplit(PC pc) 1010971522cSBarry Smith { 1020971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1030971522cSBarry Smith PetscErrorCode ierr; 1040971522cSBarry Smith PC_FieldSplitLink link = jac->head; 1050971522cSBarry Smith 1060971522cSBarry Smith PetscFunctionBegin; 1070971522cSBarry Smith while (link) { 1080971522cSBarry Smith ierr = PCDestroy(link->pc);CHKERRQ(ierr); 1090971522cSBarry Smith ierr = PetscFree2(link,link->fields);CHKERRQ(ierr); 1100971522cSBarry Smith } 1110971522cSBarry Smith if (jac->x) {ierr = PetscFree2(jac->x,jac->y);CHKERRQ(ierr);} 1120971522cSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 1130971522cSBarry Smith PetscFunctionReturn(0); 1140971522cSBarry Smith } 1150971522cSBarry Smith 1160971522cSBarry Smith #undef __FUNCT__ 1170971522cSBarry Smith #define __FUNCT__ "PCSetFromOptions_FieldSplit" 1180971522cSBarry Smith static PetscErrorCode PCSetFromOptions_FieldSplit(PC pc) 1190971522cSBarry Smith { 1200971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1210971522cSBarry Smith PetscErrorCode ierr; 1220971522cSBarry Smith PC_FieldSplitLink link = jac->head; 1230971522cSBarry Smith 1240971522cSBarry Smith PetscFunctionBegin; 1250971522cSBarry Smith if (!link) { /* user never set fields, so set them ourselves */ 1260971522cSBarry Smith link = jac->head; 1270971522cSBarry Smith } 1280971522cSBarry Smith while (link) { 1290971522cSBarry Smith ierr = PCApply(link->pc,link->x,link->y);CHKERRQ(ierr); 1300971522cSBarry Smith link = link->next; 1310971522cSBarry Smith } 1320971522cSBarry Smith PetscFunctionReturn(0); 1330971522cSBarry Smith } 1340971522cSBarry Smith 1350971522cSBarry Smith /*------------------------------------------------------------------------------------*/ 1360971522cSBarry Smith 1370971522cSBarry Smith EXTERN_C_BEGIN 1380971522cSBarry Smith #undef __FUNCT__ 1390971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields_FieldSplit" 1400971522cSBarry Smith PetscErrorCode PCFieldSplitSetFields_FieldSplit(PC pc,PetscInt n,PetscInt *fields) 1410971522cSBarry Smith { 1420971522cSBarry Smith PC_FieldSplit *jac; 1430971522cSBarry Smith PetscErrorCode ierr; 1440971522cSBarry Smith PetscInt *myfields; 1450971522cSBarry Smith PC_FieldSplitLink link,next = jac->head; 1460971522cSBarry Smith 1470971522cSBarry Smith PetscFunctionBegin; 1480971522cSBarry Smith if (n <= 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Negative number of fields requested"); 1490971522cSBarry Smith ierr = PetscMalloc2(1,struct _PC_FieldSplitLink,&link,n,PetscInt,&myfields);CHKERRQ(ierr); 1500971522cSBarry Smith ierr = PetscMemcpy(myfields,fields,n*sizeof(PetscInt));CHKERRQ(ierr); 1510971522cSBarry Smith link->fields = myfields; 1520971522cSBarry Smith link->nfields = n; 1530971522cSBarry Smith 1540971522cSBarry Smith if (!next) { 1550971522cSBarry Smith jac->head = link; 1560971522cSBarry Smith } else { 1570971522cSBarry Smith while (next->next) { 1580971522cSBarry Smith next = next->next; 1590971522cSBarry Smith } 1600971522cSBarry Smith next->next = link; 1610971522cSBarry Smith } 1620971522cSBarry Smith jac->nsplits++; 1630971522cSBarry Smith PetscFunctionReturn(0); 1640971522cSBarry Smith } 1650971522cSBarry Smith EXTERN_C_END 1660971522cSBarry Smith 1670971522cSBarry Smith 1680971522cSBarry Smith EXTERN_C_BEGIN 1690971522cSBarry Smith #undef __FUNCT__ 1700971522cSBarry Smith #define __FUNCT__ "PCFieldSplitGetSubPC_FieldSplit" 1710971522cSBarry Smith PetscErrorCode PCFieldSplitGetSubPC_FieldSplit(PC pc,PetscInt *n,PC **subpc) 1720971522cSBarry Smith { 1730971522cSBarry Smith PC_FieldSplit *jac = (PC_FieldSplit*)pc->data; 1740971522cSBarry Smith PetscErrorCode ierr; 1750971522cSBarry Smith PetscInt cnt = 0; 1760971522cSBarry Smith PC_FieldSplitLink link; 1770971522cSBarry Smith 1780971522cSBarry Smith PetscFunctionBegin; 1790971522cSBarry Smith ierr = PetscMalloc(jac->nsplits*sizeof(PC*),subpc);CHKERRQ(ierr); 1800971522cSBarry Smith while (link) { 1810971522cSBarry Smith (*subpc)[cnt++] = link->pc; 1820971522cSBarry Smith link = link->next; 1830971522cSBarry Smith } 1840971522cSBarry Smith if (cnt != jac->nsplits) SETERRQ2(PETSC_ERR_PLIB,"Corrupt PCFIELDSPLIT object: number splits in linked list %D in object %D",cnt,jac->nsplits); 1850971522cSBarry Smith *n = jac->nsplits; 1860971522cSBarry Smith PetscFunctionReturn(0); 1870971522cSBarry Smith } 1880971522cSBarry Smith EXTERN_C_END 1890971522cSBarry Smith 1900971522cSBarry Smith #undef __FUNCT__ 1910971522cSBarry Smith #define __FUNCT__ "PCFieldSplitSetFields" 1920971522cSBarry Smith /*@ 1930971522cSBarry Smith PCFieldSplitSetFields - Sets the fields for one particular split in the field split preconditioner 1940971522cSBarry Smith 1950971522cSBarry Smith Collective on PC 1960971522cSBarry Smith 1970971522cSBarry Smith Input Parameters: 1980971522cSBarry Smith + pc - the preconditioner context 1990971522cSBarry Smith . n - the number of fields in this split 2000971522cSBarry Smith . fields - the fields in this split 2010971522cSBarry Smith 2020971522cSBarry Smith Level: intermediate 2030971522cSBarry Smith 204*3d30b1ffSBarry Smith .seealso: PCFieldSplitGetSubPC(), PCFIELDSPLIT 2050971522cSBarry Smith 2060971522cSBarry Smith @*/ 2070971522cSBarry Smith PetscErrorCode PCFieldSplitSetFields(PC pc,PetscInt n, PetscInt *fields) 2080971522cSBarry Smith { 2090971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt *); 2100971522cSBarry Smith 2110971522cSBarry Smith PetscFunctionBegin; 2120971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2130971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitSetFields_C",(void (**)(void))&f);CHKERRQ(ierr); 2140971522cSBarry Smith if (f) { 2150971522cSBarry Smith ierr = (*f)(pc,n,fields);CHKERRQ(ierr); 2160971522cSBarry Smith } 2170971522cSBarry Smith PetscFunctionReturn(0); 2180971522cSBarry Smith } 2190971522cSBarry Smith 2200971522cSBarry Smith #undef __FUNCT__ 2210971522cSBarry Smith #define __FUNCT__ "PCFieldSplitGetSubPC" 2220971522cSBarry Smith /*@C 2230971522cSBarry Smith PCFieldSplitGetSubPC - Gets the PC contexts for all splits 2240971522cSBarry Smith 2250971522cSBarry Smith Collective on PC 2260971522cSBarry Smith 2270971522cSBarry Smith Input Parameter: 2280971522cSBarry Smith . pc - the preconditioner context 2290971522cSBarry Smith 2300971522cSBarry Smith Output Parameters: 2310971522cSBarry Smith + n - the number of split 2320971522cSBarry Smith - pc - the array of PC contexts 2330971522cSBarry Smith 2340971522cSBarry Smith Note: 2350971522cSBarry Smith After PCFieldSplitGetSubPC() the array of PCs IS to be freed 2360971522cSBarry Smith 2370971522cSBarry Smith You must call KSPSetUp() before calling PCFieldSplitGetSubPC(). 2380971522cSBarry Smith 2390971522cSBarry Smith Level: advanced 2400971522cSBarry Smith 2410971522cSBarry Smith .keywords: PC, 2420971522cSBarry Smith 2430971522cSBarry Smith .seealso: PCFIELDSPLIT 2440971522cSBarry Smith @*/ 2450971522cSBarry Smith PetscErrorCode PCFieldSplitGetSubPC(PC pc,PetscInt *n,PC *subpc[]) 2460971522cSBarry Smith { 2470971522cSBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt*,PC **); 2480971522cSBarry Smith 2490971522cSBarry Smith PetscFunctionBegin; 2500971522cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2510971522cSBarry Smith PetscValidIntPointer(n,2); 2520971522cSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFieldSplitGetSubPC_C",(void (**)(void))&f);CHKERRQ(ierr); 2530971522cSBarry Smith if (f) { 2540971522cSBarry Smith ierr = (*f)(pc,n,subpc);CHKERRQ(ierr); 2550971522cSBarry Smith } else { 2560971522cSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get subpc for this type of PC"); 2570971522cSBarry Smith } 2580971522cSBarry Smith PetscFunctionReturn(0); 2590971522cSBarry Smith } 2600971522cSBarry Smith 2610971522cSBarry Smith /* -------------------------------------------------------------------------------------*/ 2620971522cSBarry Smith /*MC 2630971522cSBarry Smith PCFieldSplit - Preconditioner created by combining seperate preconditioners for individual 2640971522cSBarry Smith fields or groups of fields 2650971522cSBarry Smith 2660971522cSBarry Smith 2670971522cSBarry Smith To set options on the solvers for each block append -sub_ to all the PC 2680971522cSBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_ilu_levels 1 2690971522cSBarry Smith 2700971522cSBarry Smith To set the options on the solvers seperate for each block call PCFieldSplitGetSubPC() 2710971522cSBarry Smith and set the options directly on the resulting PC object 2720971522cSBarry Smith 2730971522cSBarry Smith Level: intermediate 2740971522cSBarry Smith 2750971522cSBarry Smith Concepts: physics based preconditioners 2760971522cSBarry Smith 2770971522cSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 2780971522cSBarry Smith PCFieldSplitGetSubPC(), PCFieldSplitSetFields() 2790971522cSBarry Smith M*/ 2800971522cSBarry Smith 2810971522cSBarry Smith EXTERN_C_BEGIN 2820971522cSBarry Smith #undef __FUNCT__ 2830971522cSBarry Smith #define __FUNCT__ "PCCreate_FieldSplit" 2840971522cSBarry Smith PetscErrorCode PCCreate_FieldSplit(PC pc) 2850971522cSBarry Smith { 2860971522cSBarry Smith PetscErrorCode ierr; 2870971522cSBarry Smith PC_FieldSplit *jac; 2880971522cSBarry Smith 2890971522cSBarry Smith PetscFunctionBegin; 2900971522cSBarry Smith ierr = PetscNew(PC_FieldSplit,&jac);CHKERRQ(ierr); 2910971522cSBarry Smith PetscLogObjectMemory(pc,sizeof(PC_FieldSplit)); 2920971522cSBarry Smith ierr = PetscMemzero(jac,sizeof(PC_FieldSplit));CHKERRQ(ierr); 2930971522cSBarry Smith jac->bs = -1; 2940971522cSBarry Smith jac->nsplits = 0; 2950971522cSBarry Smith pc->data = (void*)jac; 2960971522cSBarry Smith 2970971522cSBarry Smith pc->ops->apply = PCApply_FieldSplit; 2980971522cSBarry Smith pc->ops->setup = PCSetUp_FieldSplit; 2990971522cSBarry Smith pc->ops->destroy = PCDestroy_FieldSplit; 3000971522cSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_FieldSplit; 3010971522cSBarry Smith pc->ops->view = PCView_FieldSplit; 3020971522cSBarry Smith pc->ops->applyrichardson = 0; 3030971522cSBarry Smith 3040971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitGetSubPC_C","PCFieldSplitGetSubPC_FieldSplit", 3050971522cSBarry Smith PCFieldSplitGetSubPC_FieldSplit);CHKERRQ(ierr); 3060971522cSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFieldSplitSetFields_C","PCFieldSplitSetFields_FieldSplit", 3070971522cSBarry Smith PCFieldSplitSetFields_FieldSplit);CHKERRQ(ierr); 3080971522cSBarry Smith PetscFunctionReturn(0); 3090971522cSBarry Smith } 3100971522cSBarry Smith EXTERN_C_END 3110971522cSBarry Smith 3120971522cSBarry Smith 313