1e884886fSBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> 3c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h> /*I "petscmat.h" I*/ 4e884886fSBarry Smith 5f4259b30SLisandro Dalcin PetscFunctionList MatMFFDList = NULL; 6ace3abfcSBarry Smith PetscBool MatMFFDRegisterAllCalled = PETSC_FALSE; 7e884886fSBarry Smith 87087cfbeSBarry Smith PetscClassId MATMFFD_CLASSID; 9166c7f25SBarry Smith PetscLogEvent MATMFFD_Mult; 10e884886fSBarry Smith 11ace3abfcSBarry Smith static PetscBool MatMFFDPackageInitialized = PETSC_FALSE; 12b022a5c1SBarry Smith /*@C 132390153bSJed Brown MatMFFDFinalizePackage - This function destroys everything in the MatMFFD package. It is 14b022a5c1SBarry Smith called from PetscFinalize(). 15b022a5c1SBarry Smith 16b022a5c1SBarry Smith Level: developer 17b022a5c1SBarry Smith 18db781477SPatrick Sanan .seealso: `PetscFinalize()`, `MatCreateMFFD()`, `MatCreateSNESMF()` 19b022a5c1SBarry Smith @*/ 207087cfbeSBarry Smith PetscErrorCode MatMFFDFinalizePackage(void) 21b022a5c1SBarry Smith { 22b022a5c1SBarry Smith PetscFunctionBegin; 239566063dSJacob Faibussowitsch PetscCall(PetscFunctionListDestroy(&MatMFFDList)); 24b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_FALSE; 25b022a5c1SBarry Smith MatMFFDRegisterAllCalled = PETSC_FALSE; 26b022a5c1SBarry Smith PetscFunctionReturn(0); 27b022a5c1SBarry Smith } 28b022a5c1SBarry Smith 293ec795f1SBarry Smith /*@C 303ec795f1SBarry Smith MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called 318a690491SBarry Smith from MatInitializePackage(). 323ec795f1SBarry Smith 333ec795f1SBarry Smith Level: developer 343ec795f1SBarry Smith 35db781477SPatrick Sanan .seealso: `PetscInitialize()` 363ec795f1SBarry Smith @*/ 37607a6623SBarry Smith PetscErrorCode MatMFFDInitializePackage(void) 383ec795f1SBarry Smith { 393ec795f1SBarry Smith char logList[256]; 408e81d068SLisandro Dalcin PetscBool opt,pkg; 413ec795f1SBarry Smith 423ec795f1SBarry Smith PetscFunctionBegin; 43b022a5c1SBarry Smith if (MatMFFDPackageInitialized) PetscFunctionReturn(0); 44b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_TRUE; 453ec795f1SBarry Smith /* Register Classes */ 469566063dSJacob Faibussowitsch PetscCall(PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID)); 473ec795f1SBarry Smith /* Register Constructors */ 489566063dSJacob Faibussowitsch PetscCall(MatMFFDRegisterAll()); 493ec795f1SBarry Smith /* Register Events */ 509566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister("MatMult MF",MATMFFD_CLASSID,&MATMFFD_Mult)); 51e94e781bSJacob Faibussowitsch /* Process Info */ 52e94e781bSJacob Faibussowitsch { 53e94e781bSJacob Faibussowitsch PetscClassId classids[1]; 54e94e781bSJacob Faibussowitsch 55e94e781bSJacob Faibussowitsch classids[0] = MATMFFD_CLASSID; 569566063dSJacob Faibussowitsch PetscCall(PetscInfoProcessClass("matmffd", 1, classids)); 573ec795f1SBarry Smith } 583ec795f1SBarry Smith /* Process summary exclusions */ 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt)); 603ec795f1SBarry Smith if (opt) { 619566063dSJacob Faibussowitsch PetscCall(PetscStrInList("matmffd",logList,',',&pkg)); 629566063dSJacob Faibussowitsch if (pkg) PetscCall(PetscLogEventExcludeClass(MATMFFD_CLASSID)); 633ec795f1SBarry Smith } 648e81d068SLisandro Dalcin /* Register package finalizer */ 659566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(MatMFFDFinalizePackage)); 663ec795f1SBarry Smith PetscFunctionReturn(0); 673ec795f1SBarry Smith } 683ec795f1SBarry Smith 69789d8953SBarry Smith static PetscErrorCode MatMFFDSetType_MFFD(Mat mat,MatMFFDType ftype) 70789d8953SBarry Smith { 71789d8953SBarry Smith MatMFFD ctx; 72789d8953SBarry Smith PetscBool match; 735f80ce2aSJacob Faibussowitsch PetscErrorCode (*r)(MatMFFD); 74789d8953SBarry Smith 75789d8953SBarry Smith PetscFunctionBegin; 76789d8953SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 77789d8953SBarry Smith PetscValidCharPointer(ftype,2); 789566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ctx)); 79789d8953SBarry Smith 80789d8953SBarry Smith /* already set, so just return */ 819566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)ctx,ftype,&match)); 82789d8953SBarry Smith if (match) PetscFunctionReturn(0); 83789d8953SBarry Smith 84789d8953SBarry Smith /* destroy the old one if it exists */ 859566063dSJacob Faibussowitsch if (ctx->ops->destroy) PetscCall((*ctx->ops->destroy)(ctx)); 86789d8953SBarry Smith 879566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(MatMFFDList,ftype,&r)); 885f80ce2aSJacob Faibussowitsch PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype); 899566063dSJacob Faibussowitsch PetscCall((*r)(ctx)); 909566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)ctx,ftype)); 91789d8953SBarry Smith PetscFunctionReturn(0); 92789d8953SBarry Smith } 93789d8953SBarry Smith 94e884886fSBarry Smith /*@C 95e884886fSBarry Smith MatMFFDSetType - Sets the method that is used to compute the 96e884886fSBarry Smith differencing parameter for finite differene matrix-free formulations. 97e884886fSBarry Smith 98e884886fSBarry Smith Input Parameters: 99e884886fSBarry Smith + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD() 100e884886fSBarry Smith or MatSetType(mat,MATMFFD); 101e884886fSBarry Smith - ftype - the type requested, either MATMFFD_WP or MATMFFD_DS 102e884886fSBarry Smith 103e884886fSBarry Smith Level: advanced 104e884886fSBarry Smith 105e884886fSBarry Smith Notes: 106e884886fSBarry Smith For example, such routines can compute h for use in 107e884886fSBarry Smith Jacobian-vector products of the form 108e884886fSBarry Smith 109e884886fSBarry Smith F(x+ha) - F(x) 110e884886fSBarry Smith F'(u)a ~= ---------------- 111e884886fSBarry Smith h 112e884886fSBarry Smith 113db781477SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDRegister()`, `MatMFFDSetFunction()`, `MatCreateMFFD()` 114e884886fSBarry Smith @*/ 11519fd82e9SBarry Smith PetscErrorCode MatMFFDSetType(Mat mat,MatMFFDType ftype) 116e884886fSBarry Smith { 117e884886fSBarry Smith PetscFunctionBegin; 1180700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 119e884886fSBarry Smith PetscValidCharPointer(ftype,2); 120cac4c232SBarry Smith PetscTryMethod(mat,"MatMFFDSetType_C",(Mat,MatMFFDType),(mat,ftype)); 121e884886fSBarry Smith PetscFunctionReturn(0); 122e884886fSBarry Smith } 123e884886fSBarry Smith 1245d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat,Vec); 1255d21e1f8SStefano Zampini 126e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/ 12740244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func) 128e884886fSBarry Smith { 129789d8953SBarry Smith MatMFFD ctx; 130e884886fSBarry Smith 131e884886fSBarry Smith PetscFunctionBegin; 1329566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ctx)); 133e884886fSBarry Smith ctx->funcisetbase = func; 134e884886fSBarry Smith PetscFunctionReturn(0); 135e884886fSBarry Smith } 136e884886fSBarry Smith 137e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 13840244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci) 139e884886fSBarry Smith { 140789d8953SBarry Smith MatMFFD ctx; 141e884886fSBarry Smith 142e884886fSBarry Smith PetscFunctionBegin; 1439566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ctx)); 144e884886fSBarry Smith ctx->funci = funci; 1459566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(mat,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_MFFD)); 146789d8953SBarry Smith PetscFunctionReturn(0); 1475d21e1f8SStefano Zampini } 148789d8953SBarry Smith 149789d8953SBarry Smith static PetscErrorCode MatMFFDGetH_MFFD(Mat mat,PetscScalar *h) 150789d8953SBarry Smith { 151789d8953SBarry Smith MatMFFD ctx; 152789d8953SBarry Smith 153789d8953SBarry Smith PetscFunctionBegin; 1549566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ctx)); 155789d8953SBarry Smith *h = ctx->currenth; 156e884886fSBarry Smith PetscFunctionReturn(0); 157e884886fSBarry Smith } 158e884886fSBarry Smith 15940244768SBarry Smith static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) 160bc0f08ceSBarry Smith { 161789d8953SBarry Smith MatMFFD ctx; 162bc0f08ceSBarry Smith 163bc0f08ceSBarry Smith PetscFunctionBegin; 1649566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J,&ctx)); 165bc0f08ceSBarry Smith ctx->ncurrenth = 0; 166bc0f08ceSBarry Smith PetscFunctionReturn(0); 167bc0f08ceSBarry Smith } 168e884886fSBarry Smith 1691c84c290SBarry Smith /*@C 1701c84c290SBarry Smith MatMFFDRegister - Adds a method to the MatMFFD registry. 1711c84c290SBarry Smith 1721c84c290SBarry Smith Not Collective 1731c84c290SBarry Smith 1741c84c290SBarry Smith Input Parameters: 1751c84c290SBarry Smith + name_solver - name of a new user-defined compute-h module 1761c84c290SBarry Smith - routine_create - routine to create method context 1771c84c290SBarry Smith 1781c84c290SBarry Smith Level: developer 1791c84c290SBarry Smith 1801c84c290SBarry Smith Notes: 1811c84c290SBarry Smith MatMFFDRegister() may be called multiple times to add several user-defined solvers. 1821c84c290SBarry Smith 1831c84c290SBarry Smith Sample usage: 1841c84c290SBarry Smith .vb 185bdf89e91SBarry Smith MatMFFDRegister("my_h",MyHCreate); 1861c84c290SBarry Smith .ve 1871c84c290SBarry Smith 1881c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 1891c84c290SBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1901c84c290SBarry Smith or at runtime via the option 191be7a6d03SBarry Smith $ -mat_mffd_type my_h 1921c84c290SBarry Smith 193db781477SPatrick Sanan .seealso: `MatMFFDRegisterAll()`, `MatMFFDRegisterDestroy()` 1941c84c290SBarry Smith @*/ 195bdf89e91SBarry Smith PetscErrorCode MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD)) 196e884886fSBarry Smith { 197e884886fSBarry Smith PetscFunctionBegin; 1989566063dSJacob Faibussowitsch PetscCall(MatInitializePackage()); 1999566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&MatMFFDList,sname,function)); 200e884886fSBarry Smith PetscFunctionReturn(0); 201e884886fSBarry Smith } 202e884886fSBarry Smith 203e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/ 20440244768SBarry Smith static PetscErrorCode MatDestroy_MFFD(Mat mat) 205e884886fSBarry Smith { 206789d8953SBarry Smith MatMFFD ctx; 207e884886fSBarry Smith 208e884886fSBarry Smith PetscFunctionBegin; 2099566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ctx)); 2109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->w)); 2119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->current_u)); 212cfe22f5eSBarry Smith if (ctx->current_f_allocated) { 2139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx->current_f)); 214cfe22f5eSBarry Smith } 2159566063dSJacob Faibussowitsch if (ctx->ops->destroy) PetscCall((*ctx->ops->destroy)(ctx)); 2169566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(&ctx)); 217e884886fSBarry Smith 2189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL)); 2199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL)); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL)); 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL)); 2239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL)); 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL)); 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL)); 2269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetHHistory_C",NULL)); 2279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetType_C",NULL)); 2289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatMFFDGetH_C",NULL)); 229*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFSetReuseBase_C",NULL)); 230*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)mat,"MatSNESMFGetReuseBase_C",NULL)); 231e884886fSBarry Smith PetscFunctionReturn(0); 232e884886fSBarry Smith } 233e884886fSBarry Smith 234e884886fSBarry Smith /* 235e884886fSBarry Smith MatMFFDView_MFFD - Views matrix-free parameters. 236e884886fSBarry Smith 237e884886fSBarry Smith */ 23840244768SBarry Smith static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer) 239e884886fSBarry Smith { 240789d8953SBarry Smith MatMFFD ctx; 241a283635eSDmitry Karpeev PetscBool iascii, viewbase, viewfunction; 242a283635eSDmitry Karpeev const char *prefix; 243e884886fSBarry Smith 244e884886fSBarry Smith PetscFunctionBegin; 2459566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J,&ctx)); 2469566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 247e884886fSBarry Smith if (iascii) { 2489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n")); 2499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel)); 2517adad957SLisandro Dalcin if (!((PetscObject)ctx)->type_name) { 2529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n")); 253e884886fSBarry Smith } else { 2549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name)); 255e884886fSBarry Smith } 256c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 257c92e3469SBarry Smith if (ctx->usecomplex) { 2589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Using Lyness complex number trick to compute the matrix-vector product\n")); 259c92e3469SBarry Smith } 260c92e3469SBarry Smith #endif 261e884886fSBarry Smith if (ctx->ops->view) { 2629566063dSJacob Faibussowitsch PetscCall((*ctx->ops->view)(ctx,viewer)); 263e884886fSBarry Smith } 2649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)J, &prefix)); 265a283635eSDmitry Karpeev 2669566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase)); 267a283635eSDmitry Karpeev if (viewbase) { 2689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Base:\n")); 2699566063dSJacob Faibussowitsch PetscCall(VecView(ctx->current_u, viewer)); 270a283635eSDmitry Karpeev } 2719566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction)); 272a283635eSDmitry Karpeev if (viewfunction) { 2739566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Function:\n")); 2749566063dSJacob Faibussowitsch PetscCall(VecView(ctx->current_f, viewer)); 275a283635eSDmitry Karpeev } 2769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 27711aeaf0aSBarry Smith } 278e884886fSBarry Smith PetscFunctionReturn(0); 279e884886fSBarry Smith } 280e884886fSBarry Smith 281e884886fSBarry Smith /* 282e884886fSBarry Smith MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 283e884886fSBarry Smith allows the user to indicate the beginning of a new linear solve by calling 284e884886fSBarry Smith MatAssemblyXXX() on the matrix free matrix. This then allows the 2851d0fab5eSBarry Smith MatCreateMFFD_WP() to properly compute ||U|| only the first time 286e884886fSBarry Smith in the linear solver rather than every time. 2875a576424SJed Brown 288cc2e6a90SBarry Smith This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence 289cc2e6a90SBarry Smith it must be labeled as PETSC_EXTERN 290e884886fSBarry Smith */ 2915a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt) 292e884886fSBarry Smith { 293789d8953SBarry Smith MatMFFD j; 294e884886fSBarry Smith 295e884886fSBarry Smith PetscFunctionBegin; 2969566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J,&j)); 2979566063dSJacob Faibussowitsch PetscCall(MatMFFDResetHHistory(J)); 298e884886fSBarry Smith PetscFunctionReturn(0); 299e884886fSBarry Smith } 300e884886fSBarry Smith 301e884886fSBarry Smith /* 302e884886fSBarry Smith MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 303e884886fSBarry Smith 304e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 305e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 306e884886fSBarry Smith u = current iterate 307e884886fSBarry Smith h = difference interval 308e884886fSBarry Smith */ 30940244768SBarry Smith static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y) 310e884886fSBarry Smith { 311789d8953SBarry Smith MatMFFD ctx; 312e884886fSBarry Smith PetscScalar h; 313e884886fSBarry Smith Vec w,U,F; 314ace3abfcSBarry Smith PetscBool zeroa; 315e884886fSBarry Smith 316e884886fSBarry Smith PetscFunctionBegin; 3179566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ctx)); 31828b400f6SJacob Faibussowitsch PetscCheck(ctx->current_u,PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONGSTATE,"MatMFFDSetBase() has not been called, this is often caused by forgetting to call \n\t\tMatAssemblyBegin/End on the first Mat in the SNES compute function"); 319e884886fSBarry Smith /* We log matrix-free matrix-vector products separately, so that we can 320e884886fSBarry Smith separate the performance monitoring from the cases that use conventional 321e884886fSBarry Smith storage. We may eventually modify event logging to associate events 322e884886fSBarry Smith with particular objects, hence alleviating the more general problem. */ 3239566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(MATMFFD_Mult,a,y,0,0)); 324e884886fSBarry Smith 325e884886fSBarry Smith w = ctx->w; 326e884886fSBarry Smith U = ctx->current_u; 3273ec795f1SBarry Smith F = ctx->current_f; 328e884886fSBarry Smith /* 329e884886fSBarry Smith Compute differencing parameter 330e884886fSBarry Smith */ 3314a2f8832SBarry Smith if (!((PetscObject)ctx)->type_name) { 3329566063dSJacob Faibussowitsch PetscCall(MatMFFDSetType(mat,MATMFFD_WP)); 3339566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(mat)); 334e884886fSBarry Smith } 3359566063dSJacob Faibussowitsch PetscCall((*ctx->ops->compute)(ctx,U,a,&h,&zeroa)); 336e884886fSBarry Smith if (zeroa) { 3379566063dSJacob Faibussowitsch PetscCall(VecSet(y,0.0)); 33885d14658SLisandro Dalcin PetscCall(PetscLogEventEnd(MATMFFD_Mult,a,y,0,0)); 339e884886fSBarry Smith PetscFunctionReturn(0); 340e884886fSBarry Smith } 341e884886fSBarry Smith 34208401ef6SPierre Jolivet PetscCheck(!mat->erroriffailure || !PetscIsInfOrNanScalar(h),PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h"); 343e884886fSBarry Smith if (ctx->checkh) { 3449566063dSJacob Faibussowitsch PetscCall((*ctx->checkh)(ctx->checkhctx,U,a,&h)); 345e884886fSBarry Smith } 346e884886fSBarry Smith 347e884886fSBarry Smith /* keep a record of the current differencing parameter h */ 348e884886fSBarry Smith ctx->currenth = h; 349e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX) 3509566063dSJacob Faibussowitsch PetscCall(PetscInfo(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h))); 351e884886fSBarry Smith #else 3529566063dSJacob Faibussowitsch PetscCall(PetscInfo(mat,"Current differencing parameter: %15.12e\n",(double)PetscRealPart(h))); 353e884886fSBarry Smith #endif 354e884886fSBarry Smith if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 355e884886fSBarry Smith ctx->historyh[ctx->ncurrenth] = h; 356e884886fSBarry Smith } 357e884886fSBarry Smith ctx->ncurrenth++; 358e884886fSBarry Smith 359c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 360c92e3469SBarry Smith if (ctx->usecomplex) h = PETSC_i*h; 361c92e3469SBarry Smith #endif 362c92e3469SBarry Smith 363e884886fSBarry Smith /* w = u + ha */ 3649566063dSJacob Faibussowitsch PetscCall(VecWAXPY(w,h,a,U)); 365e884886fSBarry Smith 3661b797266SDmitry Karpeev /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 3671b797266SDmitry Karpeev if (ctx->ncurrenth == 1 && ctx->current_f_allocated) { 3689566063dSJacob Faibussowitsch PetscCall((*ctx->func)(ctx->funcctx,U,F)); 36952121784SBarry Smith } 3709566063dSJacob Faibussowitsch PetscCall((*ctx->func)(ctx->funcctx,w,y)); 371e884886fSBarry Smith 372c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 373c92e3469SBarry Smith if (ctx->usecomplex) { 3749566063dSJacob Faibussowitsch PetscCall(VecImaginaryPart(y)); 375c92e3469SBarry Smith h = PetscImaginaryPart(h); 376c92e3469SBarry Smith } else { 3779566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,F)); 378c92e3469SBarry Smith } 379c92e3469SBarry Smith #else 3809566063dSJacob Faibussowitsch PetscCall(VecAXPY(y,-1.0,F)); 381c92e3469SBarry Smith #endif 3829566063dSJacob Faibussowitsch PetscCall(VecScale(y,1.0/h)); 3839566063dSJacob Faibussowitsch if (mat->nullsp) PetscCall(MatNullSpaceRemove(mat->nullsp,y)); 384e884886fSBarry Smith 3859566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(MATMFFD_Mult,a,y,0,0)); 386e884886fSBarry Smith PetscFunctionReturn(0); 387e884886fSBarry Smith } 388e884886fSBarry Smith 389e884886fSBarry Smith /* 390e884886fSBarry Smith MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 391e884886fSBarry Smith 392e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 393e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 394e884886fSBarry Smith u = current iterate 395e884886fSBarry Smith h = difference interval 396e884886fSBarry Smith */ 3975d21e1f8SStefano Zampini PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 398e884886fSBarry Smith { 399789d8953SBarry Smith MatMFFD ctx; 400e884886fSBarry Smith PetscScalar h,*aa,*ww,v; 401e884886fSBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 402e884886fSBarry Smith Vec w,U; 403e884886fSBarry Smith PetscInt i,rstart,rend; 404e884886fSBarry Smith 405e884886fSBarry Smith PetscFunctionBegin; 4069566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ctx)); 40728b400f6SJacob Faibussowitsch PetscCheck(ctx->func,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunction() first"); 40828b400f6SJacob Faibussowitsch PetscCheck(ctx->funci,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioni() first"); 409e884886fSBarry Smith w = ctx->w; 410e884886fSBarry Smith U = ctx->current_u; 4119566063dSJacob Faibussowitsch PetscCall((*ctx->func)(ctx->funcctx,U,a)); 412789d8953SBarry Smith if (ctx->funcisetbase) { 4139566063dSJacob Faibussowitsch PetscCall((*ctx->funcisetbase)(ctx->funcctx,U)); 414789d8953SBarry Smith } 4159566063dSJacob Faibussowitsch PetscCall(VecCopy(U,w)); 416e884886fSBarry Smith 4179566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(a,&rstart,&rend)); 4189566063dSJacob Faibussowitsch PetscCall(VecGetArray(a,&aa)); 419e884886fSBarry Smith for (i=rstart; i<rend; i++) { 4209566063dSJacob Faibussowitsch PetscCall(VecGetArray(w,&ww)); 421e884886fSBarry Smith h = ww[i-rstart]; 422e884886fSBarry Smith if (h == 0.0) h = 1.0; 423e884886fSBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 424e884886fSBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 425e884886fSBarry Smith h *= epsilon; 426e884886fSBarry Smith 427e884886fSBarry Smith ww[i-rstart] += h; 4289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w,&ww)); 4299566063dSJacob Faibussowitsch PetscCall((*ctx->funci)(ctx->funcctx,i,w,&v)); 430e884886fSBarry Smith aa[i-rstart] = (v - aa[i-rstart])/h; 431e884886fSBarry Smith 4329566063dSJacob Faibussowitsch PetscCall(VecGetArray(w,&ww)); 433e884886fSBarry Smith ww[i-rstart] -= h; 4349566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(w,&ww)); 435e884886fSBarry Smith } 4369566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a,&aa)); 437e884886fSBarry Smith PetscFunctionReturn(0); 438e884886fSBarry Smith } 439e884886fSBarry Smith 440d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F) 441e884886fSBarry Smith { 442789d8953SBarry Smith MatMFFD ctx; 443e884886fSBarry Smith 444e884886fSBarry Smith PetscFunctionBegin; 4459566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J,&ctx)); 4469566063dSJacob Faibussowitsch PetscCall(MatMFFDResetHHistory(J)); 447c51ad4d4SStefano Zampini if (!ctx->current_u) { 4489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U,&ctx->current_u)); 4499566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(ctx->current_u)); 450c51ad4d4SStefano Zampini } 4519566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(ctx->current_u)); 4529566063dSJacob Faibussowitsch PetscCall(VecCopy(U,ctx->current_u)); 4539566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(ctx->current_u)); 45452121784SBarry Smith if (F) { 4559566063dSJacob Faibussowitsch if (ctx->current_f_allocated) PetscCall(VecDestroy(&ctx->current_f)); 4563ec795f1SBarry Smith ctx->current_f = F; 457cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_FALSE; 458cfe22f5eSBarry Smith } else if (!ctx->current_f_allocated) { 4599566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(J,NULL,&ctx->current_f)); 460cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_TRUE; 46152121784SBarry Smith } 462e884886fSBarry Smith if (!ctx->w) { 4639566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ctx->current_u,&ctx->w)); 464e884886fSBarry Smith } 465e884886fSBarry Smith J->assembled = PETSC_TRUE; 466e884886fSBarry Smith PetscFunctionReturn(0); 467e884886fSBarry Smith } 4685a576424SJed Brown 469e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 470b2573a8aSBarry Smith 47140244768SBarry Smith static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx) 472e884886fSBarry Smith { 473789d8953SBarry Smith MatMFFD ctx; 474e884886fSBarry Smith 475e884886fSBarry Smith PetscFunctionBegin; 4769566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J,&ctx)); 477e884886fSBarry Smith ctx->checkh = fun; 478e884886fSBarry Smith ctx->checkhctx = ectx; 479e884886fSBarry Smith PetscFunctionReturn(0); 480e884886fSBarry Smith } 481e884886fSBarry Smith 4826aa9148fSLisandro Dalcin /*@C 4836aa9148fSLisandro Dalcin MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 4846aa9148fSLisandro Dalcin MatMFFD options in the database. 4856aa9148fSLisandro Dalcin 4866aa9148fSLisandro Dalcin Collective on Mat 4876aa9148fSLisandro Dalcin 488d8d19677SJose E. Roman Input Parameters: 4896aa9148fSLisandro Dalcin + A - the Mat context 4906aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names 4916aa9148fSLisandro Dalcin 4926aa9148fSLisandro Dalcin Notes: 4936aa9148fSLisandro Dalcin A hyphen (-) must NOT be given at the beginning of the prefix name. 4946aa9148fSLisandro Dalcin The first character of all runtime options is AUTOMATICALLY the hyphen. 4956aa9148fSLisandro Dalcin 4966aa9148fSLisandro Dalcin Level: advanced 4976aa9148fSLisandro Dalcin 498db781477SPatrick Sanan .seealso: `MatSetFromOptions()`, `MatCreateSNESMF()`, `MatCreateMFFD()` 4996aa9148fSLisandro Dalcin @*/ 5007087cfbeSBarry Smith PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[]) 5016aa9148fSLisandro Dalcin { 502789d8953SBarry Smith MatMFFD mfctx; 5035fd66863SKarl Rupp 5046aa9148fSLisandro Dalcin PetscFunctionBegin; 5050700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5069566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&mfctx)); 5070700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 5089566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix)); 5096aa9148fSLisandro Dalcin PetscFunctionReturn(0); 5106aa9148fSLisandro Dalcin } 5116aa9148fSLisandro Dalcin 51240244768SBarry Smith static PetscErrorCode MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat) 513e884886fSBarry Smith { 514789d8953SBarry Smith MatMFFD mfctx; 515ace3abfcSBarry Smith PetscBool flg; 516e884886fSBarry Smith char ftype[256]; 517e884886fSBarry Smith 518e884886fSBarry Smith PetscFunctionBegin; 519064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(mat,MAT_CLASSID,2); 5209566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&mfctx)); 521064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,2); 522d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)mfctx); 5239566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg)); 524e884886fSBarry Smith if (flg) { 5259566063dSJacob Faibussowitsch PetscCall(MatMFFDSetType(mat,ftype)); 526e884886fSBarry Smith } 527e884886fSBarry Smith 5289566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,NULL)); 5299566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,NULL)); 530e884886fSBarry Smith 53190d69ab7SBarry Smith flg = PETSC_FALSE; 5329566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL)); 533e884886fSBarry Smith if (flg) { 5349566063dSJacob Faibussowitsch PetscCall(MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,NULL)); 535e884886fSBarry Smith } 536c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 5379566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_mffd_complex","Use Lyness complex number trick to compute the matrix-vector product","None",mfctx->usecomplex,&mfctx->usecomplex,NULL)); 538c92e3469SBarry Smith #endif 539e884886fSBarry Smith if (mfctx->ops->setfromoptions) { 5409566063dSJacob Faibussowitsch PetscCall((*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx)); 541e884886fSBarry Smith } 542d0609cedSBarry Smith PetscOptionsEnd(); 543e884886fSBarry Smith PetscFunctionReturn(0); 544e884886fSBarry Smith } 545e884886fSBarry Smith 54640244768SBarry Smith static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period) 547bc0f08ceSBarry Smith { 548789d8953SBarry Smith MatMFFD ctx; 549bc0f08ceSBarry Smith 550bc0f08ceSBarry Smith PetscFunctionBegin; 5519566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ctx)); 552bc0f08ceSBarry Smith ctx->recomputeperiod = period; 553bc0f08ceSBarry Smith PetscFunctionReturn(0); 554bc0f08ceSBarry Smith } 555bc0f08ceSBarry Smith 55640244768SBarry Smith static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 557bc0f08ceSBarry Smith { 558789d8953SBarry Smith MatMFFD ctx; 559bc0f08ceSBarry Smith 560bc0f08ceSBarry Smith PetscFunctionBegin; 5619566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ctx)); 562bc0f08ceSBarry Smith ctx->func = func; 563bc0f08ceSBarry Smith ctx->funcctx = funcctx; 564bc0f08ceSBarry Smith PetscFunctionReturn(0); 565bc0f08ceSBarry Smith } 566bc0f08ceSBarry Smith 56740244768SBarry Smith static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error) 568bc0f08ceSBarry Smith { 569789d8953SBarry Smith MatMFFD ctx; 570bc0f08ceSBarry Smith 571bc0f08ceSBarry Smith PetscFunctionBegin; 5729566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&ctx)); 573bc0f08ceSBarry Smith if (error != PETSC_DEFAULT) ctx->error_rel = error; 574bc0f08ceSBarry Smith PetscFunctionReturn(0); 575bc0f08ceSBarry Smith } 576bc0f08ceSBarry Smith 577789d8953SBarry Smith PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J,PetscScalar history[],PetscInt nhistory) 5783b49f96aSBarry Smith { 579789d8953SBarry Smith MatMFFD ctx; 580789d8953SBarry Smith 5813b49f96aSBarry Smith PetscFunctionBegin; 5829566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(J,&ctx)); 583789d8953SBarry Smith ctx->historyh = history; 584789d8953SBarry Smith ctx->maxcurrenth = nhistory; 585789d8953SBarry Smith ctx->currenth = 0.; 5863b49f96aSBarry Smith PetscFunctionReturn(0); 5873b49f96aSBarry Smith } 5883b49f96aSBarry Smith 589e884886fSBarry Smith /*MC 590e884886fSBarry Smith MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 591e884886fSBarry Smith 592e884886fSBarry Smith Level: advanced 593e884886fSBarry Smith 594789d8953SBarry Smith Developers Note: This is implemented on top of MATSHELL to get support for scaling and shifting without requiring duplicate code 595789d8953SBarry Smith 596db781477SPatrick Sanan .seealso: `MatCreateMFFD()`, `MatCreateSNESMF()`, `MatMFFDSetFunction()`, `MatMFFDSetType()`, 597db781477SPatrick Sanan `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()` 598db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`, 599db781477SPatrick Sanan `MatMFFDGetH()`, 600e884886fSBarry Smith M*/ 6018cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 602e884886fSBarry Smith { 603e884886fSBarry Smith MatMFFD mfctx; 604e884886fSBarry Smith 605e884886fSBarry Smith PetscFunctionBegin; 6069566063dSJacob Faibussowitsch PetscCall(MatMFFDInitializePackage()); 607e884886fSBarry Smith 6089566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),NULL,NULL)); 6092205254eSKarl Rupp 610e884886fSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 611e884886fSBarry Smith mfctx->recomputeperiod = 1; 612e884886fSBarry Smith mfctx->count = 0; 613e884886fSBarry Smith mfctx->currenth = 0.0; 6140298fd71SBarry Smith mfctx->historyh = NULL; 615e884886fSBarry Smith mfctx->ncurrenth = 0; 616e884886fSBarry Smith mfctx->maxcurrenth = 0; 617f4259b30SLisandro Dalcin ((PetscObject)mfctx)->type_name = NULL; 618e884886fSBarry Smith 619e884886fSBarry Smith /* 620e884886fSBarry Smith Create the empty data structure to contain compute-h routines. 621e884886fSBarry Smith These will be filled in below from the command line options or 622e884886fSBarry Smith a later call with MatMFFDSetType() or if that is not called 623e884886fSBarry Smith then it will default in the first use of MatMult_MFFD() 624e884886fSBarry Smith */ 625f4259b30SLisandro Dalcin mfctx->ops->compute = NULL; 626f4259b30SLisandro Dalcin mfctx->ops->destroy = NULL; 627f4259b30SLisandro Dalcin mfctx->ops->view = NULL; 628f4259b30SLisandro Dalcin mfctx->ops->setfromoptions = NULL; 629f4259b30SLisandro Dalcin mfctx->hctx = NULL; 630e884886fSBarry Smith 631f4259b30SLisandro Dalcin mfctx->func = NULL; 632f4259b30SLisandro Dalcin mfctx->funcctx = NULL; 6330298fd71SBarry Smith mfctx->w = NULL; 634789d8953SBarry Smith mfctx->mat = A; 635e884886fSBarry Smith 6369566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSHELL)); 6379566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(A,mfctx)); 6389566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A,MATOP_MULT,(void (*)(void))MatMult_MFFD)); 6399566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A,MATOP_DESTROY,(void (*)(void))MatDestroy_MFFD)); 6409566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A,MATOP_VIEW,(void (*)(void))MatView_MFFD)); 6419566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A,MATOP_ASSEMBLY_END,(void (*)(void))MatAssemblyEnd_MFFD)); 6429566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(A,MATOP_SET_FROM_OPTIONS,(void (*)(void))MatSetFromOptions_MFFD)); 643e884886fSBarry Smith 6449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD)); 6459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD)); 6469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD)); 6479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD)); 6489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD)); 6499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD)); 6509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD)); 6519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD)); 6529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetHHistory_C",MatMFFDSetHHistory_MFFD)); 6539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetType_C",MatMFFDSetType_MFFD)); 6549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)A,"MatMFFDGetH_C",MatMFFDGetH_MFFD)); 6559566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)A,MATMFFD)); 656e884886fSBarry Smith PetscFunctionReturn(0); 657e884886fSBarry Smith } 658e884886fSBarry Smith 659e884886fSBarry Smith /*@ 660e884886fSBarry Smith MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 661e884886fSBarry Smith 662e884886fSBarry Smith Collective on Vec 663e884886fSBarry Smith 664e884886fSBarry Smith Input Parameters: 665fef1beadSBarry Smith + comm - MPI communicator 666fef1beadSBarry Smith . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 667fef1beadSBarry Smith This value should be the same as the local size used in creating the 668fef1beadSBarry Smith y vector for the matrix-vector product y = Ax. 669fef1beadSBarry Smith . n - This value should be the same as the local size used in creating the 670fef1beadSBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 671fef1beadSBarry Smith calculated if N is given) For square matrices n is almost always m. 672fef1beadSBarry Smith . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 673fef1beadSBarry Smith - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 674fef1beadSBarry Smith 675e884886fSBarry Smith Output Parameter: 676e884886fSBarry Smith . J - the matrix-free matrix 677e884886fSBarry Smith 6789c6ac3b3SBarry Smith Options Database Keys: call MatSetFromOptions() to trigger these 6799c6ac3b3SBarry Smith + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 680c92e3469SBarry Smith . -mat_mffd_err - square root of estimated relative error in function evaluation 681c92e3469SBarry Smith . -mat_mffd_period - how often h is recomputed, defaults to 1, every time 682c92e3469SBarry Smith . -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values 683c92e3469SBarry Smith - -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing 684c92e3469SBarry Smith (requires real valued functions but that PETSc be configured for complex numbers) 6859c6ac3b3SBarry Smith 686e884886fSBarry Smith Level: advanced 687e884886fSBarry Smith 688e884886fSBarry Smith Notes: 689e884886fSBarry Smith The matrix-free matrix context merely contains the function pointers 690e884886fSBarry Smith and work space for performing finite difference approximations of 691e884886fSBarry Smith Jacobian-vector products, F'(u)*a, 692e884886fSBarry Smith 693e884886fSBarry Smith The default code uses the following approach to compute h 694e884886fSBarry Smith 695e884886fSBarry Smith .vb 696e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 697e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 698e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 699e884886fSBarry Smith where 700e884886fSBarry Smith error_rel = square root of relative error in function evaluation 701e884886fSBarry Smith umin = minimum iterate parameter 702e884886fSBarry Smith .ve 703e884886fSBarry Smith 704ca93e954SBarry Smith You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different 705ca93e954SBarry Smith preconditioner matrix 706ca93e954SBarry Smith 707e884886fSBarry Smith The user can set the error_rel via MatMFFDSetFunctionError() and 708a7f22e61SSatish Balay umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details. 709e884886fSBarry Smith 710e884886fSBarry Smith The user should call MatDestroy() when finished with the matrix-free 711e884886fSBarry Smith matrix context. 712e884886fSBarry Smith 713e884886fSBarry Smith Options Database Keys: 714e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 7151e2ae407SBarry Smith . -mat_mffd_umin <umin> - Sets umin (for default PETSc routine that computes h only) 71667b8a455SSatish Balay - -mat_mffd_check_positivity - check positivity 717e884886fSBarry Smith 718db781477SPatrick Sanan .seealso: `MatDestroy()`, `MatMFFDSetFunctionError()`, `MatMFFDDSSetUmin()`, `MatMFFDSetFunction()` 719db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `MatCreateSNESMF()`, 720db781477SPatrick Sanan `MatMFFDGetH()`, `MatMFFDRegister()`, `MatMFFDComputeJacobian()` 721e884886fSBarry Smith 722e884886fSBarry Smith @*/ 7237087cfbeSBarry Smith PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 724e884886fSBarry Smith { 725e884886fSBarry Smith PetscFunctionBegin; 7269566063dSJacob Faibussowitsch PetscCall(MatCreate(comm,J)); 7279566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J,m,n,M,N)); 7289566063dSJacob Faibussowitsch PetscCall(MatSetType(*J,MATMFFD)); 7299566063dSJacob Faibussowitsch PetscCall(MatSetUp(*J)); 730e884886fSBarry Smith PetscFunctionReturn(0); 731e884886fSBarry Smith } 732e884886fSBarry Smith 733e884886fSBarry Smith /*@ 734e884886fSBarry Smith MatMFFDGetH - Gets the last value that was used as the differencing 735e884886fSBarry Smith parameter. 736e884886fSBarry Smith 737e884886fSBarry Smith Not Collective 738e884886fSBarry Smith 739e884886fSBarry Smith Input Parameters: 740e884886fSBarry Smith . mat - the matrix obtained with MatCreateSNESMF() 741e884886fSBarry Smith 742fd292e60Sprj- Output Parameter: 743e884886fSBarry Smith . h - the differencing step size 744e884886fSBarry Smith 745e884886fSBarry Smith Level: advanced 746e884886fSBarry Smith 747c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDSetHHistory()`, `MatCreateMFFD()`, `MATMFFD`, `MatMFFDResetHHistory()` 748e884886fSBarry Smith @*/ 7497087cfbeSBarry Smith PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 750e884886fSBarry Smith { 751e884886fSBarry Smith PetscFunctionBegin; 75288b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 753dadcf809SJacob Faibussowitsch PetscValidScalarPointer(h,2); 754cac4c232SBarry Smith PetscUseMethod(mat,"MatMFFDGetH_C",(Mat,PetscScalar*),(mat,h)); 755e884886fSBarry Smith PetscFunctionReturn(0); 756e884886fSBarry Smith } 757e884886fSBarry Smith 758e884886fSBarry Smith /*@C 759e884886fSBarry Smith MatMFFDSetFunction - Sets the function used in applying the matrix free. 760e884886fSBarry Smith 7613f9fe445SBarry Smith Logically Collective on Mat 762e884886fSBarry Smith 763e884886fSBarry Smith Input Parameters: 76414f633e2SBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD() 765e884886fSBarry Smith . func - the function to use 766e884886fSBarry Smith - funcctx - optional function context passed to function 767e884886fSBarry Smith 76814f633e2SBarry Smith Calling Sequence of func: 76914f633e2SBarry Smith $ func (void *funcctx, Vec x, Vec f) 77014f633e2SBarry Smith 77114f633e2SBarry Smith + funcctx - user provided context 77214f633e2SBarry Smith . x - input vector 77314f633e2SBarry Smith - f - computed output function 77414f633e2SBarry Smith 775e884886fSBarry Smith Level: advanced 776e884886fSBarry Smith 777e884886fSBarry Smith Notes: 778e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 779e884886fSBarry Smith matrix inside your compute Jacobian routine 780e884886fSBarry Smith 7813ec795f1SBarry Smith If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 782e884886fSBarry Smith 783c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD`, 784db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()` 785e884886fSBarry Smith @*/ 7867087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 787e884886fSBarry Smith { 788e884886fSBarry Smith PetscFunctionBegin; 78988b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 790cac4c232SBarry Smith PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx)); 791e884886fSBarry Smith PetscFunctionReturn(0); 792e884886fSBarry Smith } 793e884886fSBarry Smith 794e884886fSBarry Smith /*@C 795e884886fSBarry Smith MatMFFDSetFunctioni - Sets the function for a single component 796e884886fSBarry Smith 7973f9fe445SBarry Smith Logically Collective on Mat 798e884886fSBarry Smith 799e884886fSBarry Smith Input Parameters: 800e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 801e884886fSBarry Smith - funci - the function to use 802e884886fSBarry Smith 803e884886fSBarry Smith Level: advanced 804e884886fSBarry Smith 805e884886fSBarry Smith Notes: 806e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 80794eb4328SStefano Zampini matrix inside your compute Jacobian routine. 80894eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 809486afcceSStefano Zampini funci must not contain any MPI call as it is called inside a loop on the local portion of the vector. 810e884886fSBarry Smith 811c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()` 8121d0fab5eSBarry Smith 813e884886fSBarry Smith @*/ 8147087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 815e884886fSBarry Smith { 816e884886fSBarry Smith PetscFunctionBegin; 8170700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 818cac4c232SBarry Smith PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci)); 819e884886fSBarry Smith PetscFunctionReturn(0); 820e884886fSBarry Smith } 821e884886fSBarry Smith 822e884886fSBarry Smith /*@C 823e884886fSBarry Smith MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 824e884886fSBarry Smith 8253f9fe445SBarry Smith Logically Collective on Mat 826e884886fSBarry Smith 827e884886fSBarry Smith Input Parameters: 828e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 829e884886fSBarry Smith - func - the function to use 830e884886fSBarry Smith 831e884886fSBarry Smith Level: advanced 832e884886fSBarry Smith 833e884886fSBarry Smith Notes: 834e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 83594eb4328SStefano Zampini matrix inside your compute Jacobian routine. 83694eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 837e884886fSBarry Smith 838c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD` 839db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()`, `SNESetFunction()`, `MatGetDiagonal()` 840e884886fSBarry Smith @*/ 8417087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 842e884886fSBarry Smith { 843e884886fSBarry Smith PetscFunctionBegin; 8440700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 845cac4c232SBarry Smith PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func)); 846e884886fSBarry Smith PetscFunctionReturn(0); 847e884886fSBarry Smith } 848e884886fSBarry Smith 849e884886fSBarry Smith /*@ 850e884886fSBarry Smith MatMFFDSetPeriod - Sets how often h is recomputed, by default it is every time 851e884886fSBarry Smith 8523f9fe445SBarry Smith Logically Collective on Mat 853e884886fSBarry Smith 854e884886fSBarry Smith Input Parameters: 855e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 856e884886fSBarry Smith - period - 1 for every time, 2 for every second etc 857e884886fSBarry Smith 858e884886fSBarry Smith Options Database Keys: 85967b8a455SSatish Balay . -mat_mffd_period <period> - Sets how often h is recomputed 860e884886fSBarry Smith 861e884886fSBarry Smith Level: advanced 862e884886fSBarry Smith 863c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDGetH()`, 864db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()` 865e884886fSBarry Smith @*/ 8667087cfbeSBarry Smith PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 867e884886fSBarry Smith { 868e884886fSBarry Smith PetscFunctionBegin; 86988b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 87088b4c220SStefano Zampini PetscValidLogicalCollectiveInt(mat,period,2); 871cac4c232SBarry Smith PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period)); 872e884886fSBarry Smith PetscFunctionReturn(0); 873e884886fSBarry Smith } 874e884886fSBarry Smith 875e884886fSBarry Smith /*@ 876e884886fSBarry Smith MatMFFDSetFunctionError - Sets the error_rel for the approximation of 877e884886fSBarry Smith matrix-vector products using finite differences. 878e884886fSBarry Smith 8793f9fe445SBarry Smith Logically Collective on Mat 880e884886fSBarry Smith 881e884886fSBarry Smith Input Parameters: 882e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 883e884886fSBarry Smith - error_rel - relative error (should be set to the square root of 884e884886fSBarry Smith the relative error in the function evaluations) 885e884886fSBarry Smith 886e884886fSBarry Smith Options Database Keys: 887a2b725a8SWilliam Gropp . -mat_mffd_err <error_rel> - Sets error_rel 888e884886fSBarry Smith 889e884886fSBarry Smith Level: advanced 890e884886fSBarry Smith 891e884886fSBarry Smith Notes: 892e884886fSBarry Smith The default matrix-free matrix-vector product routine computes 893e884886fSBarry Smith .vb 894e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 895e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 896e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 897e884886fSBarry Smith .ve 898e884886fSBarry Smith 899c2e3fba1SPatrick Sanan .seealso: `MatCreateSNESMF()`, `MatMFFDGetH()`, `MatCreateMFFD()`, `MATMFFD` 900db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDResetHHistory()` 901e884886fSBarry Smith @*/ 9027087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 903e884886fSBarry Smith { 904e884886fSBarry Smith PetscFunctionBegin; 90588b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 90688b4c220SStefano Zampini PetscValidLogicalCollectiveReal(mat,error,2); 907cac4c232SBarry Smith PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error)); 908e884886fSBarry Smith PetscFunctionReturn(0); 909e884886fSBarry Smith } 910e884886fSBarry Smith 911e884886fSBarry Smith /*@ 912e884886fSBarry Smith MatMFFDSetHHistory - Sets an array to collect a history of the 913e884886fSBarry Smith differencing values (h) computed for the matrix-free product. 914e884886fSBarry Smith 9153f9fe445SBarry Smith Logically Collective on Mat 916e884886fSBarry Smith 917e884886fSBarry Smith Input Parameters: 918e884886fSBarry Smith + J - the matrix-free matrix context 919a5b23f4aSJose E. Roman . history - space to hold the history 920e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than 921e884886fSBarry Smith nhistory, then the later ones are discarded 922e884886fSBarry Smith 923e884886fSBarry Smith Level: advanced 924e884886fSBarry Smith 925e884886fSBarry Smith Notes: 926e884886fSBarry Smith Use MatMFFDResetHHistory() to reset the history counter and collect 927e884886fSBarry Smith a new batch of differencing parameters, h. 928e884886fSBarry Smith 929db781477SPatrick Sanan .seealso: `MatMFFDGetH()`, `MatCreateSNESMF()`, 930db781477SPatrick Sanan `MatMFFDResetHHistory()`, `MatMFFDSetFunctionError()` 931e884886fSBarry Smith 932e884886fSBarry Smith @*/ 9337087cfbeSBarry Smith PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 934e884886fSBarry Smith { 935e884886fSBarry Smith PetscFunctionBegin; 93688b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 937dadcf809SJacob Faibussowitsch if (history) PetscValidScalarPointer(history,2); 93888b4c220SStefano Zampini PetscValidLogicalCollectiveInt(J,nhistory,3); 939cac4c232SBarry Smith PetscUseMethod(J,"MatMFFDSetHHistory_C",(Mat,PetscScalar[],PetscInt),(J,history,nhistory)); 940e884886fSBarry Smith PetscFunctionReturn(0); 941e884886fSBarry Smith } 942e884886fSBarry Smith 943e884886fSBarry Smith /*@ 944e884886fSBarry Smith MatMFFDResetHHistory - Resets the counter to zero to begin 945e884886fSBarry Smith collecting a new set of differencing histories. 946e884886fSBarry Smith 9473f9fe445SBarry Smith Logically Collective on Mat 948e884886fSBarry Smith 949e884886fSBarry Smith Input Parameters: 950e884886fSBarry Smith . J - the matrix-free matrix context 951e884886fSBarry Smith 952e884886fSBarry Smith Level: advanced 953e884886fSBarry Smith 954e884886fSBarry Smith Notes: 955e884886fSBarry Smith Use MatMFFDSetHHistory() to create the original history counter. 956e884886fSBarry Smith 957db781477SPatrick Sanan .seealso: `MatMFFDGetH()`, `MatCreateSNESMF()`, 958db781477SPatrick Sanan `MatMFFDSetHHistory()`, `MatMFFDSetFunctionError()` 959e884886fSBarry Smith 960e884886fSBarry Smith @*/ 9617087cfbeSBarry Smith PetscErrorCode MatMFFDResetHHistory(Mat J) 962e884886fSBarry Smith { 963e884886fSBarry Smith PetscFunctionBegin; 96488b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 965cac4c232SBarry Smith PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J)); 966e884886fSBarry Smith PetscFunctionReturn(0); 967e884886fSBarry Smith } 968e884886fSBarry Smith 969487a658cSBarry Smith /*@ 970e884886fSBarry Smith MatMFFDSetBase - Sets the vector U at which matrix vector products of the 971e884886fSBarry Smith Jacobian are computed 972e884886fSBarry Smith 9733f9fe445SBarry Smith Logically Collective on Mat 974e884886fSBarry Smith 975e884886fSBarry Smith Input Parameters: 976e884886fSBarry Smith + J - the MatMFFD matrix 9773ec795f1SBarry Smith . U - the vector 978bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed 979e884886fSBarry Smith 98095452b02SPatrick Sanan Notes: 98195452b02SPatrick Sanan This is rarely used directly 982e884886fSBarry Smith 9838af5ae88SBarry Smith If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 9848af5ae88SBarry Smith point during the first MatMult() after each call to MatMFFDSetBase(). 985dff2f722SBarry Smith 986e884886fSBarry Smith Level: advanced 987e884886fSBarry Smith 988e884886fSBarry Smith @*/ 9897087cfbeSBarry Smith PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 990e884886fSBarry Smith { 991e884886fSBarry Smith PetscFunctionBegin; 9920700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 9930700a824SBarry Smith PetscValidHeaderSpecific(U,VEC_CLASSID,2); 9940700a824SBarry Smith if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 995cac4c232SBarry Smith PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F)); 996e884886fSBarry Smith PetscFunctionReturn(0); 997e884886fSBarry Smith } 998e884886fSBarry Smith 999e884886fSBarry Smith /*@C 1000e884886fSBarry Smith MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 1001e884886fSBarry Smith it to satisfy some criteria 1002e884886fSBarry Smith 10033f9fe445SBarry Smith Logically Collective on Mat 1004e884886fSBarry Smith 1005e884886fSBarry Smith Input Parameters: 1006e884886fSBarry Smith + J - the MatMFFD matrix 1007e884886fSBarry Smith . fun - the function that checks h 1008e884886fSBarry Smith - ctx - any context needed by the function 1009e884886fSBarry Smith 1010e884886fSBarry Smith Options Database Keys: 101167b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is non-negative 1012e884886fSBarry Smith 1013e884886fSBarry Smith Level: advanced 1014e884886fSBarry Smith 101595452b02SPatrick Sanan Notes: 101695452b02SPatrick Sanan For example, MatMFFDCheckPositivity() insures that all entries 1017e884886fSBarry Smith of U + h*a are non-negative 1018e884886fSBarry Smith 1019755b7f64SBarry Smith The function you provide is called after the default h has been computed and allows you to 1020755b7f64SBarry Smith modify it. 1021755b7f64SBarry Smith 1022db781477SPatrick Sanan .seealso: `MatMFFDCheckPositivity()` 1023e884886fSBarry Smith @*/ 10247087cfbeSBarry Smith PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx) 1025e884886fSBarry Smith { 1026e884886fSBarry Smith PetscFunctionBegin; 10270700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1028cac4c232SBarry Smith PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx)); 1029e884886fSBarry Smith PetscFunctionReturn(0); 1030e884886fSBarry Smith } 1031e884886fSBarry Smith 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 10363f9fe445SBarry Smith Logically 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: 104567b8a455SSatish Balay . -mat_mffd_check_positivity <bool> - Insure that U + h*a is nonnegative 1046e884886fSBarry Smith 1047e884886fSBarry Smith Level: advanced 1048e884886fSBarry Smith 104995452b02SPatrick Sanan Notes: 105095452b02SPatrick Sanan This is rarely used directly, rather it is passed as an argument to 1051e884886fSBarry Smith MatMFFDSetCheckh() 1052e884886fSBarry Smith 1053db781477SPatrick Sanan .seealso: `MatMFFDSetCheckh()` 1054e884886fSBarry Smith @*/ 10557087cfbeSBarry Smith PetscErrorCode MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h) 1056e884886fSBarry Smith { 1057e884886fSBarry Smith PetscReal val, minval; 1058e884886fSBarry Smith PetscScalar *u_vec, *a_vec; 1059e884886fSBarry Smith PetscInt i,n; 1060e884886fSBarry Smith MPI_Comm comm; 1061e884886fSBarry Smith 1062e884886fSBarry Smith PetscFunctionBegin; 106388b4c220SStefano Zampini PetscValidHeaderSpecific(U,VEC_CLASSID,2); 106488b4c220SStefano Zampini PetscValidHeaderSpecific(a,VEC_CLASSID,3); 1065dadcf809SJacob Faibussowitsch PetscValidScalarPointer(h,4); 10669566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)U,&comm)); 10679566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&u_vec)); 10689566063dSJacob Faibussowitsch PetscCall(VecGetArray(a,&a_vec)); 10699566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(U,&n)); 1070c068d9bbSLisandro Dalcin minval = PetscAbsScalar(*h)*PetscRealConstant(1.01); 1071e884886fSBarry Smith for (i=0; i<n; i++) { 1072e884886fSBarry Smith if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1073e884886fSBarry Smith val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1074e884886fSBarry Smith if (val < minval) minval = val; 1075e884886fSBarry Smith } 1076e884886fSBarry Smith } 10779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&u_vec)); 10789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(a,&a_vec)); 10791c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm)); 1080e884886fSBarry Smith if (val <= PetscAbsScalar(*h)) { 10819566063dSJacob Faibussowitsch PetscCall(PetscInfo(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val))); 1082e884886fSBarry Smith if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1083e884886fSBarry Smith else *h = -0.99*val; 1084e884886fSBarry Smith } 1085e884886fSBarry Smith PetscFunctionReturn(0); 1086e884886fSBarry Smith } 1087