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