1 /* 2 The basic KSP routines, Create, View etc. are here. 3 */ 4 #include <petsc/private/kspimpl.h> /*I "petscksp.h" I*/ 5 6 /* Logging support */ 7 PetscClassId KSP_CLASSID; 8 PetscClassId DMKSP_CLASSID; 9 PetscClassId KSPGUESS_CLASSID; 10 PetscLogEvent KSP_GMRESOrthogonalization, KSP_SetUp, KSP_Solve, KSP_SolveTranspose, KSP_MatSolve, KSP_MatSolveTranspose; 11 12 /* 13 Contains the list of registered KSP routines 14 */ 15 PetscFunctionList KSPList = NULL; 16 PetscBool KSPRegisterAllCalled = PETSC_FALSE; 17 18 /* 19 Contains the list of registered KSP monitors 20 */ 21 PetscFunctionList KSPMonitorList = NULL; 22 PetscFunctionList KSPMonitorCreateList = NULL; 23 PetscFunctionList KSPMonitorDestroyList = NULL; 24 PetscBool KSPMonitorRegisterAllCalled = PETSC_FALSE; 25 26 /*@ 27 KSPLoad - Loads a `KSP` that has been stored in a `PETSCVIEWERBINARY` with `KSPView()`. 28 29 Collective 30 31 Input Parameters: 32 + newdm - the newly loaded `KSP`, this needs to have been created with `KSPCreate()` or 33 some related function before a call to `KSPLoad()`. 34 - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` 35 36 Level: intermediate 37 38 Note: 39 The type is determined by the data in the file, any type set into the `KSP` before this call is ignored. 40 41 .seealso: [](ch_ksp), `KSP`, `PetscViewerBinaryOpen()`, `KSPView()`, `MatLoad()`, `VecLoad()` 42 @*/ 43 PetscErrorCode KSPLoad(KSP newdm, PetscViewer viewer) 44 { 45 PetscBool isbinary; 46 PetscInt classid; 47 char type[256]; 48 PC pc; 49 50 PetscFunctionBegin; 51 PetscValidHeaderSpecific(newdm, KSP_CLASSID, 1); 52 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 53 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 54 PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 55 56 PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT)); 57 PetscCheck(classid == KSP_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not KSP next in file"); 58 PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR)); 59 PetscCall(KSPSetType(newdm, type)); 60 PetscTryTypeMethod(newdm, load, viewer); 61 PetscCall(KSPGetPC(newdm, &pc)); 62 PetscCall(PCLoad(pc, viewer)); 63 PetscFunctionReturn(PETSC_SUCCESS); 64 } 65 66 #include <petscdraw.h> 67 #if defined(PETSC_HAVE_SAWS) 68 #include <petscviewersaws.h> 69 #endif 70 /*@ 71 KSPView - Prints the various parameters currently set in the `KSP` object. For example, the convergence tolerances and `KSPType`. 72 Also views the `PC` and `Mat` contained by the `KSP` with `PCView()` and `MatView()`. 73 74 Collective 75 76 Input Parameters: 77 + ksp - the Krylov space context 78 - viewer - visualization context 79 80 Options Database Key: 81 . -ksp_view - print the `KSP` data structure at the end of each `KSPSolve()` call 82 83 Level: beginner 84 85 Notes: 86 The available visualization contexts include 87 + `PETSC_VIEWER_STDOUT_SELF` - standard output (default) 88 - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard 89 output where only the first processor opens 90 the file. All other processors send their 91 data to the first processor to print. 92 93 The available formats include 94 + `PETSC_VIEWER_DEFAULT` - standard output (default) 95 - `PETSC_VIEWER_ASCII_INFO_DETAIL` - more verbose output for `PCBJACOBI` and `PCASM` 96 97 The user can open an alternative visualization context with 98 `PetscViewerASCIIOpen()` - output to a specified file. 99 100 Use `KSPViewFromOptions()` to allow the user to select many different `PetscViewerType` and formats from the options database. 101 102 In the debugger you can do call `KSPView(ksp,0)` to display the `KSP`. (The same holds for any PETSc object viewer). 103 104 .seealso: [](ch_ksp), `KSP`, `PetscViewer`, `PCView()`, `PetscViewerASCIIOpen()`, `KSPViewFromOptions()` 105 @*/ 106 PetscErrorCode KSPView(KSP ksp, PetscViewer viewer) 107 { 108 PetscBool isascii, isbinary, isdraw, isstring; 109 #if defined(PETSC_HAVE_SAWS) 110 PetscBool issaws; 111 #endif 112 113 PetscFunctionBegin; 114 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 115 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)ksp), &viewer)); 116 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 117 PetscCheckSameComm(ksp, 1, viewer, 2); 118 119 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 120 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 121 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 122 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 123 #if defined(PETSC_HAVE_SAWS) 124 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws)); 125 #endif 126 if (isascii) { 127 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ksp, viewer)); 128 PetscCall(PetscViewerASCIIPushTab(viewer)); 129 PetscTryTypeMethod(ksp, view, viewer); 130 PetscCall(PetscViewerASCIIPopTab(viewer)); 131 if (ksp->guess_zero) { 132 PetscCall(PetscViewerASCIIPrintf(viewer, " maximum iterations=%" PetscInt_FMT ", initial guess is zero\n", ksp->max_it)); 133 } else { 134 PetscCall(PetscViewerASCIIPrintf(viewer, " maximum iterations=%" PetscInt_FMT ", nonzero initial guess\n", ksp->max_it)); 135 } 136 if (ksp->min_it) PetscCall(PetscViewerASCIIPrintf(viewer, " minimum iterations=%" PetscInt_FMT "\n", ksp->min_it)); 137 if (ksp->guess_knoll) PetscCall(PetscViewerASCIIPrintf(viewer, " using preconditioner applied to right-hand side for initial guess\n")); 138 PetscCall(PetscViewerASCIIPrintf(viewer, " tolerances: relative=%g, absolute=%g, divergence=%g\n", (double)ksp->rtol, (double)ksp->abstol, (double)ksp->divtol)); 139 if (ksp->pc_side == PC_RIGHT) { 140 PetscCall(PetscViewerASCIIPrintf(viewer, " right preconditioning\n")); 141 } else if (ksp->pc_side == PC_SYMMETRIC) { 142 PetscCall(PetscViewerASCIIPrintf(viewer, " symmetric preconditioning\n")); 143 } else { 144 PetscCall(PetscViewerASCIIPrintf(viewer, " left preconditioning\n")); 145 } 146 if (ksp->guess) { 147 PetscCall(PetscViewerASCIIPushTab(viewer)); 148 PetscCall(KSPGuessView(ksp->guess, viewer)); 149 PetscCall(PetscViewerASCIIPopTab(viewer)); 150 } 151 if (ksp->dscale) PetscCall(PetscViewerASCIIPrintf(viewer, " diagonally scaled system\n")); 152 PetscCall(PetscViewerASCIIPrintf(viewer, " using %s norm type for convergence test\n", KSPNormTypes[ksp->normtype])); 153 } else if (isbinary) { 154 PetscInt classid = KSP_FILE_CLASSID; 155 MPI_Comm comm; 156 PetscMPIInt rank; 157 char type[256]; 158 159 PetscCall(PetscObjectGetComm((PetscObject)ksp, &comm)); 160 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 161 if (rank == 0) { 162 PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT)); 163 PetscCall(PetscStrncpy(type, ((PetscObject)ksp)->type_name, 256)); 164 PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR)); 165 } 166 PetscTryTypeMethod(ksp, view, viewer); 167 } else if (isstring) { 168 const char *type; 169 PetscCall(KSPGetType(ksp, &type)); 170 PetscCall(PetscViewerStringSPrintf(viewer, " KSPType: %-7.7s", type)); 171 PetscTryTypeMethod(ksp, view, viewer); 172 } else if (isdraw) { 173 PetscDraw draw; 174 char str[36]; 175 PetscReal x, y, bottom, h; 176 PetscBool flg; 177 178 PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 179 PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 180 PetscCall(PetscObjectTypeCompare((PetscObject)ksp, KSPPREONLY, &flg)); 181 if (!flg) { 182 PetscCall(PetscStrncpy(str, "KSP: ", sizeof(str))); 183 PetscCall(PetscStrlcat(str, ((PetscObject)ksp)->type_name, sizeof(str))); 184 PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h)); 185 bottom = y - h; 186 } else { 187 bottom = y; 188 } 189 PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 190 #if defined(PETSC_HAVE_SAWS) 191 } else if (issaws) { 192 PetscMPIInt rank; 193 const char *name; 194 195 PetscCall(PetscObjectGetName((PetscObject)ksp, &name)); 196 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 197 if (!((PetscObject)ksp)->amsmem && rank == 0) { 198 char dir[1024]; 199 200 PetscCall(PetscObjectViewSAWs((PetscObject)ksp, viewer)); 201 PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Objects/%s/its", name)); 202 PetscCallSAWs(SAWs_Register, (dir, &ksp->its, 1, SAWs_READ, SAWs_INT)); 203 if (!ksp->res_hist) PetscCall(KSPSetResidualHistory(ksp, NULL, PETSC_DECIDE, PETSC_TRUE)); 204 PetscCall(PetscSNPrintf(dir, 1024, "/PETSc/Objects/%s/res_hist", name)); 205 PetscCallSAWs(SAWs_Register, (dir, ksp->res_hist, 10, SAWs_READ, SAWs_DOUBLE)); 206 } 207 #endif 208 } else PetscTryTypeMethod(ksp, view, viewer); 209 if (ksp->pc) PetscCall(PCView(ksp->pc, viewer)); 210 if (isdraw) { 211 PetscDraw draw; 212 PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 213 PetscCall(PetscDrawPopCurrentPoint(draw)); 214 } 215 PetscFunctionReturn(PETSC_SUCCESS); 216 } 217 218 /*@ 219 KSPViewFromOptions - View (print) a `KSP` object based on values in the options database. Also views the `PC` and `Mat` contained by the `KSP` 220 with `PCView()` and `MatView()`. 221 222 Collective 223 224 Input Parameters: 225 + A - Krylov solver context 226 . obj - Optional object that provides the options prefix used to query the options database 227 - name - command line option 228 229 Level: intermediate 230 231 .seealso: [](ch_ksp), `KSP`, `KSPView()`, `PetscObjectViewFromOptions()`, `KSPCreate()` 232 @*/ 233 PetscErrorCode KSPViewFromOptions(KSP A, PetscObject obj, const char name[]) 234 { 235 PetscFunctionBegin; 236 PetscValidHeaderSpecific(A, KSP_CLASSID, 1); 237 PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 238 PetscFunctionReturn(PETSC_SUCCESS); 239 } 240 241 /*@ 242 KSPSetNormType - Sets the type of residual norm that is used for convergence testing in `KSPSolve()` for the given `KSP` context 243 244 Logically Collective 245 246 Input Parameters: 247 + ksp - Krylov solver context 248 - normtype - one of 249 .vb 250 KSP_NORM_NONE - skips computing the norm, this should generally only be used if you are using 251 the Krylov method as a smoother with a fixed small number of iterations. 252 Implicitly sets `KSPConvergedSkip()` as the `KSP` convergence test. 253 Note that certain algorithms such as `KSPGMRES` ALWAYS require the norm calculation, 254 for these methods the norms are still computed, they are just not used in 255 the convergence test. 256 KSP_NORM_PRECONDITIONED - the default for left-preconditioned solves, uses the 2-norm 257 of the preconditioned residual $B^{-1}(b - A x)$. 258 KSP_NORM_UNPRECONDITIONED - uses the 2-norm of the true $b - Ax$ residual. 259 KSP_NORM_NATURAL - uses the $A$ norm of the true $b - Ax$ residual; supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS` 260 .ve 261 262 Options Database Key: 263 . -ksp_norm_type <none,preconditioned,unpreconditioned,natural> - set `KSP` norm type 264 265 Level: advanced 266 267 Notes: 268 The norm is always of the equations residual $\| b - A x^n \|$ (or an approximation to that norm), they are never a norm of the error in the equation. 269 270 Not all combinations of preconditioner side (see `KSPSetPCSide()`) and norm types are supported by all Krylov methods. 271 If only one is set, PETSc tries to automatically change the other to find a compatible pair. If no such combination 272 is supported, PETSc will generate an error. 273 274 Developer Note: 275 Supported combinations of norm and preconditioner side are set using `KSPSetSupportedNorm()` for each `KSPType`. 276 277 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetCheckNormIteration()`, `KSPSetPCSide()`, `KSPGetPCSide()`, `KSPNormType` 278 @*/ 279 PetscErrorCode KSPSetNormType(KSP ksp, KSPNormType normtype) 280 { 281 PetscFunctionBegin; 282 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 283 PetscValidLogicalCollectiveEnum(ksp, normtype, 2); 284 ksp->normtype = ksp->normtype_set = normtype; 285 PetscFunctionReturn(PETSC_SUCCESS); 286 } 287 288 /*@ 289 KSPSetCheckNormIteration - Sets the first iteration at which the norm of the residual will be 290 computed and used in the convergence test of `KSPSolve()` for the given `KSP` context 291 292 Logically Collective 293 294 Input Parameters: 295 + ksp - Krylov solver context 296 - it - use -1 to check at all iterations 297 298 Level: advanced 299 300 Notes: 301 Currently only works with `KSPCG`, `KSPBCGS` and `KSPIBCGS` 302 303 Use `KSPSetNormType`(ksp,`KSP_NORM_NONE`) to never check the norm 304 305 On steps where the norm is not computed, the previous norm is still in the variable, so if you run with, for example, 306 `-ksp_monitor` the residual norm will appear to be unchanged for several iterations (though it is not really unchanged). 307 308 Certain methods such as `KSPGMRES` always compute the residual norm, this routine will not change that computation, but it will 309 prevent the computed norm from being checked. 310 311 .seealso: [](ch_ksp), `KSP`, `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetNormType()`, `KSPSetLagNorm()` 312 @*/ 313 PetscErrorCode KSPSetCheckNormIteration(KSP ksp, PetscInt it) 314 { 315 PetscFunctionBegin; 316 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 317 PetscValidLogicalCollectiveInt(ksp, it, 2); 318 ksp->chknorm = it; 319 PetscFunctionReturn(PETSC_SUCCESS); 320 } 321 322 /*@ 323 KSPSetLagNorm - Lags the residual norm calculation so that it is computed as part of the `MPI_Allreduce()` used for 324 computing the inner products needed for the next iteration. 325 326 Logically Collective 327 328 Input Parameters: 329 + ksp - Krylov solver context 330 - flg - `PETSC_TRUE` or `PETSC_FALSE` 331 332 Options Database Key: 333 . -ksp_lag_norm - lag the calculated residual norm 334 335 Level: advanced 336 337 Notes: 338 Currently only works with `KSPIBCGS`. 339 340 This can reduce communication costs at the expense of doing 341 one additional iteration because the norm used in the convergence test of `KSPSolve()` is one iteration behind the actual 342 current residual norm (which has not yet been computed due to the lag). 343 344 Use `KSPSetNormType`(ksp,`KSP_NORM_NONE`) to never check the norm 345 346 If you lag the norm and run with, for example, `-ksp_monitor`, the residual norm reported will be the lagged one. 347 348 `KSPSetCheckNormIteration()` is an alternative way of avoiding the expense of computing the residual norm at each iteration. 349 350 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSPConvergedSkip()`, `KSPSetNormType()`, `KSPSetCheckNormIteration()` 351 @*/ 352 PetscErrorCode KSPSetLagNorm(KSP ksp, PetscBool flg) 353 { 354 PetscFunctionBegin; 355 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 356 PetscValidLogicalCollectiveBool(ksp, flg, 2); 357 ksp->lagnorm = flg; 358 PetscFunctionReturn(PETSC_SUCCESS); 359 } 360 361 /*@ 362 KSPSetSupportedNorm - Sets a norm and preconditioner side supported by a `KSPType` 363 364 Logically Collective 365 366 Input Parameters: 367 + ksp - Krylov method 368 . normtype - supported norm type of the type `KSPNormType` 369 . pcside - preconditioner side, of the type `PCSide` that can be used with this `KSPNormType` 370 - priority - positive integer preference for this combination; larger values have higher priority 371 372 Level: developer 373 374 Notes: 375 This function should be called from the implementation files `KSPCreate_XXX()` to declare 376 which norms and preconditioner sides are supported. Users should not call this 377 function. 378 379 This function can be called multiple times for each combination of `KSPNormType` and `PCSide` 380 the `KSPType` supports 381 382 .seealso: [](ch_ksp), `KSP`, `KSPNormType`, `PCSide`, `KSPSetNormType()`, `KSPSetPCSide()` 383 @*/ 384 PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType normtype, PCSide pcside, PetscInt priority) 385 { 386 PetscFunctionBegin; 387 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 388 ksp->normsupporttable[normtype][pcside] = priority; 389 PetscFunctionReturn(PETSC_SUCCESS); 390 } 391 392 static PetscErrorCode KSPNormSupportTableReset_Private(KSP ksp) 393 { 394 PetscFunctionBegin; 395 PetscCall(PetscMemzero(ksp->normsupporttable, sizeof(ksp->normsupporttable))); 396 ksp->pc_side = ksp->pc_side_set; 397 ksp->normtype = ksp->normtype_set; 398 PetscFunctionReturn(PETSC_SUCCESS); 399 } 400 401 PetscErrorCode KSPSetUpNorms_Private(KSP ksp, PetscBool errorifnotsupported, KSPNormType *normtype, PCSide *pcside) 402 { 403 PetscInt i, j, best, ibest = 0, jbest = 0; 404 405 PetscFunctionBegin; 406 best = 0; 407 for (i = 0; i < KSP_NORM_MAX; i++) { 408 for (j = 0; j < PC_SIDE_MAX; j++) { 409 if ((ksp->normtype == KSP_NORM_DEFAULT || ksp->normtype == i) && (ksp->pc_side == PC_SIDE_DEFAULT || ksp->pc_side == j) && ksp->normsupporttable[i][j] > best) { 410 best = ksp->normsupporttable[i][j]; 411 ibest = i; 412 jbest = j; 413 } 414 } 415 } 416 if (best < 1 && errorifnotsupported) { 417 PetscCheck(ksp->normtype != KSP_NORM_DEFAULT || ksp->pc_side != PC_SIDE_DEFAULT, PetscObjectComm((PetscObject)ksp), PETSC_ERR_PLIB, "The %s KSP implementation did not call KSPSetSupportedNorm()", ((PetscObject)ksp)->type_name); 418 PetscCheck(ksp->normtype != KSP_NORM_DEFAULT, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSP %s does not support preconditioner side %s", ((PetscObject)ksp)->type_name, PCSides[ksp->pc_side]); 419 PetscCheck(ksp->pc_side != PC_SIDE_DEFAULT, PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSP %s does not support norm type %s", ((PetscObject)ksp)->type_name, KSPNormTypes[ksp->normtype]); 420 SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_SUP, "KSP %s does not support norm type %s with preconditioner side %s", ((PetscObject)ksp)->type_name, KSPNormTypes[ksp->normtype], PCSides[ksp->pc_side]); 421 } 422 if (normtype) *normtype = (KSPNormType)ibest; 423 if (pcside) *pcside = (PCSide)jbest; 424 PetscFunctionReturn(PETSC_SUCCESS); 425 } 426 427 /*@ 428 KSPGetNormType - Gets the `KSPNormType` that is used for convergence testing during `KSPSolve()` for this `KSP` context 429 430 Not Collective 431 432 Input Parameter: 433 . ksp - Krylov solver context 434 435 Output Parameter: 436 . normtype - the `KSPNormType` that is used for convergence testing 437 438 Level: advanced 439 440 .seealso: [](ch_ksp), `KSPNormType`, `KSPSetNormType()`, `KSPConvergedSkip()` 441 @*/ 442 PetscErrorCode KSPGetNormType(KSP ksp, KSPNormType *normtype) 443 { 444 PetscFunctionBegin; 445 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 446 PetscAssertPointer(normtype, 2); 447 PetscCall(KSPSetUpNorms_Private(ksp, PETSC_TRUE, &ksp->normtype, &ksp->pc_side)); 448 *normtype = ksp->normtype; 449 PetscFunctionReturn(PETSC_SUCCESS); 450 } 451 452 #if defined(PETSC_HAVE_SAWS) 453 #include <petscviewersaws.h> 454 #endif 455 456 /*@ 457 KSPSetOperators - Sets the matrix associated with the linear system 458 and a (possibly) different one from which the preconditioner will be built into the `KSP` context. The matrix will then be used during `KSPSolve()` 459 460 Collective 461 462 Input Parameters: 463 + ksp - the `KSP` context 464 . Amat - the matrix that defines the linear system 465 - Pmat - the matrix to be used in constructing the preconditioner, usually the same as `Amat`. 466 467 Level: beginner 468 469 Notes: 470 .vb 471 KSPSetOperators(ksp, Amat, Pmat); 472 .ve 473 is the same as 474 .vb 475 KSPGetPC(ksp, &pc); 476 PCSetOperators(pc, Amat, Pmat); 477 .ve 478 and is equivalent to 479 .vb 480 PCCreate(PetscObjectComm((PetscObject)ksp), &pc); 481 PCSetOperators(pc, Amat, Pmat); 482 KSPSetPC(ksp, pc); 483 .ve 484 485 If you know the operator `Amat` has a null space you can use `MatSetNullSpace()` and `MatSetTransposeNullSpace()` to supply the null 486 space to `Amat` and the `KSP` solvers will automatically use that null space as needed during the solution process. 487 488 All future calls to `KSPSetOperators()` must use the same size matrices, unless `KSPReset()` is called! 489 490 Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently being used from the `KSP` context. 491 492 If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then 493 first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()` 494 on it and then pass it back in your call to `KSPSetOperators()`. 495 496 Developer Notes: 497 If the operators have NOT been set with `KSPSetOperators()` then the operators 498 are created in the `PC` and returned to the user. In this case, if both operators 499 mat and pmat are requested, two DIFFERENT operators will be returned. If 500 only one is requested both operators in the `PC` will be the same (i.e. as 501 if one had called `KSPSetOperators()` with the same argument for both `Mat`s). 502 The user must set the sizes of the returned matrices and their type etc just 503 as if the user created them with `MatCreate()`. For example, 504 505 .vb 506 KSPGetOperators(ksp/pc,&mat,NULL); is equivalent to 507 set size, type, etc of mat 508 509 MatCreate(comm,&mat); 510 KSP/PCSetOperators(ksp/pc,mat,mat); 511 PetscObjectDereference((PetscObject)mat); 512 set size, type, etc of mat 513 514 and 515 516 KSP/PCGetOperators(ksp/pc,&mat,&pmat); is equivalent to 517 set size, type, etc of mat and pmat 518 519 MatCreate(comm,&mat); 520 MatCreate(comm,&pmat); 521 KSP/PCSetOperators(ksp/pc,mat,pmat); 522 PetscObjectDereference((PetscObject)mat); 523 PetscObjectDereference((PetscObject)pmat); 524 set size, type, etc of mat and pmat 525 .ve 526 527 The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy 528 of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely 529 managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look 530 at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to 531 the `SNES` object (the `SNES` object manages it for you). Similarly when you create a `KSP` 532 you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you). 533 Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when 534 it can be created for you? 535 536 .seealso: [](ch_ksp), `KSP`, `Mat`, `KSPSolve()`, `KSPGetPC()`, `PCGetOperators()`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetComputeOperators()`, `KSPSetComputeInitialGuess()`, `KSPSetComputeRHS()` 537 @*/ 538 PetscErrorCode KSPSetOperators(KSP ksp, Mat Amat, Mat Pmat) 539 { 540 PetscFunctionBegin; 541 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 542 if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 543 if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3); 544 if (Amat) PetscCheckSameComm(ksp, 1, Amat, 2); 545 if (Pmat) PetscCheckSameComm(ksp, 1, Pmat, 3); 546 if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc)); 547 PetscCall(PCSetOperators(ksp->pc, Amat, Pmat)); 548 if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX; /* so that next solve call will call PCSetUp() on new matrix */ 549 PetscFunctionReturn(PETSC_SUCCESS); 550 } 551 552 /*@ 553 KSPGetOperators - Gets the matrix associated with the linear system 554 and a (possibly) different one used to construct the preconditioner from the `KSP` context 555 556 Collective 557 558 Input Parameter: 559 . ksp - the `KSP` context 560 561 Output Parameters: 562 + Amat - the matrix that defines the linear system 563 - Pmat - the matrix to be used in constructing the preconditioner, usually the same as `Amat`. 564 565 Level: intermediate 566 567 Notes: 568 If `KSPSetOperators()` has not been called then the `KSP` object will attempt to automatically create the matrix `Amat` and return it 569 570 Use `KSPGetOperatorsSet()` to determine if matrices have been provided. 571 572 DOES NOT increase the reference counts of the matrix, so you should NOT destroy them. 573 574 .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetPC()`, `PCSetOperators()`, `KSPSetOperators()`, `KSPGetOperatorsSet()` 575 @*/ 576 PetscErrorCode KSPGetOperators(KSP ksp, Mat *Amat, Mat *Pmat) 577 { 578 PetscFunctionBegin; 579 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 580 if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc)); 581 PetscCall(PCGetOperators(ksp->pc, Amat, Pmat)); 582 PetscFunctionReturn(PETSC_SUCCESS); 583 } 584 585 /*@ 586 KSPGetOperatorsSet - Determines if the matrix associated with the linear system and 587 possibly a different one from which the preconditioner will be built have been set in the `KSP` with `KSPSetOperators()` 588 589 Not Collective, though the results on all processes will be the same 590 591 Input Parameter: 592 . ksp - the `KSP` context 593 594 Output Parameters: 595 + mat - the matrix associated with the linear system was set 596 - pmat - matrix from which the preconditioner will be built, usually the same as `mat` was set 597 598 Level: intermediate 599 600 Note: 601 This routine exists because if you call `KSPGetOperators()` on a `KSP` that does not yet have operators they are 602 automatically created in the call. 603 604 .seealso: [](ch_ksp), `KSP`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()`, `PCGetOperatorsSet()` 605 @*/ 606 PetscErrorCode KSPGetOperatorsSet(KSP ksp, PetscBool *mat, PetscBool *pmat) 607 { 608 PetscFunctionBegin; 609 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 610 if (!ksp->pc) PetscCall(KSPGetPC(ksp, &ksp->pc)); 611 PetscCall(PCGetOperatorsSet(ksp->pc, mat, pmat)); 612 PetscFunctionReturn(PETSC_SUCCESS); 613 } 614 615 /*@C 616 KSPSetPreSolve - Sets a function that is called at the beginning of each `KSPSolve()`. Used in conjunction with `KSPSetPostSolve()`. 617 618 Logically Collective 619 620 Input Parameters: 621 + ksp - the solver object 622 . presolve - the function to call before the solve, see` KSPPSolveFn` 623 - ctx - an optional context needed by the function 624 625 Level: developer 626 627 Notes: 628 The function provided here `presolve` is used to modify the right hand side, and possibly the matrix, of the linear system to be solved. 629 The function provided with `KSPSetPostSolve()` then modifies the resulting solution of that linear system to obtain the correct solution 630 to the initial linear system. 631 632 The functions `PCPreSolve()` and `PCPostSolve()` provide a similar functionality and are used, for example with `PCEISENSTAT`. 633 634 .seealso: [](ch_ksp), `KSPPSolveFn`, `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPSetPostSolve()`, `PCEISENSTAT`, `PCPreSolve()`, `PCPostSolve()` 635 @*/ 636 PetscErrorCode KSPSetPreSolve(KSP ksp, KSPPSolveFn *presolve, void *ctx) 637 { 638 PetscFunctionBegin; 639 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 640 ksp->presolve = presolve; 641 ksp->prectx = ctx; 642 PetscFunctionReturn(PETSC_SUCCESS); 643 } 644 645 /*@C 646 KSPSetPostSolve - Sets a function that is called at the end of each `KSPSolve()` (whether it converges or not). Used in conjunction with `KSPSetPreSolve()`. 647 648 Logically Collective 649 650 Input Parameters: 651 + ksp - the solver object 652 . postsolve - the function to call after the solve, see` KSPPSolveFn` 653 - ctx - an optional context needed by the function 654 655 Level: developer 656 657 .seealso: [](ch_ksp), `KSPPSolveFn`, `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPSetPreSolve()`, `PCEISENSTAT` 658 @*/ 659 PetscErrorCode KSPSetPostSolve(KSP ksp, KSPPSolveFn *postsolve, void *ctx) 660 { 661 PetscFunctionBegin; 662 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 663 ksp->postsolve = postsolve; 664 ksp->postctx = ctx; 665 PetscFunctionReturn(PETSC_SUCCESS); 666 } 667 668 /*@ 669 KSPSetNestLevel - sets the amount of nesting the `KSP` has. That is the number of levels of `KSP` above this `KSP` in a linear solve. 670 671 Collective 672 673 Input Parameters: 674 + ksp - the `KSP` 675 - level - the nest level 676 677 Level: developer 678 679 Note: 680 For example, the `KSP` in each block of a `KSPBJACOBI` has a level of 1, while the outer `KSP` has a level of 0. 681 682 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCSetKSPNestLevel()`, `PCGetKSPNestLevel()` 683 @*/ 684 PetscErrorCode KSPSetNestLevel(KSP ksp, PetscInt level) 685 { 686 PetscFunctionBegin; 687 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 688 PetscValidLogicalCollectiveInt(ksp, level, 2); 689 ksp->nestlevel = level; 690 PetscFunctionReturn(PETSC_SUCCESS); 691 } 692 693 /*@ 694 KSPGetNestLevel - gets the amount of nesting the `KSP` has 695 696 Not Collective 697 698 Input Parameter: 699 . ksp - the `KSP` 700 701 Output Parameter: 702 . level - the nest level 703 704 Level: developer 705 706 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `PCGetKSPNestLevel()` 707 @*/ 708 PetscErrorCode KSPGetNestLevel(KSP ksp, PetscInt *level) 709 { 710 PetscFunctionBegin; 711 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 712 PetscAssertPointer(level, 2); 713 *level = ksp->nestlevel; 714 PetscFunctionReturn(PETSC_SUCCESS); 715 } 716 717 /*@ 718 KSPCreate - Creates the `KSP` context. This `KSP` context is used in PETSc to solve linear systems with `KSPSolve()` 719 720 Collective 721 722 Input Parameter: 723 . comm - MPI communicator 724 725 Output Parameter: 726 . inksp - location to put the `KSP` context 727 728 Level: beginner 729 730 Note: 731 The default `KSPType` is `KSPGMRES` with a restart of 30, using modified Gram-Schmidt orthogonalization. The `KSPType` may be 732 changed with `KSPSetType()` 733 734 .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetType()` 735 @*/ 736 PetscErrorCode KSPCreate(MPI_Comm comm, KSP *inksp) 737 { 738 KSP ksp; 739 void *ctx; 740 741 PetscFunctionBegin; 742 PetscAssertPointer(inksp, 2); 743 PetscCall(KSPInitializePackage()); 744 745 PetscCall(PetscHeaderCreate(ksp, KSP_CLASSID, "KSP", "Krylov Method", "KSP", comm, KSPDestroy, KSPView)); 746 ksp->default_max_it = ksp->max_it = 10000; 747 ksp->pc_side = ksp->pc_side_set = PC_SIDE_DEFAULT; 748 749 ksp->default_rtol = ksp->rtol = 1.e-5; 750 ksp->default_abstol = ksp->abstol = PetscDefined(USE_REAL_SINGLE) ? 1.e-25 : 1.e-50; 751 ksp->default_divtol = ksp->divtol = 1.e4; 752 753 ksp->chknorm = -1; 754 ksp->normtype = ksp->normtype_set = KSP_NORM_DEFAULT; 755 ksp->rnorm = 0.0; 756 ksp->its = 0; 757 ksp->guess_zero = PETSC_TRUE; 758 ksp->calc_sings = PETSC_FALSE; 759 ksp->res_hist = NULL; 760 ksp->res_hist_alloc = NULL; 761 ksp->res_hist_len = 0; 762 ksp->res_hist_max = 0; 763 ksp->res_hist_reset = PETSC_TRUE; 764 ksp->err_hist = NULL; 765 ksp->err_hist_alloc = NULL; 766 ksp->err_hist_len = 0; 767 ksp->err_hist_max = 0; 768 ksp->err_hist_reset = PETSC_TRUE; 769 ksp->numbermonitors = 0; 770 ksp->numberreasonviews = 0; 771 ksp->setfromoptionscalled = 0; 772 ksp->nmax = PETSC_DECIDE; 773 774 PetscCall(KSPConvergedDefaultCreate(&ctx)); 775 PetscCall(KSPSetConvergenceTest(ksp, KSPConvergedDefault, ctx, KSPConvergedDefaultDestroy)); 776 ksp->ops->buildsolution = KSPBuildSolutionDefault; 777 ksp->ops->buildresidual = KSPBuildResidualDefault; 778 779 ksp->vec_sol = NULL; 780 ksp->vec_rhs = NULL; 781 ksp->pc = NULL; 782 ksp->data = NULL; 783 ksp->nwork = 0; 784 ksp->work = NULL; 785 ksp->reason = KSP_CONVERGED_ITERATING; 786 ksp->setupstage = KSP_SETUP_NEW; 787 788 PetscCall(KSPNormSupportTableReset_Private(ksp)); 789 790 *inksp = ksp; 791 PetscFunctionReturn(PETSC_SUCCESS); 792 } 793 794 /*@ 795 KSPSetType - Sets the algorithm/method to be used to solve the linear system with the given `KSP` 796 797 Logically Collective 798 799 Input Parameters: 800 + ksp - the Krylov space context 801 - type - a known method 802 803 Options Database Key: 804 . -ksp_type <method> - Sets the method; see `KSPType` or use `-help` for a list of available methods (for instance, cg or gmres) 805 806 Level: intermediate 807 808 Notes: 809 See `KSPType` for available methods (for instance, `KSPCG` or `KSPGMRES`). 810 811 Normally, it is best to use the `KSPSetFromOptions()` command and 812 then set the `KSP` type from the options database rather than by using 813 this routine. Using the options database provides the user with 814 maximum flexibility in evaluating the many different Krylov methods. 815 The `KSPSetType()` routine is provided for those situations where it 816 is necessary to set the iterative solver independently of the command 817 line or options database. This might be the case, for example, when 818 the choice of iterative solver changes during the execution of the 819 program, and the user's application is taking responsibility for 820 choosing the appropriate method. In other words, this routine is 821 not for beginners. 822 823 Developer Note: 824 `KSPRegister()` is used to add Krylov types to `KSPList` from which they are accessed by `KSPSetType()`. 825 826 .seealso: [](ch_ksp), `PCSetType()`, `KSPType`, `KSPRegister()`, `KSPCreate()`, `KSP` 827 @*/ 828 PetscErrorCode KSPSetType(KSP ksp, KSPType type) 829 { 830 PetscBool match; 831 PetscErrorCode (*r)(KSP); 832 833 PetscFunctionBegin; 834 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 835 PetscAssertPointer(type, 2); 836 837 PetscCall(PetscObjectTypeCompare((PetscObject)ksp, type, &match)); 838 if (match) PetscFunctionReturn(PETSC_SUCCESS); 839 840 PetscCall(PetscFunctionListFind(KSPList, type, &r)); 841 PetscCheck(r, PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unable to find requested KSP type %s", type); 842 /* Destroy the previous private KSP context */ 843 PetscTryTypeMethod(ksp, destroy); 844 845 /* Reinitialize function pointers in KSPOps structure */ 846 PetscCall(PetscMemzero(ksp->ops, sizeof(struct _KSPOps))); 847 ksp->ops->buildsolution = KSPBuildSolutionDefault; 848 ksp->ops->buildresidual = KSPBuildResidualDefault; 849 PetscCall(KSPNormSupportTableReset_Private(ksp)); 850 ksp->converged_neg_curve = PETSC_FALSE; // restore default 851 ksp->setupnewmatrix = PETSC_FALSE; // restore default (setup not called in case of new matrix) 852 /* Call the KSPCreate_XXX routine for this particular Krylov solver */ 853 ksp->setupstage = KSP_SETUP_NEW; 854 ksp->guess_not_read = PETSC_FALSE; // restore default 855 PetscCall((*r)(ksp)); 856 PetscCall(PetscObjectChangeTypeName((PetscObject)ksp, type)); 857 PetscFunctionReturn(PETSC_SUCCESS); 858 } 859 860 /*@ 861 KSPGetType - Gets the `KSP` type as a string from the `KSP` object. 862 863 Not Collective 864 865 Input Parameter: 866 . ksp - Krylov context 867 868 Output Parameter: 869 . type - name of the `KSP` method 870 871 Level: intermediate 872 873 .seealso: [](ch_ksp), `KSPType`, `KSP`, `KSPSetType()` 874 @*/ 875 PetscErrorCode KSPGetType(KSP ksp, KSPType *type) 876 { 877 PetscFunctionBegin; 878 PetscValidHeaderSpecific(ksp, KSP_CLASSID, 1); 879 PetscAssertPointer(type, 2); 880 *type = ((PetscObject)ksp)->type_name; 881 PetscFunctionReturn(PETSC_SUCCESS); 882 } 883 884 /*@C 885 KSPRegister - Adds a method, `KSPType`, to the Krylov subspace solver package. 886 887 Not Collective, No Fortran Support 888 889 Input Parameters: 890 + sname - name of a new user-defined solver 891 - function - routine to create method 892 893 Level: advanced 894 895 Note: 896 `KSPRegister()` may be called multiple times to add several user-defined solvers. 897 898 Example Usage: 899 .vb 900 KSPRegister("my_solver", MySolverCreate); 901 .ve 902 903 Then, your solver can be chosen with the procedural interface via 904 .vb 905 KSPSetType(ksp, "my_solver") 906 .ve 907 or at runtime via the option `-ksp_type my_solver` 908 909 .seealso: [](ch_ksp), `KSP`, `KSPType`, `KSPSetType`, `KSPRegisterAll()` 910 @*/ 911 PetscErrorCode KSPRegister(const char sname[], PetscErrorCode (*function)(KSP)) 912 { 913 PetscFunctionBegin; 914 PetscCall(KSPInitializePackage()); 915 PetscCall(PetscFunctionListAdd(&KSPList, sname, function)); 916 PetscFunctionReturn(PETSC_SUCCESS); 917 } 918 919 PetscErrorCode KSPMonitorMakeKey_Internal(const char name[], PetscViewerType vtype, PetscViewerFormat format, char key[]) 920 { 921 PetscFunctionBegin; 922 PetscCall(PetscStrncpy(key, name, PETSC_MAX_PATH_LEN)); 923 PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN)); 924 PetscCall(PetscStrlcat(key, vtype, PETSC_MAX_PATH_LEN)); 925 PetscCall(PetscStrlcat(key, ":", PETSC_MAX_PATH_LEN)); 926 PetscCall(PetscStrlcat(key, PetscViewerFormats[format], PETSC_MAX_PATH_LEN)); 927 PetscFunctionReturn(PETSC_SUCCESS); 928 } 929 930 /*@C 931 KSPMonitorRegister - Registers a Krylov subspace solver monitor routine that may be accessed with `KSPMonitorSetFromOptions()` 932 933 Not Collective 934 935 Input Parameters: 936 + name - name of a new monitor type 937 . vtype - A `PetscViewerType` for the output 938 . format - A `PetscViewerFormat` for the output 939 . monitor - Monitor routine, see `KSPMonitorRegisterFn` 940 . create - Creation routine, or `NULL` 941 - destroy - Destruction routine, or `NULL` 942 943 Level: advanced 944 945 Notes: 946 `KSPMonitorRegister()` may be called multiple times to add several user-defined monitors. 947 948 The calling sequence for the given function matches the calling sequence used by `KSPMonitorFn` functions passed to `KSPMonitorSet()` with the additional 949 requirement that its final argument be a `PetscViewerAndFormat`. 950 951 Example Usage: 952 .vb 953 KSPMonitorRegister("my_monitor", PETSCVIEWERASCII, PETSC_VIEWER_ASCII_INFO_DETAIL, MyMonitor, NULL, NULL); 954 .ve 955 956 Then, your monitor can be chosen with the procedural interface via 957 .vb 958 KSPMonitorSetFromOptions(ksp, "-ksp_monitor_my_monitor", "my_monitor", NULL) 959 .ve 960 or at runtime via the option `-ksp_monitor_my_monitor` 961 962 .seealso: [](ch_ksp), `KSP`, `KSPMonitorSet()`, `KSPMonitorRegisterAll()`, `KSPMonitorSetFromOptions()` 963 @*/ 964 PetscErrorCode KSPMonitorRegister(const char name[], PetscViewerType vtype, PetscViewerFormat format, KSPMonitorRegisterFn *monitor, KSPMonitorRegisterCreateFn *create, KSPMonitorRegisterDestroyFn *destroy) 965 { 966 char key[PETSC_MAX_PATH_LEN]; 967 968 PetscFunctionBegin; 969 PetscCall(KSPInitializePackage()); 970 PetscCall(KSPMonitorMakeKey_Internal(name, vtype, format, key)); 971 PetscCall(PetscFunctionListAdd(&KSPMonitorList, key, monitor)); 972 if (create) PetscCall(PetscFunctionListAdd(&KSPMonitorCreateList, key, create)); 973 if (destroy) PetscCall(PetscFunctionListAdd(&KSPMonitorDestroyList, key, destroy)); 974 PetscFunctionReturn(PETSC_SUCCESS); 975 } 976