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