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