1 #define PETSCSNES_DLL 2 3 #include "include/private/matimpl.h" 4 #include "src/snes/mf/snesmfj.h" /*I "petscsnes.h" I*/ 5 6 #undef __FUNCT__ 7 #define __FUNCT__ "MatMFFDComputeJacobian" 8 /*@ 9 MatMFFDComputeJacobian - Tells the matrix-free Jacobian object the new location at which 10 Jacobian matrix vector products will be computed at, i.e. J(x) * a. 11 12 Collective on SNES 13 14 Input Parameters: 15 + snes - the nonlinear solver context 16 . x - the point at which the Jacobian vector products will be performed 17 . jac - the matrix-free Jacobian object 18 . B - either the same as jac or another matrix type (ignored) 19 . flag - not relevent for matrix-free form 20 - dummy - the user context (ignored) 21 22 Level: developer 23 24 Notes: 25 This can be passed into SNESSetJacobian() when using a completely matrix-free solver, 26 that is the B matrix is also the same matrix operator. This is used when you select 27 -mat_mffd but rarely used directly by users. 28 29 .seealso: MatMFFDGetH(), MatCreateSNESMF(), 30 MatMFFDSetHHistory(), 31 MatMFFDKSPMonitor(), MatMFFDSetFunctionError(), MatMFFDCreate(), SNESSetJacobian() 32 33 @*/ 34 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDComputeJacobian(SNES *snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy) 35 { 36 PetscErrorCode ierr; 37 PetscFunctionBegin; 38 ierr = MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 39 ierr = MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 40 PetscFunctionReturn(0); 41 } 42 43 PetscFList MatSNESMPetscFList = 0; 44 PetscTruth MatSNESMFRegisterAllCalled = PETSC_FALSE; 45 46 PetscCookie PETSCSNES_DLLEXPORT MATSNESMFCTX_COOKIE = 0; 47 PetscEvent MATSNESMF_Mult = 0; 48 49 #undef __FUNCT__ 50 #define __FUNCT__ "MatSNESMFSetType" 51 /*@C 52 MatSNESMFSetType - Sets the method that is used to compute the 53 differencing parameter for finite differene matrix-free formulations. 54 55 Input Parameters: 56 + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMF() 57 or MatSetType(mat,MATMFFD); 58 - ftype - the type requested, either MATSNESMF_WP or MATSNESMF_DS 59 60 Level: advanced 61 62 Notes: 63 For example, such routines can compute h for use in 64 Jacobian-vector products of the form 65 66 F(x+ha) - F(x) 67 F'(u)a ~= ---------------- 68 h 69 70 .seealso: MatCreateSNESMF(), MatSNESMFRegisterDynamic) 71 @*/ 72 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetType(Mat mat,MatSNESMFType ftype) 73 { 74 PetscErrorCode ierr,(*r)(MatSNESMF); 75 MatSNESMF ctx = (MatSNESMF)mat->data; 76 PetscTruth match; 77 78 PetscFunctionBegin; 79 PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 80 PetscValidCharPointer(ftype,2); 81 82 /* already set, so just return */ 83 ierr = PetscTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr); 84 if (match) PetscFunctionReturn(0); 85 86 /* destroy the old one if it exists */ 87 if (ctx->ops->destroy) { 88 ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr); 89 } 90 91 ierr = PetscFListFind(ctx->comm,MatSNESMPetscFList,ftype,(void (**)(void)) &r);CHKERRQ(ierr); 92 if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatSNESMF type %s given",ftype); 93 ierr = (*r)(ctx);CHKERRQ(ierr); 94 ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr); 95 PetscFunctionReturn(0); 96 } 97 98 typedef PetscErrorCode (*FCN1)(Vec,void*); /* force argument to next function to not be extern C*/ 99 EXTERN_C_BEGIN 100 #undef __FUNCT__ 101 #define __FUNCT__ "MatSNESMFSetFunctioniBase_FD" 102 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetFunctioniBase_FD(Mat mat,FCN1 func) 103 { 104 MatSNESMF ctx = (MatSNESMF)mat->data; 105 106 PetscFunctionBegin; 107 ctx->funcisetbase = func; 108 PetscFunctionReturn(0); 109 } 110 EXTERN_C_END 111 112 typedef PetscErrorCode (*FCN2)(PetscInt,Vec,PetscScalar*,void*); /* force argument to next function to not be extern C*/ 113 EXTERN_C_BEGIN 114 #undef __FUNCT__ 115 #define __FUNCT__ "MatSNESMFSetFunctioni_FD" 116 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetFunctioni_FD(Mat mat,FCN2 funci) 117 { 118 MatSNESMF ctx = (MatSNESMF)mat->data; 119 120 PetscFunctionBegin; 121 ctx->funci = funci; 122 PetscFunctionReturn(0); 123 } 124 EXTERN_C_END 125 126 127 #undef __FUNCT__ 128 #define __FUNCT__ "MatSNESMFRegister" 129 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(MatSNESMF)) 130 { 131 PetscErrorCode ierr; 132 char fullname[PETSC_MAX_PATH_LEN]; 133 134 PetscFunctionBegin; 135 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 136 ierr = PetscFListAdd(&MatSNESMPetscFList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 137 PetscFunctionReturn(0); 138 } 139 140 141 #undef __FUNCT__ 142 #define __FUNCT__ "MatSNESMFRegisterDestroy" 143 /*@C 144 MatSNESMFRegisterDestroy - Frees the list of MatSNESMF methods that were 145 registered by MatSNESMFRegisterDynamic). 146 147 Not Collective 148 149 Level: developer 150 151 .keywords: MatSNESMF, register, destroy 152 153 .seealso: MatSNESMFRegisterDynamic), MatSNESMFRegisterAll() 154 @*/ 155 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFRegisterDestroy(void) 156 { 157 PetscErrorCode ierr; 158 159 PetscFunctionBegin; 160 ierr = PetscFListDestroy(&MatSNESMPetscFList);CHKERRQ(ierr); 161 MatSNESMFRegisterAllCalled = PETSC_FALSE; 162 PetscFunctionReturn(0); 163 } 164 165 /* ----------------------------------------------------------------------------------------*/ 166 #undef __FUNCT__ 167 #define __FUNCT__ "MatDestroy_MFFD" 168 PetscErrorCode MatDestroy_MFFD(Mat mat) 169 { 170 PetscErrorCode ierr; 171 MatSNESMF ctx = (MatSNESMF)mat->data; 172 173 PetscFunctionBegin; 174 if (ctx->w) { 175 ierr = VecDestroy(ctx->w);CHKERRQ(ierr); 176 } 177 if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);} 178 if (ctx->sp) {ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);} 179 ierr = PetscHeaderDestroy(ctx);CHKERRQ(ierr); 180 181 ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetBase_C","",PETSC_NULL);CHKERRQ(ierr); 182 ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetFunctioniBase_C","",PETSC_NULL);CHKERRQ(ierr); 183 ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetFunctioni_C","",PETSC_NULL);CHKERRQ(ierr); 184 ierr = PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetCheckh_C","",PETSC_NULL);CHKERRQ(ierr); 185 186 PetscFunctionReturn(0); 187 } 188 189 #undef __FUNCT__ 190 #define __FUNCT__ "MatView_MFFD" 191 /* 192 MatSNESMFView_MFFD - Views matrix-free parameters. 193 194 */ 195 PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer) 196 { 197 PetscErrorCode ierr; 198 MatSNESMF ctx = (MatSNESMF)J->data; 199 PetscTruth iascii; 200 201 PetscFunctionBegin; 202 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 203 if (iascii) { 204 ierr = PetscViewerASCIIPrintf(viewer," SNES matrix-free approximation:\n");CHKERRQ(ierr); 205 ierr = PetscViewerASCIIPrintf(viewer," err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr); 206 if (!ctx->type_name) { 207 ierr = PetscViewerASCIIPrintf(viewer," The compute h routine has not yet been set\n");CHKERRQ(ierr); 208 } else { 209 ierr = PetscViewerASCIIPrintf(viewer," Using %s compute h routine\n",ctx->type_name);CHKERRQ(ierr); 210 } 211 if (ctx->ops->view) { 212 ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr); 213 } 214 } else { 215 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for SNES matrix free matrix",((PetscObject)viewer)->type_name); 216 } 217 PetscFunctionReturn(0); 218 } 219 220 #undef __FUNCT__ 221 #define __FUNCT__ "MatAssemblyEnd_MFFD" 222 /* 223 MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 224 allows the user to indicate the beginning of a new linear solve by calling 225 MatAssemblyXXX() on the matrix free matrix. This then allows the 226 MatSNESMFCreate_WP() to properly compute ||U|| only the first time 227 in the linear solver rather than every time. 228 */ 229 PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt) 230 { 231 PetscErrorCode ierr; 232 MatSNESMF j = (MatSNESMF)J->data; 233 234 PetscFunctionBegin; 235 ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr); 236 if (j->usesnes) { 237 ierr = SNESGetSolution(j->snes,&j->current_u);CHKERRQ(ierr); 238 ierr = SNESGetFunction(j->snes,&j->current_f,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 239 if (!j->w) { 240 ierr = VecDuplicate(j->current_u, &j->w);CHKERRQ(ierr); 241 } 242 } 243 j->vshift = 0.0; 244 j->vscale = 1.0; 245 PetscFunctionReturn(0); 246 } 247 248 #undef __FUNCT__ 249 #define __FUNCT__ "MatMult_MFFD" 250 /* 251 MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 252 253 y ~= (F(u + ha) - F(u))/h, 254 where F = nonlinear function, as set by SNESSetFunction() 255 u = current iterate 256 h = difference interval 257 */ 258 PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y) 259 { 260 MatSNESMF ctx = (MatSNESMF)mat->data; 261 SNES snes; 262 PetscScalar h; 263 Vec w,U,F; 264 PetscErrorCode ierr,(*eval_fct)(SNES,Vec,Vec)=0; 265 PetscTruth zeroa; 266 267 PetscFunctionBegin; 268 /* We log matrix-free matrix-vector products separately, so that we can 269 separate the performance monitoring from the cases that use conventional 270 storage. We may eventually modify event logging to associate events 271 with particular objects, hence alleviating the more general problem. */ 272 ierr = PetscLogEventBegin(MATSNESMF_Mult,a,y,0,0);CHKERRQ(ierr); 273 274 snes = ctx->snes; 275 w = ctx->w; 276 U = ctx->current_u; 277 278 /* 279 Compute differencing parameter 280 */ 281 if (!ctx->ops->compute) { 282 ierr = MatSNESMFSetType(mat,MATSNESMF_WP);CHKERRQ(ierr); 283 ierr = MatSNESMFSetFromOptions(mat);CHKERRQ(ierr); 284 } 285 ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr); 286 if (zeroa) { 287 ierr = VecSet(y,0.0);CHKERRQ(ierr); 288 PetscFunctionReturn(0); 289 } 290 291 if (ctx->checkh) { 292 ierr = (*ctx->checkh)(U,a,&h,ctx->checkhctx);CHKERRQ(ierr); 293 } 294 295 /* keep a record of the current differencing parameter h */ 296 ctx->currenth = h; 297 #if defined(PETSC_USE_COMPLEX) 298 ierr = PetscInfo2(mat,"Current differencing parameter: %G + %G i\n",PetscRealPart(h),PetscImaginaryPart(h));CHKERRQ(ierr); 299 #else 300 ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr); 301 #endif 302 if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 303 ctx->historyh[ctx->ncurrenth] = h; 304 } 305 ctx->ncurrenth++; 306 307 /* w = u + ha */ 308 ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr); 309 310 if (ctx->usesnes) { 311 eval_fct = SNESComputeFunction; 312 F = ctx->current_f; 313 if (!F) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"You must call MatAssembly() even on matrix-free matrices"); 314 ierr = (*eval_fct)(snes,w,y);CHKERRQ(ierr); 315 } else { 316 F = ctx->funcvec; 317 /* compute func(U) as base for differencing */ 318 if (ctx->ncurrenth == 1) { 319 ierr = (*ctx->func)(snes,U,F,ctx->funcctx);CHKERRQ(ierr); 320 } 321 ierr = (*ctx->func)(snes,w,y,ctx->funcctx);CHKERRQ(ierr); 322 } 323 324 ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 325 ierr = VecScale(y,1.0/h);CHKERRQ(ierr); 326 327 ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr); 328 329 if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);} 330 331 ierr = PetscLogEventEnd(MATSNESMF_Mult,a,y,0,0);CHKERRQ(ierr); 332 PetscFunctionReturn(0); 333 } 334 335 #undef __FUNCT__ 336 #define __FUNCT__ "MatGetDiagonal_MFFD" 337 /* 338 MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 339 340 y ~= (F(u + ha) - F(u))/h, 341 where F = nonlinear function, as set by SNESSetFunction() 342 u = current iterate 343 h = difference interval 344 */ 345 PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 346 { 347 MatSNESMF ctx = (MatSNESMF)mat->data; 348 PetscScalar h,*aa,*ww,v; 349 PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 350 Vec w,U; 351 PetscErrorCode ierr; 352 PetscInt i,rstart,rend; 353 354 PetscFunctionBegin; 355 if (!ctx->funci) { 356 SETERRQ(PETSC_ERR_ORDER,"Requires calling MatSNESMFSetFunctioni() first"); 357 } 358 359 w = ctx->w; 360 U = ctx->current_u; 361 ierr = (*ctx->func)(0,U,a,ctx->funcctx);CHKERRQ(ierr); 362 ierr = (*ctx->funcisetbase)(U,ctx->funcctx);CHKERRQ(ierr); 363 ierr = VecCopy(U,w);CHKERRQ(ierr); 364 365 ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr); 366 ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 367 for (i=rstart; i<rend; i++) { 368 ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 369 h = ww[i-rstart]; 370 if (h == 0.0) h = 1.0; 371 #if !defined(PETSC_USE_COMPLEX) 372 if (h < umin && h >= 0.0) h = umin; 373 else if (h < 0.0 && h > -umin) h = -umin; 374 #else 375 if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 376 else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 377 #endif 378 h *= epsilon; 379 380 ww[i-rstart] += h; 381 ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 382 ierr = (*ctx->funci)(i,w,&v,ctx->funcctx);CHKERRQ(ierr); 383 aa[i-rstart] = (v - aa[i-rstart])/h; 384 385 /* possibly shift and scale result */ 386 aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart]; 387 388 ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 389 ww[i-rstart] -= h; 390 ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 391 } 392 ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 393 PetscFunctionReturn(0); 394 } 395 396 #undef __FUNCT__ 397 #define __FUNCT__ "MatShift_MFFD" 398 PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a) 399 { 400 MatSNESMF shell = (MatSNESMF)Y->data; 401 PetscFunctionBegin; 402 shell->vshift += a; 403 PetscFunctionReturn(0); 404 } 405 406 #undef __FUNCT__ 407 #define __FUNCT__ "MatScale_MFFD" 408 PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a) 409 { 410 MatSNESMF shell = (MatSNESMF)Y->data; 411 PetscFunctionBegin; 412 shell->vscale *= a; 413 PetscFunctionReturn(0); 414 } 415 416 417 #undef __FUNCT__ 418 #define __FUNCT__ "MatCreateSNESMF" 419 /*@ 420 MatCreateSNESMF - Creates a matrix-free matrix context for use with 421 a SNES solver. This matrix can be used as the Jacobian argument for 422 the routine SNESSetJacobian(). 423 424 Collective on SNES and Vec 425 426 Input Parameters: 427 + snes - the SNES context 428 - x - vector where SNES solution is to be stored. 429 430 Output Parameter: 431 . J - the matrix-free matrix 432 433 Level: advanced 434 435 Notes: 436 The matrix-free matrix context merely contains the function pointers 437 and work space for performing finite difference approximations of 438 Jacobian-vector products, F'(u)*a, 439 440 The default code uses the following approach to compute h 441 442 .vb 443 F'(u)*a = [F(u+h*a) - F(u)]/h where 444 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 445 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 446 where 447 error_rel = square root of relative error in function evaluation 448 umin = minimum iterate parameter 449 .ve 450 (see MATSNESMF_WP or MATSNESMF_DS) 451 452 The user can set the error_rel via MatSNESMFSetFunctionError() and 453 umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter 454 of the users manual for details. 455 456 The user should call MatDestroy() when finished with the matrix-free 457 matrix context. 458 459 Options Database Keys: 460 + -snes_mf_err <error_rel> - Sets error_rel 461 + -snes_mf_type - wp or ds (see MATSNESMF_WP or MATSNESMF_DS) 462 . -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only) 463 - -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h 464 465 .keywords: SNES, default, matrix-free, create, matrix 466 467 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin() 468 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateMF(), 469 MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic), MatSNESMFComputeJacobian() 470 471 @*/ 472 PetscErrorCode PETSCSNES_DLLEXPORT MatCreateSNESMF(SNES snes,Vec x,Mat *J) 473 { 474 MatSNESMF mfctx; 475 PetscErrorCode ierr; 476 477 PetscFunctionBegin; 478 ierr = MatCreateMF(x,J);CHKERRQ(ierr); 479 480 mfctx = (MatSNESMF)(*J)->data; 481 mfctx->snes = snes; 482 mfctx->usesnes = PETSC_TRUE; 483 ierr = PetscLogObjectParent(snes,*J);CHKERRQ(ierr); 484 PetscFunctionReturn(0); 485 } 486 487 EXTERN_C_BEGIN 488 #undef __FUNCT__ 489 #define __FUNCT__ "MatSNESMFSetBase_FD" 490 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetBase_FD(Mat J,Vec U) 491 { 492 PetscErrorCode ierr; 493 MatSNESMF ctx = (MatSNESMF)J->data; 494 495 PetscFunctionBegin; 496 ierr = MatSNESMFResetHHistory(J);CHKERRQ(ierr); 497 ctx->current_u = U; 498 ctx->usesnes = PETSC_FALSE; 499 if (!ctx->w) { 500 ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr); 501 } 502 J->assembled = PETSC_TRUE; 503 PetscFunctionReturn(0); 504 } 505 EXTERN_C_END 506 typedef PetscErrorCode (*FCN3)(Vec,Vec,PetscScalar*,void*); /* force argument to next function to not be extern C*/ 507 EXTERN_C_BEGIN 508 #undef __FUNCT__ 509 #define __FUNCT__ "MatSNESMFSetCheckh_FD" 510 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetCheckh_FD(Mat J,FCN3 fun,void*ectx) 511 { 512 MatSNESMF ctx = (MatSNESMF)J->data; 513 514 PetscFunctionBegin; 515 ctx->checkh = fun; 516 ctx->checkhctx = ectx; 517 PetscFunctionReturn(0); 518 } 519 EXTERN_C_END 520 521 #undef __FUNCT__ 522 #define __FUNCT__ "MatSNESMFSetFromOptions" 523 /*@ 524 MatSNESMFSetFromOptions - Sets the MatSNESMF options from the command line 525 parameter. 526 527 Collective on Mat 528 529 Input Parameters: 530 . mat - the matrix obtained with MatCreateSNESMF() 531 532 Options Database Keys: 533 + -snes_mf_type - wp or ds (see MATSNESMF_WP or MATSNESMF_DS) 534 - -snes_mf_err - square root of estimated relative error in function evaluation 535 - -snes_mf_period - how often h is recomputed, defaults to 1, everytime 536 537 Level: advanced 538 539 .keywords: SNES, matrix-free, parameters 540 541 .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(), 542 MatSNESMFResetHHistory(), MatSNESMFKSPMonitor() 543 @*/ 544 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetFromOptions(Mat mat) 545 { 546 MatSNESMF mfctx = (MatSNESMF)mat->data; 547 PetscErrorCode ierr; 548 PetscTruth flg; 549 char ftype[256]; 550 551 PetscFunctionBegin; 552 if (!MatSNESMFRegisterAllCalled) {ierr = MatSNESMFRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 553 554 ierr = PetscOptionsBegin(mfctx->comm,mfctx->prefix,"Set matrix free computation parameters","MatSNESMF");CHKERRQ(ierr); 555 ierr = PetscOptionsList("-snes_mf_type","Matrix free type","MatSNESMFSetType",MatSNESMPetscFList,mfctx->type_name,ftype,256,&flg);CHKERRQ(ierr); 556 if (flg) { 557 ierr = MatSNESMFSetType(mat,ftype);CHKERRQ(ierr); 558 } 559 560 ierr = PetscOptionsReal("-snes_mf_err","set sqrt relative error in function","MatSNESMFSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 561 ierr = PetscOptionsInt("-snes_mf_period","how often h is recomputed","MatSNESMFSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 562 if (mfctx->snes) { 563 ierr = PetscOptionsName("-snes_mf_ksp_monitor","Monitor matrix-free parameters","MatSNESMFKSPMonitor",&flg);CHKERRQ(ierr); 564 if (flg) { 565 KSP ksp; 566 ierr = SNESGetKSP(mfctx->snes,&ksp);CHKERRQ(ierr); 567 ierr = KSPMonitorSet(ksp,MatSNESMFKSPMonitor,PETSC_NULL,0);CHKERRQ(ierr); 568 } 569 } 570 ierr = PetscOptionsName("-snes_mf_check_positivity","Insure that U + h*a is nonnegative","MatSNESMFSetCheckh",&flg);CHKERRQ(ierr); 571 if (flg) { 572 ierr = MatSNESMFSetCheckh(mat,MatSNESMFCheckPositivity,0);CHKERRQ(ierr); 573 } 574 if (mfctx->ops->setfromoptions) { 575 ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr); 576 } 577 ierr = PetscOptionsEnd();CHKERRQ(ierr); 578 PetscFunctionReturn(0); 579 } 580 581 /*MC 582 MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 583 584 Level: advanced 585 586 .seealso: MatCreateMF(), MatCreateSNESMF() 587 M*/ 588 EXTERN_C_BEGIN 589 #undef __FUNCT__ 590 #define __FUNCT__ "MatCreate_MFFD" 591 PetscErrorCode PETSCSNES_DLLEXPORT MatCreate_MFFD(Mat A) 592 { 593 MatSNESMF mfctx; 594 PetscErrorCode ierr; 595 596 PetscFunctionBegin; 597 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 598 ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 599 #endif 600 601 ierr = PetscHeaderCreate(mfctx,_p_MatSNESMF,struct _MFOps,MATSNESMFCTX_COOKIE,0,"SNESMF",A->comm,MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr); 602 mfctx->sp = 0; 603 mfctx->snes = 0; 604 mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 605 mfctx->recomputeperiod = 1; 606 mfctx->count = 0; 607 mfctx->currenth = 0.0; 608 mfctx->historyh = PETSC_NULL; 609 mfctx->ncurrenth = 0; 610 mfctx->maxcurrenth = 0; 611 mfctx->type_name = 0; 612 mfctx->usesnes = PETSC_FALSE; 613 614 mfctx->vshift = 0.0; 615 mfctx->vscale = 1.0; 616 617 /* 618 Create the empty data structure to contain compute-h routines. 619 These will be filled in below from the command line options or 620 a later call with MatSNESMFSetType() or if that is not called 621 then it will default in the first use of MatMult_MFFD() 622 */ 623 mfctx->ops->compute = 0; 624 mfctx->ops->destroy = 0; 625 mfctx->ops->view = 0; 626 mfctx->ops->setfromoptions = 0; 627 mfctx->hctx = 0; 628 629 mfctx->func = 0; 630 mfctx->funcctx = 0; 631 mfctx->funcvec = 0; 632 mfctx->w = PETSC_NULL; 633 634 A->data = mfctx; 635 636 A->ops->mult = MatMult_MFFD; 637 A->ops->destroy = MatDestroy_MFFD; 638 A->ops->view = MatView_MFFD; 639 A->ops->assemblyend = MatAssemblyEnd_MFFD; 640 A->ops->getdiagonal = MatGetDiagonal_MFFD; 641 A->ops->scale = MatScale_MFFD; 642 A->ops->shift = MatShift_MFFD; 643 A->ops->setfromoptions = MatSNESMFSetFromOptions; 644 A->assembled = PETSC_TRUE; 645 646 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetBase_C","MatSNESMFSetBase_FD",MatSNESMFSetBase_FD);CHKERRQ(ierr); 647 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioniBase_C","MatSNESMFSetFunctioniBase_FD",MatSNESMFSetFunctioniBase_FD);CHKERRQ(ierr); 648 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetFunctioni_C","MatSNESMFSetFunctioni_FD",MatSNESMFSetFunctioni_FD);CHKERRQ(ierr); 649 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatSNESMFSetCheckh_C","MatSNESMFSetCheckh_FD",MatSNESMFSetCheckh_FD);CHKERRQ(ierr); 650 mfctx->mat = A; 651 ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr); 652 PetscFunctionReturn(0); 653 } 654 EXTERN_C_END 655 656 #undef __FUNCT__ 657 #define __FUNCT__ "MatCreateMF" 658 /*@ 659 MatCreateMF - Creates a matrix-free matrix. See also MatCreateSNESMF() 660 661 Collective on Vec 662 663 Input Parameters: 664 . x - vector that defines layout of the vectors and matrices 665 666 Output Parameter: 667 . J - the matrix-free matrix 668 669 Level: advanced 670 671 Notes: 672 The matrix-free matrix context merely contains the function pointers 673 and work space for performing finite difference approximations of 674 Jacobian-vector products, F'(u)*a, 675 676 The default code uses the following approach to compute h 677 678 .vb 679 F'(u)*a = [F(u+h*a) - F(u)]/h where 680 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 681 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 682 where 683 error_rel = square root of relative error in function evaluation 684 umin = minimum iterate parameter 685 .ve 686 687 The user can set the error_rel via MatSNESMFSetFunctionError() and 688 umin via MatSNESMFDefaultSetUmin(); see the nonlinear solvers chapter 689 of the users manual for details. 690 691 The user should call MatDestroy() when finished with the matrix-free 692 matrix context. 693 694 Options Database Keys: 695 + -snes_mf_err <error_rel> - Sets error_rel 696 . -snes_mf_unim <umin> - Sets umin (for default PETSc routine that computes h only) 697 . -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h 698 - -snes_mf_check_positivity 699 700 .keywords: default, matrix-free, create, matrix 701 702 .seealso: MatDestroy(), MatSNESMFSetFunctionError(), MatSNESMFDefaultSetUmin(), MatSNESMFSetFunction() 703 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), MatCreateSNESMF(), 704 MatSNESMFGetH(),MatSNESMFKSPMonitor(), MatSNESMFRegisterDynamic),, MatSNESMFComputeJacobian() 705 706 @*/ 707 PetscErrorCode PETSCSNES_DLLEXPORT MatCreateMF(Vec x,Mat *J) 708 { 709 MPI_Comm comm; 710 PetscErrorCode ierr; 711 PetscInt n,nloc; 712 713 PetscFunctionBegin; 714 ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr); 715 ierr = VecGetSize(x,&n);CHKERRQ(ierr); 716 ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr); 717 ierr = MatCreate(comm,J);CHKERRQ(ierr); 718 ierr = MatSetSizes(*J,nloc,nloc,n,n);CHKERRQ(ierr); 719 ierr = MatRegisterDynamic(MATMFFD,0,"MatCreate_MFFD",MatCreate_MFFD);CHKERRQ(ierr); 720 ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 721 PetscFunctionReturn(0); 722 } 723 724 725 #undef __FUNCT__ 726 #define __FUNCT__ "MatSNESMFGetH" 727 /*@ 728 MatSNESMFGetH - Gets the last value that was used as the differencing 729 parameter. 730 731 Not Collective 732 733 Input Parameters: 734 . mat - the matrix obtained with MatCreateSNESMF() 735 736 Output Paramter: 737 . h - the differencing step size 738 739 Level: advanced 740 741 .keywords: SNES, matrix-free, parameters 742 743 .seealso: MatCreateSNESMF(),MatSNESMFSetHHistory(), 744 MatSNESMFResetHHistory(),MatSNESMFKSPMonitor() 745 @*/ 746 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFGetH(Mat mat,PetscScalar *h) 747 { 748 MatSNESMF ctx = (MatSNESMF)mat->data; 749 750 PetscFunctionBegin; 751 *h = ctx->currenth; 752 PetscFunctionReturn(0); 753 } 754 755 #undef __FUNCT__ 756 #define __FUNCT__ "MatSNESMFKSPMonitor" 757 /*@C 758 MatSNESMFKSPMonitor - A KSP monitor for use with the PETSc 759 SNES matrix free routines. Prints the differencing parameter used at 760 each step. 761 762 Collective on KSP 763 764 Input Parameters: 765 + ksp - the Krylov solver object 766 . n - the current iteration 767 . rnorm - the current residual norm (may be preconditioned or not depending on solver and options 768 _ dummy - unused argument 769 770 Options Database: 771 . -snes_mf_ksp_monitor - turn this monitor on 772 773 Notes: Use KSPMonitorSet(ksp,MatSNESMFKSPMonitor,PETSC_NULL,PETSC_NULL); 774 775 .seealso: KSPMonitorSet() 776 777 @*/ 778 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFKSPMonitor(KSP ksp,PetscInt n,PetscReal rnorm,void *dummy) 779 { 780 MatSNESMF ctx; 781 PetscErrorCode ierr; 782 Mat mat; 783 MPI_Comm comm; 784 PetscTruth nonzeroinitialguess; 785 786 PetscFunctionBegin; 787 ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr); 788 ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr); 789 ierr = KSPGetOperators(ksp,&mat,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 790 ctx = (MatSNESMF)mat->data; 791 792 if (n > 0 || nonzeroinitialguess) { 793 #if defined(PETSC_USE_COMPLEX) 794 ierr = PetscPrintf(comm,"%D KSP Residual norm %14.12e h %G + %G i\n",n,rnorm, 795 PetscRealPart(ctx->currenth),PetscImaginaryPart(ctx->currenth));CHKERRQ(ierr); 796 #else 797 ierr = PetscPrintf(comm,"%D KSP Residual norm %14.12e h %G \n",n,rnorm,ctx->currenth);CHKERRQ(ierr); 798 #endif 799 } else { 800 ierr = PetscPrintf(comm,"%D KSP Residual norm %14.12e\n",n,rnorm);CHKERRQ(ierr); 801 } 802 PetscFunctionReturn(0); 803 } 804 805 #undef __FUNCT__ 806 #define __FUNCT__ "MatSNESMFSetFunction" 807 /*@C 808 MatSNESMFSetFunction - Sets the function used in applying the matrix free. 809 810 Collective on Mat 811 812 Input Parameters: 813 + mat - the matrix free matrix created via MatCreateSNESMF() 814 . v - workspace vector 815 . func - the function to use 816 - funcctx - optional function context passed to function 817 818 Level: advanced 819 820 Notes: 821 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 822 matrix inside your compute Jacobian routine 823 824 If this is not set then it will use the function set with SNESSetFunction() 825 826 .keywords: SNES, matrix-free, function 827 828 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 829 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 830 MatSNESMFKSPMonitor(), SNESetFunction() 831 @*/ 832 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetFunction(Mat mat,Vec v,PetscErrorCode (*func)(SNES,Vec,Vec,void *),void *funcctx) 833 { 834 MatSNESMF ctx = (MatSNESMF)mat->data; 835 836 PetscFunctionBegin; 837 ctx->func = func; 838 ctx->funcctx = funcctx; 839 ctx->funcvec = v; 840 PetscFunctionReturn(0); 841 } 842 843 #undef __FUNCT__ 844 #define __FUNCT__ "MatSNESMFSetFunctioni" 845 /*@C 846 MatSNESMFSetFunctioni - Sets the function for a single component 847 848 Collective on Mat 849 850 Input Parameters: 851 + mat - the matrix free matrix created via MatCreateSNESMF() 852 - funci - the function to use 853 854 Level: advanced 855 856 Notes: 857 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 858 matrix inside your compute Jacobian routine 859 860 861 .keywords: SNES, matrix-free, function 862 863 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 864 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 865 MatSNESMFKSPMonitor(), SNESetFunction() 866 @*/ 867 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetFunctioni(Mat mat,PetscErrorCode (*funci)(PetscInt,Vec,PetscScalar*,void *)) 868 { 869 PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(PetscInt,Vec,PetscScalar*,void *)); 870 871 PetscFunctionBegin; 872 PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 873 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioni_C",(void (**)(void))&f);CHKERRQ(ierr); 874 if (f) { 875 ierr = (*f)(mat,funci);CHKERRQ(ierr); 876 } 877 PetscFunctionReturn(0); 878 } 879 880 881 #undef __FUNCT__ 882 #define __FUNCT__ "MatSNESMFSetFunctioniBase" 883 /*@C 884 MatSNESMFSetFunctioniBase - Sets the base vector for a single component function evaluation 885 886 Collective on Mat 887 888 Input Parameters: 889 + mat - the matrix free matrix created via MatCreateSNESMF() 890 - func - the function to use 891 892 Level: advanced 893 894 Notes: 895 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 896 matrix inside your compute Jacobian routine 897 898 899 .keywords: SNES, matrix-free, function 900 901 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 902 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 903 MatSNESMFKSPMonitor(), SNESetFunction() 904 @*/ 905 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetFunctioniBase(Mat mat,PetscErrorCode (*func)(Vec,void *)) 906 { 907 PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(Vec,void *)); 908 909 PetscFunctionBegin; 910 PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 911 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatSNESMFSetFunctioniBase_C",(void (**)(void))&f);CHKERRQ(ierr); 912 if (f) { 913 ierr = (*f)(mat,func);CHKERRQ(ierr); 914 } 915 PetscFunctionReturn(0); 916 } 917 918 919 #undef __FUNCT__ 920 #define __FUNCT__ "MatSNESMFSetPeriod" 921 /*@ 922 MatSNESMFSetPeriod - Sets how often h is recomputed, by default it is everytime 923 924 Collective on Mat 925 926 Input Parameters: 927 + mat - the matrix free matrix created via MatCreateSNESMF() 928 - period - 1 for everytime, 2 for every second etc 929 930 Options Database Keys: 931 + -snes_mf_period <period> 932 933 Level: advanced 934 935 936 .keywords: SNES, matrix-free, parameters 937 938 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 939 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 940 MatSNESMFKSPMonitor() 941 @*/ 942 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetPeriod(Mat mat,PetscInt period) 943 { 944 MatSNESMF ctx = (MatSNESMF)mat->data; 945 946 PetscFunctionBegin; 947 ctx->recomputeperiod = period; 948 PetscFunctionReturn(0); 949 } 950 951 #undef __FUNCT__ 952 #define __FUNCT__ "MatSNESMFSetFunctionError" 953 /*@ 954 MatSNESMFSetFunctionError - Sets the error_rel for the approximation of 955 matrix-vector products using finite differences. 956 957 Collective on Mat 958 959 Input Parameters: 960 + mat - the matrix free matrix created via MatCreateSNESMF() 961 - error_rel - relative error (should be set to the square root of 962 the relative error in the function evaluations) 963 964 Options Database Keys: 965 + -snes_mf_err <error_rel> - Sets error_rel 966 967 Level: advanced 968 969 Notes: 970 The default matrix-free matrix-vector product routine computes 971 .vb 972 F'(u)*a = [F(u+h*a) - F(u)]/h where 973 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 974 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 975 .ve 976 977 .keywords: SNES, matrix-free, parameters 978 979 .seealso: MatCreateSNESMF(),MatSNESMFGetH(), 980 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 981 MatSNESMFKSPMonitor() 982 @*/ 983 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetFunctionError(Mat mat,PetscReal error) 984 { 985 MatSNESMF ctx = (MatSNESMF)mat->data; 986 987 PetscFunctionBegin; 988 if (error != PETSC_DEFAULT) ctx->error_rel = error; 989 PetscFunctionReturn(0); 990 } 991 992 #undef __FUNCT__ 993 #define __FUNCT__ "MatSNESMFAddNullSpace" 994 /*@ 995 MatSNESMFAddNullSpace - Provides a null space that an operator is 996 supposed to have. Since roundoff will create a small component in 997 the null space, if you know the null space you may have it 998 automatically removed. 999 1000 Collective on Mat 1001 1002 Input Parameters: 1003 + J - the matrix-free matrix context 1004 - nullsp - object created with MatNullSpaceCreate() 1005 1006 Level: advanced 1007 1008 .keywords: SNES, matrix-free, null space 1009 1010 .seealso: MatNullSpaceCreate(), MatSNESMFGetH(), MatCreateSNESMF(), 1011 MatSNESMFSetHHistory(), MatSNESMFResetHHistory(), 1012 MatSNESMFKSPMonitor(), MatSNESMFErrorRel() 1013 @*/ 1014 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFAddNullSpace(Mat J,MatNullSpace nullsp) 1015 { 1016 PetscErrorCode ierr; 1017 MatSNESMF ctx = (MatSNESMF)J->data; 1018 1019 PetscFunctionBegin; 1020 ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr); 1021 if (ctx->sp) { ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr); } 1022 ctx->sp = nullsp; 1023 PetscFunctionReturn(0); 1024 } 1025 1026 #undef __FUNCT__ 1027 #define __FUNCT__ "MatSNESMFSetHHistory" 1028 /*@ 1029 MatSNESMFSetHHistory - Sets an array to collect a history of the 1030 differencing values (h) computed for the matrix-free product. 1031 1032 Collective on Mat 1033 1034 Input Parameters: 1035 + J - the matrix-free matrix context 1036 . histroy - space to hold the history 1037 - nhistory - number of entries in history, if more entries are generated than 1038 nhistory, then the later ones are discarded 1039 1040 Level: advanced 1041 1042 Notes: 1043 Use MatSNESMFResetHHistory() to reset the history counter and collect 1044 a new batch of differencing parameters, h. 1045 1046 .keywords: SNES, matrix-free, h history, differencing history 1047 1048 .seealso: MatSNESMFGetH(), MatCreateSNESMF(), 1049 MatSNESMFResetHHistory(), 1050 MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError() 1051 1052 @*/ 1053 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 1054 { 1055 MatSNESMF ctx = (MatSNESMF)J->data; 1056 1057 PetscFunctionBegin; 1058 ctx->historyh = history; 1059 ctx->maxcurrenth = nhistory; 1060 ctx->currenth = 0; 1061 PetscFunctionReturn(0); 1062 } 1063 1064 #undef __FUNCT__ 1065 #define __FUNCT__ "MatSNESMFResetHHistory" 1066 /*@ 1067 MatSNESMFResetHHistory - Resets the counter to zero to begin 1068 collecting a new set of differencing histories. 1069 1070 Collective on Mat 1071 1072 Input Parameters: 1073 . J - the matrix-free matrix context 1074 1075 Level: advanced 1076 1077 Notes: 1078 Use MatSNESMFSetHHistory() to create the original history counter. 1079 1080 .keywords: SNES, matrix-free, h history, differencing history 1081 1082 .seealso: MatSNESMFGetH(), MatCreateSNESMF(), 1083 MatSNESMFSetHHistory(), 1084 MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError() 1085 1086 @*/ 1087 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFResetHHistory(Mat J) 1088 { 1089 MatSNESMF ctx = (MatSNESMF)J->data; 1090 1091 PetscFunctionBegin; 1092 ctx->ncurrenth = 0; 1093 PetscFunctionReturn(0); 1094 } 1095 1096 #undef __FUNCT__ 1097 #define __FUNCT__ "MatSNESMFComputeJacobian" 1098 /*@ 1099 MatSNESMFComputeJacobian - Tells the matrix-free Jacobian object the new location at which 1100 Jacobian matrix vector products will be computed at, i.e. J(x) * a. 1101 1102 Collective on SNES 1103 1104 Input Parameters: 1105 + snes - the nonlinear solver context 1106 . x - the point at which the Jacobian vector products will be performed 1107 . jac - the matrix-free Jacobian object 1108 . B - either the same as jac or another matrix type (ignored) 1109 . flag - not relevent for matrix-free form 1110 - dummy - the user context (ignored) 1111 1112 Level: developer 1113 1114 Notes: 1115 This can be passed into SNESSetJacobian() when using a completely matrix-free solver, 1116 that is the B matrix is also the same matrix operator. This is used when you select 1117 -snes_mf but rarely used directly by users. 1118 1119 .seealso: MatSNESMFGetH(), MatCreateSNESMF(), 1120 MatSNESMFSetHHistory(), 1121 MatSNESMFKSPMonitor(), MatSNESMFSetFunctionError(), MatSNESMFCreate(), SNESSetJacobian() 1122 1123 @*/ 1124 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFComputeJacobian(SNES snes,Vec x,Mat *jac,Mat *B,MatStructure *flag,void *dummy) 1125 { 1126 PetscErrorCode ierr; 1127 PetscFunctionBegin; 1128 ierr = MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1129 ierr = MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1130 PetscFunctionReturn(0); 1131 } 1132 1133 #undef __FUNCT__ 1134 #define __FUNCT__ "MatSNESMFSetBase" 1135 /*@ 1136 MatSNESMFSetBase - Sets the vector U at which matrix vector products of the 1137 Jacobian are computed 1138 1139 Collective on Mat 1140 1141 Input Parameters: 1142 + J - the MatSNESMF matrix 1143 - U - the vector 1144 1145 Notes: This is rarely used directly 1146 1147 Level: advanced 1148 1149 @*/ 1150 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetBase(Mat J,Vec U) 1151 { 1152 PetscErrorCode ierr,(*f)(Mat,Vec); 1153 1154 PetscFunctionBegin; 1155 PetscValidHeaderSpecific(J,MAT_COOKIE,1); 1156 PetscValidHeaderSpecific(U,VEC_COOKIE,2); 1157 ierr = PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetBase_C",(void (**)(void))&f);CHKERRQ(ierr); 1158 if (f) { 1159 ierr = (*f)(J,U);CHKERRQ(ierr); 1160 } 1161 PetscFunctionReturn(0); 1162 } 1163 1164 #undef __FUNCT__ 1165 #define __FUNCT__ "MatSNESMFSetCheckh" 1166 /*@C 1167 MatSNESMFSetCheckh - Sets a function that checks the computed h and adjusts 1168 it to satisfy some criteria 1169 1170 Collective on Mat 1171 1172 Input Parameters: 1173 + J - the MatSNESMF matrix 1174 . fun - the function that checks h 1175 - ctx - any context needed by the function 1176 1177 Options Database Keys: 1178 . -snes_mf_check_positivity 1179 1180 Level: advanced 1181 1182 Notes: For example, MatSNESMFSetCheckPositivity() insures that all entries 1183 of U + h*a are non-negative 1184 1185 .seealso: MatSNESMFSetCheckPositivity() 1186 @*/ 1187 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFSetCheckh(Mat J,PetscErrorCode (*fun)(Vec,Vec,PetscScalar*,void*),void* ctx) 1188 { 1189 PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(Vec,Vec,PetscScalar*,void*),void*); 1190 1191 PetscFunctionBegin; 1192 PetscValidHeaderSpecific(J,MAT_COOKIE,1); 1193 ierr = PetscObjectQueryFunction((PetscObject)J,"MatSNESMFSetCheckh_C",(void (**)(void))&f);CHKERRQ(ierr); 1194 if (f) { 1195 ierr = (*f)(J,fun,ctx);CHKERRQ(ierr); 1196 } 1197 PetscFunctionReturn(0); 1198 } 1199 1200 #undef __FUNCT__ 1201 #define __FUNCT__ "MatSNESMFSetCheckPositivity" 1202 /*@ 1203 MatSNESMFCheckPositivity - Checks that all entries in U + h*a are positive or 1204 zero, decreases h until this is satisfied. 1205 1206 Collective on Vec 1207 1208 Input Parameters: 1209 + U - base vector that is added to 1210 . a - vector that is added 1211 . h - scaling factor on a 1212 - dummy - context variable (unused) 1213 1214 Options Database Keys: 1215 . -snes_mf_check_positivity 1216 1217 Level: advanced 1218 1219 Notes: This is rarely used directly, rather it is passed as an argument to 1220 MatSNESMFSetCheckh() 1221 1222 .seealso: MatSNESMFSetCheckh() 1223 @*/ 1224 PetscErrorCode PETSCSNES_DLLEXPORT MatSNESMFCheckPositivity(Vec U,Vec a,PetscScalar *h,void *dummy) 1225 { 1226 PetscReal val, minval; 1227 PetscScalar *u_vec, *a_vec; 1228 PetscErrorCode ierr; 1229 PetscInt i,n; 1230 MPI_Comm comm; 1231 1232 PetscFunctionBegin; 1233 ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); 1234 ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); 1235 ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); 1236 ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); 1237 minval = PetscAbsScalar(*h*1.01); 1238 for(i=0;i<n;i++) { 1239 if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1240 val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1241 if (val < minval) minval = val; 1242 } 1243 } 1244 ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); 1245 ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); 1246 ierr = PetscGlobalMin(&minval,&val,comm);CHKERRQ(ierr); 1247 if (val <= PetscAbsScalar(*h)) { 1248 ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr); 1249 if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1250 else *h = -0.99*val; 1251 } 1252 PetscFunctionReturn(0); 1253 } 1254 1255 1256 1257 1258 1259 1260 1261