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