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