1*7e14e8a7SBarry Smith #ifdef PETSC_RCS_HEADER 2*7e14e8a7SBarry Smith static char vcid[] = "$Id: itcreate.c,v 1.162 1999/05/04 20:34:35 balay Exp $"; 3*7e14e8a7SBarry Smith #endif 4*7e14e8a7SBarry Smith /* 5*7e14e8a7SBarry Smith The basic KSP routines, Create, View etc. are here. 6*7e14e8a7SBarry Smith */ 7*7e14e8a7SBarry Smith #include "petsc.h" 8*7e14e8a7SBarry Smith #include "src/sles/ksp/kspimpl.h" /*I "ksp.h" I*/ 9*7e14e8a7SBarry Smith #include "sys.h" 10*7e14e8a7SBarry Smith #include "viewer.h" /*I "viewer.h" I*/ 11*7e14e8a7SBarry Smith 12*7e14e8a7SBarry Smith int KSPRegisterAllCalled = 0; 13*7e14e8a7SBarry Smith 14*7e14e8a7SBarry Smith #undef __FUNC__ 15*7e14e8a7SBarry Smith #define __FUNC__ "KSPView" 16*7e14e8a7SBarry Smith /*@ 17*7e14e8a7SBarry Smith KSPView - Prints the KSP data structure. 18*7e14e8a7SBarry Smith 19*7e14e8a7SBarry Smith Collective on KSP unless Viewer is VIEWER_STDOUT_SELF 20*7e14e8a7SBarry Smith 21*7e14e8a7SBarry Smith Input Parameters: 22*7e14e8a7SBarry Smith + ksp - the Krylov space context 23*7e14e8a7SBarry Smith - viewer - visualization context 24*7e14e8a7SBarry Smith 25*7e14e8a7SBarry Smith Note: 26*7e14e8a7SBarry Smith The available visualization contexts include 27*7e14e8a7SBarry Smith + VIEWER_STDOUT_SELF - standard output (default) 28*7e14e8a7SBarry Smith - VIEWER_STDOUT_WORLD - synchronized standard 29*7e14e8a7SBarry Smith output where only the first processor opens 30*7e14e8a7SBarry Smith the file. All other processors send their 31*7e14e8a7SBarry Smith data to the first processor to print. 32*7e14e8a7SBarry Smith 33*7e14e8a7SBarry Smith The user can open an alternative visualization context with 34*7e14e8a7SBarry Smith ViewerASCIIOpen() - output to a specified file. 35*7e14e8a7SBarry Smith 36*7e14e8a7SBarry Smith Level: developer 37*7e14e8a7SBarry Smith 38*7e14e8a7SBarry Smith .keywords: KSP, view 39*7e14e8a7SBarry Smith 40*7e14e8a7SBarry Smith .seealso: PCView(), ViewerASCIIOpen() 41*7e14e8a7SBarry Smith @*/ 42*7e14e8a7SBarry Smith int KSPView(KSP ksp,Viewer viewer) 43*7e14e8a7SBarry Smith { 44*7e14e8a7SBarry Smith char *method; 45*7e14e8a7SBarry Smith int ierr; 46*7e14e8a7SBarry Smith ViewerType vtype; 47*7e14e8a7SBarry Smith 48*7e14e8a7SBarry Smith PetscFunctionBegin; 49*7e14e8a7SBarry Smith ierr = ViewerGetType(viewer,&vtype);CHKERRQ(ierr); 50*7e14e8a7SBarry Smith if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 51*7e14e8a7SBarry Smith ierr = KSPGetType(ksp,&method);CHKERRQ(ierr); 52*7e14e8a7SBarry Smith ierr = ViewerASCIIPrintf(viewer,"KSP Object:\n");CHKERRQ(ierr); 53*7e14e8a7SBarry Smith ierr = ViewerASCIIPrintf(viewer," method: %s\n",method);CHKERRQ(ierr); 54*7e14e8a7SBarry Smith if (ksp->ops->view) { 55*7e14e8a7SBarry Smith ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 56*7e14e8a7SBarry Smith ierr = (*ksp->ops->view)(ksp,viewer);CHKERRQ(ierr); 57*7e14e8a7SBarry Smith ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 58*7e14e8a7SBarry Smith } 59*7e14e8a7SBarry Smith if (ksp->guess_zero) {ierr = ViewerASCIIPrintf(viewer," maximum iterations=%d, initial guess is zero\n",ksp->max_it);CHKERRQ(ierr);} 60*7e14e8a7SBarry Smith else {ierr = ViewerASCIIPrintf(viewer," maximum iterations=%d\n", ksp->max_it);CHKERRQ(ierr);} 61*7e14e8a7SBarry Smith ierr = ViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, divergence=%g\n",ksp->rtol, ksp->atol, ksp->divtol);CHKERRQ(ierr); 62*7e14e8a7SBarry Smith if (ksp->pc_side == PC_RIGHT) {ierr = ViewerASCIIPrintf(viewer," right preconditioning\n");CHKERRQ(ierr);} 63*7e14e8a7SBarry Smith else if (ksp->pc_side == PC_SYMMETRIC) {ierr = ViewerASCIIPrintf(viewer," symmetric preconditioning\n");CHKERRQ(ierr);} 64*7e14e8a7SBarry Smith else {ierr = ViewerASCIIPrintf(viewer," left preconditioning\n");CHKERRQ(ierr);} 65*7e14e8a7SBarry Smith } 66*7e14e8a7SBarry Smith PetscFunctionReturn(0); 67*7e14e8a7SBarry Smith } 68*7e14e8a7SBarry Smith 69*7e14e8a7SBarry Smith /* 70*7e14e8a7SBarry Smith Contains the list of registered KSP routines 71*7e14e8a7SBarry Smith */ 72*7e14e8a7SBarry Smith FList KSPList = 0; 73*7e14e8a7SBarry Smith 74*7e14e8a7SBarry Smith #undef __FUNC__ 75*7e14e8a7SBarry Smith #define __FUNC__ "KSPSetAvoidNorms" 76*7e14e8a7SBarry Smith /*@C 77*7e14e8a7SBarry Smith KSPSetAvoidNorms - Sets the KSP solver to avoid computing the residual norm 78*7e14e8a7SBarry Smith when possible. This, for example, reduces the number of collective operations 79*7e14e8a7SBarry Smith when using the Krylov method as a smoother. 80*7e14e8a7SBarry Smith 81*7e14e8a7SBarry Smith Collective on KSP 82*7e14e8a7SBarry Smith 83*7e14e8a7SBarry Smith Input Parameter: 84*7e14e8a7SBarry Smith . ksp - Krylov solver context 85*7e14e8a7SBarry Smith 86*7e14e8a7SBarry Smith Notes: 87*7e14e8a7SBarry Smith One cannot use the default convergence test routines when this option is 88*7e14e8a7SBarry Smith set, since these are based on decreases in the residual norms. Thus, this 89*7e14e8a7SBarry Smith option automatically switches to activate the KSPSkipConverged() test function. 90*7e14e8a7SBarry Smith 91*7e14e8a7SBarry Smith Currently only works with the CG, Richardson, Bi-CG-stab, CR, and CGS methods. 92*7e14e8a7SBarry Smith 93*7e14e8a7SBarry Smith Level: advanced 94*7e14e8a7SBarry Smith 95*7e14e8a7SBarry Smith .keywords: KSP, create, context, norms 96*7e14e8a7SBarry Smith 97*7e14e8a7SBarry Smith .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPSkipConverged() 98*7e14e8a7SBarry Smith @*/ 99*7e14e8a7SBarry Smith int KSPSetAvoidNorms(KSP ksp) 100*7e14e8a7SBarry Smith { 101*7e14e8a7SBarry Smith int ierr; 102*7e14e8a7SBarry Smith 103*7e14e8a7SBarry Smith PetscFunctionBegin; 104*7e14e8a7SBarry Smith PetscValidHeaderSpecific(ksp,KSP_COOKIE); 105*7e14e8a7SBarry Smith ksp->avoidnorms = PETSC_TRUE; 106*7e14e8a7SBarry Smith ierr = KSPSetConvergenceTest(ksp,KSPSkipConverged,PETSC_NULL);CHKERRQ(ierr); 107*7e14e8a7SBarry Smith PetscFunctionReturn(0); 108*7e14e8a7SBarry Smith } 109*7e14e8a7SBarry Smith 110*7e14e8a7SBarry Smith #undef __FUNC__ 111*7e14e8a7SBarry Smith #define __FUNC__ "KSPPublish_Petsc" 112*7e14e8a7SBarry Smith static int KSPPublish_Petsc(PetscObject object) 113*7e14e8a7SBarry Smith { 114*7e14e8a7SBarry Smith #if defined(HAVE_AMS) 115*7e14e8a7SBarry Smith KSP v = (KSP) object; 116*7e14e8a7SBarry Smith int ierr; 117*7e14e8a7SBarry Smith 118*7e14e8a7SBarry Smith PetscFunctionBegin; 119*7e14e8a7SBarry Smith 120*7e14e8a7SBarry Smith /* if it is already published then return */ 121*7e14e8a7SBarry Smith if (v->amem >=0 ) PetscFunctionReturn(0); 122*7e14e8a7SBarry Smith 123*7e14e8a7SBarry Smith ierr = PetscObjectPublishBaseBegin(object);CHKERRQ(ierr); 124*7e14e8a7SBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->its,1,AMS_INT,AMS_READ, 125*7e14e8a7SBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 126*7e14e8a7SBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->rnorm,1,AMS_DOUBLE,AMS_READ, 127*7e14e8a7SBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 128*7e14e8a7SBarry Smith ierr = PetscObjectPublishBaseEnd(object);CHKERRQ(ierr); 129*7e14e8a7SBarry Smith #else 130*7e14e8a7SBarry Smith PetscFunctionBegin; 131*7e14e8a7SBarry Smith #endif 132*7e14e8a7SBarry Smith 133*7e14e8a7SBarry Smith PetscFunctionReturn(0); 134*7e14e8a7SBarry Smith } 135*7e14e8a7SBarry Smith 136*7e14e8a7SBarry Smith 137*7e14e8a7SBarry Smith #undef __FUNC__ 138*7e14e8a7SBarry Smith #define __FUNC__ "KSPCreate" 139*7e14e8a7SBarry Smith /*@C 140*7e14e8a7SBarry Smith KSPCreate - Creates the default KSP context. 141*7e14e8a7SBarry Smith 142*7e14e8a7SBarry Smith Collective on MPI_Comm 143*7e14e8a7SBarry Smith 144*7e14e8a7SBarry Smith Input Parameter: 145*7e14e8a7SBarry Smith . comm - MPI communicator 146*7e14e8a7SBarry Smith 147*7e14e8a7SBarry Smith Output Parameter: 148*7e14e8a7SBarry Smith . ksp - location to put the KSP context 149*7e14e8a7SBarry Smith 150*7e14e8a7SBarry Smith Notes: 151*7e14e8a7SBarry Smith The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt 152*7e14e8a7SBarry Smith orthogonalization. 153*7e14e8a7SBarry Smith 154*7e14e8a7SBarry Smith Level: developer 155*7e14e8a7SBarry Smith 156*7e14e8a7SBarry Smith .keywords: KSP, create, context 157*7e14e8a7SBarry Smith 158*7e14e8a7SBarry Smith .seealso: KSPSetUp(), KSPSolve(), KSPDestroy() 159*7e14e8a7SBarry Smith @*/ 160*7e14e8a7SBarry Smith int KSPCreate(MPI_Comm comm,KSP *inksp) 161*7e14e8a7SBarry Smith { 162*7e14e8a7SBarry Smith KSP ksp; 163*7e14e8a7SBarry Smith 164*7e14e8a7SBarry Smith PetscFunctionBegin; 165*7e14e8a7SBarry Smith *inksp = 0; 166*7e14e8a7SBarry Smith PetscHeaderCreate(ksp,_p_KSP,struct _KSPOps,KSP_COOKIE,-1,"KSP",comm,KSPDestroy,KSPView); 167*7e14e8a7SBarry Smith PLogObjectCreate(ksp); 168*7e14e8a7SBarry Smith *inksp = ksp; 169*7e14e8a7SBarry Smith ksp->bops->publish = KSPPublish_Petsc; 170*7e14e8a7SBarry Smith 171*7e14e8a7SBarry Smith ksp->type = -1; 172*7e14e8a7SBarry Smith ksp->max_it = 10000; 173*7e14e8a7SBarry Smith ksp->pc_side = PC_LEFT; 174*7e14e8a7SBarry Smith ksp->use_pres = 0; 175*7e14e8a7SBarry Smith ksp->rtol = 1.e-5; 176*7e14e8a7SBarry Smith ksp->atol = 1.e-50; 177*7e14e8a7SBarry Smith ksp->divtol = 1.e4; 178*7e14e8a7SBarry Smith ksp->avoidnorms = PETSC_FALSE; 179*7e14e8a7SBarry Smith 180*7e14e8a7SBarry Smith ksp->rnorm = 0.0; 181*7e14e8a7SBarry Smith ksp->its = 0; 182*7e14e8a7SBarry Smith ksp->guess_zero = 1; 183*7e14e8a7SBarry Smith ksp->calc_sings = 0; 184*7e14e8a7SBarry Smith ksp->calc_res = 0; 185*7e14e8a7SBarry Smith ksp->res_hist = PETSC_NULL; 186*7e14e8a7SBarry Smith ksp->res_hist_len = 0; 187*7e14e8a7SBarry Smith ksp->res_hist_max = 0; 188*7e14e8a7SBarry Smith ksp->res_hist_reset = PETSC_TRUE; 189*7e14e8a7SBarry Smith ksp->numbermonitors = 0; 190*7e14e8a7SBarry Smith ksp->converged = KSPDefaultConverged; 191*7e14e8a7SBarry Smith ksp->ops->buildsolution = KSPDefaultBuildSolution; 192*7e14e8a7SBarry Smith ksp->ops->buildresidual = KSPDefaultBuildResidual; 193*7e14e8a7SBarry Smith 194*7e14e8a7SBarry Smith ksp->ops->setfromoptions = 0; 195*7e14e8a7SBarry Smith ksp->ops->printhelp = 0; 196*7e14e8a7SBarry Smith 197*7e14e8a7SBarry Smith ksp->vec_sol = 0; 198*7e14e8a7SBarry Smith ksp->vec_rhs = 0; 199*7e14e8a7SBarry Smith ksp->B = 0; 200*7e14e8a7SBarry Smith 201*7e14e8a7SBarry Smith ksp->ops->solve = 0; 202*7e14e8a7SBarry Smith ksp->ops->solvetrans = 0; 203*7e14e8a7SBarry Smith ksp->ops->setup = 0; 204*7e14e8a7SBarry Smith ksp->ops->destroy = 0; 205*7e14e8a7SBarry Smith 206*7e14e8a7SBarry Smith ksp->data = 0; 207*7e14e8a7SBarry Smith ksp->nwork = 0; 208*7e14e8a7SBarry Smith ksp->work = 0; 209*7e14e8a7SBarry Smith 210*7e14e8a7SBarry Smith ksp->cnvP = 0; 211*7e14e8a7SBarry Smith 212*7e14e8a7SBarry Smith ksp->setupcalled = 0; 213*7e14e8a7SBarry Smith PetscPublishAll(ksp); 214*7e14e8a7SBarry Smith PetscFunctionReturn(0); 215*7e14e8a7SBarry Smith } 216*7e14e8a7SBarry Smith 217*7e14e8a7SBarry Smith #undef __FUNC__ 218*7e14e8a7SBarry Smith #define __FUNC__ "KSPSetType" 219*7e14e8a7SBarry Smith /*@C 220*7e14e8a7SBarry Smith KSPSetType - Builds KSP for a particular solver. 221*7e14e8a7SBarry Smith 222*7e14e8a7SBarry Smith Collective on KSP 223*7e14e8a7SBarry Smith 224*7e14e8a7SBarry Smith Input Parameters: 225*7e14e8a7SBarry Smith . ksp - the Krylov space context 226*7e14e8a7SBarry Smith . itmethod - a known method 227*7e14e8a7SBarry Smith 228*7e14e8a7SBarry Smith Options Database Key: 229*7e14e8a7SBarry Smith . -ksp_type <method> - Sets the method; use -help for a list 230*7e14e8a7SBarry Smith of available methods (for instance, cg or gmres) 231*7e14e8a7SBarry Smith 232*7e14e8a7SBarry Smith Notes: 233*7e14e8a7SBarry Smith See "petsc/include/ksp.h" for available methods (for instance, 234*7e14e8a7SBarry Smith KSPCG or KSPGMRES). 235*7e14e8a7SBarry Smith 236*7e14e8a7SBarry Smith Normally, it is best to use the SLESSetFromOptions() command and 237*7e14e8a7SBarry Smith then set the KSP type from the options database rather than by using 238*7e14e8a7SBarry Smith this routine. Using the options database provides the user with 239*7e14e8a7SBarry Smith maximum flexibility in evaluating the many different Krylov methods. 240*7e14e8a7SBarry Smith The KSPSetType() routine is provided for those situations where it 241*7e14e8a7SBarry Smith is necessary to set the iterative solver independently of the command 242*7e14e8a7SBarry Smith line or options database. This might be the case, for example, when 243*7e14e8a7SBarry Smith the choice of iterative solver changes during the execution of the 244*7e14e8a7SBarry Smith program, and the user's application is taking responsibility for 245*7e14e8a7SBarry Smith choosing the appropriate method. In other words, this routine is 246*7e14e8a7SBarry Smith not for beginners. 247*7e14e8a7SBarry Smith 248*7e14e8a7SBarry Smith Level: intermediate 249*7e14e8a7SBarry Smith 250*7e14e8a7SBarry Smith .keywords: KSP, set, method 251*7e14e8a7SBarry Smith 252*7e14e8a7SBarry Smith .seealso: PCSetType() 253*7e14e8a7SBarry Smith @*/ 254*7e14e8a7SBarry Smith int KSPSetType(KSP ksp,KSPType itmethod) 255*7e14e8a7SBarry Smith { 256*7e14e8a7SBarry Smith int ierr,(*r)(KSP); 257*7e14e8a7SBarry Smith 258*7e14e8a7SBarry Smith PetscFunctionBegin; 259*7e14e8a7SBarry Smith PetscValidHeaderSpecific(ksp,KSP_COOKIE); 260*7e14e8a7SBarry Smith 261*7e14e8a7SBarry Smith if (PetscTypeCompare(ksp->type_name,itmethod)) PetscFunctionReturn(0); 262*7e14e8a7SBarry Smith 263*7e14e8a7SBarry Smith if (ksp->setupcalled) { 264*7e14e8a7SBarry Smith /* destroy the old private KSP context */ 265*7e14e8a7SBarry Smith ierr = (*ksp->ops->destroy)(ksp);CHKERRQ(ierr); 266*7e14e8a7SBarry Smith ksp->data = 0; 267*7e14e8a7SBarry Smith } 268*7e14e8a7SBarry Smith /* Get the function pointers for the iterative method requested */ 269*7e14e8a7SBarry Smith if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 270*7e14e8a7SBarry Smith 271*7e14e8a7SBarry Smith ierr = FListFind(ksp->comm, KSPList, itmethod,(int (**)(void *)) &r );CHKERRQ(ierr); 272*7e14e8a7SBarry Smith 273*7e14e8a7SBarry Smith if (!r) SETERRQ1(1,1,"Unknown KSP type given: %s",itmethod); 274*7e14e8a7SBarry Smith 275*7e14e8a7SBarry Smith if (ksp->data) PetscFree(ksp->data); 276*7e14e8a7SBarry Smith ksp->data = 0; 277*7e14e8a7SBarry Smith ksp->setupcalled = 0; 278*7e14e8a7SBarry Smith ierr = (*r)(ksp);CHKERRQ(ierr); 279*7e14e8a7SBarry Smith 280*7e14e8a7SBarry Smith if (ksp->type_name) PetscFree(ksp->type_name); 281*7e14e8a7SBarry Smith ksp->type_name = (char *) PetscMalloc((PetscStrlen(itmethod)+1)*sizeof(char));CHKPTRQ(ksp->type_name); 282*7e14e8a7SBarry Smith ierr = PetscStrcpy(ksp->type_name,itmethod);CHKERRQ(ierr); 283*7e14e8a7SBarry Smith PetscFunctionReturn(0); 284*7e14e8a7SBarry Smith } 285*7e14e8a7SBarry Smith 286*7e14e8a7SBarry Smith #undef __FUNC__ 287*7e14e8a7SBarry Smith #define __FUNC__ "KSPRegisterDestroy" 288*7e14e8a7SBarry Smith /*@C 289*7e14e8a7SBarry Smith KSPRegisterDestroy - Frees the list of KSP methods that were 290*7e14e8a7SBarry Smith registered by KSPRegister(). 291*7e14e8a7SBarry Smith 292*7e14e8a7SBarry Smith Not Collective 293*7e14e8a7SBarry Smith 294*7e14e8a7SBarry Smith Level: advanced 295*7e14e8a7SBarry Smith 296*7e14e8a7SBarry Smith .keywords: KSP, register, destroy 297*7e14e8a7SBarry Smith 298*7e14e8a7SBarry Smith .seealso: KSPRegister(), KSPRegisterAll() 299*7e14e8a7SBarry Smith @*/ 300*7e14e8a7SBarry Smith int KSPRegisterDestroy(void) 301*7e14e8a7SBarry Smith { 302*7e14e8a7SBarry Smith int ierr; 303*7e14e8a7SBarry Smith 304*7e14e8a7SBarry Smith PetscFunctionBegin; 305*7e14e8a7SBarry Smith if (KSPList) { 306*7e14e8a7SBarry Smith ierr = FListDestroy( KSPList );CHKERRQ(ierr); 307*7e14e8a7SBarry Smith KSPList = 0; 308*7e14e8a7SBarry Smith } 309*7e14e8a7SBarry Smith KSPRegisterAllCalled = 0; 310*7e14e8a7SBarry Smith PetscFunctionReturn(0); 311*7e14e8a7SBarry Smith } 312*7e14e8a7SBarry Smith 313*7e14e8a7SBarry Smith #undef __FUNC__ 314*7e14e8a7SBarry Smith #define __FUNC__ "KSPGetType" 315*7e14e8a7SBarry Smith /*@C 316*7e14e8a7SBarry Smith KSPGetType - Gets the KSP type as a string from the KSP object. 317*7e14e8a7SBarry Smith 318*7e14e8a7SBarry Smith Not Collective 319*7e14e8a7SBarry Smith 320*7e14e8a7SBarry Smith Input Parameter: 321*7e14e8a7SBarry Smith . ksp - Krylov context 322*7e14e8a7SBarry Smith 323*7e14e8a7SBarry Smith Output Parameter: 324*7e14e8a7SBarry Smith . name - name of KSP method 325*7e14e8a7SBarry Smith 326*7e14e8a7SBarry Smith Level: intermediate 327*7e14e8a7SBarry Smith 328*7e14e8a7SBarry Smith .keywords: KSP, get, method, name 329*7e14e8a7SBarry Smith 330*7e14e8a7SBarry Smith .seealso: KSPSetType() 331*7e14e8a7SBarry Smith @*/ 332*7e14e8a7SBarry Smith int KSPGetType(KSP ksp,KSPType *type) 333*7e14e8a7SBarry Smith { 334*7e14e8a7SBarry Smith PetscFunctionBegin; 335*7e14e8a7SBarry Smith *type = ksp->type_name; 336*7e14e8a7SBarry Smith PetscFunctionReturn(0); 337*7e14e8a7SBarry Smith } 338*7e14e8a7SBarry Smith 339*7e14e8a7SBarry Smith #undef __FUNC__ 340*7e14e8a7SBarry Smith #define __FUNC__ "KSPPrintHelp" 341*7e14e8a7SBarry Smith /*@ 342*7e14e8a7SBarry Smith KSPPrintHelp - Prints all options for the KSP component. 343*7e14e8a7SBarry Smith 344*7e14e8a7SBarry Smith Collective on KSP 345*7e14e8a7SBarry Smith 346*7e14e8a7SBarry Smith Input Parameter: 347*7e14e8a7SBarry Smith . ksp - the KSP context 348*7e14e8a7SBarry Smith 349*7e14e8a7SBarry Smith Options Database Keys: 350*7e14e8a7SBarry Smith + -help - Prints KSP options 351*7e14e8a7SBarry Smith - -h - Prints KSP options 352*7e14e8a7SBarry Smith 353*7e14e8a7SBarry Smith Level: developer 354*7e14e8a7SBarry Smith 355*7e14e8a7SBarry Smith .keywords: KSP, help 356*7e14e8a7SBarry Smith 357*7e14e8a7SBarry Smith .seealso: KSPSetFromOptions() 358*7e14e8a7SBarry Smith @*/ 359*7e14e8a7SBarry Smith int KSPPrintHelp(KSP ksp) 360*7e14e8a7SBarry Smith { 361*7e14e8a7SBarry Smith char p[64]; 362*7e14e8a7SBarry Smith int ierr; 363*7e14e8a7SBarry Smith 364*7e14e8a7SBarry Smith PetscFunctionBegin; 365*7e14e8a7SBarry Smith PetscValidHeaderSpecific(ksp,KSP_COOKIE); 366*7e14e8a7SBarry Smith ierr = PetscStrcpy(p,"-");CHKERRQ(ierr); 367*7e14e8a7SBarry Smith if (ksp->prefix) PetscStrcat(p,ksp->prefix); 368*7e14e8a7SBarry Smith 369*7e14e8a7SBarry Smith if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 370*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm,"KSP options -------------------------------------------------\n");CHKERRQ(ierr); 371*7e14e8a7SBarry Smith ierr = FListPrintTypes(ksp->comm,stdout,ksp->prefix,"ksp_type",KSPList);CHKERRQ(ierr);CHKERRQ(ierr); 372*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_rtol <tol>: relative tolerance, defaults to %g\n", 373*7e14e8a7SBarry Smith p,ksp->rtol);CHKERRQ(ierr); 374*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_atol <tol>: absolute tolerance, defaults to %g\n", 375*7e14e8a7SBarry Smith p,ksp->atol);CHKERRQ(ierr); 376*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_divtol <tol>: divergence tolerance, defaults to %g\n", 377*7e14e8a7SBarry Smith p,ksp->divtol);CHKERRQ(ierr); 378*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_max_it <maxit>: maximum iterations, defaults to %d\n", 379*7e14e8a7SBarry Smith p,ksp->max_it);CHKERRQ(ierr); 380*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_preres: use preconditioned residual norm in convergence test\n",p);CHKERRQ(ierr); 381*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_right_pc: use right preconditioner instead of left\n",p);CHKERRQ(ierr); 382*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_avoid_norms: do not compute residual norms (CG and Bi-CG-stab)\n",p);CHKERRQ(ierr); 383*7e14e8a7SBarry Smith 384*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," KSP Monitoring Options: Choose any of the following\n");CHKERRQ(ierr); 385*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_cancelmonitors: cancel all monitors hardwired in code\n",p);CHKERRQ(ierr); 386*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_monitor: at each iteration print (usually preconditioned) \n\ 387*7e14e8a7SBarry Smith residual norm to stdout\n",p);CHKERRQ(ierr); 388*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_smonitor: same as the above, but prints fewer digits of the\n\ 389*7e14e8a7SBarry Smith residual norm for small residual norms. This is useful to conceal\n\ 390*7e14e8a7SBarry Smith meaningless digits that may be different on different machines.\n",p);CHKERRQ(ierr); 391*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_xmonitor [x,y,w,h]: use X graphics monitor of (usually \n\ 392*7e14e8a7SBarry Smith preconditioned) residual norm\n",p);CHKERRQ(ierr); 393*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_truemonitor: at each iteration print true and preconditioned\n",p);CHKERRQ(ierr); 394*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," residual norms to stdout\n");CHKERRQ(ierr); 395*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_xtruemonitor [x,y,w,h]: use X graphics monitor of true\n",p);CHKERRQ(ierr); 396*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," residual norm\n");CHKERRQ(ierr); 397*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_singmonitor: calculate singular values during linear solve\n",p);CHKERRQ(ierr); 398*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," (only for CG and GMRES)\n");CHKERRQ(ierr); 399*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_plot_eigenvalues_explicitly\n",p);CHKERRQ(ierr); 400*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_plot_eigenvalues\n",p);CHKERRQ(ierr); 401*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_compute_eigenvalues\n",p);CHKERRQ(ierr); 402*7e14e8a7SBarry Smith ierr = (*PetscHelpPrintf)(ksp->comm," %sksp_compute_singularvalues\n",p);CHKERRQ(ierr); 403*7e14e8a7SBarry Smith 404*7e14e8a7SBarry Smith if (ksp->ops->printhelp) { 405*7e14e8a7SBarry Smith ierr = (*ksp->ops->printhelp)(ksp,p);CHKERRQ(ierr); 406*7e14e8a7SBarry Smith } 407*7e14e8a7SBarry Smith PetscFunctionReturn(0); 408*7e14e8a7SBarry Smith } 409*7e14e8a7SBarry Smith 410*7e14e8a7SBarry Smith #define MAXSETFROMOPTIONS 5 411*7e14e8a7SBarry Smith extern int numberofsetfromoptions; 412*7e14e8a7SBarry Smith extern int (*othersetfromoptions[MAXSETFROMOPTIONS])(KSP); 413*7e14e8a7SBarry Smith 414*7e14e8a7SBarry Smith #undef __FUNC__ 415*7e14e8a7SBarry Smith #define __FUNC__ "KSPSetTypeFromOptions" 416*7e14e8a7SBarry Smith /*@ 417*7e14e8a7SBarry Smith KSPSetTypeFromOptions - Sets KSP type from the options database, if not 418*7e14e8a7SBarry Smith given then sets default. 419*7e14e8a7SBarry Smith 420*7e14e8a7SBarry Smith Collective on KSP 421*7e14e8a7SBarry Smith 422*7e14e8a7SBarry Smith Input Parameters: 423*7e14e8a7SBarry Smith . ksp - the Krylov space context 424*7e14e8a7SBarry Smith 425*7e14e8a7SBarry Smith Level: developer 426*7e14e8a7SBarry Smith 427*7e14e8a7SBarry Smith .keywords: KSP, set, from, options, database 428*7e14e8a7SBarry Smith 429*7e14e8a7SBarry Smith .seealso: KSPPrintHelp(), KSPSetFromOptions(), SLESSetFromOptions(), 430*7e14e8a7SBarry Smith SLESSetTypeFromOptions() 431*7e14e8a7SBarry Smith @*/ 432*7e14e8a7SBarry Smith int KSPSetTypeFromOptions(KSP ksp) 433*7e14e8a7SBarry Smith { 434*7e14e8a7SBarry Smith int flg, ierr; 435*7e14e8a7SBarry Smith char method[256]; 436*7e14e8a7SBarry Smith 437*7e14e8a7SBarry Smith PetscFunctionBegin; 438*7e14e8a7SBarry Smith PetscValidHeaderSpecific(ksp,KSP_COOKIE); 439*7e14e8a7SBarry Smith 440*7e14e8a7SBarry Smith ierr = OptionsGetString(ksp->prefix,"-ksp_type",method,256,&flg); 441*7e14e8a7SBarry Smith if (flg) { 442*7e14e8a7SBarry Smith ierr = KSPSetType(ksp,method);CHKERRQ(ierr); 443*7e14e8a7SBarry Smith } 444*7e14e8a7SBarry Smith /* 445*7e14e8a7SBarry Smith Set the type if it was never set. 446*7e14e8a7SBarry Smith */ 447*7e14e8a7SBarry Smith if (!ksp->type_name) { 448*7e14e8a7SBarry Smith ierr = KSPSetType(ksp,KSPGMRES);CHKERRQ(ierr); 449*7e14e8a7SBarry Smith } 450*7e14e8a7SBarry Smith PetscFunctionReturn(0); 451*7e14e8a7SBarry Smith } 452*7e14e8a7SBarry Smith 453*7e14e8a7SBarry Smith #undef __FUNC__ 454*7e14e8a7SBarry Smith #define __FUNC__ "KSPSetFromOptions" 455*7e14e8a7SBarry Smith /*@ 456*7e14e8a7SBarry Smith KSPSetFromOptions - Sets KSP options from the options database. 457*7e14e8a7SBarry Smith This routine must be called before KSPSetUp() if the user is to be 458*7e14e8a7SBarry Smith allowed to set the Krylov type. 459*7e14e8a7SBarry Smith 460*7e14e8a7SBarry Smith Collective on KSP 461*7e14e8a7SBarry Smith 462*7e14e8a7SBarry Smith Input Parameters: 463*7e14e8a7SBarry Smith . ksp - the Krylov space context 464*7e14e8a7SBarry Smith 465*7e14e8a7SBarry Smith Options Database Keys: 466*7e14e8a7SBarry Smith + -ksp_max_it - maximum number of linear iterations 467*7e14e8a7SBarry Smith . -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e. 468*7e14e8a7SBarry Smith if residual norm decreases by this factor than convergence is declared 469*7e14e8a7SBarry Smith . -ksp_atol atol - absolute tolerance used in default convergence test, i.e. if residual 470*7e14e8a7SBarry Smith norm is less than this then convergence is declared 471*7e14e8a7SBarry Smith . -ksp_divtol tol - if residual norm increases by this factor than divergence is declared 472*7e14e8a7SBarry Smith . -ksp_avoid_norms - skip norms used in convergence tests (useful only when not using 473*7e14e8a7SBarry Smith convergence test (say you always want to run with 5 iterations) to 474*7e14e8a7SBarry Smith save on communication overhead 475*7e14e8a7SBarry Smith . -ksp_cancelmonitors - cancel all previous convergene monitor routines set 476*7e14e8a7SBarry Smith . -ksp_monitor - print residual norm at each iteration 477*7e14e8a7SBarry Smith . -ksp_xmonitor - plot residual norm at each iteration 478*7e14e8a7SBarry Smith . -ksp_vecmonitor - plot solution at each iteration 479*7e14e8a7SBarry Smith - -ksp_singmonitor - monitor extremem singular values at each iteration 480*7e14e8a7SBarry Smith 481*7e14e8a7SBarry Smith Notes: 482*7e14e8a7SBarry Smith To see all options, run your program with the -help option 483*7e14e8a7SBarry Smith or consult the users manual. 484*7e14e8a7SBarry Smith 485*7e14e8a7SBarry Smith Level: developer 486*7e14e8a7SBarry Smith 487*7e14e8a7SBarry Smith .keywords: KSP, set, from, options, database 488*7e14e8a7SBarry Smith 489*7e14e8a7SBarry Smith .seealso: KSPPrintHelp() 490*7e14e8a7SBarry Smith @*/ 491*7e14e8a7SBarry Smith int KSPSetFromOptions(KSP ksp) 492*7e14e8a7SBarry Smith { 493*7e14e8a7SBarry Smith int flg, ierr,loc[4], nmax = 4,i; 494*7e14e8a7SBarry Smith 495*7e14e8a7SBarry Smith PetscFunctionBegin; 496*7e14e8a7SBarry Smith PetscValidHeaderSpecific(ksp,KSP_COOKIE); 497*7e14e8a7SBarry Smith 498*7e14e8a7SBarry Smith ierr = KSPSetTypeFromOptions(ksp);CHKERRQ(ierr); 499*7e14e8a7SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 500*7e14e8a7SBarry Smith 501*7e14e8a7SBarry Smith ierr = OptionsGetInt(ksp->prefix,"-ksp_max_it",&ksp->max_it, &flg);CHKERRQ(ierr); 502*7e14e8a7SBarry Smith ierr = OptionsGetDouble(ksp->prefix,"-ksp_rtol",&ksp->rtol, &flg);CHKERRQ(ierr); 503*7e14e8a7SBarry Smith ierr = OptionsGetDouble(ksp->prefix,"-ksp_atol",&ksp->atol, &flg);CHKERRQ(ierr); 504*7e14e8a7SBarry Smith ierr = OptionsGetDouble(ksp->prefix,"-ksp_divtol",&ksp->divtol, &flg);CHKERRQ(ierr); 505*7e14e8a7SBarry Smith 506*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_avoid_norms", &flg);CHKERRQ(ierr); 507*7e14e8a7SBarry Smith if (flg) { 508*7e14e8a7SBarry Smith ierr = KSPSetAvoidNorms(ksp);CHKERRQ(ierr); 509*7e14e8a7SBarry Smith } 510*7e14e8a7SBarry Smith 511*7e14e8a7SBarry Smith /* -----------------------------------------------------------------------*/ 512*7e14e8a7SBarry Smith /* 513*7e14e8a7SBarry Smith Cancels all monitors hardwired into code before call to KSPSetFromOptions() 514*7e14e8a7SBarry Smith */ 515*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_cancelmonitors",&flg);CHKERRQ(ierr); 516*7e14e8a7SBarry Smith if (flg) { 517*7e14e8a7SBarry Smith ierr = KSPClearMonitor(ksp);CHKERRQ(ierr); 518*7e14e8a7SBarry Smith } 519*7e14e8a7SBarry Smith /* 520*7e14e8a7SBarry Smith Prints preconditioned residual norm at each iteration 521*7e14e8a7SBarry Smith */ 522*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_monitor",&flg);CHKERRQ(ierr); 523*7e14e8a7SBarry Smith if (flg) { 524*7e14e8a7SBarry Smith int rank = 0; 525*7e14e8a7SBarry Smith ierr = MPI_Comm_rank(ksp->comm,&rank);CHKERRQ(ierr); 526*7e14e8a7SBarry Smith if (!rank) { 527*7e14e8a7SBarry Smith ierr = KSPSetMonitor(ksp,KSPDefaultMonitor,0,0);CHKERRQ(ierr); 528*7e14e8a7SBarry Smith } 529*7e14e8a7SBarry Smith } 530*7e14e8a7SBarry Smith /* 531*7e14e8a7SBarry Smith Plots the vector solution 532*7e14e8a7SBarry Smith */ 533*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_vecmonitor",&flg);CHKERRQ(ierr); 534*7e14e8a7SBarry Smith if (flg) { 535*7e14e8a7SBarry Smith ierr = KSPSetMonitor(ksp,KSPVecViewMonitor,0,0);CHKERRQ(ierr); 536*7e14e8a7SBarry Smith } 537*7e14e8a7SBarry Smith /* 538*7e14e8a7SBarry Smith Prints preconditioned and true residual norm at each iteration 539*7e14e8a7SBarry Smith */ 540*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_truemonitor",&flg);CHKERRQ(ierr); 541*7e14e8a7SBarry Smith if (flg) { 542*7e14e8a7SBarry Smith ierr = KSPSetMonitor(ksp,KSPTrueMonitor,0,0);CHKERRQ(ierr); 543*7e14e8a7SBarry Smith } 544*7e14e8a7SBarry Smith /* 545*7e14e8a7SBarry Smith Prints extreme eigenvalue estimates at each iteration 546*7e14e8a7SBarry Smith */ 547*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_singmonitor",&flg);CHKERRQ(ierr); 548*7e14e8a7SBarry Smith if (flg) { 549*7e14e8a7SBarry Smith ierr = KSPSetComputeSingularValues(ksp);CHKERRQ(ierr); 550*7e14e8a7SBarry Smith ierr = KSPSetMonitor(ksp,KSPSingularValueMonitor,0,0);CHKERRQ(ierr); 551*7e14e8a7SBarry Smith } 552*7e14e8a7SBarry Smith /* 553*7e14e8a7SBarry Smith Prints preconditioned residual norm with fewer digits 554*7e14e8a7SBarry Smith */ 555*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_smonitor",&flg);CHKERRQ(ierr); 556*7e14e8a7SBarry Smith if (flg) { 557*7e14e8a7SBarry Smith int rank = 0; 558*7e14e8a7SBarry Smith ierr = MPI_Comm_rank(ksp->comm,&rank);CHKERRQ(ierr); 559*7e14e8a7SBarry Smith if (!rank) { 560*7e14e8a7SBarry Smith ierr = KSPSetMonitor(ksp,KSPDefaultSMonitor,0,0);CHKERRQ(ierr); 561*7e14e8a7SBarry Smith } 562*7e14e8a7SBarry Smith } 563*7e14e8a7SBarry Smith /* 564*7e14e8a7SBarry Smith Graphically plots preconditioned residual norm 565*7e14e8a7SBarry Smith */ 566*7e14e8a7SBarry Smith nmax = 4; 567*7e14e8a7SBarry Smith ierr = OptionsGetIntArray(ksp->prefix,"-ksp_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 568*7e14e8a7SBarry Smith if (flg) { 569*7e14e8a7SBarry Smith int rank = 0; 570*7e14e8a7SBarry Smith DrawLG lg; 571*7e14e8a7SBarry Smith ierr = MPI_Comm_rank(ksp->comm,&rank);CHKERRQ(ierr); 572*7e14e8a7SBarry Smith if (!rank) { 573*7e14e8a7SBarry Smith ierr = KSPLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg);CHKERRQ(ierr); 574*7e14e8a7SBarry Smith PLogObjectParent(ksp,(PetscObject) lg); 575*7e14e8a7SBarry Smith ierr = KSPSetMonitor(ksp,KSPLGMonitor,lg,(int (*)(void*))KSPLGMonitorDestroy);CHKERRQ(ierr); 576*7e14e8a7SBarry Smith } 577*7e14e8a7SBarry Smith } 578*7e14e8a7SBarry Smith /* 579*7e14e8a7SBarry Smith Graphically plots preconditioned and true residual norm 580*7e14e8a7SBarry Smith */ 581*7e14e8a7SBarry Smith nmax = 4; 582*7e14e8a7SBarry Smith ierr = OptionsGetIntArray(ksp->prefix,"-ksp_xtruemonitor",loc,&nmax,&flg);CHKERRQ(ierr); 583*7e14e8a7SBarry Smith if (flg){ 584*7e14e8a7SBarry Smith int rank = 0; 585*7e14e8a7SBarry Smith DrawLG lg; 586*7e14e8a7SBarry Smith ierr = MPI_Comm_rank(ksp->comm,&rank);CHKERRQ(ierr); 587*7e14e8a7SBarry Smith if (!rank) { 588*7e14e8a7SBarry Smith ierr = KSPLGTrueMonitorCreate(ksp->comm,0,0,loc[0],loc[1],loc[2],loc[3],&lg);CHKERRQ(ierr); 589*7e14e8a7SBarry Smith PLogObjectParent(ksp,(PetscObject) lg); 590*7e14e8a7SBarry Smith ierr = KSPSetMonitor(ksp,KSPLGTrueMonitor,lg,(int (*)(void*))KSPLGMonitorDestroy);CHKERRQ(ierr); 591*7e14e8a7SBarry Smith } 592*7e14e8a7SBarry Smith } 593*7e14e8a7SBarry Smith /* -----------------------------------------------------------------------*/ 594*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_preres",&flg);CHKERRQ(ierr); 595*7e14e8a7SBarry Smith if (flg) { ierr = KSPSetUsePreconditionedResidual(ksp);CHKERRQ(ierr);} 596*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_left_pc",&flg);CHKERRQ(ierr); 597*7e14e8a7SBarry Smith if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_LEFT);CHKERRQ(ierr); } 598*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_right_pc",&flg);CHKERRQ(ierr); 599*7e14e8a7SBarry Smith if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_RIGHT);CHKERRQ(ierr);} 600*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_symmetric_pc",&flg);CHKERRQ(ierr); 601*7e14e8a7SBarry Smith if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_SYMMETRIC);CHKERRQ(ierr);} 602*7e14e8a7SBarry Smith 603*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_compute_singularvalues",&flg);CHKERRQ(ierr); 604*7e14e8a7SBarry Smith if (flg) { ierr = KSPSetComputeSingularValues(ksp);CHKERRQ(ierr); } 605*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_compute_eigenvalues",&flg);CHKERRQ(ierr); 606*7e14e8a7SBarry Smith if (flg) { ierr = KSPSetComputeSingularValues(ksp);CHKERRQ(ierr); } 607*7e14e8a7SBarry Smith ierr = OptionsHasName(ksp->prefix,"-ksp_plot_eigenvalues",&flg);CHKERRQ(ierr); 608*7e14e8a7SBarry Smith if (flg) { ierr = KSPSetComputeSingularValues(ksp);CHKERRQ(ierr); } 609*7e14e8a7SBarry Smith 610*7e14e8a7SBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 611*7e14e8a7SBarry Smith ierr = (*othersetfromoptions[i])(ksp);CHKERRQ(ierr); 612*7e14e8a7SBarry Smith } 613*7e14e8a7SBarry Smith 614*7e14e8a7SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr); 615*7e14e8a7SBarry Smith if (flg) { ierr = KSPPrintHelp(ksp);CHKERRQ(ierr); } 616*7e14e8a7SBarry Smith 617*7e14e8a7SBarry Smith if (ksp->ops->setfromoptions) { 618*7e14e8a7SBarry Smith ierr = (*ksp->ops->setfromoptions)(ksp);CHKERRQ(ierr); 619*7e14e8a7SBarry Smith } 620*7e14e8a7SBarry Smith 621*7e14e8a7SBarry Smith PetscFunctionReturn(0); 622*7e14e8a7SBarry Smith } 623*7e14e8a7SBarry Smith 624*7e14e8a7SBarry Smith /*MC 625*7e14e8a7SBarry Smith KSPRegister - Adds a method to the Krylov subspace solver package. 626*7e14e8a7SBarry Smith 627*7e14e8a7SBarry Smith Synopsis: 628*7e14e8a7SBarry Smith KSPRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(KSP)) 629*7e14e8a7SBarry Smith 630*7e14e8a7SBarry Smith Not Collective 631*7e14e8a7SBarry Smith 632*7e14e8a7SBarry Smith Input Parameters: 633*7e14e8a7SBarry Smith + name_solver - name of a new user-defined solver 634*7e14e8a7SBarry Smith . path - path (either absolute or relative) the library containing this solver 635*7e14e8a7SBarry Smith . name_create - name of routine to create method context 636*7e14e8a7SBarry Smith - routine_create - routine to create method context 637*7e14e8a7SBarry Smith 638*7e14e8a7SBarry Smith Notes: 639*7e14e8a7SBarry Smith KSPRegister() may be called multiple times to add several user-defined solvers. 640*7e14e8a7SBarry Smith 641*7e14e8a7SBarry Smith If dynamic libraries are used, then the fourth input argument (routine_create) 642*7e14e8a7SBarry Smith is ignored. 643*7e14e8a7SBarry Smith 644*7e14e8a7SBarry Smith Sample usage: 645*7e14e8a7SBarry Smith .vb 646*7e14e8a7SBarry Smith KSPRegister("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a, 647*7e14e8a7SBarry Smith "MySolverCreate",MySolverCreate); 648*7e14e8a7SBarry Smith .ve 649*7e14e8a7SBarry Smith 650*7e14e8a7SBarry Smith Then, your solver can be chosen with the procedural interface via 651*7e14e8a7SBarry Smith $ KSPSetType(ksp,"my_solver") 652*7e14e8a7SBarry Smith or at runtime via the option 653*7e14e8a7SBarry Smith $ -ksp_type my_solver 654*7e14e8a7SBarry Smith 655*7e14e8a7SBarry Smith Level: advanced 656*7e14e8a7SBarry Smith 657*7e14e8a7SBarry Smith $PETSC_ARCH and $BOPT occuring in pathname will be replaced with appropriate values. 658*7e14e8a7SBarry Smith 659*7e14e8a7SBarry Smith .keywords: KSP, register 660*7e14e8a7SBarry Smith 661*7e14e8a7SBarry Smith .seealso: KSPRegisterAll(), KSPRegisterDestroy() 662*7e14e8a7SBarry Smith 663*7e14e8a7SBarry Smith M*/ 664*7e14e8a7SBarry Smith 665*7e14e8a7SBarry Smith #undef __FUNC__ 666*7e14e8a7SBarry Smith #define __FUNC__ "KSPRegister_Private" 667*7e14e8a7SBarry Smith int KSPRegister_Private(char *sname,char *path,char *name,int (*function)(KSP)) 668*7e14e8a7SBarry Smith { 669*7e14e8a7SBarry Smith int ierr; 670*7e14e8a7SBarry Smith char fullname[256]; 671*7e14e8a7SBarry Smith 672*7e14e8a7SBarry Smith PetscFunctionBegin; 673*7e14e8a7SBarry Smith ierr = PetscStrcpy(fullname,path);CHKERRQ(ierr); 674*7e14e8a7SBarry Smith PetscStrcat(fullname,":");PetscStrcat(fullname,name); 675*7e14e8a7SBarry Smith ierr = FListAdd_Private(&KSPList,sname,fullname,(int (*)(void*))function);CHKERRQ(ierr); 676*7e14e8a7SBarry Smith PetscFunctionReturn(0); 677*7e14e8a7SBarry Smith } 678