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