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