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