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