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