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