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