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