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