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