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