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