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