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