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