1e884886fSBarry Smith 2af0996ceSBarry Smith #include <petsc/private/matimpl.h> 3c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h> /*I "petscmat.h" I*/ 4e884886fSBarry Smith 5140e18c1SBarry Smith PetscFunctionList MatMFFDList = 0; 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 #undef __FUNCT__ 13b022a5c1SBarry Smith #define __FUNCT__ "MatMFFDFinalizePackage" 14b022a5c1SBarry Smith /*@C 152390153bSJed Brown MatMFFDFinalizePackage - This function destroys everything in the MatMFFD package. It is 16b022a5c1SBarry Smith called from PetscFinalize(). 17b022a5c1SBarry Smith 18b022a5c1SBarry Smith Level: developer 19b022a5c1SBarry Smith 202390153bSJed Brown .keywords: Petsc, destroy, package 21b022a5c1SBarry Smith .seealso: PetscFinalize() 22b022a5c1SBarry Smith @*/ 237087cfbeSBarry Smith PetscErrorCode MatMFFDFinalizePackage(void) 24b022a5c1SBarry Smith { 25a703d84cSPatrick Lacasse PetscErrorCode ierr; 26a703d84cSPatrick Lacasse 27b022a5c1SBarry Smith PetscFunctionBegin; 2837e93019SBarry Smith ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr); 29b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_FALSE; 30b022a5c1SBarry Smith MatMFFDRegisterAllCalled = PETSC_FALSE; 31b022a5c1SBarry Smith PetscFunctionReturn(0); 32b022a5c1SBarry Smith } 33b022a5c1SBarry Smith 34e884886fSBarry Smith #undef __FUNCT__ 353ec795f1SBarry Smith #define __FUNCT__ "MatMFFDInitializePackage" 363ec795f1SBarry Smith /*@C 373ec795f1SBarry Smith MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called 383ec795f1SBarry Smith from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to MatCreate_MFFD() 393ec795f1SBarry Smith when using static libraries. 403ec795f1SBarry Smith 413ec795f1SBarry Smith Level: developer 423ec795f1SBarry Smith 433ec795f1SBarry Smith .keywords: Vec, initialize, package 443ec795f1SBarry Smith .seealso: PetscInitialize() 453ec795f1SBarry Smith @*/ 46607a6623SBarry Smith PetscErrorCode MatMFFDInitializePackage(void) 473ec795f1SBarry Smith { 483ec795f1SBarry Smith char logList[256]; 493ec795f1SBarry Smith char *className; 50ace3abfcSBarry Smith PetscBool opt; 513ec795f1SBarry Smith PetscErrorCode ierr; 523ec795f1SBarry Smith 533ec795f1SBarry Smith PetscFunctionBegin; 54b022a5c1SBarry Smith if (MatMFFDPackageInitialized) PetscFunctionReturn(0); 55b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_TRUE; 563ec795f1SBarry Smith /* Register Classes */ 570700a824SBarry Smith ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr); 583ec795f1SBarry Smith /* Register Constructors */ 59607a6623SBarry Smith ierr = MatMFFDRegisterAll();CHKERRQ(ierr); 603ec795f1SBarry Smith /* Register Events */ 610700a824SBarry Smith ierr = PetscLogEventRegister("MatMult MF", MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr); 623ec795f1SBarry Smith 633ec795f1SBarry Smith /* Process info exclusions */ 64c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL, "-info_exclude", logList, 256, &opt);CHKERRQ(ierr); 653ec795f1SBarry Smith if (opt) { 663ec795f1SBarry Smith ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr); 673ec795f1SBarry Smith if (className) { 680700a824SBarry Smith ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr); 693ec795f1SBarry Smith } 703ec795f1SBarry Smith } 713ec795f1SBarry Smith /* Process summary exclusions */ 72c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL, "-log_summary_exclude", logList, 256, &opt);CHKERRQ(ierr); 733ec795f1SBarry Smith if (opt) { 743ec795f1SBarry Smith ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr); 753ec795f1SBarry Smith if (className) { 760700a824SBarry Smith ierr = PetscLogEventDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr); 773ec795f1SBarry Smith } 783ec795f1SBarry Smith } 79b022a5c1SBarry Smith ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr); 803ec795f1SBarry Smith PetscFunctionReturn(0); 813ec795f1SBarry Smith } 823ec795f1SBarry Smith 833ec795f1SBarry Smith #undef __FUNCT__ 84e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetType" 85e884886fSBarry Smith /*@C 86e884886fSBarry Smith MatMFFDSetType - Sets the method that is used to compute the 87e884886fSBarry Smith differencing parameter for finite differene matrix-free formulations. 88e884886fSBarry Smith 89e884886fSBarry Smith Input Parameters: 90e884886fSBarry Smith + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD() 91e884886fSBarry Smith or MatSetType(mat,MATMFFD); 92e884886fSBarry Smith - ftype - the type requested, either MATMFFD_WP or MATMFFD_DS 93e884886fSBarry Smith 94e884886fSBarry Smith Level: advanced 95e884886fSBarry Smith 96e884886fSBarry Smith Notes: 97e884886fSBarry Smith For example, such routines can compute h for use in 98e884886fSBarry Smith Jacobian-vector products of the form 99e884886fSBarry Smith 100e884886fSBarry Smith F(x+ha) - F(x) 101e884886fSBarry Smith F'(u)a ~= ---------------- 102e884886fSBarry Smith h 103e884886fSBarry Smith 1041c84c290SBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction() 105e884886fSBarry Smith @*/ 10619fd82e9SBarry Smith PetscErrorCode MatMFFDSetType(Mat mat,MatMFFDType ftype) 107e884886fSBarry Smith { 108e884886fSBarry Smith PetscErrorCode ierr,(*r)(MatMFFD); 109e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 110ace3abfcSBarry Smith PetscBool match; 111e884886fSBarry Smith 112e884886fSBarry Smith PetscFunctionBegin; 1130700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 114e884886fSBarry Smith PetscValidCharPointer(ftype,2); 115e884886fSBarry Smith 116251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 1179a13eae7SJed Brown if (!match) PetscFunctionReturn(0); 1189a13eae7SJed Brown 119e884886fSBarry Smith /* already set, so just return */ 120251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr); 121e884886fSBarry Smith if (match) PetscFunctionReturn(0); 122e884886fSBarry Smith 123e884886fSBarry Smith /* destroy the old one if it exists */ 124e884886fSBarry Smith if (ctx->ops->destroy) { 125e884886fSBarry Smith ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr); 126e884886fSBarry Smith } 127e884886fSBarry Smith 1281c9cd337SJed Brown ierr = PetscFunctionListFind(MatMFFDList,ftype,&r);CHKERRQ(ierr); 129e32f2f54SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype); 130e884886fSBarry Smith ierr = (*r)(ctx);CHKERRQ(ierr); 131e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr); 132e884886fSBarry Smith PetscFunctionReturn(0); 133e884886fSBarry Smith } 134e884886fSBarry Smith 135e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/ 136e884886fSBarry Smith #undef __FUNCT__ 137c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase_MFFD" 1387087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func) 139e884886fSBarry Smith { 140e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 141e884886fSBarry Smith 142e884886fSBarry Smith PetscFunctionBegin; 143e884886fSBarry Smith ctx->funcisetbase = func; 144e884886fSBarry Smith PetscFunctionReturn(0); 145e884886fSBarry Smith } 146e884886fSBarry Smith 147e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 148e884886fSBarry Smith #undef __FUNCT__ 149c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni_MFFD" 1507087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci) 151e884886fSBarry Smith { 152e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 153e884886fSBarry Smith 154e884886fSBarry Smith PetscFunctionBegin; 155e884886fSBarry Smith ctx->funci = funci; 156e884886fSBarry Smith PetscFunctionReturn(0); 157e884886fSBarry Smith } 158e884886fSBarry Smith 159bc0f08ceSBarry Smith #undef __FUNCT__ 160bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory_MFFD" 161bc0f08ceSBarry Smith PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) 162bc0f08ceSBarry Smith { 163bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 164bc0f08ceSBarry Smith 165bc0f08ceSBarry Smith PetscFunctionBegin; 166bc0f08ceSBarry Smith ctx->ncurrenth = 0; 167bc0f08ceSBarry Smith PetscFunctionReturn(0); 168bc0f08ceSBarry Smith } 169e884886fSBarry Smith 170e884886fSBarry Smith #undef __FUNCT__ 171e884886fSBarry Smith #define __FUNCT__ "MatMFFDRegister" 1721c84c290SBarry Smith /*@C 1731c84c290SBarry Smith MatMFFDRegister - Adds a method to the MatMFFD registry. 1741c84c290SBarry Smith 1751c84c290SBarry Smith Not Collective 1761c84c290SBarry Smith 1771c84c290SBarry Smith Input Parameters: 1781c84c290SBarry Smith + name_solver - name of a new user-defined compute-h module 1791c84c290SBarry Smith - routine_create - routine to create method context 1801c84c290SBarry Smith 1811c84c290SBarry Smith Level: developer 1821c84c290SBarry Smith 1831c84c290SBarry Smith Notes: 1841c84c290SBarry Smith MatMFFDRegister() may be called multiple times to add several user-defined solvers. 1851c84c290SBarry Smith 1861c84c290SBarry Smith Sample usage: 1871c84c290SBarry Smith .vb 188bdf89e91SBarry Smith MatMFFDRegister("my_h",MyHCreate); 1891c84c290SBarry Smith .ve 1901c84c290SBarry Smith 1911c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 1921c84c290SBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1931c84c290SBarry Smith or at runtime via the option 194be7a6d03SBarry Smith $ -mat_mffd_type my_h 1951c84c290SBarry Smith 1961c84c290SBarry Smith .keywords: MatMFFD, register 1971c84c290SBarry Smith 1981c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1991c84c290SBarry Smith @*/ 200bdf89e91SBarry Smith PetscErrorCode MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD)) 201e884886fSBarry Smith { 202e884886fSBarry Smith PetscErrorCode ierr; 203e884886fSBarry Smith 204e884886fSBarry Smith PetscFunctionBegin; 205a240a19fSJed Brown ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr); 206e884886fSBarry Smith PetscFunctionReturn(0); 207e884886fSBarry Smith } 208e884886fSBarry Smith 209e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/ 210e884886fSBarry Smith #undef __FUNCT__ 211e884886fSBarry Smith #define __FUNCT__ "MatDestroy_MFFD" 212e884886fSBarry Smith PetscErrorCode MatDestroy_MFFD(Mat mat) 213e884886fSBarry Smith { 214e884886fSBarry Smith PetscErrorCode ierr; 215e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 216e884886fSBarry Smith 217e884886fSBarry Smith PetscFunctionBegin; 2186bf464f9SBarry Smith ierr = VecDestroy(&ctx->w);CHKERRQ(ierr); 2196bf464f9SBarry Smith ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr); 2206bf464f9SBarry Smith ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr); 2216bf464f9SBarry Smith ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr); 222cfe22f5eSBarry Smith if (ctx->current_f_allocated) { 2236bf464f9SBarry Smith ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr); 224cfe22f5eSBarry Smith } 225e884886fSBarry Smith if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);} 2266bf464f9SBarry Smith ierr = PetscHeaderDestroy(&ctx);CHKERRQ(ierr); 227bf0cc555SLisandro Dalcin mat->data = 0; 228e884886fSBarry Smith 229bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr); 230bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr); 231bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr); 232bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr); 233bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr); 234bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr); 235bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr); 236bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr); 237e884886fSBarry Smith PetscFunctionReturn(0); 238e884886fSBarry Smith } 239e884886fSBarry Smith 240e884886fSBarry Smith #undef __FUNCT__ 241e884886fSBarry Smith #define __FUNCT__ "MatView_MFFD" 242e884886fSBarry Smith /* 243e884886fSBarry Smith MatMFFDView_MFFD - Views matrix-free parameters. 244e884886fSBarry Smith 245e884886fSBarry Smith */ 246e884886fSBarry Smith PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer) 247e884886fSBarry Smith { 248e884886fSBarry Smith PetscErrorCode ierr; 249e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 250a283635eSDmitry Karpeev PetscBool iascii, viewbase, viewfunction; 251a283635eSDmitry Karpeev const char *prefix; 252e884886fSBarry Smith 253e884886fSBarry Smith PetscFunctionBegin; 254251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 255e884886fSBarry Smith if (iascii) { 256a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr); 257a283635eSDmitry Karpeev ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 25857622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr); 2597adad957SLisandro Dalcin if (!((PetscObject)ctx)->type_name) { 260e884886fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr); 261e884886fSBarry Smith } else { 2627adad957SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr); 263e884886fSBarry Smith } 264e884886fSBarry Smith if (ctx->ops->view) { 265e884886fSBarry Smith ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr); 266e884886fSBarry Smith } 267a283635eSDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr); 268a283635eSDmitry Karpeev 269c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr); 270a283635eSDmitry Karpeev if (viewbase) { 271a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr); 272a283635eSDmitry Karpeev ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr); 273a283635eSDmitry Karpeev } 274c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr); 275a283635eSDmitry Karpeev if (viewfunction) { 276a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr); 277a283635eSDmitry Karpeev ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr); 278a283635eSDmitry Karpeev } 279a283635eSDmitry Karpeev ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 28011aeaf0aSBarry Smith } 281e884886fSBarry Smith PetscFunctionReturn(0); 282e884886fSBarry Smith } 283e884886fSBarry Smith 284e884886fSBarry Smith #undef __FUNCT__ 285e884886fSBarry Smith #define __FUNCT__ "MatAssemblyEnd_MFFD" 286e884886fSBarry Smith /* 287e884886fSBarry Smith MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 288e884886fSBarry Smith allows the user to indicate the beginning of a new linear solve by calling 289e884886fSBarry Smith MatAssemblyXXX() on the matrix free matrix. This then allows the 2901d0fab5eSBarry Smith MatCreateMFFD_WP() to properly compute ||U|| only the first time 291e884886fSBarry Smith in the linear solver rather than every time. 2925a576424SJed Brown 2935a576424SJed Brown This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library. 294e884886fSBarry Smith */ 2955a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt) 296e884886fSBarry Smith { 297e884886fSBarry Smith PetscErrorCode ierr; 298e884886fSBarry Smith MatMFFD j = (MatMFFD)J->data; 299e884886fSBarry Smith 300e884886fSBarry Smith PetscFunctionBegin; 301e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 302e884886fSBarry Smith j->vshift = 0.0; 303e884886fSBarry Smith j->vscale = 1.0; 304e884886fSBarry Smith PetscFunctionReturn(0); 305e884886fSBarry Smith } 306e884886fSBarry Smith 307e884886fSBarry Smith #undef __FUNCT__ 308e884886fSBarry Smith #define __FUNCT__ "MatMult_MFFD" 309e884886fSBarry Smith /* 310e884886fSBarry Smith MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 311e884886fSBarry Smith 312e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 313e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 314e884886fSBarry Smith u = current iterate 315e884886fSBarry Smith h = difference interval 316e884886fSBarry Smith */ 317e884886fSBarry Smith PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y) 318e884886fSBarry Smith { 319e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 320e884886fSBarry Smith PetscScalar h; 321e884886fSBarry Smith Vec w,U,F; 322e884886fSBarry Smith PetscErrorCode ierr; 323ace3abfcSBarry Smith PetscBool zeroa; 324e884886fSBarry Smith 325e884886fSBarry Smith PetscFunctionBegin; 326ce94432eSBarry Smith if (!ctx->current_u) SETERRQ(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"); 327e884886fSBarry Smith /* We log matrix-free matrix-vector products separately, so that we can 328e884886fSBarry Smith separate the performance monitoring from the cases that use conventional 329e884886fSBarry Smith storage. We may eventually modify event logging to associate events 330e884886fSBarry Smith with particular objects, hence alleviating the more general problem. */ 331e884886fSBarry Smith ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 332e884886fSBarry Smith 333e884886fSBarry Smith w = ctx->w; 334e884886fSBarry Smith U = ctx->current_u; 3353ec795f1SBarry Smith F = ctx->current_f; 336e884886fSBarry Smith /* 337e884886fSBarry Smith Compute differencing parameter 338e884886fSBarry Smith */ 339e884886fSBarry Smith if (!ctx->ops->compute) { 340e884886fSBarry Smith ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr); 3419c6ac3b3SBarry Smith ierr = MatSetFromOptions(mat);CHKERRQ(ierr); 342e884886fSBarry Smith } 343e884886fSBarry Smith ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr); 344e884886fSBarry Smith if (zeroa) { 345e884886fSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 346e884886fSBarry Smith PetscFunctionReturn(0); 347e884886fSBarry Smith } 348e884886fSBarry Smith 349*84d44b13SHong Zhang if (mat->erroriffailure && PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h"); 350e884886fSBarry Smith if (ctx->checkh) { 351e884886fSBarry Smith ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr); 352e884886fSBarry Smith } 353e884886fSBarry Smith 354e884886fSBarry Smith /* keep a record of the current differencing parameter h */ 355e884886fSBarry Smith ctx->currenth = h; 356e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX) 35757622a8eSBarry Smith ierr = PetscInfo2(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h));CHKERRQ(ierr); 358e884886fSBarry Smith #else 359e884886fSBarry Smith ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr); 360e884886fSBarry Smith #endif 361e884886fSBarry Smith if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 362e884886fSBarry Smith ctx->historyh[ctx->ncurrenth] = h; 363e884886fSBarry Smith } 364e884886fSBarry Smith ctx->ncurrenth++; 365e884886fSBarry Smith 366e884886fSBarry Smith /* w = u + ha */ 367c3bb7e23SBarry Smith if (ctx->drscale) { 368c3bb7e23SBarry Smith ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr); 369c3bb7e23SBarry Smith ierr = VecAYPX(U,h,w);CHKERRQ(ierr); 370c3bb7e23SBarry Smith } else { 371e884886fSBarry Smith ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr); 372c3bb7e23SBarry Smith } 373e884886fSBarry Smith 3741b797266SDmitry Karpeev /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 3751b797266SDmitry Karpeev if (ctx->ncurrenth == 1 && ctx->current_f_allocated) { 376e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr); 37752121784SBarry Smith } 378e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr); 379e884886fSBarry Smith 380e884886fSBarry Smith ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 381e884886fSBarry Smith ierr = VecScale(y,1.0/h);CHKERRQ(ierr); 382e884886fSBarry Smith 383c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 384e884886fSBarry Smith ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr); 385c3bb7e23SBarry Smith } 386c3bb7e23SBarry Smith if (ctx->dlscale) { 387c3bb7e23SBarry Smith ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr); 388c3bb7e23SBarry Smith } 389c3bb7e23SBarry Smith if (ctx->dshift) { 390c3bb7e23SBarry Smith ierr = VecPointwiseMult(ctx->dshift,a,U);CHKERRQ(ierr); 391c3bb7e23SBarry Smith ierr = VecAXPY(y,1.0,U);CHKERRQ(ierr); 392c3bb7e23SBarry Smith } 393e884886fSBarry Smith 39439601f4eSBarry Smith if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);} 395e884886fSBarry Smith 396e884886fSBarry Smith ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 397e884886fSBarry Smith PetscFunctionReturn(0); 398e884886fSBarry Smith } 399e884886fSBarry Smith 400e884886fSBarry Smith #undef __FUNCT__ 401e884886fSBarry Smith #define __FUNCT__ "MatGetDiagonal_MFFD" 402e884886fSBarry Smith /* 403e884886fSBarry Smith MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 404e884886fSBarry Smith 405e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 406e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 407e884886fSBarry Smith u = current iterate 408e884886fSBarry Smith h = difference interval 409e884886fSBarry Smith */ 410e884886fSBarry Smith PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 411e884886fSBarry Smith { 412e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 413e884886fSBarry Smith PetscScalar h,*aa,*ww,v; 414e884886fSBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 415e884886fSBarry Smith Vec w,U; 416e884886fSBarry Smith PetscErrorCode ierr; 417e884886fSBarry Smith PetscInt i,rstart,rend; 418e884886fSBarry Smith 419e884886fSBarry Smith PetscFunctionBegin; 420e7e72b3dSBarry Smith if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first"); 421e884886fSBarry Smith 422e884886fSBarry Smith w = ctx->w; 423e884886fSBarry Smith U = ctx->current_u; 424e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr); 425e884886fSBarry Smith ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr); 426e884886fSBarry Smith ierr = VecCopy(U,w);CHKERRQ(ierr); 427e884886fSBarry Smith 428e884886fSBarry Smith ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr); 429e884886fSBarry Smith ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 430e884886fSBarry Smith for (i=rstart; i<rend; i++) { 431e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 432e884886fSBarry Smith h = ww[i-rstart]; 433e884886fSBarry Smith if (h == 0.0) h = 1.0; 434e884886fSBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 435e884886fSBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 436e884886fSBarry Smith h *= epsilon; 437e884886fSBarry Smith 438e884886fSBarry Smith ww[i-rstart] += h; 439e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 440e884886fSBarry Smith ierr = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr); 441e884886fSBarry Smith aa[i-rstart] = (v - aa[i-rstart])/h; 442e884886fSBarry Smith 443e884886fSBarry Smith /* possibly shift and scale result */ 444c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 445e884886fSBarry Smith aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart]; 446c3bb7e23SBarry Smith } 447e884886fSBarry Smith 448e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 449e884886fSBarry Smith ww[i-rstart] -= h; 450e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 451e884886fSBarry Smith } 452e884886fSBarry Smith ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 453e884886fSBarry Smith PetscFunctionReturn(0); 454e884886fSBarry Smith } 455e884886fSBarry Smith 456e884886fSBarry Smith #undef __FUNCT__ 457c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalScale_MFFD" 458c3bb7e23SBarry Smith PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr) 459c3bb7e23SBarry Smith { 460c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 461c3bb7e23SBarry Smith PetscErrorCode ierr; 462c3bb7e23SBarry Smith 463c3bb7e23SBarry Smith PetscFunctionBegin; 464c3bb7e23SBarry Smith if (ll && !aij->dlscale) { 465c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr); 466c3bb7e23SBarry Smith } 467c3bb7e23SBarry Smith if (rr && !aij->drscale) { 468c3bb7e23SBarry Smith ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr); 469c3bb7e23SBarry Smith } 470c3bb7e23SBarry Smith if (ll) { 471c3bb7e23SBarry Smith ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr); 472c3bb7e23SBarry Smith } 473c3bb7e23SBarry Smith if (rr) { 474c3bb7e23SBarry Smith ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr); 475c3bb7e23SBarry Smith } 476c3bb7e23SBarry Smith PetscFunctionReturn(0); 477c3bb7e23SBarry Smith } 478c3bb7e23SBarry Smith 479c3bb7e23SBarry Smith #undef __FUNCT__ 480c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalSet_MFFD" 481c3bb7e23SBarry Smith PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode) 482c3bb7e23SBarry Smith { 483c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 484c3bb7e23SBarry Smith PetscErrorCode ierr; 485c3bb7e23SBarry Smith 486c3bb7e23SBarry Smith PetscFunctionBegin; 487ce94432eSBarry Smith if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES"); 488c3bb7e23SBarry Smith if (!aij->dshift) { 489c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr); 490c3bb7e23SBarry Smith } 491c3bb7e23SBarry Smith ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr); 492c3bb7e23SBarry Smith PetscFunctionReturn(0); 493c3bb7e23SBarry Smith } 494c3bb7e23SBarry Smith 495c3bb7e23SBarry Smith #undef __FUNCT__ 496e884886fSBarry Smith #define __FUNCT__ "MatShift_MFFD" 497e884886fSBarry Smith PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a) 498e884886fSBarry Smith { 499e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 5005fd66863SKarl Rupp 501e884886fSBarry Smith PetscFunctionBegin; 502e884886fSBarry Smith shell->vshift += a; 503e884886fSBarry Smith PetscFunctionReturn(0); 504e884886fSBarry Smith } 505e884886fSBarry Smith 506e884886fSBarry Smith #undef __FUNCT__ 507e884886fSBarry Smith #define __FUNCT__ "MatScale_MFFD" 508e884886fSBarry Smith PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a) 509e884886fSBarry Smith { 510e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 5115fd66863SKarl Rupp 512e884886fSBarry Smith PetscFunctionBegin; 513e884886fSBarry Smith shell->vscale *= a; 514e884886fSBarry Smith PetscFunctionReturn(0); 515e884886fSBarry Smith } 516e884886fSBarry Smith 517e884886fSBarry Smith #undef __FUNCT__ 518c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetBase_MFFD" 519d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F) 520e884886fSBarry Smith { 521e884886fSBarry Smith PetscErrorCode ierr; 522e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 523e884886fSBarry Smith 524e884886fSBarry Smith PetscFunctionBegin; 525e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 5262205254eSKarl Rupp 527e884886fSBarry Smith ctx->current_u = U; 52852121784SBarry Smith if (F) { 5296bf464f9SBarry Smith if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);} 5303ec795f1SBarry Smith ctx->current_f = F; 531cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_FALSE; 532cfe22f5eSBarry Smith } else if (!ctx->current_f_allocated) { 53306c4b6baSBarry Smith ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr); 5342205254eSKarl Rupp 535cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_TRUE; 53652121784SBarry Smith } 537e884886fSBarry Smith if (!ctx->w) { 538e884886fSBarry Smith ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr); 539e884886fSBarry Smith } 540e884886fSBarry Smith J->assembled = PETSC_TRUE; 541e884886fSBarry Smith PetscFunctionReturn(0); 542e884886fSBarry Smith } 5435a576424SJed Brown 544e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 545b2573a8aSBarry Smith 546e884886fSBarry Smith #undef __FUNCT__ 547c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh_MFFD" 5487087cfbeSBarry Smith PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx) 549e884886fSBarry Smith { 550e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 551e884886fSBarry Smith 552e884886fSBarry Smith PetscFunctionBegin; 553e884886fSBarry Smith ctx->checkh = fun; 554e884886fSBarry Smith ctx->checkhctx = ectx; 555e884886fSBarry Smith PetscFunctionReturn(0); 556e884886fSBarry Smith } 557e884886fSBarry Smith 558e884886fSBarry Smith #undef __FUNCT__ 5596aa9148fSLisandro Dalcin #define __FUNCT__ "MatMFFDSetOptionsPrefix" 5606aa9148fSLisandro Dalcin /*@C 5616aa9148fSLisandro Dalcin MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 5626aa9148fSLisandro Dalcin MatMFFD options in the database. 5636aa9148fSLisandro Dalcin 5646aa9148fSLisandro Dalcin Collective on Mat 5656aa9148fSLisandro Dalcin 5666aa9148fSLisandro Dalcin Input Parameter: 5676aa9148fSLisandro Dalcin + A - the Mat context 5686aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names 5696aa9148fSLisandro Dalcin 5706aa9148fSLisandro Dalcin Notes: 5716aa9148fSLisandro Dalcin A hyphen (-) must NOT be given at the beginning of the prefix name. 5726aa9148fSLisandro Dalcin The first character of all runtime options is AUTOMATICALLY the hyphen. 5736aa9148fSLisandro Dalcin 5746aa9148fSLisandro Dalcin Level: advanced 5756aa9148fSLisandro Dalcin 5766aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters 5776aa9148fSLisandro Dalcin 5789c6ac3b3SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF() 5796aa9148fSLisandro Dalcin @*/ 5807087cfbeSBarry Smith PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[]) 5816aa9148fSLisandro Dalcin 5826aa9148fSLisandro Dalcin { 5830298fd71SBarry Smith MatMFFD mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL; 5846aa9148fSLisandro Dalcin PetscErrorCode ierr; 5855fd66863SKarl Rupp 5866aa9148fSLisandro Dalcin PetscFunctionBegin; 5870700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5880700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 5896aa9148fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr); 5906aa9148fSLisandro Dalcin PetscFunctionReturn(0); 5916aa9148fSLisandro Dalcin } 5926aa9148fSLisandro Dalcin 5936aa9148fSLisandro Dalcin #undef __FUNCT__ 5949c6ac3b3SBarry Smith #define __FUNCT__ "MatSetFromOptions_MFFD" 5954416b707SBarry Smith PetscErrorCode MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat) 596e884886fSBarry Smith { 59771f08665SBarry Smith MatMFFD mfctx = (MatMFFD)mat->data; 598e884886fSBarry Smith PetscErrorCode ierr; 599ace3abfcSBarry Smith PetscBool flg; 600e884886fSBarry Smith char ftype[256]; 601e884886fSBarry Smith 602e884886fSBarry Smith PetscFunctionBegin; 6030700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 6040700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 6053194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr); 606a264d7a6SBarry Smith ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr); 607e884886fSBarry Smith if (flg) { 608e884886fSBarry Smith ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr); 609e884886fSBarry Smith } 610e884886fSBarry Smith 611e884886fSBarry Smith ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 612e884886fSBarry Smith ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 613e884886fSBarry Smith 61490d69ab7SBarry Smith flg = PETSC_FALSE; 6150298fd71SBarry Smith ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr); 616e884886fSBarry Smith if (flg) { 617e884886fSBarry Smith ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr); 618e884886fSBarry Smith } 619e884886fSBarry Smith if (mfctx->ops->setfromoptions) { 620e55864a3SBarry Smith ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr); 621e884886fSBarry Smith } 622e884886fSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 623e884886fSBarry Smith PetscFunctionReturn(0); 624e884886fSBarry Smith } 625e884886fSBarry Smith 626bc0f08ceSBarry Smith #undef __FUNCT__ 627bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod_MFFD" 628bc0f08ceSBarry Smith PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period) 629bc0f08ceSBarry Smith { 630bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 631bc0f08ceSBarry Smith 632bc0f08ceSBarry Smith PetscFunctionBegin; 633bc0f08ceSBarry Smith PetscValidLogicalCollectiveInt(mat,period,2); 634bc0f08ceSBarry Smith ctx->recomputeperiod = period; 635bc0f08ceSBarry Smith PetscFunctionReturn(0); 636bc0f08ceSBarry Smith } 637bc0f08ceSBarry Smith 638bc0f08ceSBarry Smith #undef __FUNCT__ 639bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunction_MFFD" 640bc0f08ceSBarry Smith PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 641bc0f08ceSBarry Smith { 642bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 643bc0f08ceSBarry Smith 644bc0f08ceSBarry Smith PetscFunctionBegin; 645bc0f08ceSBarry Smith ctx->func = func; 646bc0f08ceSBarry Smith ctx->funcctx = funcctx; 647bc0f08ceSBarry Smith PetscFunctionReturn(0); 648bc0f08ceSBarry Smith } 649bc0f08ceSBarry Smith 650bc0f08ceSBarry Smith #undef __FUNCT__ 651bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError_MFFD" 652bc0f08ceSBarry Smith PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error) 653bc0f08ceSBarry Smith { 654bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 655bc0f08ceSBarry Smith 656bc0f08ceSBarry Smith PetscFunctionBegin; 657bc0f08ceSBarry Smith PetscValidLogicalCollectiveReal(mat,error,2); 658bc0f08ceSBarry Smith if (error != PETSC_DEFAULT) ctx->error_rel = error; 659bc0f08ceSBarry Smith PetscFunctionReturn(0); 660bc0f08ceSBarry Smith } 661bc0f08ceSBarry Smith 6623b49f96aSBarry Smith #undef __FUNCT__ 6633b49f96aSBarry Smith #define __FUNCT__ "MatMissingDiagonal_MFFD" 6643b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MFFD(Mat A,PetscBool *missing,PetscInt *d) 6653b49f96aSBarry Smith { 6663b49f96aSBarry Smith PetscFunctionBegin; 6673b49f96aSBarry Smith *missing = PETSC_FALSE; 6683b49f96aSBarry Smith PetscFunctionReturn(0); 6693b49f96aSBarry Smith } 6703b49f96aSBarry Smith 671e884886fSBarry Smith /*MC 672e884886fSBarry Smith MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 673e884886fSBarry Smith 674e884886fSBarry Smith Level: advanced 675e884886fSBarry Smith 676d2d6cebeSBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction() 677e884886fSBarry Smith M*/ 678e884886fSBarry Smith #undef __FUNCT__ 679e884886fSBarry Smith #define __FUNCT__ "MatCreate_MFFD" 6808cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 681e884886fSBarry Smith { 682e884886fSBarry Smith MatMFFD mfctx; 683e884886fSBarry Smith PetscErrorCode ierr; 684e884886fSBarry Smith 685e884886fSBarry Smith PetscFunctionBegin; 686607a6623SBarry Smith ierr = MatMFFDInitializePackage();CHKERRQ(ierr); 687e884886fSBarry Smith 68873107ff1SLisandro Dalcin ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr); 6892205254eSKarl Rupp 690e884886fSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 691e884886fSBarry Smith mfctx->recomputeperiod = 1; 692e884886fSBarry Smith mfctx->count = 0; 693e884886fSBarry Smith mfctx->currenth = 0.0; 6940298fd71SBarry Smith mfctx->historyh = NULL; 695e884886fSBarry Smith mfctx->ncurrenth = 0; 696e884886fSBarry Smith mfctx->maxcurrenth = 0; 6977adad957SLisandro Dalcin ((PetscObject)mfctx)->type_name = 0; 698e884886fSBarry Smith 699e884886fSBarry Smith mfctx->vshift = 0.0; 700e884886fSBarry Smith mfctx->vscale = 1.0; 701e884886fSBarry Smith 702e884886fSBarry Smith /* 703e884886fSBarry Smith Create the empty data structure to contain compute-h routines. 704e884886fSBarry Smith These will be filled in below from the command line options or 705e884886fSBarry Smith a later call with MatMFFDSetType() or if that is not called 706e884886fSBarry Smith then it will default in the first use of MatMult_MFFD() 707e884886fSBarry Smith */ 708e884886fSBarry Smith mfctx->ops->compute = 0; 709e884886fSBarry Smith mfctx->ops->destroy = 0; 710e884886fSBarry Smith mfctx->ops->view = 0; 711e884886fSBarry Smith mfctx->ops->setfromoptions = 0; 712e884886fSBarry Smith mfctx->hctx = 0; 713e884886fSBarry Smith 714e884886fSBarry Smith mfctx->func = 0; 715e884886fSBarry Smith mfctx->funcctx = 0; 7160298fd71SBarry Smith mfctx->w = NULL; 717e884886fSBarry Smith 718e884886fSBarry Smith A->data = mfctx; 719e884886fSBarry Smith 720e884886fSBarry Smith A->ops->mult = MatMult_MFFD; 721e884886fSBarry Smith A->ops->destroy = MatDestroy_MFFD; 722e884886fSBarry Smith A->ops->view = MatView_MFFD; 723e884886fSBarry Smith A->ops->assemblyend = MatAssemblyEnd_MFFD; 724e884886fSBarry Smith A->ops->getdiagonal = MatGetDiagonal_MFFD; 725e884886fSBarry Smith A->ops->scale = MatScale_MFFD; 726e884886fSBarry Smith A->ops->shift = MatShift_MFFD; 727c3bb7e23SBarry Smith A->ops->diagonalscale = MatDiagonalScale_MFFD; 728c3bb7e23SBarry Smith A->ops->diagonalset = MatDiagonalSet_MFFD; 7299c6ac3b3SBarry Smith A->ops->setfromoptions = MatSetFromOptions_MFFD; 7303b49f96aSBarry Smith A->ops->missingdiagonal = MatMissingDiagonal_MFFD; 731e884886fSBarry Smith A->assembled = PETSC_TRUE; 732e884886fSBarry Smith 73326283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 73426283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 735ee1cef2cSJed Brown 736bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr); 737bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr); 738bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr); 739bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr); 740bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr); 741bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr); 742bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr); 743bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr); 7442205254eSKarl Rupp 745e884886fSBarry Smith mfctx->mat = A; 7462205254eSKarl Rupp 747e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr); 748e884886fSBarry Smith PetscFunctionReturn(0); 749e884886fSBarry Smith } 750e884886fSBarry Smith 751e884886fSBarry Smith #undef __FUNCT__ 752e884886fSBarry Smith #define __FUNCT__ "MatCreateMFFD" 753e884886fSBarry Smith /*@ 754e884886fSBarry Smith MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 755e884886fSBarry Smith 756e884886fSBarry Smith Collective on Vec 757e884886fSBarry Smith 758e884886fSBarry Smith Input Parameters: 759fef1beadSBarry Smith + comm - MPI communicator 760fef1beadSBarry Smith . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 761fef1beadSBarry Smith This value should be the same as the local size used in creating the 762fef1beadSBarry Smith y vector for the matrix-vector product y = Ax. 763fef1beadSBarry Smith . n - This value should be the same as the local size used in creating the 764fef1beadSBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 765fef1beadSBarry Smith calculated if N is given) For square matrices n is almost always m. 766fef1beadSBarry Smith . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 767fef1beadSBarry Smith - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 768fef1beadSBarry Smith 769e884886fSBarry Smith 770e884886fSBarry Smith Output Parameter: 771e884886fSBarry Smith . J - the matrix-free matrix 772e884886fSBarry Smith 7739c6ac3b3SBarry Smith Options Database Keys: call MatSetFromOptions() to trigger these 7749c6ac3b3SBarry Smith + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 7759c6ac3b3SBarry Smith - -mat_mffd_err - square root of estimated relative error in function evaluation 7769c6ac3b3SBarry Smith - -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 7779c6ac3b3SBarry Smith 7789c6ac3b3SBarry Smith 779e884886fSBarry Smith Level: advanced 780e884886fSBarry Smith 781e884886fSBarry Smith Notes: 782e884886fSBarry Smith The matrix-free matrix context merely contains the function pointers 783e884886fSBarry Smith and work space for performing finite difference approximations of 784e884886fSBarry Smith Jacobian-vector products, F'(u)*a, 785e884886fSBarry Smith 786e884886fSBarry Smith The default code uses the following approach to compute h 787e884886fSBarry Smith 788e884886fSBarry Smith .vb 789e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 790e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 791e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 792e884886fSBarry Smith where 793e884886fSBarry Smith error_rel = square root of relative error in function evaluation 794e884886fSBarry Smith umin = minimum iterate parameter 795e884886fSBarry Smith .ve 796e884886fSBarry Smith 797ca93e954SBarry Smith You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different 798ca93e954SBarry Smith preconditioner matrix 799ca93e954SBarry Smith 800e884886fSBarry Smith The user can set the error_rel via MatMFFDSetFunctionError() and 801a7f22e61SSatish Balay umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details. 802e884886fSBarry Smith 803e884886fSBarry Smith The user should call MatDestroy() when finished with the matrix-free 804e884886fSBarry Smith matrix context. 805e884886fSBarry Smith 806e884886fSBarry Smith Options Database Keys: 807e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 808e884886fSBarry Smith . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only) 809e884886fSBarry Smith - -mat_mffd_check_positivity 810e884886fSBarry Smith 811e884886fSBarry Smith .keywords: default, matrix-free, create, matrix 812e884886fSBarry Smith 8131d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 814e884886fSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 81581242352SJed Brown MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian() 816e884886fSBarry Smith 817e884886fSBarry Smith @*/ 8187087cfbeSBarry Smith PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 819e884886fSBarry Smith { 820e884886fSBarry Smith PetscErrorCode ierr; 821e884886fSBarry Smith 822e884886fSBarry Smith PetscFunctionBegin; 823e884886fSBarry Smith ierr = MatCreate(comm,J);CHKERRQ(ierr); 824fef1beadSBarry Smith ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr); 825e884886fSBarry Smith ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 826be7c243fSJed Brown ierr = MatSetUp(*J);CHKERRQ(ierr); 827e884886fSBarry Smith PetscFunctionReturn(0); 828e884886fSBarry Smith } 829e884886fSBarry Smith 830e884886fSBarry Smith 831e884886fSBarry Smith #undef __FUNCT__ 832e884886fSBarry Smith #define __FUNCT__ "MatMFFDGetH" 833e884886fSBarry Smith /*@ 834e884886fSBarry Smith MatMFFDGetH - Gets the last value that was used as the differencing 835e884886fSBarry Smith parameter. 836e884886fSBarry Smith 837e884886fSBarry Smith Not Collective 838e884886fSBarry Smith 839e884886fSBarry Smith Input Parameters: 840e884886fSBarry Smith . mat - the matrix obtained with MatCreateSNESMF() 841e884886fSBarry Smith 842e884886fSBarry Smith Output Paramter: 843e884886fSBarry Smith . h - the differencing step size 844e884886fSBarry Smith 845e884886fSBarry Smith Level: advanced 846e884886fSBarry Smith 847e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 848e884886fSBarry Smith 8491d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory() 850e884886fSBarry Smith @*/ 8517087cfbeSBarry Smith PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 852e884886fSBarry Smith { 853e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 854bc0f08ceSBarry Smith PetscErrorCode ierr; 855bc0f08ceSBarry Smith PetscBool match; 856e884886fSBarry Smith 857e884886fSBarry Smith PetscFunctionBegin; 858251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 859ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 860bc0f08ceSBarry Smith 861e884886fSBarry Smith *h = ctx->currenth; 862e884886fSBarry Smith PetscFunctionReturn(0); 863e884886fSBarry Smith } 864e884886fSBarry Smith 865e884886fSBarry Smith #undef __FUNCT__ 866e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunction" 867e884886fSBarry Smith /*@C 868e884886fSBarry Smith MatMFFDSetFunction - Sets the function used in applying the matrix free. 869e884886fSBarry Smith 8703f9fe445SBarry Smith Logically Collective on Mat 871e884886fSBarry Smith 872e884886fSBarry Smith Input Parameters: 87314f633e2SBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD() 874e884886fSBarry Smith . func - the function to use 875e884886fSBarry Smith - funcctx - optional function context passed to function 876e884886fSBarry Smith 87714f633e2SBarry Smith Calling Sequence of func: 87814f633e2SBarry Smith $ func (void *funcctx, Vec x, Vec f) 87914f633e2SBarry Smith 88014f633e2SBarry Smith + funcctx - user provided context 88114f633e2SBarry Smith . x - input vector 88214f633e2SBarry Smith - f - computed output function 88314f633e2SBarry Smith 884e884886fSBarry Smith Level: advanced 885e884886fSBarry Smith 886e884886fSBarry Smith Notes: 887e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 888e884886fSBarry Smith matrix inside your compute Jacobian routine 889e884886fSBarry Smith 8903ec795f1SBarry Smith If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 891e884886fSBarry Smith 892e884886fSBarry Smith .keywords: SNES, matrix-free, function 893e884886fSBarry Smith 8941d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD, 8951d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 896e884886fSBarry Smith @*/ 8977087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 898e884886fSBarry Smith { 899bc0f08ceSBarry Smith PetscErrorCode ierr; 900e884886fSBarry Smith 901e884886fSBarry Smith PetscFunctionBegin; 902bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr); 903e884886fSBarry Smith PetscFunctionReturn(0); 904e884886fSBarry Smith } 905e884886fSBarry Smith 906e884886fSBarry Smith #undef __FUNCT__ 907e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni" 908e884886fSBarry Smith /*@C 909e884886fSBarry Smith MatMFFDSetFunctioni - Sets the function for a single component 910e884886fSBarry Smith 9113f9fe445SBarry Smith Logically Collective on Mat 912e884886fSBarry Smith 913e884886fSBarry Smith Input Parameters: 914e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 915e884886fSBarry Smith - funci - the function to use 916e884886fSBarry Smith 917e884886fSBarry Smith Level: advanced 918e884886fSBarry Smith 919e884886fSBarry Smith Notes: 920e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 921e884886fSBarry Smith matrix inside your compute Jacobian routine 922e884886fSBarry Smith 923e884886fSBarry Smith 924e884886fSBarry Smith .keywords: SNES, matrix-free, function 925e884886fSBarry Smith 9261d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 9271d0fab5eSBarry Smith 928e884886fSBarry Smith @*/ 9297087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 930e884886fSBarry Smith { 9314ac538c5SBarry Smith PetscErrorCode ierr; 932e884886fSBarry Smith 933e884886fSBarry Smith PetscFunctionBegin; 9340700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 9354ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr); 936e884886fSBarry Smith PetscFunctionReturn(0); 937e884886fSBarry Smith } 938e884886fSBarry Smith 939e884886fSBarry Smith 940e884886fSBarry Smith #undef __FUNCT__ 941e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase" 942e884886fSBarry Smith /*@C 943e884886fSBarry Smith MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 944e884886fSBarry Smith 9453f9fe445SBarry Smith Logically Collective on Mat 946e884886fSBarry Smith 947e884886fSBarry Smith Input Parameters: 948e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 949e884886fSBarry Smith - func - the function to use 950e884886fSBarry Smith 951e884886fSBarry Smith Level: advanced 952e884886fSBarry Smith 953e884886fSBarry Smith Notes: 954e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 955e884886fSBarry Smith matrix inside your compute Jacobian routine 956e884886fSBarry Smith 957e884886fSBarry Smith 958e884886fSBarry Smith .keywords: SNES, matrix-free, function 959e884886fSBarry Smith 960e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 9611d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 962e884886fSBarry Smith @*/ 9637087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 964e884886fSBarry Smith { 9654ac538c5SBarry Smith PetscErrorCode ierr; 966e884886fSBarry Smith 967e884886fSBarry Smith PetscFunctionBegin; 9680700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 9694ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr); 970e884886fSBarry Smith PetscFunctionReturn(0); 971e884886fSBarry Smith } 972e884886fSBarry Smith 973e884886fSBarry Smith #undef __FUNCT__ 974e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod" 975e884886fSBarry Smith /*@ 976e884886fSBarry Smith MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 977e884886fSBarry Smith 9783f9fe445SBarry Smith Logically Collective on Mat 979e884886fSBarry Smith 980e884886fSBarry Smith Input Parameters: 981e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 982e884886fSBarry Smith - period - 1 for everytime, 2 for every second etc 983e884886fSBarry Smith 984e884886fSBarry Smith Options Database Keys: 985e884886fSBarry Smith + -mat_mffd_period <period> 986e884886fSBarry Smith 987e884886fSBarry Smith Level: advanced 988e884886fSBarry Smith 989e884886fSBarry Smith 990e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 991e884886fSBarry Smith 992e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), 9931d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 994e884886fSBarry Smith @*/ 9957087cfbeSBarry Smith PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 996e884886fSBarry Smith { 997bc0f08ceSBarry Smith PetscErrorCode ierr; 998e884886fSBarry Smith 999e884886fSBarry Smith PetscFunctionBegin; 1000bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr); 1001e884886fSBarry Smith PetscFunctionReturn(0); 1002e884886fSBarry Smith } 1003e884886fSBarry Smith 1004e884886fSBarry Smith #undef __FUNCT__ 1005e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError" 1006e884886fSBarry Smith /*@ 1007e884886fSBarry Smith MatMFFDSetFunctionError - Sets the error_rel for the approximation of 1008e884886fSBarry Smith matrix-vector products using finite differences. 1009e884886fSBarry Smith 10103f9fe445SBarry Smith Logically Collective on Mat 1011e884886fSBarry Smith 1012e884886fSBarry Smith Input Parameters: 1013e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 1014e884886fSBarry Smith - error_rel - relative error (should be set to the square root of 1015e884886fSBarry Smith the relative error in the function evaluations) 1016e884886fSBarry Smith 1017e884886fSBarry Smith Options Database Keys: 1018e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 1019e884886fSBarry Smith 1020e884886fSBarry Smith Level: advanced 1021e884886fSBarry Smith 1022e884886fSBarry Smith Notes: 1023e884886fSBarry Smith The default matrix-free matrix-vector product routine computes 1024e884886fSBarry Smith .vb 1025e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 1026e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 1027e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 1028e884886fSBarry Smith .ve 1029e884886fSBarry Smith 1030e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 1031e884886fSBarry Smith 1032e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 10331d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 1034e884886fSBarry Smith @*/ 10357087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 1036e884886fSBarry Smith { 1037bc0f08ceSBarry Smith PetscErrorCode ierr; 1038e884886fSBarry Smith 1039e884886fSBarry Smith PetscFunctionBegin; 1040bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr); 1041e884886fSBarry Smith PetscFunctionReturn(0); 1042e884886fSBarry Smith } 1043e884886fSBarry Smith 1044e884886fSBarry Smith #undef __FUNCT__ 1045e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetHHistory" 1046e884886fSBarry Smith /*@ 1047e884886fSBarry Smith MatMFFDSetHHistory - Sets an array to collect a history of the 1048e884886fSBarry Smith differencing values (h) computed for the matrix-free product. 1049e884886fSBarry Smith 10503f9fe445SBarry Smith Logically Collective on Mat 1051e884886fSBarry Smith 1052e884886fSBarry Smith Input Parameters: 1053e884886fSBarry Smith + J - the matrix-free matrix context 1054e884886fSBarry Smith . histroy - space to hold the history 1055e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than 1056e884886fSBarry Smith nhistory, then the later ones are discarded 1057e884886fSBarry Smith 1058e884886fSBarry Smith Level: advanced 1059e884886fSBarry Smith 1060e884886fSBarry Smith Notes: 1061e884886fSBarry Smith Use MatMFFDResetHHistory() to reset the history counter and collect 1062e884886fSBarry Smith a new batch of differencing parameters, h. 1063e884886fSBarry Smith 1064e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1065e884886fSBarry Smith 1066e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 10671d0fab5eSBarry Smith MatMFFDResetHHistory(), MatMFFDSetFunctionError() 1068e884886fSBarry Smith 1069e884886fSBarry Smith @*/ 10707087cfbeSBarry Smith PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 1071e884886fSBarry Smith { 1072e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 1073bc0f08ceSBarry Smith PetscErrorCode ierr; 1074bc0f08ceSBarry Smith PetscBool match; 1075e884886fSBarry Smith 1076e884886fSBarry Smith PetscFunctionBegin; 1077251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr); 1078ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 1079e884886fSBarry Smith ctx->historyh = history; 1080e884886fSBarry Smith ctx->maxcurrenth = nhistory; 108175567043SBarry Smith ctx->currenth = 0.; 1082e884886fSBarry Smith PetscFunctionReturn(0); 1083e884886fSBarry Smith } 1084e884886fSBarry Smith 1085bc0f08ceSBarry Smith 1086e884886fSBarry Smith #undef __FUNCT__ 1087e884886fSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory" 1088e884886fSBarry Smith /*@ 1089e884886fSBarry Smith MatMFFDResetHHistory - Resets the counter to zero to begin 1090e884886fSBarry Smith collecting a new set of differencing histories. 1091e884886fSBarry Smith 10923f9fe445SBarry Smith Logically Collective on Mat 1093e884886fSBarry Smith 1094e884886fSBarry Smith Input Parameters: 1095e884886fSBarry Smith . J - the matrix-free matrix context 1096e884886fSBarry Smith 1097e884886fSBarry Smith Level: advanced 1098e884886fSBarry Smith 1099e884886fSBarry Smith Notes: 1100e884886fSBarry Smith Use MatMFFDSetHHistory() to create the original history counter. 1101e884886fSBarry Smith 1102e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1103e884886fSBarry Smith 1104e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 11051d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDSetFunctionError() 1106e884886fSBarry Smith 1107e884886fSBarry Smith @*/ 11087087cfbeSBarry Smith PetscErrorCode MatMFFDResetHHistory(Mat J) 1109e884886fSBarry Smith { 1110bc0f08ceSBarry Smith PetscErrorCode ierr; 1111e884886fSBarry Smith 1112e884886fSBarry Smith PetscFunctionBegin; 1113bc0f08ceSBarry Smith ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr); 1114e884886fSBarry Smith PetscFunctionReturn(0); 1115e884886fSBarry Smith } 1116e884886fSBarry Smith 1117e884886fSBarry Smith 1118e884886fSBarry Smith #undef __FUNCT__ 1119e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetBase" 1120e884886fSBarry Smith /*@ 1121e884886fSBarry Smith MatMFFDSetBase - Sets the vector U at which matrix vector products of the 1122e884886fSBarry Smith Jacobian are computed 1123e884886fSBarry Smith 11243f9fe445SBarry Smith Logically Collective on Mat 1125e884886fSBarry Smith 1126e884886fSBarry Smith Input Parameters: 1127e884886fSBarry Smith + J - the MatMFFD matrix 11283ec795f1SBarry Smith . U - the vector 1129bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed 1130e884886fSBarry Smith 1131e884886fSBarry Smith Notes: This is rarely used directly 1132e884886fSBarry Smith 11338af5ae88SBarry Smith If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 11348af5ae88SBarry Smith point during the first MatMult() after each call to MatMFFDSetBase(). 1135dff2f722SBarry Smith 1136e884886fSBarry Smith Level: advanced 1137e884886fSBarry Smith 1138e884886fSBarry Smith @*/ 11397087cfbeSBarry Smith PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 1140e884886fSBarry Smith { 11414ac538c5SBarry Smith PetscErrorCode ierr; 1142e884886fSBarry Smith 1143e884886fSBarry Smith PetscFunctionBegin; 11440700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 11450700a824SBarry Smith PetscValidHeaderSpecific(U,VEC_CLASSID,2); 11460700a824SBarry Smith if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 11474ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr); 1148e884886fSBarry Smith PetscFunctionReturn(0); 1149e884886fSBarry Smith } 1150e884886fSBarry Smith 1151e884886fSBarry Smith #undef __FUNCT__ 1152e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh" 1153e884886fSBarry Smith /*@C 1154e884886fSBarry Smith MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 1155e884886fSBarry Smith it to satisfy some criteria 1156e884886fSBarry Smith 11573f9fe445SBarry Smith Logically Collective on Mat 1158e884886fSBarry Smith 1159e884886fSBarry Smith Input Parameters: 1160e884886fSBarry Smith + J - the MatMFFD matrix 1161e884886fSBarry Smith . fun - the function that checks h 1162e884886fSBarry Smith - ctx - any context needed by the function 1163e884886fSBarry Smith 1164e884886fSBarry Smith Options Database Keys: 1165e884886fSBarry Smith . -mat_mffd_check_positivity 1166e884886fSBarry Smith 1167e884886fSBarry Smith Level: advanced 1168e884886fSBarry Smith 1169e884886fSBarry Smith Notes: For example, MatMFFDSetCheckPositivity() insures that all entries 1170e884886fSBarry Smith of U + h*a are non-negative 1171e884886fSBarry Smith 1172e884886fSBarry Smith .seealso: MatMFFDSetCheckPositivity() 1173e884886fSBarry Smith @*/ 11747087cfbeSBarry Smith PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx) 1175e884886fSBarry Smith { 11764ac538c5SBarry Smith PetscErrorCode ierr; 1177e884886fSBarry Smith 1178e884886fSBarry Smith PetscFunctionBegin; 11790700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 11804ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr); 1181e884886fSBarry Smith PetscFunctionReturn(0); 1182e884886fSBarry Smith } 1183e884886fSBarry Smith 1184e884886fSBarry Smith #undef __FUNCT__ 1185e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckPositivity" 1186e884886fSBarry Smith /*@ 1187e884886fSBarry Smith MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1188e884886fSBarry Smith zero, decreases h until this is satisfied. 1189e884886fSBarry Smith 11903f9fe445SBarry Smith Logically Collective on Vec 1191e884886fSBarry Smith 1192e884886fSBarry Smith Input Parameters: 1193e884886fSBarry Smith + U - base vector that is added to 1194e884886fSBarry Smith . a - vector that is added 1195e884886fSBarry Smith . h - scaling factor on a 1196e884886fSBarry Smith - dummy - context variable (unused) 1197e884886fSBarry Smith 1198e884886fSBarry Smith Options Database Keys: 1199e884886fSBarry Smith . -mat_mffd_check_positivity 1200e884886fSBarry Smith 1201e884886fSBarry Smith Level: advanced 1202e884886fSBarry Smith 1203e884886fSBarry Smith Notes: This is rarely used directly, rather it is passed as an argument to 1204e884886fSBarry Smith MatMFFDSetCheckh() 1205e884886fSBarry Smith 1206e884886fSBarry Smith .seealso: MatMFFDSetCheckh() 1207e884886fSBarry Smith @*/ 12087087cfbeSBarry Smith PetscErrorCode MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h) 1209e884886fSBarry Smith { 1210e884886fSBarry Smith PetscReal val, minval; 1211e884886fSBarry Smith PetscScalar *u_vec, *a_vec; 1212e884886fSBarry Smith PetscErrorCode ierr; 1213e884886fSBarry Smith PetscInt i,n; 1214e884886fSBarry Smith MPI_Comm comm; 1215e884886fSBarry Smith 1216e884886fSBarry Smith PetscFunctionBegin; 1217e884886fSBarry Smith ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); 1218e884886fSBarry Smith ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); 1219e884886fSBarry Smith ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); 1220e884886fSBarry Smith ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); 1221e884886fSBarry Smith minval = PetscAbsScalar(*h*1.01); 1222e884886fSBarry Smith for (i=0; i<n; i++) { 1223e884886fSBarry Smith if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1224e884886fSBarry Smith val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1225e884886fSBarry Smith if (val < minval) minval = val; 1226e884886fSBarry Smith } 1227e884886fSBarry Smith } 1228e884886fSBarry Smith ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); 1229e884886fSBarry Smith ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); 1230b2566f29SBarry Smith ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr); 1231e884886fSBarry Smith if (val <= PetscAbsScalar(*h)) { 12326712e2f1SBarry Smith ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr); 1233e884886fSBarry Smith if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1234e884886fSBarry Smith else *h = -0.99*val; 1235e884886fSBarry Smith } 1236e884886fSBarry Smith PetscFunctionReturn(0); 1237e884886fSBarry Smith } 1238e884886fSBarry Smith 1239e884886fSBarry Smith 1240e884886fSBarry Smith 1241e884886fSBarry Smith 1242e884886fSBarry Smith 1243e884886fSBarry Smith 1244e884886fSBarry Smith 1245