1 2 #include <petsc/private/matimpl.h> 3 #include <../src/mat/impls/mffd/mffdimpl.h> /*I "petscmat.h" I*/ 4 5 PetscFunctionList MatMFFDList = NULL; 6 PetscBool MatMFFDRegisterAllCalled = PETSC_FALSE; 7 8 PetscClassId MATMFFD_CLASSID; 9 PetscLogEvent MATMFFD_Mult; 10 11 static PetscBool MatMFFDPackageInitialized = PETSC_FALSE; 12 /*@C 13 MatMFFDFinalizePackage - This function destroys everything in the MatMFFD package. It is 14 called from PetscFinalize(). 15 16 Level: developer 17 18 .seealso: PetscFinalize(), MatCreateMFFD(), MatCreateSNESMF() 19 @*/ 20 PetscErrorCode MatMFFDFinalizePackage(void) 21 { 22 PetscFunctionBegin; 23 CHKERRQ(PetscFunctionListDestroy(&MatMFFDList)); 24 MatMFFDPackageInitialized = PETSC_FALSE; 25 MatMFFDRegisterAllCalled = PETSC_FALSE; 26 PetscFunctionReturn(0); 27 } 28 29 /*@C 30 MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called 31 from MatInitializePackage(). 32 33 Level: developer 34 35 .seealso: PetscInitialize() 36 @*/ 37 PetscErrorCode MatMFFDInitializePackage(void) 38 { 39 char logList[256]; 40 PetscBool opt,pkg; 41 42 PetscFunctionBegin; 43 if (MatMFFDPackageInitialized) PetscFunctionReturn(0); 44 MatMFFDPackageInitialized = PETSC_TRUE; 45 /* Register Classes */ 46 CHKERRQ(PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID)); 47 /* Register Constructors */ 48 CHKERRQ(MatMFFDRegisterAll()); 49 /* Register Events */ 50 CHKERRQ(PetscLogEventRegister("MatMult MF",MATMFFD_CLASSID,&MATMFFD_Mult)); 51 /* Process Info */ 52 { 53 PetscClassId classids[1]; 54 55 classids[0] = MATMFFD_CLASSID; 56 CHKERRQ(PetscInfoProcessClass("matmffd", 1, classids)); 57 } 58 /* Process summary exclusions */ 59 CHKERRQ(PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt)); 60 if (opt) { 61 CHKERRQ(PetscStrInList("matmffd",logList,',',&pkg)); 62 if (pkg) CHKERRQ(PetscLogEventExcludeClass(MATMFFD_CLASSID)); 63 } 64 /* Register package finalizer */ 65 CHKERRQ(PetscRegisterFinalize(MatMFFDFinalizePackage)); 66 PetscFunctionReturn(0); 67 } 68 69 static PetscErrorCode MatMFFDSetType_MFFD(Mat mat,MatMFFDType ftype) 70 { 71 MatMFFD ctx; 72 PetscBool match; 73 PetscErrorCode (*r)(MatMFFD); 74 75 PetscFunctionBegin; 76 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 77 PetscValidCharPointer(ftype,2); 78 CHKERRQ(MatShellGetContext(mat,&ctx)); 79 80 /* already set, so just return */ 81 CHKERRQ(PetscObjectTypeCompare((PetscObject)ctx,ftype,&match)); 82 if (match) PetscFunctionReturn(0); 83 84 /* destroy the old one if it exists */ 85 if (ctx->ops->destroy) CHKERRQ((*ctx->ops->destroy)(ctx)); 86 87 CHKERRQ(PetscFunctionListFind(MatMFFDList,ftype,&r)); 88 PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype); 89 CHKERRQ((*r)(ctx)); 90 CHKERRQ(PetscObjectChangeTypeName((PetscObject)ctx,ftype)); 91 PetscFunctionReturn(0); 92 } 93 94 /*@C 95 MatMFFDSetType - Sets the method that is used to compute the 96 differencing parameter for finite differene matrix-free formulations. 97 98 Input Parameters: 99 + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD() 100 or MatSetType(mat,MATMFFD); 101 - ftype - the type requested, either MATMFFD_WP or MATMFFD_DS 102 103 Level: advanced 104 105 Notes: 106 For example, such routines can compute h for use in 107 Jacobian-vector products of the form 108 109 F(x+ha) - F(x) 110 F'(u)a ~= ---------------- 111 h 112 113 .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction(), MatCreateMFFD() 114 @*/ 115 PetscErrorCode MatMFFDSetType(Mat mat,MatMFFDType ftype) 116 { 117 PetscFunctionBegin; 118 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 119 PetscValidCharPointer(ftype,2); 120 CHKERRQ(PetscTryMethod(mat,"MatMFFDSetType_C",(Mat,MatMFFDType),(mat,ftype))); 121 PetscFunctionReturn(0); 122 } 123 124 static PetscErrorCode MatGetDiagonal_MFFD(Mat,Vec); 125 126 typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/ 127 static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func) 128 { 129 MatMFFD ctx; 130 131 PetscFunctionBegin; 132 CHKERRQ(MatShellGetContext(mat,&ctx)); 133 ctx->funcisetbase = func; 134 PetscFunctionReturn(0); 135 } 136 137 typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 138 static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci) 139 { 140 MatMFFD ctx; 141 142 PetscFunctionBegin; 143 CHKERRQ(MatShellGetContext(mat,&ctx)); 144 ctx->funci = funci; 145 CHKERRQ(MatShellSetOperation(mat,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_MFFD)); 146 PetscFunctionReturn(0); 147 } 148 149 static PetscErrorCode MatMFFDGetH_MFFD(Mat mat,PetscScalar *h) 150 { 151 MatMFFD ctx; 152 153 PetscFunctionBegin; 154 CHKERRQ(MatShellGetContext(mat,&ctx)); 155 *h = ctx->currenth; 156 PetscFunctionReturn(0); 157 } 158 159 static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) 160 { 161 MatMFFD ctx; 162 163 PetscFunctionBegin; 164 CHKERRQ(MatShellGetContext(J,&ctx)); 165 ctx->ncurrenth = 0; 166 PetscFunctionReturn(0); 167 } 168 169 /*@C 170 MatMFFDRegister - Adds a method to the MatMFFD registry. 171 172 Not Collective 173 174 Input Parameters: 175 + name_solver - name of a new user-defined compute-h module 176 - routine_create - routine to create method context 177 178 Level: developer 179 180 Notes: 181 MatMFFDRegister() may be called multiple times to add several user-defined solvers. 182 183 Sample usage: 184 .vb 185 MatMFFDRegister("my_h",MyHCreate); 186 .ve 187 188 Then, your solver can be chosen with the procedural interface via 189 $ MatMFFDSetType(mfctx,"my_h") 190 or at runtime via the option 191 $ -mat_mffd_type my_h 192 193 .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 194 @*/ 195 PetscErrorCode MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD)) 196 { 197 PetscFunctionBegin; 198 CHKERRQ(MatInitializePackage()); 199 CHKERRQ(PetscFunctionListAdd(&MatMFFDList,sname,function)); 200 PetscFunctionReturn(0); 201 } 202 203 /* ----------------------------------------------------------------------------------------*/ 204 static PetscErrorCode MatDestroy_MFFD(Mat mat) 205 { 206 MatMFFD ctx; 207 208 PetscFunctionBegin; 209 CHKERRQ(MatShellGetContext(mat,&ctx)); 210 CHKERRQ(VecDestroy(&ctx->w)); 211 CHKERRQ(VecDestroy(&ctx->current_u)); 212 if (ctx->current_f_allocated) { 213 CHKERRQ(VecDestroy(&ctx->current_f)); 214 } 215 if (ctx->ops->destroy) CHKERRQ((*ctx->ops->destroy)(ctx)); 216 CHKERRQ(PetscHeaderDestroy(&ctx)); 217 218 CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL)); 219 CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL)); 220 CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL)); 221 CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL)); 222 CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL)); 223 CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL)); 224 CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL)); 225 CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL)); 226 CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetHHistory_C",NULL)); 227 CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetType_C",NULL)); 228 CHKERRQ(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDGetH_C",NULL)); 229 PetscFunctionReturn(0); 230 } 231 232 /* 233 MatMFFDView_MFFD - Views matrix-free parameters. 234 235 */ 236 static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer) 237 { 238 MatMFFD ctx; 239 PetscBool iascii, viewbase, viewfunction; 240 const char *prefix; 241 242 PetscFunctionBegin; 243 CHKERRQ(MatShellGetContext(J,&ctx)); 244 CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 245 if (iascii) { 246 CHKERRQ(PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n")); 247 CHKERRQ(PetscViewerASCIIPushTab(viewer)); 248 CHKERRQ(PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel)); 249 if (!((PetscObject)ctx)->type_name) { 250 CHKERRQ(PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n")); 251 } else { 252 CHKERRQ(PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name)); 253 } 254 #if defined(PETSC_USE_COMPLEX) 255 if (ctx->usecomplex) { 256 CHKERRQ(PetscViewerASCIIPrintf(viewer,"Using Lyness complex number trick to compute the matrix-vector product\n")); 257 } 258 #endif 259 if (ctx->ops->view) { 260 CHKERRQ((*ctx->ops->view)(ctx,viewer)); 261 } 262 CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject)J, &prefix)); 263 264 CHKERRQ(PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase)); 265 if (viewbase) { 266 CHKERRQ(PetscViewerASCIIPrintf(viewer, "Base:\n")); 267 CHKERRQ(VecView(ctx->current_u, viewer)); 268 } 269 CHKERRQ(PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction)); 270 if (viewfunction) { 271 CHKERRQ(PetscViewerASCIIPrintf(viewer, "Function:\n")); 272 CHKERRQ(VecView(ctx->current_f, viewer)); 273 } 274 CHKERRQ(PetscViewerASCIIPopTab(viewer)); 275 } 276 PetscFunctionReturn(0); 277 } 278 279 /* 280 MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 281 allows the user to indicate the beginning of a new linear solve by calling 282 MatAssemblyXXX() on the matrix free matrix. This then allows the 283 MatCreateMFFD_WP() to properly compute ||U|| only the first time 284 in the linear solver rather than every time. 285 286 This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence 287 it must be labeled as PETSC_EXTERN 288 */ 289 PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt) 290 { 291 MatMFFD j; 292 293 PetscFunctionBegin; 294 CHKERRQ(MatShellGetContext(J,&j)); 295 CHKERRQ(MatMFFDResetHHistory(J)); 296 PetscFunctionReturn(0); 297 } 298 299 /* 300 MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 301 302 y ~= (F(u + ha) - F(u))/h, 303 where F = nonlinear function, as set by SNESSetFunction() 304 u = current iterate 305 h = difference interval 306 */ 307 static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y) 308 { 309 MatMFFD ctx; 310 PetscScalar h; 311 Vec w,U,F; 312 PetscBool zeroa; 313 314 PetscFunctionBegin; 315 CHKERRQ(MatShellGetContext(mat,&ctx)); 316 PetscCheck(ctx->current_u,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MatMFFDSetBase() has not been called, this is often caused by forgetting to call \n\t\tMatAssemblyBegin/End on the first Mat in the SNES compute function"); 317 /* We log matrix-free matrix-vector products separately, so that we can 318 separate the performance monitoring from the cases that use conventional 319 storage. We may eventually modify event logging to associate events 320 with particular objects, hence alleviating the more general problem. */ 321 CHKERRQ(PetscLogEventBegin(MATMFFD_Mult,a,y,0,0)); 322 323 w = ctx->w; 324 U = ctx->current_u; 325 F = ctx->current_f; 326 /* 327 Compute differencing parameter 328 */ 329 if (!((PetscObject)ctx)->type_name) { 330 CHKERRQ(MatMFFDSetType(mat,MATMFFD_WP)); 331 CHKERRQ(MatSetFromOptions(mat)); 332 } 333 CHKERRQ((*ctx->ops->compute)(ctx,U,a,&h,&zeroa)); 334 if (zeroa) { 335 CHKERRQ(VecSet(y,0.0)); 336 PetscFunctionReturn(0); 337 } 338 339 PetscCheckFalse(mat->erroriffailure && PetscIsInfOrNanScalar(h),PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h"); 340 if (ctx->checkh) { 341 CHKERRQ((*ctx->checkh)(ctx->checkhctx,U,a,&h)); 342 } 343 344 /* keep a record of the current differencing parameter h */ 345 ctx->currenth = h; 346 #if defined(PETSC_USE_COMPLEX) 347 CHKERRQ(PetscInfo(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h))); 348 #else 349 CHKERRQ(PetscInfo(mat,"Current differencing parameter: %15.12e\n",(double)PetscRealPart(h))); 350 #endif 351 if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 352 ctx->historyh[ctx->ncurrenth] = h; 353 } 354 ctx->ncurrenth++; 355 356 #if defined(PETSC_USE_COMPLEX) 357 if (ctx->usecomplex) h = PETSC_i*h; 358 #endif 359 360 /* w = u + ha */ 361 CHKERRQ(VecWAXPY(w,h,a,U)); 362 363 /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 364 if (ctx->ncurrenth == 1 && ctx->current_f_allocated) { 365 CHKERRQ((*ctx->func)(ctx->funcctx,U,F)); 366 } 367 CHKERRQ((*ctx->func)(ctx->funcctx,w,y)); 368 369 #if defined(PETSC_USE_COMPLEX) 370 if (ctx->usecomplex) { 371 CHKERRQ(VecImaginaryPart(y)); 372 h = PetscImaginaryPart(h); 373 } else { 374 CHKERRQ(VecAXPY(y,-1.0,F)); 375 } 376 #else 377 CHKERRQ(VecAXPY(y,-1.0,F)); 378 #endif 379 CHKERRQ(VecScale(y,1.0/h)); 380 if (mat->nullsp) CHKERRQ(MatNullSpaceRemove(mat->nullsp,y)); 381 382 CHKERRQ(PetscLogEventEnd(MATMFFD_Mult,a,y,0,0)); 383 PetscFunctionReturn(0); 384 } 385 386 /* 387 MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 388 389 y ~= (F(u + ha) - F(u))/h, 390 where F = nonlinear function, as set by SNESSetFunction() 391 u = current iterate 392 h = difference interval 393 */ 394 PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 395 { 396 MatMFFD ctx; 397 PetscScalar h,*aa,*ww,v; 398 PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 399 Vec w,U; 400 PetscInt i,rstart,rend; 401 402 PetscFunctionBegin; 403 CHKERRQ(MatShellGetContext(mat,&ctx)); 404 PetscCheck(ctx->func,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunction() first"); 405 PetscCheck(ctx->funci,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioni() first"); 406 w = ctx->w; 407 U = ctx->current_u; 408 CHKERRQ((*ctx->func)(ctx->funcctx,U,a)); 409 if (ctx->funcisetbase) { 410 CHKERRQ((*ctx->funcisetbase)(ctx->funcctx,U)); 411 } 412 CHKERRQ(VecCopy(U,w)); 413 414 CHKERRQ(VecGetOwnershipRange(a,&rstart,&rend)); 415 CHKERRQ(VecGetArray(a,&aa)); 416 for (i=rstart; i<rend; i++) { 417 CHKERRQ(VecGetArray(w,&ww)); 418 h = ww[i-rstart]; 419 if (h == 0.0) h = 1.0; 420 if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 421 else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 422 h *= epsilon; 423 424 ww[i-rstart] += h; 425 CHKERRQ(VecRestoreArray(w,&ww)); 426 CHKERRQ((*ctx->funci)(ctx->funcctx,i,w,&v)); 427 aa[i-rstart] = (v - aa[i-rstart])/h; 428 429 CHKERRQ(VecGetArray(w,&ww)); 430 ww[i-rstart] -= h; 431 CHKERRQ(VecRestoreArray(w,&ww)); 432 } 433 CHKERRQ(VecRestoreArray(a,&aa)); 434 PetscFunctionReturn(0); 435 } 436 437 PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F) 438 { 439 MatMFFD ctx; 440 441 PetscFunctionBegin; 442 CHKERRQ(MatShellGetContext(J,&ctx)); 443 CHKERRQ(MatMFFDResetHHistory(J)); 444 if (!ctx->current_u) { 445 CHKERRQ(VecDuplicate(U,&ctx->current_u)); 446 CHKERRQ(VecLockReadPush(ctx->current_u)); 447 } 448 CHKERRQ(VecLockReadPop(ctx->current_u)); 449 CHKERRQ(VecCopy(U,ctx->current_u)); 450 CHKERRQ(VecLockReadPush(ctx->current_u)); 451 if (F) { 452 if (ctx->current_f_allocated) CHKERRQ(VecDestroy(&ctx->current_f)); 453 ctx->current_f = F; 454 ctx->current_f_allocated = PETSC_FALSE; 455 } else if (!ctx->current_f_allocated) { 456 CHKERRQ(MatCreateVecs(J,NULL,&ctx->current_f)); 457 ctx->current_f_allocated = PETSC_TRUE; 458 } 459 if (!ctx->w) { 460 CHKERRQ(VecDuplicate(ctx->current_u,&ctx->w)); 461 } 462 J->assembled = PETSC_TRUE; 463 PetscFunctionReturn(0); 464 } 465 466 typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 467 468 static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx) 469 { 470 MatMFFD ctx; 471 472 PetscFunctionBegin; 473 CHKERRQ(MatShellGetContext(J,&ctx)); 474 ctx->checkh = fun; 475 ctx->checkhctx = ectx; 476 PetscFunctionReturn(0); 477 } 478 479 /*@C 480 MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 481 MatMFFD options in the database. 482 483 Collective on Mat 484 485 Input Parameters: 486 + A - the Mat context 487 - prefix - the prefix to prepend to all option names 488 489 Notes: 490 A hyphen (-) must NOT be given at the beginning of the prefix name. 491 The first character of all runtime options is AUTOMATICALLY the hyphen. 492 493 Level: advanced 494 495 .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD() 496 @*/ 497 PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[]) 498 { 499 MatMFFD mfctx; 500 501 PetscFunctionBegin; 502 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 503 CHKERRQ(MatShellGetContext(mat,&mfctx)); 504 PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 505 CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix)); 506 PetscFunctionReturn(0); 507 } 508 509 static PetscErrorCode MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat) 510 { 511 MatMFFD mfctx; 512 PetscErrorCode ierr; 513 PetscBool flg; 514 char ftype[256]; 515 516 PetscFunctionBegin; 517 PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 518 CHKERRQ(MatShellGetContext(mat,&mfctx)); 519 PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,2); 520 ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr); 521 CHKERRQ(PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg)); 522 if (flg) { 523 CHKERRQ(MatMFFDSetType(mat,ftype)); 524 } 525 526 CHKERRQ(PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,NULL)); 527 CHKERRQ(PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,NULL)); 528 529 flg = PETSC_FALSE; 530 CHKERRQ(PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL)); 531 if (flg) { 532 CHKERRQ(MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,NULL)); 533 } 534 #if defined(PETSC_USE_COMPLEX) 535 CHKERRQ(PetscOptionsBool("-mat_mffd_complex","Use Lyness complex number trick to compute the matrix-vector product","None",mfctx->usecomplex,&mfctx->usecomplex,NULL)); 536 #endif 537 if (mfctx->ops->setfromoptions) { 538 CHKERRQ((*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx)); 539 } 540 ierr = PetscOptionsEnd();CHKERRQ(ierr); 541 PetscFunctionReturn(0); 542 } 543 544 static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period) 545 { 546 MatMFFD ctx; 547 548 PetscFunctionBegin; 549 CHKERRQ(MatShellGetContext(mat,&ctx)); 550 ctx->recomputeperiod = period; 551 PetscFunctionReturn(0); 552 } 553 554 static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 555 { 556 MatMFFD ctx; 557 558 PetscFunctionBegin; 559 CHKERRQ(MatShellGetContext(mat,&ctx)); 560 ctx->func = func; 561 ctx->funcctx = funcctx; 562 PetscFunctionReturn(0); 563 } 564 565 static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error) 566 { 567 MatMFFD ctx; 568 569 PetscFunctionBegin; 570 CHKERRQ(MatShellGetContext(mat,&ctx)); 571 if (error != PETSC_DEFAULT) ctx->error_rel = error; 572 PetscFunctionReturn(0); 573 } 574 575 PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J,PetscScalar history[],PetscInt nhistory) 576 { 577 MatMFFD ctx; 578 579 PetscFunctionBegin; 580 CHKERRQ(MatShellGetContext(J,&ctx)); 581 ctx->historyh = history; 582 ctx->maxcurrenth = nhistory; 583 ctx->currenth = 0.; 584 PetscFunctionReturn(0); 585 } 586 587 /*MC 588 MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 589 590 Level: advanced 591 592 Developers Note: This is implemented on top of MATSHELL to get support for scaling and shifting without requiring duplicate code 593 594 .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(), 595 MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 596 MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 597 MatMFFDGetH(), 598 M*/ 599 PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 600 { 601 MatMFFD mfctx; 602 603 PetscFunctionBegin; 604 CHKERRQ(MatMFFDInitializePackage()); 605 606 CHKERRQ(PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),NULL,NULL)); 607 608 mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 609 mfctx->recomputeperiod = 1; 610 mfctx->count = 0; 611 mfctx->currenth = 0.0; 612 mfctx->historyh = NULL; 613 mfctx->ncurrenth = 0; 614 mfctx->maxcurrenth = 0; 615 ((PetscObject)mfctx)->type_name = NULL; 616 617 /* 618 Create the empty data structure to contain compute-h routines. 619 These will be filled in below from the command line options or 620 a later call with MatMFFDSetType() or if that is not called 621 then it will default in the first use of MatMult_MFFD() 622 */ 623 mfctx->ops->compute = NULL; 624 mfctx->ops->destroy = NULL; 625 mfctx->ops->view = NULL; 626 mfctx->ops->setfromoptions = NULL; 627 mfctx->hctx = NULL; 628 629 mfctx->func = NULL; 630 mfctx->funcctx = NULL; 631 mfctx->w = NULL; 632 mfctx->mat = A; 633 634 CHKERRQ(MatSetType(A,MATSHELL)); 635 CHKERRQ(MatShellSetContext(A,mfctx)); 636 CHKERRQ(MatShellSetOperation(A,MATOP_MULT,(void (*)(void))MatMult_MFFD)); 637 CHKERRQ(MatShellSetOperation(A,MATOP_DESTROY,(void (*)(void))MatDestroy_MFFD)); 638 CHKERRQ(MatShellSetOperation(A,MATOP_VIEW,(void (*)(void))MatView_MFFD)); 639 CHKERRQ(MatShellSetOperation(A,MATOP_ASSEMBLY_END,(void (*)(void))MatAssemblyEnd_MFFD)); 640 CHKERRQ(MatShellSetOperation(A,MATOP_SET_FROM_OPTIONS,(void (*)(void))MatSetFromOptions_MFFD)); 641 642 CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD)); 643 CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD)); 644 CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD)); 645 CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD)); 646 CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD)); 647 CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD)); 648 CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD)); 649 CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD)); 650 CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetHHistory_C",MatMFFDSetHHistory_MFFD)); 651 CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetType_C",MatMFFDSetType_MFFD)); 652 CHKERRQ(PetscObjectComposeFunction((PetscObject)A,"MatMFFDGetH_C",MatMFFDGetH_MFFD)); 653 CHKERRQ(PetscObjectChangeTypeName((PetscObject)A,MATMFFD)); 654 PetscFunctionReturn(0); 655 } 656 657 /*@ 658 MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 659 660 Collective on Vec 661 662 Input Parameters: 663 + comm - MPI communicator 664 . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 665 This value should be the same as the local size used in creating the 666 y vector for the matrix-vector product y = Ax. 667 . n - This value should be the same as the local size used in creating the 668 x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 669 calculated if N is given) For square matrices n is almost always m. 670 . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 671 - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 672 673 Output Parameter: 674 . J - the matrix-free matrix 675 676 Options Database Keys: call MatSetFromOptions() to trigger these 677 + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 678 . -mat_mffd_err - square root of estimated relative error in function evaluation 679 . -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 680 . -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values 681 - -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing 682 (requires real valued functions but that PETSc be configured for complex numbers) 683 684 Level: advanced 685 686 Notes: 687 The matrix-free matrix context merely contains the function pointers 688 and work space for performing finite difference approximations of 689 Jacobian-vector products, F'(u)*a, 690 691 The default code uses the following approach to compute h 692 693 .vb 694 F'(u)*a = [F(u+h*a) - F(u)]/h where 695 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 696 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 697 where 698 error_rel = square root of relative error in function evaluation 699 umin = minimum iterate parameter 700 .ve 701 702 You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different 703 preconditioner matrix 704 705 The user can set the error_rel via MatMFFDSetFunctionError() and 706 umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details. 707 708 The user should call MatDestroy() when finished with the matrix-free 709 matrix context. 710 711 Options Database Keys: 712 + -mat_mffd_err <error_rel> - Sets error_rel 713 . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only) 714 - -mat_mffd_check_positivity - check positivity 715 716 .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 717 MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 718 MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian() 719 720 @*/ 721 PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 722 { 723 PetscFunctionBegin; 724 CHKERRQ(MatCreate(comm,J)); 725 CHKERRQ(MatSetSizes(*J,m,n,M,N)); 726 CHKERRQ(MatSetType(*J,MATMFFD)); 727 CHKERRQ(MatSetUp(*J)); 728 PetscFunctionReturn(0); 729 } 730 731 /*@ 732 MatMFFDGetH - Gets the last value that was used as the differencing 733 parameter. 734 735 Not Collective 736 737 Input Parameters: 738 . mat - the matrix obtained with MatCreateSNESMF() 739 740 Output Parameter: 741 . h - the differencing step size 742 743 Level: advanced 744 745 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory() 746 @*/ 747 PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 748 { 749 PetscFunctionBegin; 750 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 751 PetscValidScalarPointer(h,2); 752 CHKERRQ(PetscUseMethod(mat,"MatMFFDGetH_C",(Mat,PetscScalar*),(mat,h))); 753 PetscFunctionReturn(0); 754 } 755 756 /*@C 757 MatMFFDSetFunction - Sets the function used in applying the matrix free. 758 759 Logically Collective on Mat 760 761 Input Parameters: 762 + mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD() 763 . func - the function to use 764 - funcctx - optional function context passed to function 765 766 Calling Sequence of func: 767 $ func (void *funcctx, Vec x, Vec f) 768 769 + funcctx - user provided context 770 . x - input vector 771 - f - computed output function 772 773 Level: advanced 774 775 Notes: 776 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 777 matrix inside your compute Jacobian routine 778 779 If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 780 781 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD, 782 MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 783 @*/ 784 PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 785 { 786 PetscFunctionBegin; 787 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 788 CHKERRQ(PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx))); 789 PetscFunctionReturn(0); 790 } 791 792 /*@C 793 MatMFFDSetFunctioni - Sets the function for a single component 794 795 Logically Collective on Mat 796 797 Input Parameters: 798 + mat - the matrix free matrix created via MatCreateSNESMF() 799 - funci - the function to use 800 801 Level: advanced 802 803 Notes: 804 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 805 matrix inside your compute Jacobian routine. 806 This function is necessary to compute the diagonal of the matrix. 807 funci must not contain any MPI call as it is called inside a loop on the local portion of the vector. 808 809 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 810 811 @*/ 812 PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 813 { 814 PetscFunctionBegin; 815 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 816 CHKERRQ(PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci))); 817 PetscFunctionReturn(0); 818 } 819 820 /*@C 821 MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 822 823 Logically 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 This function is necessary to compute the diagonal of the matrix. 835 836 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 837 MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 838 @*/ 839 PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 840 { 841 PetscFunctionBegin; 842 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 843 CHKERRQ(PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func))); 844 PetscFunctionReturn(0); 845 } 846 847 /*@ 848 MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 849 850 Logically Collective on Mat 851 852 Input Parameters: 853 + mat - the matrix free matrix created via MatCreateSNESMF() 854 - period - 1 for everytime, 2 for every second etc 855 856 Options Database Keys: 857 . -mat_mffd_period <period> - Sets how often h is recomputed 858 859 Level: advanced 860 861 .seealso: MatCreateSNESMF(),MatMFFDGetH(), 862 MatMFFDSetHHistory(), MatMFFDResetHHistory() 863 @*/ 864 PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 865 { 866 PetscFunctionBegin; 867 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 868 PetscValidLogicalCollectiveInt(mat,period,2); 869 CHKERRQ(PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period))); 870 PetscFunctionReturn(0); 871 } 872 873 /*@ 874 MatMFFDSetFunctionError - Sets the error_rel for the approximation of 875 matrix-vector products using finite differences. 876 877 Logically Collective on Mat 878 879 Input Parameters: 880 + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 881 - error_rel - relative error (should be set to the square root of 882 the relative error in the function evaluations) 883 884 Options Database Keys: 885 . -mat_mffd_err <error_rel> - Sets error_rel 886 887 Level: advanced 888 889 Notes: 890 The default matrix-free matrix-vector product routine computes 891 .vb 892 F'(u)*a = [F(u+h*a) - F(u)]/h where 893 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 894 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 895 .ve 896 897 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 898 MatMFFDSetHHistory(), MatMFFDResetHHistory() 899 @*/ 900 PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 901 { 902 PetscFunctionBegin; 903 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 904 PetscValidLogicalCollectiveReal(mat,error,2); 905 CHKERRQ(PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error))); 906 PetscFunctionReturn(0); 907 } 908 909 /*@ 910 MatMFFDSetHHistory - Sets an array to collect a history of the 911 differencing values (h) computed for the matrix-free product. 912 913 Logically Collective on Mat 914 915 Input Parameters: 916 + J - the matrix-free matrix context 917 . history - space to hold the history 918 - nhistory - number of entries in history, if more entries are generated than 919 nhistory, then the later ones are discarded 920 921 Level: advanced 922 923 Notes: 924 Use MatMFFDResetHHistory() to reset the history counter and collect 925 a new batch of differencing parameters, h. 926 927 .seealso: MatMFFDGetH(), MatCreateSNESMF(), 928 MatMFFDResetHHistory(), MatMFFDSetFunctionError() 929 930 @*/ 931 PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 932 { 933 PetscFunctionBegin; 934 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 935 if (history) PetscValidScalarPointer(history,2); 936 PetscValidLogicalCollectiveInt(J,nhistory,3); 937 CHKERRQ(PetscUseMethod(J,"MatMFFDSetHHistory_C",(Mat,PetscScalar[],PetscInt),(J,history,nhistory))); 938 PetscFunctionReturn(0); 939 } 940 941 /*@ 942 MatMFFDResetHHistory - Resets the counter to zero to begin 943 collecting a new set of differencing histories. 944 945 Logically Collective on Mat 946 947 Input Parameters: 948 . J - the matrix-free matrix context 949 950 Level: advanced 951 952 Notes: 953 Use MatMFFDSetHHistory() to create the original history counter. 954 955 .seealso: MatMFFDGetH(), MatCreateSNESMF(), 956 MatMFFDSetHHistory(), MatMFFDSetFunctionError() 957 958 @*/ 959 PetscErrorCode MatMFFDResetHHistory(Mat J) 960 { 961 PetscFunctionBegin; 962 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 963 CHKERRQ(PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J))); 964 PetscFunctionReturn(0); 965 } 966 967 /*@ 968 MatMFFDSetBase - Sets the vector U at which matrix vector products of the 969 Jacobian are computed 970 971 Logically Collective on Mat 972 973 Input Parameters: 974 + J - the MatMFFD matrix 975 . U - the vector 976 - F - (optional) vector that contains F(u) if it has been already computed 977 978 Notes: 979 This is rarely used directly 980 981 If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 982 point during the first MatMult() after each call to MatMFFDSetBase(). 983 984 Level: advanced 985 986 @*/ 987 PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 988 { 989 PetscFunctionBegin; 990 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 991 PetscValidHeaderSpecific(U,VEC_CLASSID,2); 992 if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 993 CHKERRQ(PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F))); 994 PetscFunctionReturn(0); 995 } 996 997 /*@C 998 MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 999 it to satisfy some criteria 1000 1001 Logically Collective on Mat 1002 1003 Input Parameters: 1004 + J - the MatMFFD matrix 1005 . fun - the function that checks h 1006 - ctx - any context needed by the function 1007 1008 Options Database Keys: 1009 . -mat_mffd_check_positivity <bool> - Insure that U + h*a is non-negative 1010 1011 Level: advanced 1012 1013 Notes: 1014 For example, MatMFFDCheckPositivity() insures that all entries 1015 of U + h*a are non-negative 1016 1017 The function you provide is called after the default h has been computed and allows you to 1018 modify it. 1019 1020 .seealso: MatMFFDCheckPositivity() 1021 @*/ 1022 PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx) 1023 { 1024 PetscFunctionBegin; 1025 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1026 CHKERRQ(PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx))); 1027 PetscFunctionReturn(0); 1028 } 1029 1030 /*@ 1031 MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1032 zero, decreases h until this is satisfied. 1033 1034 Logically Collective on Vec 1035 1036 Input Parameters: 1037 + U - base vector that is added to 1038 . a - vector that is added 1039 . h - scaling factor on a 1040 - dummy - context variable (unused) 1041 1042 Options Database Keys: 1043 . -mat_mffd_check_positivity <bool> - Insure that U + h*a is nonnegative 1044 1045 Level: advanced 1046 1047 Notes: 1048 This is rarely used directly, rather it is passed as an argument to 1049 MatMFFDSetCheckh() 1050 1051 .seealso: MatMFFDSetCheckh() 1052 @*/ 1053 PetscErrorCode MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h) 1054 { 1055 PetscReal val, minval; 1056 PetscScalar *u_vec, *a_vec; 1057 PetscInt i,n; 1058 MPI_Comm comm; 1059 1060 PetscFunctionBegin; 1061 PetscValidHeaderSpecific(U,VEC_CLASSID,2); 1062 PetscValidHeaderSpecific(a,VEC_CLASSID,3); 1063 PetscValidScalarPointer(h,4); 1064 CHKERRQ(PetscObjectGetComm((PetscObject)U,&comm)); 1065 CHKERRQ(VecGetArray(U,&u_vec)); 1066 CHKERRQ(VecGetArray(a,&a_vec)); 1067 CHKERRQ(VecGetLocalSize(U,&n)); 1068 minval = PetscAbsScalar(*h)*PetscRealConstant(1.01); 1069 for (i=0; i<n; i++) { 1070 if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1071 val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1072 if (val < minval) minval = val; 1073 } 1074 } 1075 CHKERRQ(VecRestoreArray(U,&u_vec)); 1076 CHKERRQ(VecRestoreArray(a,&a_vec)); 1077 CHKERRMPI(MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm)); 1078 if (val <= PetscAbsScalar(*h)) { 1079 CHKERRQ(PetscInfo(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val))); 1080 if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1081 else *h = -0.99*val; 1082 } 1083 PetscFunctionReturn(0); 1084 } 1085