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