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