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 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for matrix-free matrix",((PetscObject)viewer)->type_name); 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 PetscFunctionBegin; 524 shell->vshift += a; 525 PetscFunctionReturn(0); 526 } 527 528 #undef __FUNCT__ 529 #define __FUNCT__ "MatScale_MFFD" 530 PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a) 531 { 532 MatMFFD shell = (MatMFFD)Y->data; 533 PetscFunctionBegin; 534 shell->vscale *= a; 535 PetscFunctionReturn(0); 536 } 537 538 EXTERN_C_BEGIN 539 #undef __FUNCT__ 540 #define __FUNCT__ "MatMFFDSetBase_MFFD" 541 PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F) 542 { 543 PetscErrorCode ierr; 544 MatMFFD ctx = (MatMFFD)J->data; 545 546 PetscFunctionBegin; 547 ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 548 ctx->current_u = U; 549 if (F) { 550 if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);} 551 ctx->current_f = F; 552 ctx->current_f_allocated = PETSC_FALSE; 553 } else if (!ctx->current_f_allocated) { 554 ierr = VecDuplicate(ctx->current_u, &ctx->current_f);CHKERRQ(ierr); 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 EXTERN_C_END 564 typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 565 EXTERN_C_BEGIN 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 EXTERN_C_END 578 579 #undef __FUNCT__ 580 #define __FUNCT__ "MatMFFDSetOptionsPrefix" 581 /*@C 582 MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 583 MatMFFD options in the database. 584 585 Collective on Mat 586 587 Input Parameter: 588 + A - the Mat context 589 - prefix - the prefix to prepend to all option names 590 591 Notes: 592 A hyphen (-) must NOT be given at the beginning of the prefix name. 593 The first character of all runtime options is AUTOMATICALLY the hyphen. 594 595 Level: advanced 596 597 .keywords: SNES, matrix-free, parameters 598 599 .seealso: MatSetFromOptions(), MatCreateSNESMF() 600 @*/ 601 PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[]) 602 603 { 604 MatMFFD mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)PETSC_NULL; 605 PetscErrorCode ierr; 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,PETSC_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 EXTERN_C_BEGIN 647 #undef __FUNCT__ 648 #define __FUNCT__ "MatMFFDSetPeriod_MFFD" 649 PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period) 650 { 651 MatMFFD ctx = (MatMFFD)mat->data; 652 653 PetscFunctionBegin; 654 PetscValidLogicalCollectiveInt(mat,period,2); 655 ctx->recomputeperiod = period; 656 PetscFunctionReturn(0); 657 } 658 EXTERN_C_END 659 660 EXTERN_C_BEGIN 661 #undef __FUNCT__ 662 #define __FUNCT__ "MatMFFDSetFunction_MFFD" 663 PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 664 { 665 MatMFFD ctx = (MatMFFD)mat->data; 666 667 PetscFunctionBegin; 668 ctx->func = func; 669 ctx->funcctx = funcctx; 670 PetscFunctionReturn(0); 671 } 672 EXTERN_C_END 673 674 EXTERN_C_BEGIN 675 #undef __FUNCT__ 676 #define __FUNCT__ "MatMFFDSetFunctionError_MFFD" 677 PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error) 678 { 679 MatMFFD ctx = (MatMFFD)mat->data; 680 681 PetscFunctionBegin; 682 PetscValidLogicalCollectiveReal(mat,error,2); 683 if (error != PETSC_DEFAULT) ctx->error_rel = error; 684 PetscFunctionReturn(0); 685 } 686 EXTERN_C_END 687 688 /*MC 689 MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 690 691 Level: advanced 692 693 .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction() 694 M*/ 695 EXTERN_C_BEGIN 696 #undef __FUNCT__ 697 #define __FUNCT__ "MatCreate_MFFD" 698 PetscErrorCode MatCreate_MFFD(Mat A) 699 { 700 MatMFFD mfctx; 701 PetscErrorCode ierr; 702 703 PetscFunctionBegin; 704 #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 705 ierr = MatMFFDInitializePackage(PETSC_NULL);CHKERRQ(ierr); 706 #endif 707 708 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); 709 mfctx->sp = 0; 710 mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 711 mfctx->recomputeperiod = 1; 712 mfctx->count = 0; 713 mfctx->currenth = 0.0; 714 mfctx->historyh = PETSC_NULL; 715 mfctx->ncurrenth = 0; 716 mfctx->maxcurrenth = 0; 717 ((PetscObject)mfctx)->type_name = 0; 718 719 mfctx->vshift = 0.0; 720 mfctx->vscale = 1.0; 721 722 /* 723 Create the empty data structure to contain compute-h routines. 724 These will be filled in below from the command line options or 725 a later call with MatMFFDSetType() or if that is not called 726 then it will default in the first use of MatMult_MFFD() 727 */ 728 mfctx->ops->compute = 0; 729 mfctx->ops->destroy = 0; 730 mfctx->ops->view = 0; 731 mfctx->ops->setfromoptions = 0; 732 mfctx->hctx = 0; 733 734 mfctx->func = 0; 735 mfctx->funcctx = 0; 736 mfctx->w = PETSC_NULL; 737 738 A->data = mfctx; 739 740 A->ops->mult = MatMult_MFFD; 741 A->ops->destroy = MatDestroy_MFFD; 742 A->ops->view = MatView_MFFD; 743 A->ops->assemblyend = MatAssemblyEnd_MFFD; 744 A->ops->getdiagonal = MatGetDiagonal_MFFD; 745 A->ops->scale = MatScale_MFFD; 746 A->ops->shift = MatShift_MFFD; 747 A->ops->diagonalscale = MatDiagonalScale_MFFD; 748 A->ops->diagonalset = MatDiagonalSet_MFFD; 749 A->ops->setfromoptions = MatSetFromOptions_MFFD; 750 A->assembled = PETSC_TRUE; 751 752 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 753 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 754 755 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetBase_C","MatMFFDSetBase_MFFD",MatMFFDSetBase_MFFD);CHKERRQ(ierr); 756 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioniBase_C","MatMFFDSetFunctioniBase_MFFD",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr); 757 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioni_C","MatMFFDSetFunctioni_MFFD",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr); 758 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunction_C","MatMFFDSetFunction_MFFD",MatMFFDSetFunction_MFFD);CHKERRQ(ierr); 759 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetCheckh_C","MatMFFDSetCheckh_MFFD",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr); 760 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetPeriod_C","MatMFFDSetPeriod_MFFD",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr); 761 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctionError_C","MatMFFDSetFunctionError_MFFD",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr); 762 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDResetHHistory_C","MatMFFDResetHHistory_MFFD",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr); 763 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDAddNullSpace_C","MatMFFDAddNullSpace_MFFD",MatMFFDAddNullSpace_MFFD);CHKERRQ(ierr); 764 mfctx->mat = A; 765 ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr); 766 PetscFunctionReturn(0); 767 } 768 EXTERN_C_END 769 770 #undef __FUNCT__ 771 #define __FUNCT__ "MatCreateMFFD" 772 /*@ 773 MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 774 775 Collective on Vec 776 777 Input Parameters: 778 + comm - MPI communicator 779 . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 780 This value should be the same as the local size used in creating the 781 y vector for the matrix-vector product y = Ax. 782 . n - This value should be the same as the local size used in creating the 783 x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 784 calculated if N is given) For square matrices n is almost always m. 785 . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 786 - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 787 788 789 Output Parameter: 790 . J - the matrix-free matrix 791 792 Options Database Keys: call MatSetFromOptions() to trigger these 793 + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 794 - -mat_mffd_err - square root of estimated relative error in function evaluation 795 - -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 796 797 798 Level: advanced 799 800 Notes: 801 The matrix-free matrix context merely contains the function pointers 802 and work space for performing finite difference approximations of 803 Jacobian-vector products, F'(u)*a, 804 805 The default code uses the following approach to compute h 806 807 .vb 808 F'(u)*a = [F(u+h*a) - F(u)]/h where 809 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 810 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 811 where 812 error_rel = square root of relative error in function evaluation 813 umin = minimum iterate parameter 814 .ve 815 816 The user can set the error_rel via MatMFFDSetFunctionError() and 817 umin via MatMFFDDSSetUmin(); see the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A> for details. 818 819 The user should call MatDestroy() when finished with the matrix-free 820 matrix context. 821 822 Options Database Keys: 823 + -mat_mffd_err <error_rel> - Sets error_rel 824 . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only) 825 - -mat_mffd_check_positivity 826 827 .keywords: default, matrix-free, create, matrix 828 829 .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 830 MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 831 MatMFFDGetH(), MatMFFDRegisterDynamic), MatMFFDComputeJacobian() 832 833 @*/ 834 PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 835 { 836 PetscErrorCode ierr; 837 838 PetscFunctionBegin; 839 ierr = MatCreate(comm,J);CHKERRQ(ierr); 840 ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr); 841 ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 842 ierr = MatSetUp(*J);CHKERRQ(ierr); 843 PetscFunctionReturn(0); 844 } 845 846 847 #undef __FUNCT__ 848 #define __FUNCT__ "MatMFFDGetH" 849 /*@ 850 MatMFFDGetH - Gets the last value that was used as the differencing 851 parameter. 852 853 Not Collective 854 855 Input Parameters: 856 . mat - the matrix obtained with MatCreateSNESMF() 857 858 Output Paramter: 859 . h - the differencing step size 860 861 Level: advanced 862 863 .keywords: SNES, matrix-free, parameters 864 865 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory() 866 @*/ 867 PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 868 { 869 MatMFFD ctx = (MatMFFD)mat->data; 870 PetscErrorCode ierr; 871 PetscBool match; 872 873 PetscFunctionBegin; 874 ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 875 if (!match) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 876 877 *h = ctx->currenth; 878 PetscFunctionReturn(0); 879 } 880 881 #undef __FUNCT__ 882 #define __FUNCT__ "MatMFFDSetFunction" 883 /*@C 884 MatMFFDSetFunction - Sets the function used in applying the matrix free. 885 886 Logically Collective on Mat 887 888 Input Parameters: 889 + mat - the matrix free matrix created via MatCreateSNESMF() 890 . func - the function to use 891 - funcctx - optional function context passed to function 892 893 Level: advanced 894 895 Notes: 896 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 897 matrix inside your compute Jacobian routine 898 899 If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 900 901 .keywords: SNES, matrix-free, function 902 903 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD, 904 MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 905 @*/ 906 PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 907 { 908 PetscErrorCode ierr; 909 910 PetscFunctionBegin; 911 ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr); 912 PetscFunctionReturn(0); 913 } 914 915 #undef __FUNCT__ 916 #define __FUNCT__ "MatMFFDSetFunctioni" 917 /*@C 918 MatMFFDSetFunctioni - Sets the function for a single component 919 920 Logically Collective on Mat 921 922 Input Parameters: 923 + mat - the matrix free matrix created via MatCreateSNESMF() 924 - funci - the function to use 925 926 Level: advanced 927 928 Notes: 929 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 930 matrix inside your compute Jacobian routine 931 932 933 .keywords: SNES, matrix-free, function 934 935 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 936 937 @*/ 938 PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 939 { 940 PetscErrorCode ierr; 941 942 PetscFunctionBegin; 943 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 944 ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr); 945 PetscFunctionReturn(0); 946 } 947 948 949 #undef __FUNCT__ 950 #define __FUNCT__ "MatMFFDSetFunctioniBase" 951 /*@C 952 MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 953 954 Logically Collective on Mat 955 956 Input Parameters: 957 + mat - the matrix free matrix created via MatCreateSNESMF() 958 - func - the function to use 959 960 Level: advanced 961 962 Notes: 963 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 964 matrix inside your compute Jacobian routine 965 966 967 .keywords: SNES, matrix-free, function 968 969 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 970 MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 971 @*/ 972 PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 973 { 974 PetscErrorCode ierr; 975 976 PetscFunctionBegin; 977 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 978 ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr); 979 PetscFunctionReturn(0); 980 } 981 982 #undef __FUNCT__ 983 #define __FUNCT__ "MatMFFDSetPeriod" 984 /*@ 985 MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 986 987 Logically Collective on Mat 988 989 Input Parameters: 990 + mat - the matrix free matrix created via MatCreateSNESMF() 991 - period - 1 for everytime, 2 for every second etc 992 993 Options Database Keys: 994 + -mat_mffd_period <period> 995 996 Level: advanced 997 998 999 .keywords: SNES, matrix-free, parameters 1000 1001 .seealso: MatCreateSNESMF(),MatMFFDGetH(), 1002 MatMFFDSetHHistory(), MatMFFDResetHHistory() 1003 @*/ 1004 PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 1005 { 1006 PetscErrorCode ierr; 1007 1008 PetscFunctionBegin; 1009 ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr); 1010 PetscFunctionReturn(0); 1011 } 1012 1013 #undef __FUNCT__ 1014 #define __FUNCT__ "MatMFFDSetFunctionError" 1015 /*@ 1016 MatMFFDSetFunctionError - Sets the error_rel for the approximation of 1017 matrix-vector products using finite differences. 1018 1019 Logically Collective on Mat 1020 1021 Input Parameters: 1022 + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 1023 - error_rel - relative error (should be set to the square root of 1024 the relative error in the function evaluations) 1025 1026 Options Database Keys: 1027 + -mat_mffd_err <error_rel> - Sets error_rel 1028 1029 Level: advanced 1030 1031 Notes: 1032 The default matrix-free matrix-vector product routine computes 1033 .vb 1034 F'(u)*a = [F(u+h*a) - F(u)]/h where 1035 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 1036 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 1037 .ve 1038 1039 .keywords: SNES, matrix-free, parameters 1040 1041 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 1042 MatMFFDSetHHistory(), MatMFFDResetHHistory() 1043 @*/ 1044 PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 1045 { 1046 PetscErrorCode ierr; 1047 1048 PetscFunctionBegin; 1049 ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr); 1050 PetscFunctionReturn(0); 1051 } 1052 1053 #undef __FUNCT__ 1054 #define __FUNCT__ "MatMFFDAddNullSpace" 1055 /*@ 1056 MatMFFDAddNullSpace - Provides a null space that an operator is 1057 supposed to have. Since roundoff will create a small component in 1058 the null space, if you know the null space you may have it 1059 automatically removed. 1060 1061 Logically Collective on Mat 1062 1063 Input Parameters: 1064 + J - the matrix-free matrix context 1065 - nullsp - object created with MatNullSpaceCreate() 1066 1067 Level: advanced 1068 1069 .keywords: SNES, matrix-free, null space 1070 1071 .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD 1072 MatMFFDSetHHistory(), MatMFFDResetHHistory() 1073 @*/ 1074 PetscErrorCode MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp) 1075 { 1076 PetscErrorCode ierr; 1077 1078 PetscFunctionBegin; 1079 ierr = PetscTryMethod(J,"MatMFFDAddNullSpace_C",(Mat,MatNullSpace),(J,nullsp));CHKERRQ(ierr); 1080 PetscFunctionReturn(0); 1081 } 1082 1083 #undef __FUNCT__ 1084 #define __FUNCT__ "MatMFFDSetHHistory" 1085 /*@ 1086 MatMFFDSetHHistory - Sets an array to collect a history of the 1087 differencing values (h) computed for the matrix-free product. 1088 1089 Logically Collective on Mat 1090 1091 Input Parameters: 1092 + J - the matrix-free matrix context 1093 . histroy - space to hold the history 1094 - nhistory - number of entries in history, if more entries are generated than 1095 nhistory, then the later ones are discarded 1096 1097 Level: advanced 1098 1099 Notes: 1100 Use MatMFFDResetHHistory() to reset the history counter and collect 1101 a new batch of differencing parameters, h. 1102 1103 .keywords: SNES, matrix-free, h history, differencing history 1104 1105 .seealso: MatMFFDGetH(), MatCreateSNESMF(), 1106 MatMFFDResetHHistory(), MatMFFDSetFunctionError() 1107 1108 @*/ 1109 PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 1110 { 1111 MatMFFD ctx = (MatMFFD)J->data; 1112 PetscErrorCode ierr; 1113 PetscBool match; 1114 1115 PetscFunctionBegin; 1116 ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr); 1117 if (!match) SETERRQ(((PetscObject)J)->comm,PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 1118 ctx->historyh = history; 1119 ctx->maxcurrenth = nhistory; 1120 ctx->currenth = 0.; 1121 PetscFunctionReturn(0); 1122 } 1123 1124 1125 #undef __FUNCT__ 1126 #define __FUNCT__ "MatMFFDResetHHistory" 1127 /*@ 1128 MatMFFDResetHHistory - Resets the counter to zero to begin 1129 collecting a new set of differencing histories. 1130 1131 Logically Collective on Mat 1132 1133 Input Parameters: 1134 . J - the matrix-free matrix context 1135 1136 Level: advanced 1137 1138 Notes: 1139 Use MatMFFDSetHHistory() to create the original history counter. 1140 1141 .keywords: SNES, matrix-free, h history, differencing history 1142 1143 .seealso: MatMFFDGetH(), MatCreateSNESMF(), 1144 MatMFFDSetHHistory(), MatMFFDSetFunctionError() 1145 1146 @*/ 1147 PetscErrorCode MatMFFDResetHHistory(Mat J) 1148 { 1149 PetscErrorCode ierr; 1150 1151 PetscFunctionBegin; 1152 ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr); 1153 PetscFunctionReturn(0); 1154 } 1155 1156 1157 #undef __FUNCT__ 1158 #define __FUNCT__ "MatMFFDSetBase" 1159 /*@ 1160 MatMFFDSetBase - Sets the vector U at which matrix vector products of the 1161 Jacobian are computed 1162 1163 Logically Collective on Mat 1164 1165 Input Parameters: 1166 + J - the MatMFFD matrix 1167 . U - the vector 1168 - F - (optional) vector that contains F(u) if it has been already computed 1169 1170 Notes: This is rarely used directly 1171 1172 If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 1173 point during the first MatMult() after each call to MatMFFDSetBase(). 1174 1175 Level: advanced 1176 1177 @*/ 1178 PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 1179 { 1180 PetscErrorCode ierr; 1181 1182 PetscFunctionBegin; 1183 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1184 PetscValidHeaderSpecific(U,VEC_CLASSID,2); 1185 if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 1186 ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr); 1187 PetscFunctionReturn(0); 1188 } 1189 1190 #undef __FUNCT__ 1191 #define __FUNCT__ "MatMFFDSetCheckh" 1192 /*@C 1193 MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 1194 it to satisfy some criteria 1195 1196 Logically Collective on Mat 1197 1198 Input Parameters: 1199 + J - the MatMFFD matrix 1200 . fun - the function that checks h 1201 - ctx - any context needed by the function 1202 1203 Options Database Keys: 1204 . -mat_mffd_check_positivity 1205 1206 Level: advanced 1207 1208 Notes: For example, MatMFFDSetCheckPositivity() insures that all entries 1209 of U + h*a are non-negative 1210 1211 .seealso: MatMFFDSetCheckPositivity() 1212 @*/ 1213 PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void* ctx) 1214 { 1215 PetscErrorCode ierr; 1216 1217 PetscFunctionBegin; 1218 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1219 ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr); 1220 PetscFunctionReturn(0); 1221 } 1222 1223 #undef __FUNCT__ 1224 #define __FUNCT__ "MatMFFDSetCheckPositivity" 1225 /*@ 1226 MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1227 zero, decreases h until this is satisfied. 1228 1229 Logically Collective on Vec 1230 1231 Input Parameters: 1232 + U - base vector that is added to 1233 . a - vector that is added 1234 . h - scaling factor on a 1235 - dummy - context variable (unused) 1236 1237 Options Database Keys: 1238 . -mat_mffd_check_positivity 1239 1240 Level: advanced 1241 1242 Notes: This is rarely used directly, rather it is passed as an argument to 1243 MatMFFDSetCheckh() 1244 1245 .seealso: MatMFFDSetCheckh() 1246 @*/ 1247 PetscErrorCode MatMFFDCheckPositivity(void* dummy,Vec U,Vec a,PetscScalar *h) 1248 { 1249 PetscReal val, minval; 1250 PetscScalar *u_vec, *a_vec; 1251 PetscErrorCode ierr; 1252 PetscInt i,n; 1253 MPI_Comm comm; 1254 1255 PetscFunctionBegin; 1256 ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); 1257 ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); 1258 ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); 1259 ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); 1260 minval = PetscAbsScalar(*h*1.01); 1261 for (i=0;i<n;i++) { 1262 if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1263 val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1264 if (val < minval) minval = val; 1265 } 1266 } 1267 ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); 1268 ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); 1269 ierr = MPI_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr); 1270 if (val <= PetscAbsScalar(*h)) { 1271 ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr); 1272 if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1273 else *h = -0.99*val; 1274 } 1275 PetscFunctionReturn(0); 1276 } 1277 1278 1279 1280 1281 1282 1283 1284