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