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