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