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