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