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