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