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 VecView(U,0); 285 /* 286 Compute differencing parameter 287 */ 288 if (!ctx->ops->compute) { 289 ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr); 290 ierr = MatMFFDSetFromOptions(mat);CHKERRQ(ierr); 291 } 292 ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr); 293 if (zeroa) { 294 ierr = VecSet(y,0.0);CHKERRQ(ierr); 295 PetscFunctionReturn(0); 296 } 297 298 if (h != h) SETERRQ(PETSC_ERR_PLIB,"Computed Nan differencing parameter h"); 299 if (ctx->checkh) { 300 ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr); 301 } 302 303 /* keep a record of the current differencing parameter h */ 304 ctx->currenth = h; 305 #if defined(PETSC_USE_COMPLEX) 306 ierr = PetscInfo2(mat,"Current differencing parameter: %G + %G i\n",PetscRealPart(h),PetscImaginaryPart(h));CHKERRQ(ierr); 307 #else 308 ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr); 309 #endif 310 if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 311 ctx->historyh[ctx->ncurrenth] = h; 312 } 313 ctx->ncurrenth++; 314 315 /* w = u + ha */ 316 ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr); 317 318 /* compute func(U) as base for differencing */ 319 /* if (ctx->ncurrenth == 1) {*/ 320 ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr); 321 /* }*/ 322 ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr); 323 324 ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 325 ierr = VecScale(y,1.0/h);CHKERRQ(ierr); 326 327 ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr); 328 329 if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);} 330 331 ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 332 PetscFunctionReturn(0); 333 } 334 335 #undef __FUNCT__ 336 #define __FUNCT__ "MatGetDiagonal_MFFD" 337 /* 338 MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 339 340 y ~= (F(u + ha) - F(u))/h, 341 where F = nonlinear function, as set by SNESSetFunction() 342 u = current iterate 343 h = difference interval 344 */ 345 PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 346 { 347 MatMFFD ctx = (MatMFFD)mat->data; 348 PetscScalar h,*aa,*ww,v; 349 PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 350 Vec w,U; 351 PetscErrorCode ierr; 352 PetscInt i,rstart,rend; 353 354 PetscFunctionBegin; 355 if (!ctx->funci) { 356 SETERRQ(PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first"); 357 } 358 359 w = ctx->w; 360 U = ctx->current_u; 361 ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr); 362 ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr); 363 ierr = VecCopy(U,w);CHKERRQ(ierr); 364 365 ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr); 366 ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 367 for (i=rstart; i<rend; i++) { 368 ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 369 h = ww[i-rstart]; 370 if (h == 0.0) h = 1.0; 371 #if !defined(PETSC_USE_COMPLEX) 372 if (h < umin && h >= 0.0) h = umin; 373 else if (h < 0.0 && h > -umin) h = -umin; 374 #else 375 if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 376 else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 377 #endif 378 h *= epsilon; 379 380 ww[i-rstart] += h; 381 ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 382 ierr = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr); 383 aa[i-rstart] = (v - aa[i-rstart])/h; 384 385 /* possibly shift and scale result */ 386 aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart]; 387 388 ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 389 ww[i-rstart] -= h; 390 ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 391 } 392 ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 393 PetscFunctionReturn(0); 394 } 395 396 #undef __FUNCT__ 397 #define __FUNCT__ "MatShift_MFFD" 398 PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a) 399 { 400 MatMFFD shell = (MatMFFD)Y->data; 401 PetscFunctionBegin; 402 shell->vshift += a; 403 PetscFunctionReturn(0); 404 } 405 406 #undef __FUNCT__ 407 #define __FUNCT__ "MatScale_MFFD" 408 PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a) 409 { 410 MatMFFD shell = (MatMFFD)Y->data; 411 PetscFunctionBegin; 412 shell->vscale *= a; 413 PetscFunctionReturn(0); 414 } 415 416 EXTERN_C_BEGIN 417 #undef __FUNCT__ 418 #define __FUNCT__ "MatMFFDSetBase_FD" 419 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase_FD(Mat J,Vec U,Vec F) 420 { 421 PetscErrorCode ierr; 422 MatMFFD ctx = (MatMFFD)J->data; 423 424 PetscFunctionBegin; 425 ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 426 ctx->current_u = U; 427 ctx->current_f = F; 428 if (!ctx->w) { 429 ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr); 430 } 431 J->assembled = PETSC_TRUE; 432 PetscFunctionReturn(0); 433 } 434 EXTERN_C_END 435 typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 436 EXTERN_C_BEGIN 437 #undef __FUNCT__ 438 #define __FUNCT__ "MatMFFDSetCheckh_FD" 439 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh_FD(Mat J,FCN3 fun,void*ectx) 440 { 441 MatMFFD ctx = (MatMFFD)J->data; 442 443 PetscFunctionBegin; 444 ctx->checkh = fun; 445 ctx->checkhctx = ectx; 446 PetscFunctionReturn(0); 447 } 448 EXTERN_C_END 449 450 #undef __FUNCT__ 451 #define __FUNCT__ "MatMFFDSetFromOptions" 452 /*@ 453 MatMFFDSetFromOptions - Sets the MatMFFD options from the command line 454 parameter. 455 456 Collective on Mat 457 458 Input Parameters: 459 . mat - the matrix obtained with MatCreateMFFD() or MatCreateSNESMF() 460 461 Options Database Keys: 462 + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 463 - -mat_mffd_err - square root of estimated relative error in function evaluation 464 - -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 465 466 Level: advanced 467 468 .keywords: SNES, matrix-free, parameters 469 470 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), 471 MatMFFDResetHHistory(), MatMFFDKSPMonitor() 472 @*/ 473 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat mat) 474 { 475 MatMFFD mfctx = (MatMFFD)mat->data; 476 PetscErrorCode ierr; 477 PetscTruth flg; 478 char ftype[256]; 479 480 PetscFunctionBegin; 481 ierr = PetscOptionsBegin(mfctx->comm,mfctx->prefix,"Set matrix free computation parameters","MatMFFD");CHKERRQ(ierr); 482 ierr = PetscOptionsList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDPetscFList,mfctx->type_name,ftype,256,&flg);CHKERRQ(ierr); 483 if (flg) { 484 ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr); 485 } 486 487 ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 488 ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 489 490 ierr = PetscOptionsName("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",&flg);CHKERRQ(ierr); 491 if (flg) { 492 ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr); 493 } 494 if (mfctx->ops->setfromoptions) { 495 ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr); 496 } 497 ierr = PetscOptionsEnd();CHKERRQ(ierr); 498 PetscFunctionReturn(0); 499 } 500 501 /*MC 502 MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 503 504 Level: advanced 505 506 .seealso: MatCreateMFFD(), MatCreateSNESMF() 507 M*/ 508 EXTERN_C_BEGIN 509 #undef __FUNCT__ 510 #define __FUNCT__ "MatCreate_MFFD" 511 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MFFD(Mat A) 512 { 513 MatMFFD mfctx; 514 PetscErrorCode ierr; 515 516 PetscFunctionBegin; 517 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 518 ierr = MatMFFDInitializePackage(PETSC_NULL);CHKERRQ(ierr); 519 #endif 520 521 ierr = PetscHeaderCreate(mfctx,_p_MatMFFD,struct _MFOps,MATMFFD_COOKIE,0,"MatMFFD",A->comm,MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr); 522 mfctx->sp = 0; 523 mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 524 mfctx->recomputeperiod = 1; 525 mfctx->count = 0; 526 mfctx->currenth = 0.0; 527 mfctx->historyh = PETSC_NULL; 528 mfctx->ncurrenth = 0; 529 mfctx->maxcurrenth = 0; 530 mfctx->type_name = 0; 531 532 mfctx->vshift = 0.0; 533 mfctx->vscale = 1.0; 534 535 /* 536 Create the empty data structure to contain compute-h routines. 537 These will be filled in below from the command line options or 538 a later call with MatMFFDSetType() or if that is not called 539 then it will default in the first use of MatMult_MFFD() 540 */ 541 mfctx->ops->compute = 0; 542 mfctx->ops->destroy = 0; 543 mfctx->ops->view = 0; 544 mfctx->ops->setfromoptions = 0; 545 mfctx->hctx = 0; 546 547 mfctx->func = 0; 548 mfctx->funcctx = 0; 549 mfctx->w = PETSC_NULL; 550 551 A->data = mfctx; 552 553 A->ops->mult = MatMult_MFFD; 554 A->ops->destroy = MatDestroy_MFFD; 555 A->ops->view = MatView_MFFD; 556 A->ops->assemblyend = MatAssemblyEnd_MFFD; 557 A->ops->getdiagonal = MatGetDiagonal_MFFD; 558 A->ops->scale = MatScale_MFFD; 559 A->ops->shift = MatShift_MFFD; 560 A->ops->setfromoptions = MatMFFDSetFromOptions; 561 A->assembled = PETSC_TRUE; 562 563 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetBase_C","MatMFFDSetBase_FD",MatMFFDSetBase_FD);CHKERRQ(ierr); 564 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioniBase_C","MatMFFDSetFunctioniBase_FD",MatMFFDSetFunctioniBase_FD);CHKERRQ(ierr); 565 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioni_C","MatMFFDSetFunctioni_FD",MatMFFDSetFunctioni_FD);CHKERRQ(ierr); 566 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetCheckh_C","MatMFFDSetCheckh_FD",MatMFFDSetCheckh_FD);CHKERRQ(ierr); 567 mfctx->mat = A; 568 ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr); 569 PetscFunctionReturn(0); 570 } 571 EXTERN_C_END 572 573 #undef __FUNCT__ 574 #define __FUNCT__ "MatCreateMFFD" 575 /*@ 576 MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 577 578 Collective on Vec 579 580 Input Parameters: 581 . x - vector that defines layout of the vectors and matrices 582 583 Output Parameter: 584 . J - the matrix-free matrix 585 586 Level: advanced 587 588 Notes: 589 The matrix-free matrix context merely contains the function pointers 590 and work space for performing finite difference approximations of 591 Jacobian-vector products, F'(u)*a, 592 593 The default code uses the following approach to compute h 594 595 .vb 596 F'(u)*a = [F(u+h*a) - F(u)]/h where 597 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 598 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 599 where 600 error_rel = square root of relative error in function evaluation 601 umin = minimum iterate parameter 602 .ve 603 604 The user can set the error_rel via MatMFFDSetFunctionError() and 605 umin via MatMFFDDefaultSetUmin(); see the nonlinear solvers chapter 606 of the users manual for details. 607 608 The user should call MatDestroy() when finished with the matrix-free 609 matrix context. 610 611 Options Database Keys: 612 + -mat_mffd_err <error_rel> - Sets error_rel 613 . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only) 614 . -mat_mffd_ksp_monitor - KSP monitor routine that prints differencing h 615 - -mat_mffd_check_positivity 616 617 .keywords: default, matrix-free, create, matrix 618 619 .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDefaultSetUmin(), MatMFFDSetFunction() 620 MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 621 MatMFFDGetH(),MatMFFDKSPMonitor(), MatMFFDRegisterDynamic),, MatMFFDComputeJacobian() 622 623 @*/ 624 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(Vec x,Mat *J) 625 { 626 MPI_Comm comm; 627 PetscErrorCode ierr; 628 PetscInt n,nloc; 629 630 PetscFunctionBegin; 631 ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr); 632 ierr = VecGetSize(x,&n);CHKERRQ(ierr); 633 ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr); 634 ierr = MatCreate(comm,J);CHKERRQ(ierr); 635 ierr = MatSetSizes(*J,nloc,nloc,n,n);CHKERRQ(ierr); 636 ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 637 PetscFunctionReturn(0); 638 } 639 640 641 #undef __FUNCT__ 642 #define __FUNCT__ "MatMFFDGetH" 643 /*@ 644 MatMFFDGetH - Gets the last value that was used as the differencing 645 parameter. 646 647 Not Collective 648 649 Input Parameters: 650 . mat - the matrix obtained with MatCreateSNESMF() 651 652 Output Paramter: 653 . h - the differencing step size 654 655 Level: advanced 656 657 .keywords: SNES, matrix-free, parameters 658 659 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD 660 MatMFFDResetHHistory(),MatMFFDKSPMonitor() 661 @*/ 662 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat mat,PetscScalar *h) 663 { 664 MatMFFD ctx = (MatMFFD)mat->data; 665 666 PetscFunctionBegin; 667 *h = ctx->currenth; 668 PetscFunctionReturn(0); 669 } 670 671 672 #undef __FUNCT__ 673 #define __FUNCT__ "MatMFFDSetFunction" 674 /*@C 675 MatMFFDSetFunction - Sets the function used in applying the matrix free. 676 677 Collective on Mat 678 679 Input Parameters: 680 + mat - the matrix free matrix created via MatCreateSNESMF() 681 . func - the function to use 682 - funcctx - optional function context passed to function 683 684 Level: advanced 685 686 Notes: 687 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 688 matrix inside your compute Jacobian routine 689 690 If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 691 692 .keywords: SNES, matrix-free, function 693 694 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 695 MatMFFDSetHHistory(), MatMFFDResetHHistory(), 696 MatMFFDKSPMonitor(), SNESetFunction() 697 @*/ 698 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 699 { 700 MatMFFD ctx = (MatMFFD)mat->data; 701 702 PetscFunctionBegin; 703 ctx->func = func; 704 ctx->funcctx = funcctx; 705 PetscFunctionReturn(0); 706 } 707 708 #undef __FUNCT__ 709 #define __FUNCT__ "MatMFFDSetFunctioni" 710 /*@C 711 MatMFFDSetFunctioni - Sets the function for a single component 712 713 Collective on Mat 714 715 Input Parameters: 716 + mat - the matrix free matrix created via MatCreateSNESMF() 717 - funci - the function to use 718 719 Level: advanced 720 721 Notes: 722 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 723 matrix inside your compute Jacobian routine 724 725 726 .keywords: SNES, matrix-free, function 727 728 .seealso: MatCreateSNESMF(),MatMFFDGetH(), 729 MatMFFDSetHHistory(), MatMFFDResetHHistory(), 730 MatMFFDKSPMonitor(), SNESetFunction() 731 @*/ 732 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 733 { 734 PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)); 735 736 PetscFunctionBegin; 737 PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 738 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",(void (**)(void))&f);CHKERRQ(ierr); 739 if (f) { 740 ierr = (*f)(mat,funci);CHKERRQ(ierr); 741 } 742 PetscFunctionReturn(0); 743 } 744 745 746 #undef __FUNCT__ 747 #define __FUNCT__ "MatMFFDSetFunctioniBase" 748 /*@C 749 MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 750 751 Collective on Mat 752 753 Input Parameters: 754 + mat - the matrix free matrix created via MatCreateSNESMF() 755 - func - the function to use 756 757 Level: advanced 758 759 Notes: 760 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 761 matrix inside your compute Jacobian routine 762 763 764 .keywords: SNES, matrix-free, function 765 766 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 767 MatMFFDSetHHistory(), MatMFFDResetHHistory(), 768 MatMFFDKSPMonitor(), SNESetFunction() 769 @*/ 770 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 771 { 772 PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(void*,Vec)); 773 774 PetscFunctionBegin; 775 PetscValidHeaderSpecific(mat,MAT_COOKIE,1); 776 ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",(void (**)(void))&f);CHKERRQ(ierr); 777 if (f) { 778 ierr = (*f)(mat,func);CHKERRQ(ierr); 779 } 780 PetscFunctionReturn(0); 781 } 782 783 784 #undef __FUNCT__ 785 #define __FUNCT__ "MatMFFDSetPeriod" 786 /*@ 787 MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 788 789 Collective on Mat 790 791 Input Parameters: 792 + mat - the matrix free matrix created via MatCreateSNESMF() 793 - period - 1 for everytime, 2 for every second etc 794 795 Options Database Keys: 796 + -mat_mffd_period <period> 797 798 Level: advanced 799 800 801 .keywords: SNES, matrix-free, parameters 802 803 .seealso: MatCreateSNESMF(),MatMFFDGetH(), 804 MatMFFDSetHHistory(), MatMFFDResetHHistory(), 805 MatMFFDKSPMonitor() 806 @*/ 807 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat mat,PetscInt period) 808 { 809 MatMFFD ctx = (MatMFFD)mat->data; 810 811 PetscFunctionBegin; 812 ctx->recomputeperiod = period; 813 PetscFunctionReturn(0); 814 } 815 816 #undef __FUNCT__ 817 #define __FUNCT__ "MatMFFDSetFunctionError" 818 /*@ 819 MatMFFDSetFunctionError - Sets the error_rel for the approximation of 820 matrix-vector products using finite differences. 821 822 Collective on Mat 823 824 Input Parameters: 825 + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 826 - error_rel - relative error (should be set to the square root of 827 the relative error in the function evaluations) 828 829 Options Database Keys: 830 + -mat_mffd_err <error_rel> - Sets error_rel 831 832 Level: advanced 833 834 Notes: 835 The default matrix-free matrix-vector product routine computes 836 .vb 837 F'(u)*a = [F(u+h*a) - F(u)]/h where 838 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 839 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 840 .ve 841 842 .keywords: SNES, matrix-free, parameters 843 844 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 845 MatMFFDSetHHistory(), MatMFFDResetHHistory(), 846 MatMFFDKSPMonitor() 847 @*/ 848 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat mat,PetscReal error) 849 { 850 MatMFFD ctx = (MatMFFD)mat->data; 851 852 PetscFunctionBegin; 853 if (error != PETSC_DEFAULT) ctx->error_rel = error; 854 PetscFunctionReturn(0); 855 } 856 857 #undef __FUNCT__ 858 #define __FUNCT__ "MatMFFDAddNullSpace" 859 /*@ 860 MatMFFDAddNullSpace - Provides a null space that an operator is 861 supposed to have. Since roundoff will create a small component in 862 the null space, if you know the null space you may have it 863 automatically removed. 864 865 Collective on Mat 866 867 Input Parameters: 868 + J - the matrix-free matrix context 869 - nullsp - object created with MatNullSpaceCreate() 870 871 Level: advanced 872 873 .keywords: SNES, matrix-free, null space 874 875 .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD 876 MatMFFDSetHHistory(), MatMFFDResetHHistory(), 877 MatMFFDKSPMonitor(), MatMFFDErrorRel() 878 @*/ 879 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp) 880 { 881 PetscErrorCode ierr; 882 MatMFFD ctx = (MatMFFD)J->data; 883 884 PetscFunctionBegin; 885 ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr); 886 if (ctx->sp) { ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr); } 887 ctx->sp = nullsp; 888 PetscFunctionReturn(0); 889 } 890 891 #undef __FUNCT__ 892 #define __FUNCT__ "MatMFFDSetHHistory" 893 /*@ 894 MatMFFDSetHHistory - Sets an array to collect a history of the 895 differencing values (h) computed for the matrix-free product. 896 897 Collective on Mat 898 899 Input Parameters: 900 + J - the matrix-free matrix context 901 . histroy - space to hold the history 902 - nhistory - number of entries in history, if more entries are generated than 903 nhistory, then the later ones are discarded 904 905 Level: advanced 906 907 Notes: 908 Use MatMFFDResetHHistory() to reset the history counter and collect 909 a new batch of differencing parameters, h. 910 911 .keywords: SNES, matrix-free, h history, differencing history 912 913 .seealso: MatMFFDGetH(), MatCreateSNESMF(), 914 MatMFFDResetHHistory(), 915 MatMFFDKSPMonitor(), MatMFFDSetFunctionError() 916 917 @*/ 918 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 919 { 920 MatMFFD ctx = (MatMFFD)J->data; 921 922 PetscFunctionBegin; 923 ctx->historyh = history; 924 ctx->maxcurrenth = nhistory; 925 ctx->currenth = 0; 926 PetscFunctionReturn(0); 927 } 928 929 #undef __FUNCT__ 930 #define __FUNCT__ "MatMFFDResetHHistory" 931 /*@ 932 MatMFFDResetHHistory - Resets the counter to zero to begin 933 collecting a new set of differencing histories. 934 935 Collective on Mat 936 937 Input Parameters: 938 . J - the matrix-free matrix context 939 940 Level: advanced 941 942 Notes: 943 Use MatMFFDSetHHistory() to create the original history counter. 944 945 .keywords: SNES, matrix-free, h history, differencing history 946 947 .seealso: MatMFFDGetH(), MatCreateSNESMF(), 948 MatMFFDSetHHistory(), 949 MatMFFDKSPMonitor(), MatMFFDSetFunctionError() 950 951 @*/ 952 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat J) 953 { 954 MatMFFD ctx = (MatMFFD)J->data; 955 956 PetscFunctionBegin; 957 ctx->ncurrenth = 0; 958 PetscFunctionReturn(0); 959 } 960 961 962 #undef __FUNCT__ 963 #define __FUNCT__ "MatMFFDSetBase" 964 /*@ 965 MatMFFDSetBase - Sets the vector U at which matrix vector products of the 966 Jacobian are computed 967 968 Collective on Mat 969 970 Input Parameters: 971 + J - the MatMFFD matrix 972 . U - the vector 973 - F - vector that contains F(u) if it has been already computed 974 975 Notes: This is rarely used directly 976 977 Level: advanced 978 979 @*/ 980 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat J,Vec U,Vec F) 981 { 982 PetscErrorCode ierr,(*f)(Mat,Vec,Vec); 983 984 PetscFunctionBegin; 985 PetscValidHeaderSpecific(J,MAT_COOKIE,1); 986 PetscValidHeaderSpecific(U,VEC_COOKIE,2); 987 ierr = PetscObjectQueryFunction((PetscObject)J,"MatMFFDSetBase_C",(void (**)(void))&f);CHKERRQ(ierr); 988 if (f) { 989 ierr = (*f)(J,U,F);CHKERRQ(ierr); 990 } 991 PetscFunctionReturn(0); 992 } 993 994 #undef __FUNCT__ 995 #define __FUNCT__ "MatMFFDSetCheckh" 996 /*@C 997 MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 998 it to satisfy some criteria 999 1000 Collective on Mat 1001 1002 Input Parameters: 1003 + J - the MatMFFD matrix 1004 . fun - the function that checks h 1005 - ctx - any context needed by the function 1006 1007 Options Database Keys: 1008 . -mat_mffd_check_positivity 1009 1010 Level: advanced 1011 1012 Notes: For example, MatMFFDSetCheckPositivity() insures that all entries 1013 of U + h*a are non-negative 1014 1015 .seealso: MatMFFDSetCheckPositivity() 1016 @*/ 1017 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void* ctx) 1018 { 1019 PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*); 1020 1021 PetscFunctionBegin; 1022 PetscValidHeaderSpecific(J,MAT_COOKIE,1); 1023 ierr = PetscObjectQueryFunction((PetscObject)J,"MatMFFDSetCheckh_C",(void (**)(void))&f);CHKERRQ(ierr); 1024 if (f) { 1025 ierr = (*f)(J,fun,ctx);CHKERRQ(ierr); 1026 } 1027 PetscFunctionReturn(0); 1028 } 1029 1030 #undef __FUNCT__ 1031 #define __FUNCT__ "MatMFFDSetCheckPositivity" 1032 /*@ 1033 MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1034 zero, decreases h until this is satisfied. 1035 1036 Collective on Vec 1037 1038 Input Parameters: 1039 + U - base vector that is added to 1040 . a - vector that is added 1041 . h - scaling factor on a 1042 - dummy - context variable (unused) 1043 1044 Options Database Keys: 1045 . -mat_mffd_check_positivity 1046 1047 Level: advanced 1048 1049 Notes: This is rarely used directly, rather it is passed as an argument to 1050 MatMFFDSetCheckh() 1051 1052 .seealso: MatMFFDSetCheckh() 1053 @*/ 1054 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void* dummy,Vec U,Vec a,PetscScalar *h) 1055 { 1056 PetscReal val, minval; 1057 PetscScalar *u_vec, *a_vec; 1058 PetscErrorCode ierr; 1059 PetscInt i,n; 1060 MPI_Comm comm; 1061 1062 PetscFunctionBegin; 1063 ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); 1064 ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); 1065 ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); 1066 ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); 1067 minval = PetscAbsScalar(*h*1.01); 1068 for(i=0;i<n;i++) { 1069 if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1070 val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1071 if (val < minval) minval = val; 1072 } 1073 } 1074 ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); 1075 ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); 1076 ierr = PetscGlobalMin(&minval,&val,comm);CHKERRQ(ierr); 1077 if (val <= PetscAbsScalar(*h)) { 1078 ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr); 1079 if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1080 else *h = -0.99*val; 1081 } 1082 PetscFunctionReturn(0); 1083 } 1084 1085 1086 1087 1088 1089 1090 1091