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