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