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