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