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