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