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