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 PetscLogObjectCreate(ksp); 286 *inksp = ksp; 287 ksp->bops->publish = KSPPublish_Petsc; 288 289 ksp->type = -1; 290 ksp->max_it = 10000; 291 ksp->pc_side = PC_LEFT; 292 ksp->rtol = 1.e-5; 293 ksp->abstol = 1.e-50; 294 ksp->divtol = 1.e4; 295 296 ksp->normtype = KSP_PRECONDITIONED_NORM; 297 ksp->rnorm = 0.0; 298 ksp->its = 0; 299 ksp->guess_zero = PETSC_TRUE; 300 ksp->calc_sings = PETSC_FALSE; 301 ksp->res_hist = PETSC_NULL; 302 ksp->res_hist_len = 0; 303 ksp->res_hist_max = 0; 304 ksp->res_hist_reset = PETSC_TRUE; 305 ksp->numbermonitors = 0; 306 ksp->converged = KSPDefaultConverged; 307 ksp->ops->buildsolution = KSPDefaultBuildSolution; 308 ksp->ops->buildresidual = KSPDefaultBuildResidual; 309 310 ksp->ops->setfromoptions = 0; 311 312 ksp->vec_sol = 0; 313 ksp->vec_rhs = 0; 314 ksp->pc = 0; 315 316 ksp->ops->solve = 0; 317 ksp->ops->setup = 0; 318 ksp->ops->destroy = 0; 319 320 ksp->data = 0; 321 ksp->nwork = 0; 322 ksp->work = 0; 323 324 ksp->cnvP = 0; 325 326 ksp->reason = KSP_CONVERGED_ITERATING; 327 328 ksp->setupcalled = 0; 329 ierr = PetscPublishAll(ksp);CHKERRQ(ierr); 330 ierr = PCCreate(comm,&ksp->pc);CHKERRQ(ierr); 331 PetscFunctionReturn(0); 332 } 333 334 #undef __FUNCT__ 335 #define __FUNCT__ "KSPSetType" 336 /*@C 337 KSPSetType - Builds KSP for a particular solver. 338 339 Collective on KSP 340 341 Input Parameters: 342 + ksp - the Krylov space context 343 - type - a known method 344 345 Options Database Key: 346 . -ksp_type <method> - Sets the method; use -help for a list 347 of available methods (for instance, cg or gmres) 348 349 Notes: 350 See "petsc/include/petscksp.h" for available methods (for instance, 351 KSPCG or KSPGMRES). 352 353 Normally, it is best to use the KSPSetFromOptions() command and 354 then set the KSP type from the options database rather than by using 355 this routine. Using the options database provides the user with 356 maximum flexibility in evaluating the many different Krylov methods. 357 The KSPSetType() routine is provided for those situations where it 358 is necessary to set the iterative solver independently of the command 359 line or options database. This might be the case, for example, when 360 the choice of iterative solver changes during the execution of the 361 program, and the user's application is taking responsibility for 362 choosing the appropriate method. In other words, this routine is 363 not for beginners. 364 365 Level: intermediate 366 367 .keywords: KSP, set, method 368 369 .seealso: PCSetType(), KSPType 370 371 @*/ 372 PetscErrorCode KSPSetType(KSP ksp,const KSPType type) 373 { 374 PetscErrorCode ierr,(*r)(KSP); 375 PetscTruth match; 376 377 PetscFunctionBegin; 378 PetscValidHeaderSpecific(ksp,KSP_COOKIE,1); 379 PetscValidCharPointer(type,2); 380 381 ierr = PetscTypeCompare((PetscObject)ksp,type,&match);CHKERRQ(ierr); 382 if (match) PetscFunctionReturn(0); 383 384 if (ksp->data) { 385 /* destroy the old private KSP context */ 386 ierr = (*ksp->ops->destroy)(ksp);CHKERRQ(ierr); 387 ksp->data = 0; 388 } 389 /* Get the function pointers for the iterative method requested */ 390 if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 391 ierr = PetscFListFind(ksp->comm,KSPList,type,(void (**)(void)) &r);CHKERRQ(ierr); 392 if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown KSP type given: %s",type); 393 ksp->setupcalled = 0; 394 ierr = (*r)(ksp);CHKERRQ(ierr); 395 ierr = PetscObjectChangeTypeName((PetscObject)ksp,type);CHKERRQ(ierr); 396 PetscFunctionReturn(0); 397 } 398 399 #undef __FUNCT__ 400 #define __FUNCT__ "KSPRegisterDestroy" 401 /*@C 402 KSPRegisterDestroy - Frees the list of KSP methods that were 403 registered by KSPRegisterDynamic(). 404 405 Not Collective 406 407 Level: advanced 408 409 .keywords: KSP, register, destroy 410 411 .seealso: KSPRegisterDynamic(), KSPRegisterAll() 412 @*/ 413 PetscErrorCode KSPRegisterDestroy(void) 414 { 415 PetscErrorCode ierr; 416 417 PetscFunctionBegin; 418 if (KSPList) { 419 ierr = PetscFListDestroy(&KSPList);CHKERRQ(ierr); 420 KSPList = 0; 421 } 422 KSPRegisterAllCalled = PETSC_FALSE; 423 PetscFunctionReturn(0); 424 } 425 426 #undef __FUNCT__ 427 #define __FUNCT__ "KSPGetType" 428 /*@C 429 KSPGetType - Gets the KSP type as a string from the KSP object. 430 431 Not Collective 432 433 Input Parameter: 434 . ksp - Krylov context 435 436 Output Parameter: 437 . name - name of KSP method 438 439 Level: intermediate 440 441 .keywords: KSP, get, method, name 442 443 .seealso: KSPSetType() 444 @*/ 445 PetscErrorCode KSPGetType(KSP ksp,KSPType *type) 446 { 447 PetscFunctionBegin; 448 PetscValidHeaderSpecific(ksp,KSP_COOKIE,1); 449 PetscValidPointer(type,2); 450 *type = ksp->type_name; 451 PetscFunctionReturn(0); 452 } 453 454 #undef __FUNCT__ 455 #define __FUNCT__ "KSPSetFromOptions" 456 /*@ 457 KSPSetFromOptions - Sets KSP options from the options database. 458 This routine must be called before KSPSetUp() if the user is to be 459 allowed to set the Krylov type. 460 461 Collective on KSP 462 463 Input Parameters: 464 . ksp - the Krylov space context 465 466 Options Database Keys: 467 + -ksp_max_it - maximum number of linear iterations 468 . -ksp_rtol rtol - relative tolerance used in default determination of convergence, i.e. 469 if residual norm decreases by this factor than convergence is declared 470 . -ksp_atol abstol - absolute tolerance used in default convergence test, i.e. if residual 471 norm is less than this then convergence is declared 472 . -ksp_divtol tol - if residual norm increases by this factor than divergence is declared 473 . -ksp_norm_type - none - skip norms used in convergence tests (useful only when not using 474 $ convergence test (say you always want to run with 5 iterations) to 475 $ save on communication overhead 476 $ preconditioned - default for left preconditioning 477 $ unpreconditioned - see KSPSetNormType() 478 $ natural - see KSPSetNormType() 479 . -ksp_constant_null_space - assume the operator (matrix) has the constant vector in its null space 480 . -ksp_test_null_space - tests the null space set with KSPSetNullSpace() to see if it truly is a null space 481 . -ksp_knoll - compute initial guess by applying the preconditioner to the right hand side 482 . -ksp_cancelmonitors - cancel all previous convergene monitor routines set 483 . -ksp_monitor - print residual norm at each iteration 484 . -ksp_xmonitor - plot residual norm at each iteration 485 . -ksp_vecmonitor - plot solution at each iteration 486 - -ksp_singmonitor - monitor extremem singular values at each iteration 487 488 Notes: 489 To see all options, run your program with the -help option 490 or consult the users manual. 491 492 Level: beginner 493 494 .keywords: KSP, set, from, options, database 495 496 .seealso: 497 @*/ 498 PetscErrorCode KSPSetFromOptions(KSP ksp) 499 { 500 PetscErrorCode ierr; 501 PetscInt indx; 502 char type[256]; 503 const char *stype[] = {"none","preconditioned","unpreconditioned","natural"}; 504 PetscTruth flg; 505 506 PetscFunctionBegin; 507 PetscValidHeaderSpecific(ksp,KSP_COOKIE,1); 508 ierr = PCSetFromOptions(ksp->pc);CHKERRQ(ierr); 509 510 if (!KSPRegisterAllCalled) {ierr = KSPRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 511 ierr = PetscOptionsBegin(ksp->comm,ksp->prefix,"Krylov Method (KSP) Options","KSP");CHKERRQ(ierr); 512 ierr = PetscOptionsList("-ksp_type","Krylov method","KSPSetType",KSPList,(char*)(ksp->type_name?ksp->type_name:KSPGMRES),type,256,&flg);CHKERRQ(ierr); 513 if (flg) { 514 ierr = KSPSetType(ksp,type);CHKERRQ(ierr); 515 } 516 /* 517 Set the type if it was never set. 518 */ 519 if (!ksp->type_name) { 520 ierr = KSPSetType(ksp,KSPGMRES);CHKERRQ(ierr); 521 } 522 523 ierr = PetscOptionsInt("-ksp_max_it","Maximum number of iterations","KSPSetTolerances",ksp->max_it,&ksp->max_it,PETSC_NULL);CHKERRQ(ierr); 524 ierr = PetscOptionsReal("-ksp_rtol","Relative decrease in residual norm","KSPSetTolerances",ksp->rtol,&ksp->rtol,PETSC_NULL);CHKERRQ(ierr); 525 ierr = PetscOptionsReal("-ksp_atol","Absolute value of residual norm","KSPSetTolerances",ksp->abstol,&ksp->abstol,PETSC_NULL);CHKERRQ(ierr); 526 ierr = PetscOptionsReal("-ksp_divtol","Residual norm increase cause divergence","KSPSetTolerances",ksp->divtol,&ksp->divtol,PETSC_NULL);CHKERRQ(ierr); 527 ierr = PetscOptionsLogical("-ksp_knoll","Use preconditioner applied to b for initial guess","KSPSetInitialGuessKnoll",ksp->guess_knoll, 528 &ksp->guess_knoll,PETSC_NULL);CHKERRQ(ierr); 529 530 ierr = PetscOptionsEList("-ksp_norm_type","KSP Norm type","KSPSetNormType",stype,4,"preconditioned",&indx,&flg);CHKERRQ(ierr); 531 if (flg) { 532 switch (indx) { 533 case 0: 534 ierr = KSPSetNormType(ksp,KSP_NO_NORM);CHKERRQ(ierr); 535 ierr = KSPSetConvergenceTest(ksp,KSPSkipConverged,0);CHKERRQ(ierr); 536 break; 537 case 1: 538 ierr = KSPSetNormType(ksp,KSP_PRECONDITIONED_NORM);CHKERRQ(ierr); 539 break; 540 case 2: 541 ierr = KSPSetNormType(ksp,KSP_UNPRECONDITIONED_NORM);CHKERRQ(ierr); 542 break; 543 case 3: 544 ierr = KSPSetNormType(ksp,KSP_NATURAL_NORM);CHKERRQ(ierr); 545 break; 546 } 547 } 548 549 ierr = PetscOptionsName("-ksp_diagonal_scale","Diagonal scale matrix before building preconditioner","KSPSetDiagonalScale",&flg);CHKERRQ(ierr); 550 if (flg) { 551 ierr = KSPSetDiagonalScale(ksp,PETSC_TRUE);CHKERRQ(ierr); 552 } 553 ierr = PetscOptionsName("-ksp_diagonal_scale_fix","Fix diagonaled scaled matrix after solve","KSPSetDiagonalScaleFix",&flg);CHKERRQ(ierr); 554 if (flg) { 555 ierr = KSPSetDiagonalScaleFix(ksp,PETSC_TRUE);CHKERRQ(ierr); 556 } 557 558 559 ierr = PetscOptionsName("-ksp_constant_null_space","Add constant null space to Krylov solver","KSPSetNullSpace",&flg);CHKERRQ(ierr); 560 if (flg) { 561 MatNullSpace nsp; 562 563 ierr = MatNullSpaceCreate(ksp->comm,PETSC_TRUE,0,0,&nsp);CHKERRQ(ierr); 564 ierr = KSPSetNullSpace(ksp,nsp);CHKERRQ(ierr); 565 ierr = MatNullSpaceDestroy(nsp);CHKERRQ(ierr); 566 } 567 568 /* option is actually checked in KSPSetUp() */ 569 if (ksp->nullsp) { 570 ierr = PetscOptionsName("-ksp_test_null_space","Is provided null space correct","None",&flg);CHKERRQ(ierr); 571 } 572 573 /* 574 Prints reason for convergence or divergence of each linear solve 575 */ 576 ierr = PetscOptionsName("-ksp_converged_reason","Print reason for converged or diverged","KSPSolve",&flg);CHKERRQ(ierr); 577 if (flg) { 578 ksp->printreason = PETSC_TRUE; 579 } 580 581 ierr = PetscOptionsName("-ksp_cancelmonitors","Remove any hardwired monitor routines","KSPClearMonitor",&flg);CHKERRQ(ierr); 582 /* -----------------------------------------------------------------------*/ 583 /* 584 Cancels all monitors hardwired into code before call to KSPSetFromOptions() 585 */ 586 if (flg) { 587 ierr = KSPClearMonitor(ksp);CHKERRQ(ierr); 588 } 589 /* 590 Prints preconditioned residual norm at each iteration 591 */ 592 ierr = PetscOptionsName("-ksp_monitor","Monitor preconditioned residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr); 593 if (flg) { 594 ierr = KSPSetMonitor(ksp,KSPDefaultMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 595 } 596 /* 597 Plots the vector solution 598 */ 599 ierr = PetscOptionsName("-ksp_vecmonitor","Monitor solution graphically","KSPSetMonitor",&flg);CHKERRQ(ierr); 600 if (flg) { 601 ierr = KSPSetMonitor(ksp,KSPVecViewMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 602 } 603 /* 604 Prints preconditioned and true residual norm at each iteration 605 */ 606 ierr = PetscOptionsName("-ksp_truemonitor","Monitor true (unpreconditioned) residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr); 607 if (flg) { 608 ierr = KSPSetMonitor(ksp,KSPTrueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 609 } 610 /* 611 Prints extreme eigenvalue estimates at each iteration 612 */ 613 ierr = PetscOptionsName("-ksp_singmonitor","Monitor singular values","KSPSetMonitor",&flg);CHKERRQ(ierr); 614 if (flg) { 615 ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); 616 ierr = KSPSetMonitor(ksp,KSPSingularValueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 617 } 618 /* 619 Prints preconditioned residual norm with fewer digits 620 */ 621 ierr = PetscOptionsName("-ksp_smonitor","Monitor preconditioned residual norm with fewer digitis","KSPSetMonitor",&flg);CHKERRQ(ierr); 622 if (flg) { 623 ierr = KSPSetMonitor(ksp,KSPDefaultSMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 624 } 625 /* 626 Graphically plots preconditioned residual norm 627 */ 628 ierr = PetscOptionsName("-ksp_xmonitor","Monitor graphically preconditioned residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr); 629 if (flg) { 630 ierr = KSPSetMonitor(ksp,KSPLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 631 } 632 /* 633 Graphically plots preconditioned and true residual norm 634 */ 635 ierr = PetscOptionsName("-ksp_xtruemonitor","Monitor graphically true residual norm","KSPSetMonitor",&flg);CHKERRQ(ierr); 636 if (flg){ 637 ierr = KSPSetMonitor(ksp,KSPLGTrueMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 638 } 639 640 /* -----------------------------------------------------------------------*/ 641 642 ierr = PetscOptionsLogicalGroupBegin("-ksp_left_pc","Use left preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr); 643 if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_LEFT);CHKERRQ(ierr); } 644 ierr = PetscOptionsLogicalGroup("-ksp_right_pc","Use right preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr); 645 if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_RIGHT);CHKERRQ(ierr);} 646 ierr = PetscOptionsLogicalGroupEnd("-ksp_symmetric_pc","Use symmetric (factorized) preconditioning","KSPSetPreconditionerSide",&flg);CHKERRQ(ierr); 647 if (flg) { ierr = KSPSetPreconditionerSide(ksp,PC_SYMMETRIC);CHKERRQ(ierr);} 648 649 ierr = PetscOptionsName("-ksp_compute_singularvalues","Compute singular values of preconditioned operator","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr); 650 if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); } 651 ierr = PetscOptionsName("-ksp_compute_eigenvalues","Compute eigenvalues of preconditioned operator","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr); 652 if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); } 653 ierr = PetscOptionsName("-ksp_plot_eigenvalues","Scatter plot extreme eigenvalues","KSPSetComputeSingularValues",&flg);CHKERRQ(ierr); 654 if (flg) { ierr = KSPSetComputeSingularValues(ksp,PETSC_TRUE);CHKERRQ(ierr); } 655 656 if (ksp->ops->setfromoptions) { 657 ierr = (*ksp->ops->setfromoptions)(ksp);CHKERRQ(ierr); 658 } 659 ierr = PetscOptionsName("-ksp_view","View linear solver parameters","KSPView",&flg);CHKERRQ(ierr); 660 ierr = PetscOptionsEnd();CHKERRQ(ierr); 661 PetscFunctionReturn(0); 662 } 663 664 #undef __FUNCT__ 665 #define __FUNCT__ "KSPRegister" 666 /*@C 667 KSPRegister - See KSPRegisterDynamic() 668 669 Level: advanced 670 @*/ 671 PetscErrorCode KSPRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(KSP)) 672 { 673 PetscErrorCode ierr; 674 char fullname[PETSC_MAX_PATH_LEN]; 675 676 PetscFunctionBegin; 677 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 678 ierr = PetscFListAdd(&KSPList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 679 PetscFunctionReturn(0); 680 } 681 682 #undef __FUNCT__ 683 #define __FUNCT__ "KSPSetNullSpace" 684 /*@C 685 KSPSetNullSpace - Sets the null space of the operator 686 687 Collective on KSP 688 689 Input Parameters: 690 + ksp - the Krylov space object 691 - nullsp - the null space of the operator 692 693 Level: advanced 694 695 .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPGetNullSpace() 696 @*/ 697 PetscErrorCode KSPSetNullSpace(KSP ksp,MatNullSpace nullsp) 698 { 699 PetscErrorCode ierr; 700 701 PetscFunctionBegin; 702 if (ksp->nullsp) { 703 ierr = MatNullSpaceDestroy(ksp->nullsp);CHKERRQ(ierr); 704 } 705 ksp->nullsp = nullsp; 706 ierr = PetscObjectReference((PetscObject)ksp->nullsp);CHKERRQ(ierr); 707 PetscFunctionReturn(0); 708 } 709 710 #undef __FUNCT__ 711 #define __FUNCT__ "KSPGetNullSpace" 712 /*@C 713 KSPGetNullSpace - Gets the null space of the operator 714 715 Collective on KSP 716 717 Input Parameters: 718 + ksp - the Krylov space object 719 - nullsp - the null space of the operator 720 721 Level: advanced 722 723 .seealso: KSPSetOperators(), MatNullSpaceCreate(), KSPSetNullSpace() 724 @*/ 725 PetscErrorCode KSPGetNullSpace(KSP ksp,MatNullSpace *nullsp) 726 { 727 PetscFunctionBegin; 728 *nullsp = ksp->nullsp; 729 PetscFunctionReturn(0); 730 } 731 732