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