1 /* 2 The basic KSP routines, Create, View etc. are here. 3 */ 4 #include "src/ksp/ksp/kspimpl.h" /*I "petscksp.h" I*/ 5 #include "petscsys.h" 6 7 /* Logging support */ 8 PetscCookie KSP_COOKIE = 0; 9 PetscEvent KSP_GMRESOrthogonalization = 0, KSP_SetUp = 0, KSP_Solve = 0; 10 11 12 PetscTruth KSPRegisterAllCalled = PETSC_FALSE; 13 14 #undef __FUNCT__ 15 #define __FUNCT__ "KSPView" 16 /*@C 17 KSPView - Prints the KSP data structure. 18 19 Collective on KSP 20 21 Input Parameters: 22 + ksp - the Krylov space context 23 - viewer - visualization context 24 25 Options Database Keys: 26 . -ksp_view - print the ksp data structure at the end of a KSPSolve call 27 28 Note: 29 The available visualization contexts include 30 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 31 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 32 output where only the first processor opens 33 the file. All other processors send their 34 data to the first processor to print. 35 36 The user can open an alternative visualization context with 37 PetscViewerASCIIOpen() - output to a specified file. 38 39 Level: beginner 40 41 .keywords: KSP, view 42 43 .seealso: PCView(), PetscViewerASCIIOpen() 44 @*/ 45 PetscErrorCode KSPView(KSP ksp,PetscViewer viewer) 46 { 47 char *type; 48 PetscErrorCode ierr; 49 PetscTruth iascii; 50 51 PetscFunctionBegin; 52 PetscValidHeaderSpecific(ksp,KSP_COOKIE,1); 53 if (!viewer) viewer = PETSC_VIEWER_STDOUT_(ksp->comm); 54 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 55 PetscCheckSameComm(ksp,1,viewer,2); 56 57 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 58 if (iascii) { 59 ierr = KSPGetType(ksp,&type);CHKERRQ(ierr); 60 if (ksp->prefix) { 61 ierr = PetscViewerASCIIPrintf(viewer,"KSP Object:(%s)\n",ksp->prefix);CHKERRQ(ierr); 62 } else { 63 ierr = PetscViewerASCIIPrintf(viewer,"KSP Object:\n");CHKERRQ(ierr); 64 } 65 if (type) { 66 ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 67 } else { 68 ierr = PetscViewerASCIIPrintf(viewer," type: not yet set\n");CHKERRQ(ierr); 69 } 70 if (ksp->ops->view) { 71 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 72 ierr = (*ksp->ops->view)(ksp,viewer);CHKERRQ(ierr); 73 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 74 } 75 if (ksp->guess_zero) {ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, initial guess is zero\n",ksp->max_it);CHKERRQ(ierr);} 76 else {ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D\n", ksp->max_it);CHKERRQ(ierr);} 77 if (ksp->guess_knoll) {ierr = PetscViewerASCIIPrintf(viewer," using preconditioner applied to right hand side for initial guess\n");CHKERRQ(ierr);} 78 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, divergence=%g\n",ksp->rtol,ksp->abstol,ksp->divtol);CHKERRQ(ierr); 79 if (ksp->pc_side == PC_RIGHT) {ierr = PetscViewerASCIIPrintf(viewer," right preconditioning\n");CHKERRQ(ierr);} 80 else if (ksp->pc_side == PC_SYMMETRIC) {ierr = PetscViewerASCIIPrintf(viewer," symmetric preconditioning\n");CHKERRQ(ierr);} 81 else {ierr = PetscViewerASCIIPrintf(viewer," left preconditioning\n");CHKERRQ(ierr);} 82 } else { 83 if (ksp->ops->view) { 84 ierr = (*ksp->ops->view)(ksp,viewer);CHKERRQ(ierr); 85 } 86 } 87 ierr = PCView(ksp->pc,viewer);CHKERRQ(ierr); 88 PetscFunctionReturn(0); 89 } 90 91 /* 92 Contains the list of registered KSP routines 93 */ 94 PetscFList KSPList = 0; 95 96 #undef __FUNCT__ 97 #define __FUNCT__ "KSPSetNormType" 98 /*@C 99 KSPSetNormType - Sets the norm that is used for convergence testing. 100 101 Collective on KSP 102 103 Input Parameter: 104 + ksp - Krylov solver context 105 - normtype - one of 106 $ KSP_NO_NORM - skips computing the norm, this should only be used if you are using 107 $ the Krylov method as a smoother with a fixed small number of iterations. 108 $ You must also call KSPSetConvergenceTest(ksp,KSPSkipConverged,PETSC_NULL); 109 $ supported only by CG, Richardson, Bi-CG-stab, CR, and CGS methods. 110 $ KSP_PRECONDITIONED_NORM - the default for left preconditioned solves, uses the l2 norm 111 $ of the preconditioned residual 112 $ KSP_UNPRECONDITIONED_NORM - uses the l2 norm of the true b - Ax residual, supported only by 113 $ CG, CHEBYCHEV, and RICHARDSON 114 $ KSP_NATURAL_NORM - supported by cg, cr, and cgs 115 116 117 Options Database Key: 118 . -ksp_norm_type <none,preconditioned,unpreconditioned,natural> 119 120 Notes: 121 Currently only works with the CG, Richardson, Bi-CG-stab, CR, and CGS methods. 122 123 Level: advanced 124 125 .keywords: KSP, create, context, norms 126 127 .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSPSkipConverged() 128 @*/ 129 PetscErrorCode KSPSetNormType(KSP ksp,KSPNormType normtype) 130 { 131 PetscFunctionBegin; 132 PetscValidHeaderSpecific(ksp,KSP_COOKIE,1); 133 ksp->normtype = normtype; 134 if (normtype == KSP_NO_NORM) { 135 PetscLogInfo(ksp,"KSPSetNormType:Warning seting KSPNormType to skip computing the norm\n\ 136 make sure you set the KSP convergence test to KSPSkipConvergence\n"); 137 } 138 PetscFunctionReturn(0); 139 } 140 141 #undef __FUNCT__ 142 #define __FUNCT__ "KSPPublish_Petsc" 143 static PetscErrorCode KSPPublish_Petsc(PetscObject obj) 144 { 145 PetscFunctionBegin; 146 PetscFunctionReturn(0); 147 } 148 149 #undef __FUNCT__ 150 #define __FUNCT__ "KSPSetOperators" 151 /*@ 152 KSPSetOperators - Sets the matrix associated with the linear system 153 and a (possibly) different one associated with the preconditioner. 154 155 Collective on KSP and Mat 156 157 Input Parameters: 158 + ksp - the KSP context 159 . Amat - the matrix associated with the linear system 160 . Pmat - the matrix to be used in constructing the preconditioner, usually the 161 same as Amat. 162 - flag - flag indicating information about the preconditioner matrix structure 163 during successive linear solves. This flag is ignored the first time a 164 linear system is solved, and thus is irrelevant when solving just one linear 165 system. 166 167 Notes: 168 The flag can be used to eliminate unnecessary work in the preconditioner 169 during the repeated solution of linear systems of the same size. The 170 available options are 171 $ SAME_PRECONDITIONER - 172 $ Pmat is identical during successive linear solves. 173 $ This option is intended for folks who are using 174 $ different Amat and Pmat matrices and want to reuse the 175 $ same preconditioner matrix. For example, this option 176 $ saves work by not recomputing incomplete factorization 177 $ for ILU/ICC preconditioners. 178 $ SAME_NONZERO_PATTERN - 179 $ Pmat has the same nonzero structure during 180 $ successive linear solves. 181 $ DIFFERENT_NONZERO_PATTERN - 182 $ Pmat does not have the same nonzero structure. 183 184 Caution: 185 If you specify SAME_NONZERO_PATTERN, PETSc believes your assertion 186 and does not check the structure of the matrix. If you erroneously 187 claim that the structure is the same when it actually is not, the new 188 preconditioner will not function correctly. Thus, use this optimization 189 feature carefully! 190 191 If in doubt about whether your preconditioner matrix has changed 192 structure or not, use the flag DIFFERENT_NONZERO_PATTERN. 193 194 Level: beginner 195 196 .keywords: KSP, set, operators, matrix, preconditioner, linear system 197 198 .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPGetOperators() 199 @*/ 200 PetscErrorCode KSPSetOperators(KSP ksp,Mat Amat,Mat Pmat,MatStructure flag) 201 { 202 PetscErrorCode ierr; 203 204 PetscFunctionBegin; 205 PetscValidHeaderSpecific(ksp,KSP_COOKIE,1); 206 if (Amat) PetscValidHeaderSpecific(Amat,MAT_COOKIE,2); 207 if (Pmat) PetscValidHeaderSpecific(Pmat,MAT_COOKIE,3); 208 ierr = PCSetOperators(ksp->pc,Amat,Pmat,flag);CHKERRQ(ierr); 209 if (ksp->setupcalled > 1) ksp->setupcalled = 1; /* so that next solve call will call setup */ 210 PetscFunctionReturn(0); 211 } 212 213 #undef __FUNCT__ 214 #define __FUNCT__ "KSPGetOperators" 215 /*@ 216 KSPGetOperators - Gets the matrix associated with the linear system 217 and a (possibly) different one associated with the preconditioner. 218 219 Collective on KSP and Mat 220 221 Input Parameter: 222 . ksp - the KSP context 223 224 Output Parameters: 225 + Amat - the matrix associated with the linear system 226 . Pmat - the matrix to be used in constructing the preconditioner, usually the 227 same as Amat. 228 - flag - flag indicating information about the preconditioner matrix structure 229 during successive linear solves. This flag is ignored the first time a 230 linear system is solved, and thus is irrelevant when solving just one linear 231 system. 232 233 Level: intermediate 234 235 .keywords: KSP, set, get, operators, matrix, preconditioner, linear system 236 237 .seealso: KSPSolve(), KSPGetPC(), PCGetOperators(), PCSetOperators(), KSPSetOperators() 238 @*/ 239 PetscErrorCode KSPGetOperators(KSP ksp,Mat *Amat,Mat *Pmat,MatStructure *flag) 240 { 241 PetscErrorCode ierr; 242 243 PetscFunctionBegin; 244 PetscValidHeaderSpecific(ksp,KSP_COOKIE,1); 245 ierr = PCGetOperators(ksp->pc,Amat,Pmat,flag);CHKERRQ(ierr); 246 PetscFunctionReturn(0); 247 } 248 249 #undef __FUNCT__ 250 #define __FUNCT__ "KSPCreate" 251 /*@C 252 KSPCreate - Creates the default KSP context. 253 254 Collective on MPI_Comm 255 256 Input Parameter: 257 . comm - MPI communicator 258 259 Output Parameter: 260 . ksp - location to put the KSP context 261 262 Notes: 263 The default KSP type is GMRES with a restart of 30, using modified Gram-Schmidt 264 orthogonalization. 265 266 Level: beginner 267 268 .keywords: KSP, create, context 269 270 .seealso: KSPSetUp(), KSPSolve(), KSPDestroy(), KSP 271 @*/ 272 PetscErrorCode KSPCreate(MPI_Comm comm,KSP *inksp) 273 { 274 KSP ksp; 275 PetscErrorCode ierr; 276 277 PetscFunctionBegin; 278 PetscValidPointer(inksp,2); 279 *inksp = 0; 280 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 281 ierr = KSPInitializePackage(PETSC_NULL);CHKERRQ(ierr); 282 #endif 283 284 PetscHeaderCreate(ksp,_p_KSP,struct _KSPOps,KSP_COOKIE,-1,"KSP",comm,KSPDestroy,KSPView); 285 *inksp = ksp; 286 ksp->bops->publish = KSPPublish_Petsc; 287 288 ksp->type = -1; 289 ksp->max_it = 10000; 290 ksp->pc_side = PC_LEFT; 291 ksp->rtol = 1.e-5; 292 ksp->abstol = 1.e-50; 293 ksp->divtol = 1.e4; 294 295 ksp->normtype = KSP_PRECONDITIONED_NORM; 296 ksp->rnorm = 0.0; 297 ksp->its = 0; 298 ksp->guess_zero = PETSC_TRUE; 299 ksp->calc_sings = PETSC_FALSE; 300 ksp->res_hist = PETSC_NULL; 301 ksp->res_hist_len = 0; 302 ksp->res_hist_max = 0; 303 ksp->res_hist_reset = PETSC_TRUE; 304 ksp->numbermonitors = 0; 305 ksp->converged = KSPDefaultConverged; 306 ksp->ops->buildsolution = KSPDefaultBuildSolution; 307 ksp->ops->buildresidual = KSPDefaultBuildResidual; 308 309 ksp->ops->setfromoptions = 0; 310 311 ksp->vec_sol = 0; 312 ksp->vec_rhs = 0; 313 ksp->pc = 0; 314 315 ksp->ops->solve = 0; 316 ksp->ops->setup = 0; 317 ksp->ops->destroy = 0; 318 319 ksp->data = 0; 320 ksp->nwork = 0; 321 ksp->work = 0; 322 323 ksp->cnvP = 0; 324 325 ksp->reason = KSP_CONVERGED_ITERATING; 326 327 ksp->setupcalled = 0; 328 ierr = PetscPublishAll(ksp);CHKERRQ(ierr); 329 ierr = PCCreate(comm,&ksp->pc);CHKERRQ(ierr); 330 PetscFunctionReturn(0); 331 } 332 333 #undef __FUNCT__ 334 #define __FUNCT__ "KSPSetType" 335 /*@C 336 KSPSetType - Builds KSP for a particular solver. 337 338 Collective on KSP 339 340 Input Parameters: 341 + ksp - the Krylov space context 342 - type - a known method 343 344 Options Database Key: 345 . -ksp_type <method> - Sets the method; use -help for a list 346 of available methods (for instance, cg or gmres) 347 348 Notes: 349 See "petsc/include/petscksp.h" for available methods (for instance, 350 KSPCG or KSPGMRES). 351 352 Normally, it is best to use the KSPSetFromOptions() command and 353 then set the KSP type from the options database rather than by using 354 this routine. Using the options database provides the user with 355 maximum flexibility in evaluating the many different Krylov methods. 356 The KSPSetType() routine is provided for those situations where it 357 is necessary to set the iterative solver independently of the command 358 line or options database. This might be the case, for example, when 359 the choice of iterative solver changes during the execution of the 360 program, and the user's application is taking responsibility for 361 choosing the appropriate method. In other words, this routine is 362 not for beginners. 363 364 Level: intermediate 365 366 .keywords: KSP, set, method 367 368 .seealso: PCSetType(), KSPType 369 370 @*/ 371 PetscErrorCode KSPSetType(KSP ksp,const KSPType type) 372 { 373 PetscErrorCode ierr,(*r)(KSP); 374 PetscTruth match; 375 376 PetscFunctionBegin; 377 PetscValidHeaderSpecific(ksp,KSP_COOKIE,1); 378 PetscValidCharPointer(type,2); 379 380 ierr = PetscTypeCompare((PetscObject)ksp,type,&match);CHKERRQ(ierr); 381 if (match) PetscFunctionReturn(0); 382 383 if (ksp->data) { 384 /* destroy the old private KSP context */ 385 ierr = (*ksp->ops->destroy)(ksp);CHKERRQ(ierr); 386 ksp->data = 0; 387 } 388 /* Get the function pointers for the iterative method requested */ 389 if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 390 ierr = PetscFListFind(ksp->comm,KSPList,type,(void (**)(void)) &r);CHKERRQ(ierr); 391 if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown KSP type given: %s",type); 392 ksp->setupcalled = 0; 393 ierr = (*r)(ksp);CHKERRQ(ierr); 394 ierr = PetscObjectChangeTypeName((PetscObject)ksp,type);CHKERRQ(ierr); 395 PetscFunctionReturn(0); 396 } 397 398 #undef __FUNCT__ 399 #define __FUNCT__ "KSPRegisterDestroy" 400 /*@C 401 KSPRegisterDestroy - Frees the list of KSP methods that were 402 registered by KSPRegisterDynamic(). 403 404 Not Collective 405 406 Level: advanced 407 408 .keywords: KSP, register, destroy 409 410 .seealso: KSPRegisterDynamic(), KSPRegisterAll() 411 @*/ 412 PetscErrorCode KSPRegisterDestroy(void) 413 { 414 PetscErrorCode ierr; 415 416 PetscFunctionBegin; 417 if (KSPList) { 418 ierr = PetscFListDestroy(&KSPList);CHKERRQ(ierr); 419 KSPList = 0; 420 } 421 KSPRegisterAllCalled = PETSC_FALSE; 422 PetscFunctionReturn(0); 423 } 424 425 #undef __FUNCT__ 426 #define __FUNCT__ "KSPGetType" 427 /*@C 428 KSPGetType - Gets the KSP type as a string from the KSP object. 429 430 Not Collective 431 432 Input Parameter: 433 . ksp - Krylov context 434 435 Output Parameter: 436 . name - name of KSP method 437 438 Level: intermediate 439 440 .keywords: KSP, get, method, name 441 442 .seealso: KSPSetType() 443 @*/ 444 PetscErrorCode KSPGetType(KSP ksp,KSPType *type) 445 { 446 PetscFunctionBegin; 447 PetscValidHeaderSpecific(ksp,KSP_COOKIE,1); 448 PetscValidPointer(type,2); 449 *type = ksp->type_name; 450 PetscFunctionReturn(0); 451 } 452 453 #undef __FUNCT__ 454 #define __FUNCT__ "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 abstol - 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_norm_type - none - 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 $ preconditioned - default for left preconditioning 476 $ unpreconditioned - see KSPSetNormType() 477 $ natural - see KSPSetNormType() 478 . -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space 479 . -ksp_test_null_space - tests the null space set with KSPSetNullSpace() to see if it truly is a null space 480 . -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side 481 . -ksp_cancelmonitors - cancel all previous convergene monitor routines set 482 . -ksp_monitor - print residual norm at each iteration 483 . -ksp_xmonitor - plot residual norm at each iteration 484 . -ksp_vecmonitor - plot solution at each iteration 485 - -ksp_singmonitor - monitor extremem singular values at each iteration 486 487 Notes: 488 To see all options, run your program with the -help option 489 or consult the users manual. 490 491 Level: beginner 492 493 .keywords: KSP, set, from, options, database 494 495 .seealso: 496 @*/ 497 PetscErrorCode KSPSetFromOptions(KSP ksp) 498 { 499 PetscErrorCode ierr; 500 PetscInt indx; 501 char type[256]; 502 const char *stype[] = {"none","preconditioned","unpreconditioned","natural"}; 503 PetscTruth flg; 504 505 PetscFunctionBegin; 506 PetscValidHeaderSpecific(ksp,KSP_COOKIE,1); 507 ierr = PCSetFromOptions(ksp->pc);CHKERRQ(ierr); 508 509 if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 510 ierr = PetscOptionsBegin(ksp->comm,ksp->prefix,"Krylov Method (KSP) Options","KSP");CHKERRQ(ierr); 511 ierr = PetscOptionsList("-ksp_type","Krylov method","KSPSetType",KSPList,(char*)(ksp->type_name?ksp->type_name:KSPGMRES),type,256,&flg);CHKERRQ(ierr); 512 if (flg) { 513 ierr = KSPSetType(ksp,type);CHKERRQ(ierr); 514 } 515 /* 516 Set the type if it was never set. 517 */ 518 if (!ksp->type_name) { 519 ierr = KSPSetType(ksp,KSPGMRES);CHKERRQ(ierr); 520 } 521 522 ierr = PetscOptionsInt("-ksp_max_it","Maximum number of iterations","KSPSetTolerances",ksp->max_it,&ksp->max_it,PETSC_NULL);CHKERRQ(ierr); 523 ierr = PetscOptionsReal("-ksp_rtol","Relative decrease in residual norm","KSPSetTolerances",ksp->rtol,&ksp->rtol,PETSC_NULL);CHKERRQ(ierr); 524 ierr = PetscOptionsReal("-ksp_atol","Absolute value of residual norm","KSPSetTolerances",ksp->abstol,&ksp->abstol,PETSC_NULL);CHKERRQ(ierr); 525 ierr = PetscOptionsReal("-ksp_divtol","Residual norm increase cause divergence","KSPSetTolerances",ksp->divtol,&ksp->divtol,PETSC_NULL);CHKERRQ(ierr); 526 ierr = PetscOptionsLogical("-ksp_knoll","Use preconditioner applied to b for initial guess","KSPSetInitialGuessKnoll",ksp->guess_knoll, 527 &ksp->guess_knoll,PETSC_NULL);CHKERRQ(ierr); 528 529 ierr = PetscOptionsEList("-ksp_norm_type","KSP Norm type","KSPSetNormType",stype,4,"preconditioned",&indx,&flg);CHKERRQ(ierr); 530 if (flg) { 531 switch (indx) { 532 case 0: 533 ierr = KSPSetNormType(ksp,KSP_NO_NORM);CHKERRQ(ierr); 534 ierr = KSPSetConvergenceTest(ksp,KSPSkipConverged,0);CHKERRQ(ierr); 535 break; 536 case 1: 537 ierr = KSPSetNormType(ksp,KSP_PRECONDITIONED_NORM);CHKERRQ(ierr); 538 break; 539 case 2: 540 ierr = KSPSetNormType(ksp,KSP_UNPRECONDITIONED_NORM);CHKERRQ(ierr); 541 break; 542 case 3: 543 ierr = KSPSetNormType(ksp,KSP_NATURAL_NORM);CHKERRQ(ierr); 544 break; 545 } 546 } 547 548 ierr = PetscOptionsName("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",&flg);CHKERRQ(ierr); 549 if (flg) { 550 ierr = KSPSetDiagonalScale(ksp,PETSC_TRUE);CHKERRQ(ierr); 551 } 552 ierr = PetscOptionsName("-ksp_diagonal_scale_fix","Fix diagonaled scaled matrix after solve","KSPSetDiagonalScaleFix",&flg);CHKERRQ(ierr); 553 if (flg) { 554 ierr = KSPSetDiagonalScaleFix(ksp,PETSC_TRUE);CHKERRQ(ierr); 555 } 556 557 558 ierr = PetscOptionsName("-ksp_constant_null_space","Add constant null space to Krylov solver","KSPSetNullSpace",&flg);CHKERRQ(ierr); 559 if (flg) { 560 MatNullSpace nsp; 561 562 ierr = MatNullSpaceCreate(ksp->comm,PETSC_TRUE,0,0,&nsp);CHKERRQ(ierr); 563 ierr = KSPSetNullSpace(ksp,nsp);CHKERRQ(ierr); 564 ierr = MatNullSpaceDestroy(nsp);CHKERRQ(ierr); 565 } 566 567 /* option is actually checked in KSPSetUp() */ 568 if (ksp->nullsp) { 569 ierr = PetscOptionsName("-ksp_test_null_space","Is provided null space correct","None",&flg);CHKERRQ(ierr); 570 } 571 572 /* 573 Prints reason for convergence or divergence of each linear solve 574 */ 575 ierr = PetscOptionsName("-ksp_converged_reason","Print reason for converged or diverged","KSPSolve",&flg);CHKERRQ(ierr); 576 if (flg) { 577 ksp->printreason = PETSC_TRUE; 578 } 579 580 ierr = PetscOptionsName("-ksp_cancelmonitors","Remove any hardwired monitor routines","KSPClearMonitor",&flg);CHKERRQ(ierr); 581 /* -----------------------------------------------------------------------*/ 582 /* 583 Cancels all monitors hardwired into code before call to KSPSetFromOptions() 584 */ 585 if (flg) { 586 ierr = KSPClearMonitor(ksp);CHKERRQ(ierr); 587 } 588 /* 589 Prints preconditioned residual norm at each iteration 590 */ 591 ierr = PetscOptionsName("-ksp_monitor","Monitor preconditioned residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr); 592 if (flg) { 593 ierr = KSPSetMonitor(ksp,KSPDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 594 } 595 /* 596 Plots the vector solution 597 */ 598 ierr = PetscOptionsName("-ksp_vecmonitor","Monitor solution graphically","KSPSetMonitor",&flg);CHKERRQ(ierr); 599 if (flg) { 600 ierr = KSPSetMonitor(ksp,KSPVecViewMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 601 } 602 /* 603 Prints preconditioned and true residual norm at each iteration 604 */ 605 ierr = PetscOptionsName("-ksp_truemonitor","Monitor true (unpreconditioned) residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr); 606 if (flg) { 607 ierr = KSPSetMonitor(ksp,KSPTrueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 608 } 609 /* 610 Prints extreme eigenvalue estimates at each iteration 611 */ 612 ierr = PetscOptionsName("-ksp_singmonitor","Monitor singular values","KSPSetMonitor",&flg);CHKERRQ(ierr); 613 if (flg) { 614 ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); 615 ierr = KSPSetMonitor(ksp,KSPSingularValueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 616 } 617 /* 618 Prints preconditioned residual norm with fewer digits 619 */ 620 ierr = PetscOptionsName("-ksp_smonitor","Monitor preconditioned residual norm with fewer digitis","KSPSetMonitor",&flg);CHKERRQ(ierr); 621 if (flg) { 622 ierr = KSPSetMonitor(ksp,KSPDefaultSMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 623 } 624 /* 625 Graphically plots preconditioned residual norm 626 */ 627 ierr = PetscOptionsName("-ksp_xmonitor","Monitor graphically preconditioned residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr); 628 if (flg) { 629 ierr = KSPSetMonitor(ksp,KSPLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 630 } 631 /* 632 Graphically plots preconditioned and true residual norm 633 */ 634 ierr = PetscOptionsName("-ksp_xtruemonitor","Monitor graphically true residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr); 635 if (flg){ 636 ierr = KSPSetMonitor(ksp,KSPLGTrueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 637 } 638 639 /* -----------------------------------------------------------------------*/ 640 641 ierr = PetscOptionsLogicalGroupBegin("-ksp_left_pc","Use left preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr); 642 if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_LEFT);CHKERRQ(ierr); } 643 ierr = PetscOptionsLogicalGroup("-ksp_right_pc","Use right preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr); 644 if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_RIGHT);CHKERRQ(ierr);} 645 ierr = PetscOptionsLogicalGroupEnd("-ksp_symmetric_pc","Use symmetric (factorized) preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr); 646 if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_SYMMETRIC);CHKERRQ(ierr);} 647 648 ierr = PetscOptionsName("-ksp_compute_singularvalues","Compute singular values of preconditioned operator","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr); 649 if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); } 650 ierr = PetscOptionsName("-ksp_compute_eigenvalues","Compute eigenvalues of preconditioned operator","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr); 651 if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); } 652 ierr = PetscOptionsName("-ksp_plot_eigenvalues","Scatter plot extreme eigenvalues","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr); 653 if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); } 654 655 if (ksp->ops->setfromoptions) { 656 ierr = (*ksp->ops->setfromoptions)(ksp);CHKERRQ(ierr); 657 } 658 ierr = PetscOptionsName("-ksp_view","View linear solver parameters","KSPView",&flg);CHKERRQ(ierr); 659 ierr = PetscOptionsEnd();CHKERRQ(ierr); 660 PetscFunctionReturn(0); 661 } 662 663 #undef __FUNCT__ 664 #define __FUNCT__ "KSPRegister" 665 /*@C 666 KSPRegister - See KSPRegisterDynamic() 667 668 Level: advanced 669 @*/ 670 PetscErrorCode KSPRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(KSP)) 671 { 672 PetscErrorCode ierr; 673 char fullname[PETSC_MAX_PATH_LEN]; 674 675 PetscFunctionBegin; 676 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 677 ierr = PetscFListAdd(&KSPList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 678 PetscFunctionReturn(0); 679 } 680 681 #undef __FUNCT__ 682 #define __FUNCT__ "KSPSetNullSpace" 683 /*@C 684 KSPSetNullSpace - Sets the null space of the operator 685 686 Collective on KSP 687 688 Input Parameters: 689 + ksp - the Krylov space object 690 - nullsp - the null space of the operator 691 692 Level: advanced 693 694 .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPGetNullSpace() 695 @*/ 696 PetscErrorCode KSPSetNullSpace(KSP ksp,MatNullSpace nullsp) 697 { 698 PetscErrorCode ierr; 699 700 PetscFunctionBegin; 701 if (ksp->nullsp) { 702 ierr = MatNullSpaceDestroy(ksp->nullsp);CHKERRQ(ierr); 703 } 704 ksp->nullsp = nullsp; 705 ierr = PetscObjectReference((PetscObject)ksp->nullsp);CHKERRQ(ierr); 706 PetscFunctionReturn(0); 707 } 708 709 #undef __FUNCT__ 710 #define __FUNCT__ "KSPGetNullSpace" 711 /*@C 712 KSPGetNullSpace - Gets the null space of the operator 713 714 Collective on KSP 715 716 Input Parameters: 717 + ksp - the Krylov space object 718 - nullsp - the null space of the operator 719 720 Level: advanced 721 722 .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPSetNullSpace() 723 @*/ 724 PetscErrorCode KSPGetNullSpace(KSP ksp,MatNullSpace *nullsp) 725 { 726 PetscFunctionBegin; 727 *nullsp = ksp->nullsp; 728 PetscFunctionReturn(0); 729 } 730 731