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 = MatSetFromOptions(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: MatSetFromOptions(), 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__ "MatSetFromOptions_MFFD" 581 PetscErrorCode MatSetFromOptions_MFFD(Mat mat) 582 { 583 MatMFFD mfctx = (MatMFFD)mat->data; 584 PetscErrorCode ierr; 585 PetscBool flg; 586 char ftype[256]; 587 588 PetscFunctionBegin; 589 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 590 PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 591 ierr = PetscOptionsBegin(((PetscObject)mfctx)->comm,((PetscObject)mfctx)->prefix,"Set matrix free computation parameters","MatMFFD");CHKERRQ(ierr); 592 ierr = PetscOptionsList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr); 593 if (flg) { 594 ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr); 595 } 596 597 ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 598 ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 599 600 flg = PETSC_FALSE; 601 ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 602 if (flg) { 603 ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr); 604 } 605 if (mfctx->ops->setfromoptions) { 606 ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr); 607 } 608 ierr = PetscOptionsEnd();CHKERRQ(ierr); 609 PetscFunctionReturn(0); 610 } 611 612 /*MC 613 MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 614 615 Level: advanced 616 617 .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction() 618 M*/ 619 EXTERN_C_BEGIN 620 #undef __FUNCT__ 621 #define __FUNCT__ "MatCreate_MFFD" 622 PetscErrorCode MatCreate_MFFD(Mat A) 623 { 624 MatMFFD mfctx; 625 PetscErrorCode ierr; 626 627 PetscFunctionBegin; 628 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 629 ierr = MatMFFDInitializePackage(PETSC_NULL);CHKERRQ(ierr); 630 #endif 631 632 ierr = PetscHeaderCreate(mfctx,_p_MatMFFD,struct _MFOps,MATMFFD_CLASSID,0,"MatMFFD",((PetscObject)A)->comm,MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr); 633 mfctx->sp = 0; 634 mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 635 mfctx->recomputeperiod = 1; 636 mfctx->count = 0; 637 mfctx->currenth = 0.0; 638 mfctx->historyh = PETSC_NULL; 639 mfctx->ncurrenth = 0; 640 mfctx->maxcurrenth = 0; 641 ((PetscObject)mfctx)->type_name = 0; 642 643 mfctx->vshift = 0.0; 644 mfctx->vscale = 1.0; 645 646 /* 647 Create the empty data structure to contain compute-h routines. 648 These will be filled in below from the command line options or 649 a later call with MatMFFDSetType() or if that is not called 650 then it will default in the first use of MatMult_MFFD() 651 */ 652 mfctx->ops->compute = 0; 653 mfctx->ops->destroy = 0; 654 mfctx->ops->view = 0; 655 mfctx->ops->setfromoptions = 0; 656 mfctx->hctx = 0; 657 658 mfctx->func = 0; 659 mfctx->funcctx = 0; 660 mfctx->w = PETSC_NULL; 661 662 A->data = mfctx; 663 664 A->ops->mult = MatMult_MFFD; 665 A->ops->destroy = MatDestroy_MFFD; 666 A->ops->view = MatView_MFFD; 667 A->ops->assemblyend = MatAssemblyEnd_MFFD; 668 A->ops->getdiagonal = MatGetDiagonal_MFFD; 669 A->ops->scale = MatScale_MFFD; 670 A->ops->shift = MatShift_MFFD; 671 A->ops->diagonalscale = MatDiagonalScale_MFFD; 672 A->ops->diagonalset = MatDiagonalSet_MFFD; 673 A->ops->setfromoptions = MatSetFromOptions_MFFD; 674 A->assembled = PETSC_TRUE; 675 676 ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr); 677 ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr); 678 ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 679 ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 680 681 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetBase_C","MatMFFDSetBase_MFFD",MatMFFDSetBase_MFFD);CHKERRQ(ierr); 682 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioniBase_C","MatMFFDSetFunctioniBase_MFFD",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr); 683 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioni_C","MatMFFDSetFunctioni_MFFD",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr); 684 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetCheckh_C","MatMFFDSetCheckh_MFFD",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr); 685 mfctx->mat = A; 686 ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr); 687 PetscFunctionReturn(0); 688 } 689 EXTERN_C_END 690 691 #undef __FUNCT__ 692 #define __FUNCT__ "MatCreateMFFD" 693 /*@ 694 MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 695 696 Collective on Vec 697 698 Input Parameters: 699 + comm - MPI communicator 700 . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 701 This value should be the same as the local size used in creating the 702 y vector for the matrix-vector product y = Ax. 703 . n - This value should be the same as the local size used in creating the 704 x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 705 calculated if N is given) For square matrices n is almost always m. 706 . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 707 - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 708 709 710 Output Parameter: 711 . J - the matrix-free matrix 712 713 Options Database Keys: call MatSetFromOptions() to trigger these 714 + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 715 - -mat_mffd_err - square root of estimated relative error in function evaluation 716 - -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 717 718 719 Level: advanced 720 721 Notes: 722 The matrix-free matrix context merely contains the function pointers 723 and work space for performing finite difference approximations of 724 Jacobian-vector products, F'(u)*a, 725 726 The default code uses the following approach to compute h 727 728 .vb 729 F'(u)*a = [F(u+h*a) - F(u)]/h where 730 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 731 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 732 where 733 error_rel = square root of relative error in function evaluation 734 umin = minimum iterate parameter 735 .ve 736 737 The user can set the error_rel via MatMFFDSetFunctionError() and 738 umin via MatMFFDDSSetUmin(); see the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A> for details. 739 740 The user should call MatDestroy() when finished with the matrix-free 741 matrix context. 742 743 Options Database Keys: 744 + -mat_mffd_err <error_rel> - Sets error_rel 745 . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only) 746 - -mat_mffd_check_positivity 747 748 .keywords: default, matrix-free, create, matrix 749 750 .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 751 MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 752 MatMFFDGetH(), MatMFFDRegisterDynamic), MatMFFDComputeJacobian() 753 754 @*/ 755 PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 756 { 757 PetscErrorCode ierr; 758 759 PetscFunctionBegin; 760 ierr = MatCreate(comm,J);CHKERRQ(ierr); 761 ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr); 762 ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 763 PetscFunctionReturn(0); 764 } 765 766 767 #undef __FUNCT__ 768 #define __FUNCT__ "MatMFFDGetH" 769 /*@ 770 MatMFFDGetH - Gets the last value that was used as the differencing 771 parameter. 772 773 Not Collective 774 775 Input Parameters: 776 . mat - the matrix obtained with MatCreateSNESMF() 777 778 Output Paramter: 779 . h - the differencing step size 780 781 Level: advanced 782 783 .keywords: SNES, matrix-free, parameters 784 785 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory() 786 @*/ 787 PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 788 { 789 MatMFFD ctx = (MatMFFD)mat->data; 790 791 PetscFunctionBegin; 792 *h = ctx->currenth; 793 PetscFunctionReturn(0); 794 } 795 796 797 #undef __FUNCT__ 798 #define __FUNCT__ "MatMFFDSetFunction" 799 /*@C 800 MatMFFDSetFunction - Sets the function used in applying the matrix free. 801 802 Logically Collective on Mat 803 804 Input Parameters: 805 + mat - the matrix free matrix created via MatCreateSNESMF() 806 . func - the function to use 807 - funcctx - optional function context passed to function 808 809 Level: advanced 810 811 Notes: 812 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 813 matrix inside your compute Jacobian routine 814 815 If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 816 817 .keywords: SNES, matrix-free, function 818 819 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD, 820 MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 821 @*/ 822 PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 823 { 824 MatMFFD ctx = (MatMFFD)mat->data; 825 826 PetscFunctionBegin; 827 ctx->func = func; 828 ctx->funcctx = funcctx; 829 PetscFunctionReturn(0); 830 } 831 832 #undef __FUNCT__ 833 #define __FUNCT__ "MatMFFDSetFunctioni" 834 /*@C 835 MatMFFDSetFunctioni - Sets the function for a single component 836 837 Logically Collective on Mat 838 839 Input Parameters: 840 + mat - the matrix free matrix created via MatCreateSNESMF() 841 - funci - the function to use 842 843 Level: advanced 844 845 Notes: 846 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 847 matrix inside your compute Jacobian routine 848 849 850 .keywords: SNES, matrix-free, function 851 852 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 853 854 @*/ 855 PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 856 { 857 PetscErrorCode ierr; 858 859 PetscFunctionBegin; 860 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 861 ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr); 862 PetscFunctionReturn(0); 863 } 864 865 866 #undef __FUNCT__ 867 #define __FUNCT__ "MatMFFDSetFunctioniBase" 868 /*@C 869 MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 870 871 Logically Collective on Mat 872 873 Input Parameters: 874 + mat - the matrix free matrix created via MatCreateSNESMF() 875 - func - the function to use 876 877 Level: advanced 878 879 Notes: 880 If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 881 matrix inside your compute Jacobian routine 882 883 884 .keywords: SNES, matrix-free, function 885 886 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 887 MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 888 @*/ 889 PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 890 { 891 PetscErrorCode ierr; 892 893 PetscFunctionBegin; 894 PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 895 ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr); 896 PetscFunctionReturn(0); 897 } 898 899 900 #undef __FUNCT__ 901 #define __FUNCT__ "MatMFFDSetPeriod" 902 /*@ 903 MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 904 905 Logically Collective on Mat 906 907 Input Parameters: 908 + mat - the matrix free matrix created via MatCreateSNESMF() 909 - period - 1 for everytime, 2 for every second etc 910 911 Options Database Keys: 912 + -mat_mffd_period <period> 913 914 Level: advanced 915 916 917 .keywords: SNES, matrix-free, parameters 918 919 .seealso: MatCreateSNESMF(),MatMFFDGetH(), 920 MatMFFDSetHHistory(), MatMFFDResetHHistory() 921 @*/ 922 PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 923 { 924 MatMFFD ctx = (MatMFFD)mat->data; 925 926 PetscFunctionBegin; 927 PetscValidLogicalCollectiveInt(mat,period,2); 928 ctx->recomputeperiod = period; 929 PetscFunctionReturn(0); 930 } 931 932 #undef __FUNCT__ 933 #define __FUNCT__ "MatMFFDSetFunctionError" 934 /*@ 935 MatMFFDSetFunctionError - Sets the error_rel for the approximation of 936 matrix-vector products using finite differences. 937 938 Logically Collective on Mat 939 940 Input Parameters: 941 + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 942 - error_rel - relative error (should be set to the square root of 943 the relative error in the function evaluations) 944 945 Options Database Keys: 946 + -mat_mffd_err <error_rel> - Sets error_rel 947 948 Level: advanced 949 950 Notes: 951 The default matrix-free matrix-vector product routine computes 952 .vb 953 F'(u)*a = [F(u+h*a) - F(u)]/h where 954 h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 955 = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 956 .ve 957 958 .keywords: SNES, matrix-free, parameters 959 960 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 961 MatMFFDSetHHistory(), MatMFFDResetHHistory() 962 @*/ 963 PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 964 { 965 MatMFFD ctx = (MatMFFD)mat->data; 966 967 PetscFunctionBegin; 968 PetscValidLogicalCollectiveReal(mat,error,2); 969 if (error != PETSC_DEFAULT) ctx->error_rel = error; 970 PetscFunctionReturn(0); 971 } 972 973 #undef __FUNCT__ 974 #define __FUNCT__ "MatMFFDAddNullSpace" 975 /*@ 976 MatMFFDAddNullSpace - Provides a null space that an operator is 977 supposed to have. Since roundoff will create a small component in 978 the null space, if you know the null space you may have it 979 automatically removed. 980 981 Logically Collective on Mat 982 983 Input Parameters: 984 + J - the matrix-free matrix context 985 - nullsp - object created with MatNullSpaceCreate() 986 987 Level: advanced 988 989 .keywords: SNES, matrix-free, null space 990 991 .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD 992 MatMFFDSetHHistory(), MatMFFDResetHHistory() 993 @*/ 994 PetscErrorCode MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp) 995 { 996 PetscErrorCode ierr; 997 MatMFFD ctx = (MatMFFD)J->data; 998 999 PetscFunctionBegin; 1000 ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr); 1001 if (ctx->sp) { ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr); } 1002 ctx->sp = nullsp; 1003 PetscFunctionReturn(0); 1004 } 1005 1006 #undef __FUNCT__ 1007 #define __FUNCT__ "MatMFFDSetHHistory" 1008 /*@ 1009 MatMFFDSetHHistory - Sets an array to collect a history of the 1010 differencing values (h) computed for the matrix-free product. 1011 1012 Logically Collective on Mat 1013 1014 Input Parameters: 1015 + J - the matrix-free matrix context 1016 . histroy - space to hold the history 1017 - nhistory - number of entries in history, if more entries are generated than 1018 nhistory, then the later ones are discarded 1019 1020 Level: advanced 1021 1022 Notes: 1023 Use MatMFFDResetHHistory() to reset the history counter and collect 1024 a new batch of differencing parameters, h. 1025 1026 .keywords: SNES, matrix-free, h history, differencing history 1027 1028 .seealso: MatMFFDGetH(), MatCreateSNESMF(), 1029 MatMFFDResetHHistory(), MatMFFDSetFunctionError() 1030 1031 @*/ 1032 PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 1033 { 1034 MatMFFD ctx = (MatMFFD)J->data; 1035 1036 PetscFunctionBegin; 1037 ctx->historyh = history; 1038 ctx->maxcurrenth = nhistory; 1039 ctx->currenth = 0.; 1040 PetscFunctionReturn(0); 1041 } 1042 1043 #undef __FUNCT__ 1044 #define __FUNCT__ "MatMFFDResetHHistory" 1045 /*@ 1046 MatMFFDResetHHistory - Resets the counter to zero to begin 1047 collecting a new set of differencing histories. 1048 1049 Logically Collective on Mat 1050 1051 Input Parameters: 1052 . J - the matrix-free matrix context 1053 1054 Level: advanced 1055 1056 Notes: 1057 Use MatMFFDSetHHistory() to create the original history counter. 1058 1059 .keywords: SNES, matrix-free, h history, differencing history 1060 1061 .seealso: MatMFFDGetH(), MatCreateSNESMF(), 1062 MatMFFDSetHHistory(), MatMFFDSetFunctionError() 1063 1064 @*/ 1065 PetscErrorCode MatMFFDResetHHistory(Mat J) 1066 { 1067 MatMFFD ctx = (MatMFFD)J->data; 1068 1069 PetscFunctionBegin; 1070 ctx->ncurrenth = 0; 1071 PetscFunctionReturn(0); 1072 } 1073 1074 1075 #undef __FUNCT__ 1076 #define __FUNCT__ "MatMFFDSetBase" 1077 /*@ 1078 MatMFFDSetBase - Sets the vector U at which matrix vector products of the 1079 Jacobian are computed 1080 1081 Logically Collective on Mat 1082 1083 Input Parameters: 1084 + J - the MatMFFD matrix 1085 . U - the vector 1086 - F - (optional) vector that contains F(u) if it has been already computed 1087 1088 Notes: This is rarely used directly 1089 1090 If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 1091 point during the first MatMult() after each call to MatMFFDSetBase(). 1092 1093 Level: advanced 1094 1095 @*/ 1096 PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 1097 { 1098 PetscErrorCode ierr; 1099 1100 PetscFunctionBegin; 1101 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1102 PetscValidHeaderSpecific(U,VEC_CLASSID,2); 1103 if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 1104 ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr); 1105 PetscFunctionReturn(0); 1106 } 1107 1108 #undef __FUNCT__ 1109 #define __FUNCT__ "MatMFFDSetCheckh" 1110 /*@C 1111 MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 1112 it to satisfy some criteria 1113 1114 Logically Collective on Mat 1115 1116 Input Parameters: 1117 + J - the MatMFFD matrix 1118 . fun - the function that checks h 1119 - ctx - any context needed by the function 1120 1121 Options Database Keys: 1122 . -mat_mffd_check_positivity 1123 1124 Level: advanced 1125 1126 Notes: For example, MatMFFDSetCheckPositivity() insures that all entries 1127 of U + h*a are non-negative 1128 1129 .seealso: MatMFFDSetCheckPositivity() 1130 @*/ 1131 PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void* ctx) 1132 { 1133 PetscErrorCode ierr; 1134 1135 PetscFunctionBegin; 1136 PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1137 ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr); 1138 PetscFunctionReturn(0); 1139 } 1140 1141 #undef __FUNCT__ 1142 #define __FUNCT__ "MatMFFDSetCheckPositivity" 1143 /*@ 1144 MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1145 zero, decreases h until this is satisfied. 1146 1147 Logically Collective on Vec 1148 1149 Input Parameters: 1150 + U - base vector that is added to 1151 . a - vector that is added 1152 . h - scaling factor on a 1153 - dummy - context variable (unused) 1154 1155 Options Database Keys: 1156 . -mat_mffd_check_positivity 1157 1158 Level: advanced 1159 1160 Notes: This is rarely used directly, rather it is passed as an argument to 1161 MatMFFDSetCheckh() 1162 1163 .seealso: MatMFFDSetCheckh() 1164 @*/ 1165 PetscErrorCode MatMFFDCheckPositivity(void* dummy,Vec U,Vec a,PetscScalar *h) 1166 { 1167 PetscReal val, minval; 1168 PetscScalar *u_vec, *a_vec; 1169 PetscErrorCode ierr; 1170 PetscInt i,n; 1171 MPI_Comm comm; 1172 1173 PetscFunctionBegin; 1174 ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); 1175 ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); 1176 ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); 1177 ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); 1178 minval = PetscAbsScalar(*h*1.01); 1179 for(i=0;i<n;i++) { 1180 if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1181 val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1182 if (val < minval) minval = val; 1183 } 1184 } 1185 ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); 1186 ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); 1187 ierr = MPI_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr); 1188 if (val <= PetscAbsScalar(*h)) { 1189 ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr); 1190 if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1191 else *h = -0.99*val; 1192 } 1193 PetscFunctionReturn(0); 1194 } 1195 1196 1197 1198 1199 1200 1201 1202