14b9ad928SBarry Smith /* 24b9ad928SBarry Smith Routines to set PC methods and options. 34b9ad928SBarry Smith */ 44b9ad928SBarry Smith 54b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 64b9ad928SBarry Smith #include "petscsys.h" 74b9ad928SBarry Smith 84b9ad928SBarry Smith PetscTruth PCRegisterAllCalled = PETSC_FALSE; 94b9ad928SBarry Smith /* 104b9ad928SBarry Smith Contains the list of registered KSP routines 114b9ad928SBarry Smith */ 124b9ad928SBarry Smith PetscFList PCList = 0; 134b9ad928SBarry Smith 144b9ad928SBarry Smith #undef __FUNCT__ 154b9ad928SBarry Smith #define __FUNCT__ "PCSetType" 164b9ad928SBarry Smith /*@C 174b9ad928SBarry Smith PCSetType - Builds PC for a particular preconditioner. 184b9ad928SBarry Smith 194b9ad928SBarry Smith Collective on PC 204b9ad928SBarry Smith 214b9ad928SBarry Smith Input Parameter: 224b9ad928SBarry Smith + pc - the preconditioner context. 234b9ad928SBarry Smith - type - a known method 244b9ad928SBarry Smith 254b9ad928SBarry Smith Options Database Key: 264b9ad928SBarry Smith . -pc_type <type> - Sets PC type 274b9ad928SBarry Smith 284b9ad928SBarry Smith Use -help for a list of available methods (for instance, 294b9ad928SBarry Smith jacobi or bjacobi) 304b9ad928SBarry Smith 314b9ad928SBarry Smith Notes: 324b9ad928SBarry Smith See "petsc/include/petscpc.h" for available methods (for instance, 334b9ad928SBarry Smith PCJACOBI, PCILU, or PCBJACOBI). 344b9ad928SBarry Smith 354b9ad928SBarry Smith Normally, it is best to use the KSPSetFromOptions() command and 364b9ad928SBarry Smith then set the PC type from the options database rather than by using 374b9ad928SBarry Smith this routine. Using the options database provides the user with 384b9ad928SBarry Smith maximum flexibility in evaluating the many different preconditioners. 394b9ad928SBarry Smith The PCSetType() routine is provided for those situations where it 404b9ad928SBarry Smith is necessary to set the preconditioner independently of the command 414b9ad928SBarry Smith line or options database. This might be the case, for example, when 424b9ad928SBarry Smith the choice of preconditioner changes during the execution of the 434b9ad928SBarry Smith program, and the user's application is taking responsibility for 444b9ad928SBarry Smith choosing the appropriate preconditioner. In other words, this 454b9ad928SBarry Smith routine is not for beginners. 464b9ad928SBarry Smith 474b9ad928SBarry Smith Level: intermediate 484b9ad928SBarry Smith 494b9ad928SBarry Smith .keywords: PC, set, method, type 504b9ad928SBarry Smith 514b9ad928SBarry Smith .seealso: KSPSetType(), PCType 524b9ad928SBarry Smith 534b9ad928SBarry Smith @*/ 54dfbe8321SBarry Smith PetscErrorCode PCSetType(PC pc,const PCType type) 554b9ad928SBarry Smith { 56dfbe8321SBarry Smith PetscErrorCode ierr,(*r)(PC); 574b9ad928SBarry Smith PetscTruth match; 584b9ad928SBarry Smith 594b9ad928SBarry Smith PetscFunctionBegin; 604482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 614482741eSBarry Smith PetscValidCharPointer(type,2); 624b9ad928SBarry Smith 634b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)pc,type,&match);CHKERRQ(ierr); 644b9ad928SBarry Smith if (match) PetscFunctionReturn(0); 654b9ad928SBarry Smith 664b9ad928SBarry Smith if (pc->ops->destroy) {ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr);} 674b9ad928SBarry Smith ierr = PetscFListDestroy(&pc->qlist);CHKERRQ(ierr); 684b9ad928SBarry Smith pc->data = 0; 694b9ad928SBarry Smith pc->setupcalled = 0; 704b9ad928SBarry Smith 714b9ad928SBarry Smith /* Get the function pointers for the method requested */ 724b9ad928SBarry Smith if (!PCRegisterAllCalled) {ierr = PCRegisterAll(0);CHKERRQ(ierr);} 734b9ad928SBarry Smith 744b9ad928SBarry Smith /* Determine the PCCreateXXX routine for a particular preconditioner */ 754b9ad928SBarry Smith ierr = PetscFListFind(pc->comm,PCList,type,(void (**)(void)) &r);CHKERRQ(ierr); 76958c9bccSBarry Smith if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested PC type %s",type); 774b9ad928SBarry Smith if (pc->data) {ierr = PetscFree(pc->data);CHKERRQ(ierr);} 784b9ad928SBarry Smith 79*6849ba73SBarry Smith pc->ops->destroy = (PetscErrorCode (*)(PC)) 0; 80*6849ba73SBarry Smith pc->ops->view = (PetscErrorCode (*)(PC,PetscViewer)) 0; 81*6849ba73SBarry Smith pc->ops->apply = (PetscErrorCode (*)(PC,Vec,Vec)) 0; 82*6849ba73SBarry Smith pc->ops->setup = (PetscErrorCode (*)(PC)) 0; 83*6849ba73SBarry Smith pc->ops->applyrichardson = (PetscErrorCode (*)(PC,Vec,Vec,Vec,PetscReal,PetscReal,PetscReal,int)) 0; 84*6849ba73SBarry Smith pc->ops->applyBA = (PetscErrorCode (*)(PC,int,Vec,Vec,Vec)) 0; 85*6849ba73SBarry Smith pc->ops->setfromoptions = (PetscErrorCode (*)(PC)) 0; 86*6849ba73SBarry Smith pc->ops->applytranspose = (PetscErrorCode (*)(PC,Vec,Vec)) 0; 87*6849ba73SBarry Smith pc->ops->applyBAtranspose = (PetscErrorCode (*)(PC,int,Vec,Vec,Vec)) 0; 88*6849ba73SBarry Smith pc->ops->presolve = (PetscErrorCode (*)(PC,KSP,Vec,Vec)) 0; 89*6849ba73SBarry Smith pc->ops->postsolve = (PetscErrorCode (*)(PC,KSP,Vec,Vec)) 0; 90*6849ba73SBarry Smith pc->ops->getfactoredmatrix = (PetscErrorCode (*)(PC,Mat*)) 0; 91*6849ba73SBarry Smith pc->ops->applysymmetricleft = (PetscErrorCode (*)(PC,Vec,Vec)) 0; 92*6849ba73SBarry Smith pc->ops->applysymmetricright = (PetscErrorCode (*)(PC,Vec,Vec)) 0; 93*6849ba73SBarry Smith pc->ops->setuponblocks = (PetscErrorCode (*)(PC)) 0; 94*6849ba73SBarry Smith pc->modifysubmatrices = (PetscErrorCode (*)(PC,int,const IS[],const IS[],Mat[],void*)) 0; 954b9ad928SBarry Smith 964b9ad928SBarry Smith /* Call the PCCreateXXX routine for this particular preconditioner */ 974b9ad928SBarry Smith ierr = (*r)(pc);CHKERRQ(ierr); 984b9ad928SBarry Smith 994b9ad928SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)pc,type);CHKERRQ(ierr); 1004b9ad928SBarry Smith PetscFunctionReturn(0); 1014b9ad928SBarry Smith } 1024b9ad928SBarry Smith 1034b9ad928SBarry Smith #undef __FUNCT__ 1044b9ad928SBarry Smith #define __FUNCT__ "PCRegisterDestroy" 1054b9ad928SBarry Smith /*@C 1064b9ad928SBarry Smith PCRegisterDestroy - Frees the list of preconditioners that were 1074b9ad928SBarry Smith registered by PCRegisterDynamic(). 1084b9ad928SBarry Smith 1094b9ad928SBarry Smith Not Collective 1104b9ad928SBarry Smith 1114b9ad928SBarry Smith Level: advanced 1124b9ad928SBarry Smith 1134b9ad928SBarry Smith .keywords: PC, register, destroy 1144b9ad928SBarry Smith 1154b9ad928SBarry Smith .seealso: PCRegisterAll(), PCRegisterAll() 1164b9ad928SBarry Smith 1174b9ad928SBarry Smith @*/ 118dfbe8321SBarry Smith PetscErrorCode PCRegisterDestroy(void) 1194b9ad928SBarry Smith { 120dfbe8321SBarry Smith PetscErrorCode ierr; 1214b9ad928SBarry Smith 1224b9ad928SBarry Smith PetscFunctionBegin; 1234b9ad928SBarry Smith if (PCList) { 1244b9ad928SBarry Smith ierr = PetscFListDestroy(&PCList);CHKERRQ(ierr); 1254b9ad928SBarry Smith PCList = 0; 1264b9ad928SBarry Smith } 1274b9ad928SBarry Smith PCRegisterAllCalled = PETSC_FALSE; 1284b9ad928SBarry Smith PetscFunctionReturn(0); 1294b9ad928SBarry Smith } 1304b9ad928SBarry Smith 1314b9ad928SBarry Smith #undef __FUNCT__ 1324b9ad928SBarry Smith #define __FUNCT__ "PCGetType" 1334b9ad928SBarry Smith /*@C 1344b9ad928SBarry Smith PCGetType - Gets the PC method type and name (as a string) from the PC 1354b9ad928SBarry Smith context. 1364b9ad928SBarry Smith 1374b9ad928SBarry Smith Not Collective 1384b9ad928SBarry Smith 1394b9ad928SBarry Smith Input Parameter: 1404b9ad928SBarry Smith . pc - the preconditioner context 1414b9ad928SBarry Smith 1424b9ad928SBarry Smith Output Parameter: 1434b9ad928SBarry Smith . name - name of preconditioner 1444b9ad928SBarry Smith 1454b9ad928SBarry Smith Level: intermediate 1464b9ad928SBarry Smith 1474b9ad928SBarry Smith .keywords: PC, get, method, name, type 1484b9ad928SBarry Smith 1494b9ad928SBarry Smith .seealso: PCSetType() 1504b9ad928SBarry Smith 1514b9ad928SBarry Smith @*/ 152dfbe8321SBarry Smith PetscErrorCode PCGetType(PC pc,PCType *meth) 1534b9ad928SBarry Smith { 1544b9ad928SBarry Smith PetscFunctionBegin; 1554b9ad928SBarry Smith *meth = (PCType) pc->type_name; 1564b9ad928SBarry Smith PetscFunctionReturn(0); 1574b9ad928SBarry Smith } 1584b9ad928SBarry Smith 159dfbe8321SBarry Smith EXTERN PetscErrorCode PCGetDefaultType_Private(PC,const char*[]); 160566e8bf2SBarry Smith 1614b9ad928SBarry Smith #undef __FUNCT__ 1624b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions" 1634b9ad928SBarry Smith /*@ 1644b9ad928SBarry Smith PCSetFromOptions - Sets PC options from the options database. 1654b9ad928SBarry Smith This routine must be called before PCSetUp() if the user is to be 1664b9ad928SBarry Smith allowed to set the preconditioner method. 1674b9ad928SBarry Smith 1684b9ad928SBarry Smith Collective on PC 1694b9ad928SBarry Smith 1704b9ad928SBarry Smith Input Parameter: 1714b9ad928SBarry Smith . pc - the preconditioner context 1724b9ad928SBarry Smith 1734b9ad928SBarry Smith Level: developer 1744b9ad928SBarry Smith 1754b9ad928SBarry Smith .keywords: PC, set, from, options, database 1764b9ad928SBarry Smith 1774b9ad928SBarry Smith .seealso: 1784b9ad928SBarry Smith 1794b9ad928SBarry Smith @*/ 180dfbe8321SBarry Smith PetscErrorCode PCSetFromOptions(PC pc) 1814b9ad928SBarry Smith { 182dfbe8321SBarry Smith PetscErrorCode ierr; 1832fc52814SBarry Smith char type[256]; 1842fc52814SBarry Smith const char *def; 1854b9ad928SBarry Smith PetscTruth flg; 1864b9ad928SBarry Smith 1874b9ad928SBarry Smith PetscFunctionBegin; 1884482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1894b9ad928SBarry Smith 1904b9ad928SBarry Smith if (!PCRegisterAllCalled) {ierr = PCRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1914b9ad928SBarry Smith ierr = PetscOptionsBegin(pc->comm,pc->prefix,"Preconditioner (PC) Options","PC");CHKERRQ(ierr); 1924b9ad928SBarry Smith if (!pc->type_name) { 193566e8bf2SBarry Smith ierr = PCGetDefaultType_Private(pc,&def);CHKERRQ(ierr); 1944b9ad928SBarry Smith } else { 1954b9ad928SBarry Smith def = pc->type_name; 1964b9ad928SBarry Smith } 1974b9ad928SBarry Smith 1984b9ad928SBarry Smith ierr = PetscOptionsList("-pc_type","Preconditioner","PCSetType",PCList,def,type,256,&flg);CHKERRQ(ierr); 1994b9ad928SBarry Smith if (flg) { 2004b9ad928SBarry Smith ierr = PCSetType(pc,type);CHKERRQ(ierr); 201566e8bf2SBarry Smith } else if (!pc->type_name){ 202566e8bf2SBarry Smith ierr = PCSetType(pc,def);CHKERRQ(ierr); 2034b9ad928SBarry Smith } 2044b9ad928SBarry Smith 2054b9ad928SBarry Smith if (pc->ops->setfromoptions) { 2064b9ad928SBarry Smith ierr = (*pc->ops->setfromoptions)(pc);CHKERRQ(ierr); 2074b9ad928SBarry Smith } 2084b9ad928SBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 2094b9ad928SBarry Smith #if defined(__cplusplus) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && defined(PETSC_HAVE_CXX_NAMESPACE) 2104b9ad928SBarry Smith ierr = PCESISetFromOptions(pc);CHKERRQ(ierr); 2114b9ad928SBarry Smith #endif 2124b9ad928SBarry Smith PetscFunctionReturn(0); 2134b9ad928SBarry Smith } 214