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 /*@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 182390153bSJed Brown .keywords: Petsc, destroy, package 19755b7f64SBarry Smith .seealso: PetscFinalize(), MatCreateMFFD(), MatCreateSNESMF() 20b022a5c1SBarry Smith @*/ 217087cfbeSBarry Smith PetscErrorCode MatMFFDFinalizePackage(void) 22b022a5c1SBarry Smith { 23a703d84cSPatrick Lacasse PetscErrorCode ierr; 24a703d84cSPatrick Lacasse 25b022a5c1SBarry Smith PetscFunctionBegin; 2637e93019SBarry Smith ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr); 27b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_FALSE; 28b022a5c1SBarry Smith MatMFFDRegisterAllCalled = PETSC_FALSE; 29b022a5c1SBarry Smith PetscFunctionReturn(0); 30b022a5c1SBarry Smith } 31b022a5c1SBarry Smith 323ec795f1SBarry Smith /*@C 333ec795f1SBarry Smith MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called 343ec795f1SBarry Smith from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to MatCreate_MFFD() 353ec795f1SBarry Smith when using static libraries. 363ec795f1SBarry Smith 373ec795f1SBarry Smith Level: developer 383ec795f1SBarry Smith 393ec795f1SBarry Smith .keywords: Vec, initialize, package 403ec795f1SBarry Smith .seealso: PetscInitialize() 413ec795f1SBarry Smith @*/ 42607a6623SBarry Smith PetscErrorCode MatMFFDInitializePackage(void) 433ec795f1SBarry Smith { 443ec795f1SBarry Smith char logList[256]; 453ec795f1SBarry Smith char *className; 46ace3abfcSBarry Smith PetscBool opt; 473ec795f1SBarry Smith PetscErrorCode ierr; 483ec795f1SBarry Smith 493ec795f1SBarry Smith PetscFunctionBegin; 50b022a5c1SBarry Smith if (MatMFFDPackageInitialized) PetscFunctionReturn(0); 51b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_TRUE; 523ec795f1SBarry Smith /* Register Classes */ 530700a824SBarry Smith ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr); 543ec795f1SBarry Smith /* Register Constructors */ 55607a6623SBarry Smith ierr = MatMFFDRegisterAll();CHKERRQ(ierr); 563ec795f1SBarry Smith /* Register Events */ 570700a824SBarry Smith ierr = PetscLogEventRegister("MatMult MF", MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr); 583ec795f1SBarry Smith 593ec795f1SBarry Smith /* Process info exclusions */ 60c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,NULL, "-info_exclude", logList, 256, &opt);CHKERRQ(ierr); 613ec795f1SBarry Smith if (opt) { 623ec795f1SBarry Smith ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr); 633ec795f1SBarry Smith if (className) { 640700a824SBarry Smith ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr); 653ec795f1SBarry Smith } 663ec795f1SBarry Smith } 673ec795f1SBarry Smith /* Process summary exclusions */ 687bf5a629SBarry Smith ierr = PetscOptionsGetString(NULL,NULL, "-log_exclude", logList, 256, &opt);CHKERRQ(ierr); 693ec795f1SBarry Smith if (opt) { 703ec795f1SBarry Smith ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr); 713ec795f1SBarry Smith if (className) { 720700a824SBarry Smith ierr = PetscLogEventDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr); 733ec795f1SBarry Smith } 743ec795f1SBarry Smith } 75b022a5c1SBarry Smith ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr); 763ec795f1SBarry Smith PetscFunctionReturn(0); 773ec795f1SBarry Smith } 783ec795f1SBarry Smith 79e884886fSBarry Smith /*@C 80e884886fSBarry Smith MatMFFDSetType - Sets the method that is used to compute the 81e884886fSBarry Smith differencing parameter for finite differene matrix-free formulations. 82e884886fSBarry Smith 83e884886fSBarry Smith Input Parameters: 84e884886fSBarry Smith + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD() 85e884886fSBarry Smith or MatSetType(mat,MATMFFD); 86e884886fSBarry Smith - ftype - the type requested, either MATMFFD_WP or MATMFFD_DS 87e884886fSBarry Smith 88e884886fSBarry Smith Level: advanced 89e884886fSBarry Smith 90e884886fSBarry Smith Notes: 91e884886fSBarry Smith For example, such routines can compute h for use in 92e884886fSBarry Smith Jacobian-vector products of the form 93e884886fSBarry Smith 94e884886fSBarry Smith F(x+ha) - F(x) 95e884886fSBarry Smith F'(u)a ~= ---------------- 96e884886fSBarry Smith h 97e884886fSBarry Smith 98755b7f64SBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction(), MatCreateMFFD() 99e884886fSBarry Smith @*/ 10019fd82e9SBarry Smith PetscErrorCode MatMFFDSetType(Mat mat,MatMFFDType ftype) 101e884886fSBarry Smith { 102e884886fSBarry Smith PetscErrorCode ierr,(*r)(MatMFFD); 103e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 104ace3abfcSBarry Smith PetscBool match; 105e884886fSBarry Smith 106e884886fSBarry Smith PetscFunctionBegin; 1070700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 108e884886fSBarry Smith PetscValidCharPointer(ftype,2); 109e884886fSBarry Smith 110251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 1119a13eae7SJed Brown if (!match) PetscFunctionReturn(0); 1129a13eae7SJed Brown 113e884886fSBarry Smith /* already set, so just return */ 114251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr); 115e884886fSBarry Smith if (match) PetscFunctionReturn(0); 116e884886fSBarry Smith 117e884886fSBarry Smith /* destroy the old one if it exists */ 118e884886fSBarry Smith if (ctx->ops->destroy) { 119e884886fSBarry Smith ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr); 120e884886fSBarry Smith } 121e884886fSBarry Smith 1221c9cd337SJed Brown ierr = PetscFunctionListFind(MatMFFDList,ftype,&r);CHKERRQ(ierr); 123e32f2f54SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype); 124e884886fSBarry Smith ierr = (*r)(ctx);CHKERRQ(ierr); 125e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr); 126e884886fSBarry Smith PetscFunctionReturn(0); 127e884886fSBarry Smith } 128e884886fSBarry Smith 129e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/ 13040244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func) 131e884886fSBarry Smith { 132e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 133e884886fSBarry Smith 134e884886fSBarry Smith PetscFunctionBegin; 135e884886fSBarry Smith ctx->funcisetbase = func; 136e884886fSBarry Smith PetscFunctionReturn(0); 137e884886fSBarry Smith } 138e884886fSBarry Smith 139e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 14040244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci) 141e884886fSBarry Smith { 142e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 143e884886fSBarry Smith 144e884886fSBarry Smith PetscFunctionBegin; 145e884886fSBarry Smith ctx->funci = funci; 146e884886fSBarry Smith PetscFunctionReturn(0); 147e884886fSBarry Smith } 148e884886fSBarry Smith 14940244768SBarry Smith static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) 150bc0f08ceSBarry Smith { 151bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 152bc0f08ceSBarry Smith 153bc0f08ceSBarry Smith PetscFunctionBegin; 154bc0f08ceSBarry Smith ctx->ncurrenth = 0; 155bc0f08ceSBarry Smith PetscFunctionReturn(0); 156bc0f08ceSBarry Smith } 157e884886fSBarry Smith 1581c84c290SBarry Smith /*@C 1591c84c290SBarry Smith MatMFFDRegister - Adds a method to the MatMFFD registry. 1601c84c290SBarry Smith 1611c84c290SBarry Smith Not Collective 1621c84c290SBarry Smith 1631c84c290SBarry Smith Input Parameters: 1641c84c290SBarry Smith + name_solver - name of a new user-defined compute-h module 1651c84c290SBarry Smith - routine_create - routine to create method context 1661c84c290SBarry Smith 1671c84c290SBarry Smith Level: developer 1681c84c290SBarry Smith 1691c84c290SBarry Smith Notes: 1701c84c290SBarry Smith MatMFFDRegister() may be called multiple times to add several user-defined solvers. 1711c84c290SBarry Smith 1721c84c290SBarry Smith Sample usage: 1731c84c290SBarry Smith .vb 174bdf89e91SBarry Smith MatMFFDRegister("my_h",MyHCreate); 1751c84c290SBarry Smith .ve 1761c84c290SBarry Smith 1771c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 1781c84c290SBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1791c84c290SBarry Smith or at runtime via the option 180be7a6d03SBarry Smith $ -mat_mffd_type my_h 1811c84c290SBarry Smith 1821c84c290SBarry Smith .keywords: MatMFFD, register 1831c84c290SBarry Smith 1841c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1851c84c290SBarry Smith @*/ 186bdf89e91SBarry Smith PetscErrorCode MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD)) 187e884886fSBarry Smith { 188e884886fSBarry Smith PetscErrorCode ierr; 189e884886fSBarry Smith 190e884886fSBarry Smith PetscFunctionBegin; 191a240a19fSJed Brown ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr); 192e884886fSBarry Smith PetscFunctionReturn(0); 193e884886fSBarry Smith } 194e884886fSBarry Smith 195e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/ 19640244768SBarry Smith static PetscErrorCode MatDestroy_MFFD(Mat mat) 197e884886fSBarry Smith { 198e884886fSBarry Smith PetscErrorCode ierr; 199e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 200e884886fSBarry Smith 201e884886fSBarry Smith PetscFunctionBegin; 2026bf464f9SBarry Smith ierr = VecDestroy(&ctx->w);CHKERRQ(ierr); 2036bf464f9SBarry Smith ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr); 2046bf464f9SBarry Smith ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr); 2056bf464f9SBarry Smith ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr); 206*c51ad4d4SStefano Zampini ierr = VecDestroy(&ctx->dshiftw);CHKERRQ(ierr); 207*c51ad4d4SStefano Zampini ierr = VecDestroy(&ctx->current_u);CHKERRQ(ierr); 208cfe22f5eSBarry Smith if (ctx->current_f_allocated) { 2096bf464f9SBarry Smith ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr); 210cfe22f5eSBarry Smith } 211e884886fSBarry Smith if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);} 2126bf464f9SBarry Smith ierr = PetscHeaderDestroy(&ctx);CHKERRQ(ierr); 213bf0cc555SLisandro Dalcin mat->data = 0; 214e884886fSBarry Smith 215bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr); 216bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr); 217bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr); 218bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr); 219bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr); 220bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr); 221bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr); 222bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr); 223e884886fSBarry Smith PetscFunctionReturn(0); 224e884886fSBarry Smith } 225e884886fSBarry Smith 22688b4c220SStefano Zampini static PetscErrorCode MatSetUp_MFFD(Mat mat) 22788b4c220SStefano Zampini { 22888b4c220SStefano Zampini PetscErrorCode ierr; 22988b4c220SStefano Zampini 23088b4c220SStefano Zampini PetscFunctionBegin; 23188b4c220SStefano Zampini ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr); 23288b4c220SStefano Zampini ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr); 23388b4c220SStefano Zampini PetscFunctionReturn(0); 23488b4c220SStefano Zampini } 23588b4c220SStefano Zampini 236e884886fSBarry Smith /* 237e884886fSBarry Smith MatMFFDView_MFFD - Views matrix-free parameters. 238e884886fSBarry Smith 239e884886fSBarry Smith */ 24040244768SBarry Smith static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer) 241e884886fSBarry Smith { 242e884886fSBarry Smith PetscErrorCode ierr; 243e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 244a283635eSDmitry Karpeev PetscBool iascii, viewbase, viewfunction; 245a283635eSDmitry Karpeev const char *prefix; 246e884886fSBarry Smith 247e884886fSBarry Smith PetscFunctionBegin; 248251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 249e884886fSBarry Smith if (iascii) { 250a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr); 251a283635eSDmitry Karpeev ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 25257622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr); 2537adad957SLisandro Dalcin if (!((PetscObject)ctx)->type_name) { 254e884886fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr); 255e884886fSBarry Smith } else { 2567adad957SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr); 257e884886fSBarry Smith } 258e884886fSBarry Smith if (ctx->ops->view) { 259e884886fSBarry Smith ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr); 260e884886fSBarry Smith } 261a283635eSDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr); 262a283635eSDmitry Karpeev 263c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr); 264a283635eSDmitry Karpeev if (viewbase) { 265a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr); 266a283635eSDmitry Karpeev ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr); 267a283635eSDmitry Karpeev } 268c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr); 269a283635eSDmitry Karpeev if (viewfunction) { 270a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr); 271a283635eSDmitry Karpeev ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr); 272a283635eSDmitry Karpeev } 273a283635eSDmitry Karpeev ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 27411aeaf0aSBarry Smith } 275e884886fSBarry Smith PetscFunctionReturn(0); 276e884886fSBarry Smith } 277e884886fSBarry Smith 278e884886fSBarry Smith /* 279e884886fSBarry Smith MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 280e884886fSBarry Smith allows the user to indicate the beginning of a new linear solve by calling 281e884886fSBarry Smith MatAssemblyXXX() on the matrix free matrix. This then allows the 2821d0fab5eSBarry Smith MatCreateMFFD_WP() to properly compute ||U|| only the first time 283e884886fSBarry Smith in the linear solver rather than every time. 2845a576424SJed Brown 285cc2e6a90SBarry Smith This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence 286cc2e6a90SBarry Smith it must be labeled as PETSC_EXTERN 287e884886fSBarry Smith */ 2885a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt) 289e884886fSBarry Smith { 290e884886fSBarry Smith PetscErrorCode ierr; 291e884886fSBarry Smith MatMFFD j = (MatMFFD)J->data; 292e884886fSBarry Smith 293e884886fSBarry Smith PetscFunctionBegin; 294e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 295e884886fSBarry Smith j->vshift = 0.0; 296e884886fSBarry Smith j->vscale = 1.0; 297e884886fSBarry Smith PetscFunctionReturn(0); 298e884886fSBarry Smith } 299e884886fSBarry Smith 300e884886fSBarry Smith /* 301e884886fSBarry Smith MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 302e884886fSBarry Smith 303e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 304e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 305e884886fSBarry Smith u = current iterate 306e884886fSBarry Smith h = difference interval 307e884886fSBarry Smith */ 30840244768SBarry Smith static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y) 309e884886fSBarry Smith { 310e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 311e884886fSBarry Smith PetscScalar h; 312e884886fSBarry Smith Vec w,U,F; 313e884886fSBarry Smith PetscErrorCode ierr; 314ace3abfcSBarry Smith PetscBool zeroa; 315e884886fSBarry Smith 316e884886fSBarry Smith PetscFunctionBegin; 317ce94432eSBarry 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"); 318e884886fSBarry Smith /* We log matrix-free matrix-vector products separately, so that we can 319e884886fSBarry Smith separate the performance monitoring from the cases that use conventional 320e884886fSBarry Smith storage. We may eventually modify event logging to associate events 321e884886fSBarry Smith with particular objects, hence alleviating the more general problem. */ 322e884886fSBarry Smith ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 323e884886fSBarry Smith 324e884886fSBarry Smith w = ctx->w; 325e884886fSBarry Smith U = ctx->current_u; 3263ec795f1SBarry Smith F = ctx->current_f; 327e884886fSBarry Smith /* 328e884886fSBarry Smith Compute differencing parameter 329e884886fSBarry Smith */ 3304a2f8832SBarry Smith if (!((PetscObject)ctx)->type_name) { 331e884886fSBarry Smith ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr); 3329c6ac3b3SBarry Smith ierr = MatSetFromOptions(mat);CHKERRQ(ierr); 333e884886fSBarry Smith } 334e884886fSBarry Smith ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr); 335e884886fSBarry Smith if (zeroa) { 336e884886fSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 337e884886fSBarry Smith PetscFunctionReturn(0); 338e884886fSBarry Smith } 339e884886fSBarry Smith 34084d44b13SHong Zhang if (mat->erroriffailure && PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h"); 341e884886fSBarry Smith if (ctx->checkh) { 342e884886fSBarry Smith ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr); 343e884886fSBarry Smith } 344e884886fSBarry Smith 345e884886fSBarry Smith /* keep a record of the current differencing parameter h */ 346e884886fSBarry Smith ctx->currenth = h; 347e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX) 34857622a8eSBarry Smith ierr = PetscInfo2(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h));CHKERRQ(ierr); 349e884886fSBarry Smith #else 350e884886fSBarry Smith ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr); 351e884886fSBarry Smith #endif 352e884886fSBarry Smith if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 353e884886fSBarry Smith ctx->historyh[ctx->ncurrenth] = h; 354e884886fSBarry Smith } 355e884886fSBarry Smith ctx->ncurrenth++; 356e884886fSBarry Smith 357e884886fSBarry Smith /* w = u + ha */ 358c3bb7e23SBarry Smith if (ctx->drscale) { 359c3bb7e23SBarry Smith ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr); 360c3bb7e23SBarry Smith ierr = VecAYPX(U,h,w);CHKERRQ(ierr); 361c3bb7e23SBarry Smith } else { 362e884886fSBarry Smith ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr); 363c3bb7e23SBarry Smith } 364e884886fSBarry Smith 3651b797266SDmitry Karpeev /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 3661b797266SDmitry Karpeev if (ctx->ncurrenth == 1 && ctx->current_f_allocated) { 367e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr); 36852121784SBarry Smith } 369e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr); 370e884886fSBarry Smith 371e884886fSBarry Smith ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 372e884886fSBarry Smith ierr = VecScale(y,1.0/h);CHKERRQ(ierr); 373e884886fSBarry Smith 374c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 375e884886fSBarry Smith ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr); 376c3bb7e23SBarry Smith } 377c3bb7e23SBarry Smith if (ctx->dlscale) { 378c3bb7e23SBarry Smith ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr); 379c3bb7e23SBarry Smith } 380c3bb7e23SBarry Smith if (ctx->dshift) { 381*c51ad4d4SStefano Zampini if (!ctx->dshiftw) { 382*c51ad4d4SStefano Zampini ierr = VecDuplicate(y,&ctx->dshiftw);CHKERRQ(ierr); 383*c51ad4d4SStefano Zampini } 384*c51ad4d4SStefano Zampini ierr = VecPointwiseMult(ctx->dshift,a,ctx->dshiftw);CHKERRQ(ierr); 385*c51ad4d4SStefano Zampini ierr = VecAXPY(y,1.0,ctx->dshiftw);CHKERRQ(ierr); 386c3bb7e23SBarry Smith } 387e884886fSBarry Smith 38839601f4eSBarry Smith if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);} 389e884886fSBarry Smith 390e884886fSBarry Smith ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 391e884886fSBarry Smith PetscFunctionReturn(0); 392e884886fSBarry Smith } 393e884886fSBarry Smith 394e884886fSBarry Smith /* 395e884886fSBarry Smith MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 396e884886fSBarry Smith 397e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 398e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 399e884886fSBarry Smith u = current iterate 400e884886fSBarry Smith h = difference interval 401e884886fSBarry Smith */ 40240244768SBarry Smith static PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 403e884886fSBarry Smith { 404e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 405e884886fSBarry Smith PetscScalar h,*aa,*ww,v; 406e884886fSBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 407e884886fSBarry Smith Vec w,U; 408e884886fSBarry Smith PetscErrorCode ierr; 409e884886fSBarry Smith PetscInt i,rstart,rend; 410e884886fSBarry Smith 411e884886fSBarry Smith PetscFunctionBegin; 412e7e72b3dSBarry Smith if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first"); 413e884886fSBarry Smith 414e884886fSBarry Smith w = ctx->w; 415e884886fSBarry Smith U = ctx->current_u; 416e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr); 417e884886fSBarry Smith ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr); 418e884886fSBarry Smith ierr = VecCopy(U,w);CHKERRQ(ierr); 419e884886fSBarry Smith 420e884886fSBarry Smith ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr); 421e884886fSBarry Smith ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 422e884886fSBarry Smith for (i=rstart; i<rend; i++) { 423e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 424e884886fSBarry Smith h = ww[i-rstart]; 425e884886fSBarry Smith if (h == 0.0) h = 1.0; 426e884886fSBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 427e884886fSBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 428e884886fSBarry Smith h *= epsilon; 429e884886fSBarry Smith 430e884886fSBarry Smith ww[i-rstart] += h; 431e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 432e884886fSBarry Smith ierr = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr); 433e884886fSBarry Smith aa[i-rstart] = (v - aa[i-rstart])/h; 434e884886fSBarry Smith 435e884886fSBarry Smith /* possibly shift and scale result */ 436c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 437e884886fSBarry Smith aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart]; 438c3bb7e23SBarry Smith } 439e884886fSBarry Smith 440e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 441e884886fSBarry Smith ww[i-rstart] -= h; 442e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 443e884886fSBarry Smith } 444e884886fSBarry Smith ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 445e884886fSBarry Smith PetscFunctionReturn(0); 446e884886fSBarry Smith } 447e884886fSBarry Smith 44840244768SBarry Smith static PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr) 449c3bb7e23SBarry Smith { 450c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 451c3bb7e23SBarry Smith PetscErrorCode ierr; 452c3bb7e23SBarry Smith 453c3bb7e23SBarry Smith PetscFunctionBegin; 454c3bb7e23SBarry Smith if (ll && !aij->dlscale) { 455c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr); 456c3bb7e23SBarry Smith } 457c3bb7e23SBarry Smith if (rr && !aij->drscale) { 458c3bb7e23SBarry Smith ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr); 459c3bb7e23SBarry Smith } 460c3bb7e23SBarry Smith if (ll) { 461c3bb7e23SBarry Smith ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr); 462c3bb7e23SBarry Smith } 463c3bb7e23SBarry Smith if (rr) { 464c3bb7e23SBarry Smith ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr); 465c3bb7e23SBarry Smith } 466c3bb7e23SBarry Smith PetscFunctionReturn(0); 467c3bb7e23SBarry Smith } 468c3bb7e23SBarry Smith 46940244768SBarry Smith static PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode) 470c3bb7e23SBarry Smith { 471c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 472c3bb7e23SBarry Smith PetscErrorCode ierr; 473c3bb7e23SBarry Smith 474c3bb7e23SBarry Smith PetscFunctionBegin; 475ce94432eSBarry Smith if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES"); 476c3bb7e23SBarry Smith if (!aij->dshift) { 477c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr); 478c3bb7e23SBarry Smith } 479c3bb7e23SBarry Smith ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr); 480c3bb7e23SBarry Smith PetscFunctionReturn(0); 481c3bb7e23SBarry Smith } 482c3bb7e23SBarry Smith 48340244768SBarry Smith static PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a) 484e884886fSBarry Smith { 485e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 4865fd66863SKarl Rupp 487e884886fSBarry Smith PetscFunctionBegin; 488e884886fSBarry Smith shell->vshift += a; 489e884886fSBarry Smith PetscFunctionReturn(0); 490e884886fSBarry Smith } 491e884886fSBarry Smith 49240244768SBarry Smith static PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a) 493e884886fSBarry Smith { 494e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 4955fd66863SKarl Rupp 496e884886fSBarry Smith PetscFunctionBegin; 497e884886fSBarry Smith shell->vscale *= a; 498e884886fSBarry Smith PetscFunctionReturn(0); 499e884886fSBarry Smith } 500e884886fSBarry Smith 501d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F) 502e884886fSBarry Smith { 503e884886fSBarry Smith PetscErrorCode ierr; 504e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 505e884886fSBarry Smith 506e884886fSBarry Smith PetscFunctionBegin; 507e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 508*c51ad4d4SStefano Zampini if (!ctx->current_u) { 509*c51ad4d4SStefano Zampini ierr = VecDuplicate(U,&ctx->current_u);CHKERRQ(ierr); 510*c51ad4d4SStefano Zampini ierr = VecLockPush(ctx->current_u);CHKERRQ(ierr); 511*c51ad4d4SStefano Zampini } 512*c51ad4d4SStefano Zampini ierr = VecLockPop(ctx->current_u);CHKERRQ(ierr); 513*c51ad4d4SStefano Zampini ierr = VecCopy(U,ctx->current_u);CHKERRQ(ierr); 514*c51ad4d4SStefano Zampini ierr = VecLockPush(ctx->current_u);CHKERRQ(ierr); 51552121784SBarry Smith if (F) { 5166bf464f9SBarry Smith if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);} 5173ec795f1SBarry Smith ctx->current_f = F; 518cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_FALSE; 519cfe22f5eSBarry Smith } else if (!ctx->current_f_allocated) { 52006c4b6baSBarry Smith ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr); 5212205254eSKarl Rupp 522cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_TRUE; 52352121784SBarry Smith } 524e884886fSBarry Smith if (!ctx->w) { 525e884886fSBarry Smith ierr = VecDuplicate(ctx->current_u,&ctx->w);CHKERRQ(ierr); 526e884886fSBarry Smith } 527e884886fSBarry Smith J->assembled = PETSC_TRUE; 528e884886fSBarry Smith PetscFunctionReturn(0); 529e884886fSBarry Smith } 5305a576424SJed Brown 531e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 532b2573a8aSBarry Smith 53340244768SBarry Smith static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx) 534e884886fSBarry Smith { 535e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 536e884886fSBarry Smith 537e884886fSBarry Smith PetscFunctionBegin; 538e884886fSBarry Smith ctx->checkh = fun; 539e884886fSBarry Smith ctx->checkhctx = ectx; 540e884886fSBarry Smith PetscFunctionReturn(0); 541e884886fSBarry Smith } 542e884886fSBarry Smith 5436aa9148fSLisandro Dalcin /*@C 5446aa9148fSLisandro Dalcin MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 5456aa9148fSLisandro Dalcin MatMFFD options in the database. 5466aa9148fSLisandro Dalcin 5476aa9148fSLisandro Dalcin Collective on Mat 5486aa9148fSLisandro Dalcin 5496aa9148fSLisandro Dalcin Input Parameter: 5506aa9148fSLisandro Dalcin + A - the Mat context 5516aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names 5526aa9148fSLisandro Dalcin 5536aa9148fSLisandro Dalcin Notes: 5546aa9148fSLisandro Dalcin A hyphen (-) must NOT be given at the beginning of the prefix name. 5556aa9148fSLisandro Dalcin The first character of all runtime options is AUTOMATICALLY the hyphen. 5566aa9148fSLisandro Dalcin 5576aa9148fSLisandro Dalcin Level: advanced 5586aa9148fSLisandro Dalcin 5596aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters 5606aa9148fSLisandro Dalcin 561755b7f64SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD() 5626aa9148fSLisandro Dalcin @*/ 5637087cfbeSBarry Smith PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[]) 5646aa9148fSLisandro Dalcin 5656aa9148fSLisandro Dalcin { 5660298fd71SBarry Smith MatMFFD mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL; 5676aa9148fSLisandro Dalcin PetscErrorCode ierr; 5685fd66863SKarl Rupp 5696aa9148fSLisandro Dalcin PetscFunctionBegin; 5700700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5710700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 5726aa9148fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr); 5736aa9148fSLisandro Dalcin PetscFunctionReturn(0); 5746aa9148fSLisandro Dalcin } 5756aa9148fSLisandro Dalcin 57640244768SBarry Smith static PetscErrorCode MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat) 577e884886fSBarry Smith { 57871f08665SBarry Smith MatMFFD mfctx = (MatMFFD)mat->data; 579e884886fSBarry Smith PetscErrorCode ierr; 580ace3abfcSBarry Smith PetscBool flg; 581e884886fSBarry Smith char ftype[256]; 582e884886fSBarry Smith 583e884886fSBarry Smith PetscFunctionBegin; 5840700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5850700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 5863194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr); 587a264d7a6SBarry Smith ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr); 588e884886fSBarry Smith if (flg) { 589e884886fSBarry Smith ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr); 590e884886fSBarry Smith } 591e884886fSBarry Smith 592e884886fSBarry Smith ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 593e884886fSBarry Smith ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 594e884886fSBarry Smith 59590d69ab7SBarry Smith flg = PETSC_FALSE; 5960298fd71SBarry Smith ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr); 597e884886fSBarry Smith if (flg) { 598e884886fSBarry Smith ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr); 599e884886fSBarry Smith } 600e884886fSBarry Smith if (mfctx->ops->setfromoptions) { 601e55864a3SBarry Smith ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr); 602e884886fSBarry Smith } 603e884886fSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 604e884886fSBarry Smith PetscFunctionReturn(0); 605e884886fSBarry Smith } 606e884886fSBarry Smith 60740244768SBarry Smith static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period) 608bc0f08ceSBarry Smith { 609bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 610bc0f08ceSBarry Smith 611bc0f08ceSBarry Smith PetscFunctionBegin; 612bc0f08ceSBarry Smith ctx->recomputeperiod = period; 613bc0f08ceSBarry Smith PetscFunctionReturn(0); 614bc0f08ceSBarry Smith } 615bc0f08ceSBarry Smith 61640244768SBarry Smith static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 617bc0f08ceSBarry Smith { 618bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 619bc0f08ceSBarry Smith 620bc0f08ceSBarry Smith PetscFunctionBegin; 621bc0f08ceSBarry Smith ctx->func = func; 622bc0f08ceSBarry Smith ctx->funcctx = funcctx; 623bc0f08ceSBarry Smith PetscFunctionReturn(0); 624bc0f08ceSBarry Smith } 625bc0f08ceSBarry Smith 62640244768SBarry Smith static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error) 627bc0f08ceSBarry Smith { 628bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 629bc0f08ceSBarry Smith 630bc0f08ceSBarry Smith PetscFunctionBegin; 631bc0f08ceSBarry Smith if (error != PETSC_DEFAULT) ctx->error_rel = error; 632bc0f08ceSBarry Smith PetscFunctionReturn(0); 633bc0f08ceSBarry Smith } 634bc0f08ceSBarry Smith 6353b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MFFD(Mat A,PetscBool *missing,PetscInt *d) 6363b49f96aSBarry Smith { 6373b49f96aSBarry Smith PetscFunctionBegin; 6383b49f96aSBarry Smith *missing = PETSC_FALSE; 6393b49f96aSBarry Smith PetscFunctionReturn(0); 6403b49f96aSBarry Smith } 6413b49f96aSBarry Smith 642e884886fSBarry Smith /*MC 643e884886fSBarry Smith MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 644e884886fSBarry Smith 645e884886fSBarry Smith Level: advanced 646e884886fSBarry Smith 647755b7f64SBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(), 648755b7f64SBarry Smith MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 649755b7f64SBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 650755b7f64SBarry Smith MatMFFDGetH(), 651e884886fSBarry Smith M*/ 6528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 653e884886fSBarry Smith { 654e884886fSBarry Smith MatMFFD mfctx; 655e884886fSBarry Smith PetscErrorCode ierr; 656e884886fSBarry Smith 657e884886fSBarry Smith PetscFunctionBegin; 658607a6623SBarry Smith ierr = MatMFFDInitializePackage();CHKERRQ(ierr); 659e884886fSBarry Smith 66073107ff1SLisandro Dalcin ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr); 6612205254eSKarl Rupp 662e884886fSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 663e884886fSBarry Smith mfctx->recomputeperiod = 1; 664e884886fSBarry Smith mfctx->count = 0; 665e884886fSBarry Smith mfctx->currenth = 0.0; 6660298fd71SBarry Smith mfctx->historyh = NULL; 667e884886fSBarry Smith mfctx->ncurrenth = 0; 668e884886fSBarry Smith mfctx->maxcurrenth = 0; 6697adad957SLisandro Dalcin ((PetscObject)mfctx)->type_name = 0; 670e884886fSBarry Smith 671e884886fSBarry Smith mfctx->vshift = 0.0; 672e884886fSBarry Smith mfctx->vscale = 1.0; 673e884886fSBarry Smith 674e884886fSBarry Smith /* 675e884886fSBarry Smith Create the empty data structure to contain compute-h routines. 676e884886fSBarry Smith These will be filled in below from the command line options or 677e884886fSBarry Smith a later call with MatMFFDSetType() or if that is not called 678e884886fSBarry Smith then it will default in the first use of MatMult_MFFD() 679e884886fSBarry Smith */ 680e884886fSBarry Smith mfctx->ops->compute = 0; 681e884886fSBarry Smith mfctx->ops->destroy = 0; 682e884886fSBarry Smith mfctx->ops->view = 0; 683e884886fSBarry Smith mfctx->ops->setfromoptions = 0; 684e884886fSBarry Smith mfctx->hctx = 0; 685e884886fSBarry Smith 686e884886fSBarry Smith mfctx->func = 0; 687e884886fSBarry Smith mfctx->funcctx = 0; 6880298fd71SBarry Smith mfctx->w = NULL; 689e884886fSBarry Smith 690e884886fSBarry Smith A->data = mfctx; 691e884886fSBarry Smith 692e884886fSBarry Smith A->ops->mult = MatMult_MFFD; 693e884886fSBarry Smith A->ops->destroy = MatDestroy_MFFD; 69488b4c220SStefano Zampini A->ops->setup = MatSetUp_MFFD; 695e884886fSBarry Smith A->ops->view = MatView_MFFD; 696e884886fSBarry Smith A->ops->assemblyend = MatAssemblyEnd_MFFD; 697e884886fSBarry Smith A->ops->getdiagonal = MatGetDiagonal_MFFD; 698e884886fSBarry Smith A->ops->scale = MatScale_MFFD; 699e884886fSBarry Smith A->ops->shift = MatShift_MFFD; 700c3bb7e23SBarry Smith A->ops->diagonalscale = MatDiagonalScale_MFFD; 701c3bb7e23SBarry Smith A->ops->diagonalset = MatDiagonalSet_MFFD; 7029c6ac3b3SBarry Smith A->ops->setfromoptions = MatSetFromOptions_MFFD; 7033b49f96aSBarry Smith A->ops->missingdiagonal = MatMissingDiagonal_MFFD; 704e884886fSBarry Smith A->assembled = PETSC_TRUE; 705e884886fSBarry Smith 706bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr); 707bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr); 708bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr); 709bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr); 710bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr); 711bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr); 712bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr); 713bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr); 7142205254eSKarl Rupp 715e884886fSBarry Smith mfctx->mat = A; 7162205254eSKarl Rupp 717e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr); 718e884886fSBarry Smith PetscFunctionReturn(0); 719e884886fSBarry Smith } 720e884886fSBarry Smith 721e884886fSBarry Smith /*@ 722e884886fSBarry Smith MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 723e884886fSBarry Smith 724e884886fSBarry Smith Collective on Vec 725e884886fSBarry Smith 726e884886fSBarry Smith Input Parameters: 727fef1beadSBarry Smith + comm - MPI communicator 728fef1beadSBarry Smith . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 729fef1beadSBarry Smith This value should be the same as the local size used in creating the 730fef1beadSBarry Smith y vector for the matrix-vector product y = Ax. 731fef1beadSBarry Smith . n - This value should be the same as the local size used in creating the 732fef1beadSBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 733fef1beadSBarry Smith calculated if N is given) For square matrices n is almost always m. 734fef1beadSBarry Smith . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 735fef1beadSBarry Smith - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 736fef1beadSBarry Smith 737e884886fSBarry Smith 738e884886fSBarry Smith Output Parameter: 739e884886fSBarry Smith . J - the matrix-free matrix 740e884886fSBarry Smith 7419c6ac3b3SBarry Smith Options Database Keys: call MatSetFromOptions() to trigger these 7429c6ac3b3SBarry Smith + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 7439c6ac3b3SBarry Smith - -mat_mffd_err - square root of estimated relative error in function evaluation 7449c6ac3b3SBarry Smith - -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 7459c6ac3b3SBarry Smith 7469c6ac3b3SBarry Smith 747e884886fSBarry Smith Level: advanced 748e884886fSBarry Smith 749e884886fSBarry Smith Notes: 750e884886fSBarry Smith The matrix-free matrix context merely contains the function pointers 751e884886fSBarry Smith and work space for performing finite difference approximations of 752e884886fSBarry Smith Jacobian-vector products, F'(u)*a, 753e884886fSBarry Smith 754e884886fSBarry Smith The default code uses the following approach to compute h 755e884886fSBarry Smith 756e884886fSBarry Smith .vb 757e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 758e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 759e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 760e884886fSBarry Smith where 761e884886fSBarry Smith error_rel = square root of relative error in function evaluation 762e884886fSBarry Smith umin = minimum iterate parameter 763e884886fSBarry Smith .ve 764e884886fSBarry Smith 765ca93e954SBarry Smith You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different 766ca93e954SBarry Smith preconditioner matrix 767ca93e954SBarry Smith 768e884886fSBarry Smith The user can set the error_rel via MatMFFDSetFunctionError() and 769a7f22e61SSatish Balay umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details. 770e884886fSBarry Smith 771e884886fSBarry Smith The user should call MatDestroy() when finished with the matrix-free 772e884886fSBarry Smith matrix context. 773e884886fSBarry Smith 774e884886fSBarry Smith Options Database Keys: 775e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 776e884886fSBarry Smith . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only) 777e884886fSBarry Smith - -mat_mffd_check_positivity 778e884886fSBarry Smith 779e884886fSBarry Smith .keywords: default, matrix-free, create, matrix 780e884886fSBarry Smith 7811d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 782e884886fSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 78381242352SJed Brown MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian() 784e884886fSBarry Smith 785e884886fSBarry Smith @*/ 7867087cfbeSBarry Smith PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 787e884886fSBarry Smith { 788e884886fSBarry Smith PetscErrorCode ierr; 789e884886fSBarry Smith 790e884886fSBarry Smith PetscFunctionBegin; 791e884886fSBarry Smith ierr = MatCreate(comm,J);CHKERRQ(ierr); 792fef1beadSBarry Smith ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr); 793e884886fSBarry Smith ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 794be7c243fSJed Brown ierr = MatSetUp(*J);CHKERRQ(ierr); 795e884886fSBarry Smith PetscFunctionReturn(0); 796e884886fSBarry Smith } 797e884886fSBarry Smith 798e884886fSBarry Smith /*@ 799e884886fSBarry Smith MatMFFDGetH - Gets the last value that was used as the differencing 800e884886fSBarry Smith parameter. 801e884886fSBarry Smith 802e884886fSBarry Smith Not Collective 803e884886fSBarry Smith 804e884886fSBarry Smith Input Parameters: 805e884886fSBarry Smith . mat - the matrix obtained with MatCreateSNESMF() 806e884886fSBarry Smith 807e884886fSBarry Smith Output Paramter: 808e884886fSBarry Smith . h - the differencing step size 809e884886fSBarry Smith 810e884886fSBarry Smith Level: advanced 811e884886fSBarry Smith 812e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 813e884886fSBarry Smith 8141d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory() 815e884886fSBarry Smith @*/ 8167087cfbeSBarry Smith PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 817e884886fSBarry Smith { 818e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 819bc0f08ceSBarry Smith PetscErrorCode ierr; 820bc0f08ceSBarry Smith PetscBool match; 821e884886fSBarry Smith 822e884886fSBarry Smith PetscFunctionBegin; 82388b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 82488b4c220SStefano Zampini PetscValidPointer(h,2); 825251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 826ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 827bc0f08ceSBarry Smith 828e884886fSBarry Smith *h = ctx->currenth; 829e884886fSBarry Smith PetscFunctionReturn(0); 830e884886fSBarry Smith } 831e884886fSBarry Smith 832e884886fSBarry Smith /*@C 833e884886fSBarry Smith MatMFFDSetFunction - Sets the function used in applying the matrix free. 834e884886fSBarry Smith 8353f9fe445SBarry Smith Logically Collective on Mat 836e884886fSBarry Smith 837e884886fSBarry Smith Input Parameters: 83814f633e2SBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD() 839e884886fSBarry Smith . func - the function to use 840e884886fSBarry Smith - funcctx - optional function context passed to function 841e884886fSBarry Smith 84214f633e2SBarry Smith Calling Sequence of func: 84314f633e2SBarry Smith $ func (void *funcctx, Vec x, Vec f) 84414f633e2SBarry Smith 84514f633e2SBarry Smith + funcctx - user provided context 84614f633e2SBarry Smith . x - input vector 84714f633e2SBarry Smith - f - computed output function 84814f633e2SBarry Smith 849e884886fSBarry Smith Level: advanced 850e884886fSBarry Smith 851e884886fSBarry Smith Notes: 852e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 853e884886fSBarry Smith matrix inside your compute Jacobian routine 854e884886fSBarry Smith 8553ec795f1SBarry Smith If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 856e884886fSBarry Smith 857e884886fSBarry Smith .keywords: SNES, matrix-free, function 858e884886fSBarry Smith 8591d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD, 8601d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 861e884886fSBarry Smith @*/ 8627087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 863e884886fSBarry Smith { 864bc0f08ceSBarry Smith PetscErrorCode ierr; 865e884886fSBarry Smith 866e884886fSBarry Smith PetscFunctionBegin; 86788b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 868bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr); 869e884886fSBarry Smith PetscFunctionReturn(0); 870e884886fSBarry Smith } 871e884886fSBarry Smith 872e884886fSBarry Smith /*@C 873e884886fSBarry Smith MatMFFDSetFunctioni - Sets the function for a single component 874e884886fSBarry Smith 8753f9fe445SBarry Smith Logically Collective on Mat 876e884886fSBarry Smith 877e884886fSBarry Smith Input Parameters: 878e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 879e884886fSBarry Smith - funci - the function to use 880e884886fSBarry Smith 881e884886fSBarry Smith Level: advanced 882e884886fSBarry Smith 883e884886fSBarry Smith Notes: 884e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 885e884886fSBarry Smith matrix inside your compute Jacobian routine 886e884886fSBarry Smith 887e884886fSBarry Smith 888e884886fSBarry Smith .keywords: SNES, matrix-free, function 889e884886fSBarry Smith 8901d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 8911d0fab5eSBarry Smith 892e884886fSBarry Smith @*/ 8937087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 894e884886fSBarry Smith { 8954ac538c5SBarry Smith PetscErrorCode ierr; 896e884886fSBarry Smith 897e884886fSBarry Smith PetscFunctionBegin; 8980700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 8994ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr); 900e884886fSBarry Smith PetscFunctionReturn(0); 901e884886fSBarry Smith } 902e884886fSBarry Smith 903e884886fSBarry Smith /*@C 904e884886fSBarry Smith MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 905e884886fSBarry Smith 9063f9fe445SBarry Smith Logically Collective on Mat 907e884886fSBarry Smith 908e884886fSBarry Smith Input Parameters: 909e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 910e884886fSBarry Smith - func - the function to use 911e884886fSBarry Smith 912e884886fSBarry Smith Level: advanced 913e884886fSBarry Smith 914e884886fSBarry Smith Notes: 915e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 916e884886fSBarry Smith matrix inside your compute Jacobian routine 917e884886fSBarry Smith 918e884886fSBarry Smith 919e884886fSBarry Smith .keywords: SNES, matrix-free, function 920e884886fSBarry Smith 921e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 9221d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 923e884886fSBarry Smith @*/ 9247087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 925e884886fSBarry Smith { 9264ac538c5SBarry Smith PetscErrorCode ierr; 927e884886fSBarry Smith 928e884886fSBarry Smith PetscFunctionBegin; 9290700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 9304ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr); 931e884886fSBarry Smith PetscFunctionReturn(0); 932e884886fSBarry Smith } 933e884886fSBarry Smith 934e884886fSBarry Smith /*@ 935e884886fSBarry Smith MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 936e884886fSBarry Smith 9373f9fe445SBarry Smith Logically Collective on Mat 938e884886fSBarry Smith 939e884886fSBarry Smith Input Parameters: 940e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 941e884886fSBarry Smith - period - 1 for everytime, 2 for every second etc 942e884886fSBarry Smith 943e884886fSBarry Smith Options Database Keys: 944e884886fSBarry Smith + -mat_mffd_period <period> 945e884886fSBarry Smith 946e884886fSBarry Smith Level: advanced 947e884886fSBarry Smith 948e884886fSBarry Smith 949e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 950e884886fSBarry Smith 951e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), 9521d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 953e884886fSBarry Smith @*/ 9547087cfbeSBarry Smith PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 955e884886fSBarry Smith { 956bc0f08ceSBarry Smith PetscErrorCode ierr; 957e884886fSBarry Smith 958e884886fSBarry Smith PetscFunctionBegin; 95988b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 96088b4c220SStefano Zampini PetscValidLogicalCollectiveInt(mat,period,2); 961bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr); 962e884886fSBarry Smith PetscFunctionReturn(0); 963e884886fSBarry Smith } 964e884886fSBarry Smith 965e884886fSBarry Smith /*@ 966e884886fSBarry Smith MatMFFDSetFunctionError - Sets the error_rel for the approximation of 967e884886fSBarry Smith matrix-vector products using finite differences. 968e884886fSBarry Smith 9693f9fe445SBarry Smith Logically Collective on Mat 970e884886fSBarry Smith 971e884886fSBarry Smith Input Parameters: 972e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 973e884886fSBarry Smith - error_rel - relative error (should be set to the square root of 974e884886fSBarry Smith the relative error in the function evaluations) 975e884886fSBarry Smith 976e884886fSBarry Smith Options Database Keys: 977e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 978e884886fSBarry Smith 979e884886fSBarry Smith Level: advanced 980e884886fSBarry Smith 981e884886fSBarry Smith Notes: 982e884886fSBarry Smith The default matrix-free matrix-vector product routine computes 983e884886fSBarry Smith .vb 984e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 985e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 986e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 987e884886fSBarry Smith .ve 988e884886fSBarry Smith 989e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 990e884886fSBarry Smith 991e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 9921d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 993e884886fSBarry Smith @*/ 9947087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 995e884886fSBarry Smith { 996bc0f08ceSBarry Smith PetscErrorCode ierr; 997e884886fSBarry Smith 998e884886fSBarry Smith PetscFunctionBegin; 99988b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 100088b4c220SStefano Zampini PetscValidLogicalCollectiveReal(mat,error,2); 1001bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr); 1002e884886fSBarry Smith PetscFunctionReturn(0); 1003e884886fSBarry Smith } 1004e884886fSBarry Smith 1005e884886fSBarry Smith /*@ 1006e884886fSBarry Smith MatMFFDSetHHistory - Sets an array to collect a history of the 1007e884886fSBarry Smith differencing values (h) computed for the matrix-free product. 1008e884886fSBarry Smith 10093f9fe445SBarry Smith Logically Collective on Mat 1010e884886fSBarry Smith 1011e884886fSBarry Smith Input Parameters: 1012e884886fSBarry Smith + J - the matrix-free matrix context 1013e884886fSBarry Smith . histroy - space to hold the history 1014e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than 1015e884886fSBarry Smith nhistory, then the later ones are discarded 1016e884886fSBarry Smith 1017e884886fSBarry Smith Level: advanced 1018e884886fSBarry Smith 1019e884886fSBarry Smith Notes: 1020e884886fSBarry Smith Use MatMFFDResetHHistory() to reset the history counter and collect 1021e884886fSBarry Smith a new batch of differencing parameters, h. 1022e884886fSBarry Smith 1023e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1024e884886fSBarry Smith 1025e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 10261d0fab5eSBarry Smith MatMFFDResetHHistory(), MatMFFDSetFunctionError() 1027e884886fSBarry Smith 1028e884886fSBarry Smith @*/ 10297087cfbeSBarry Smith PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 1030e884886fSBarry Smith { 1031e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 1032bc0f08ceSBarry Smith PetscErrorCode ierr; 1033bc0f08ceSBarry Smith PetscBool match; 1034e884886fSBarry Smith 1035e884886fSBarry Smith PetscFunctionBegin; 103688b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 103788b4c220SStefano Zampini if (history) PetscValidPointer(history,2); 103888b4c220SStefano Zampini PetscValidLogicalCollectiveInt(J,nhistory,3); 1039251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr); 1040ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 1041e884886fSBarry Smith ctx->historyh = history; 1042e884886fSBarry Smith ctx->maxcurrenth = nhistory; 104375567043SBarry Smith ctx->currenth = 0.; 1044e884886fSBarry Smith PetscFunctionReturn(0); 1045e884886fSBarry Smith } 1046e884886fSBarry Smith 1047e884886fSBarry Smith /*@ 1048e884886fSBarry Smith MatMFFDResetHHistory - Resets the counter to zero to begin 1049e884886fSBarry Smith collecting a new set of differencing histories. 1050e884886fSBarry Smith 10513f9fe445SBarry Smith Logically Collective on Mat 1052e884886fSBarry Smith 1053e884886fSBarry Smith Input Parameters: 1054e884886fSBarry Smith . J - the matrix-free matrix context 1055e884886fSBarry Smith 1056e884886fSBarry Smith Level: advanced 1057e884886fSBarry Smith 1058e884886fSBarry Smith Notes: 1059e884886fSBarry Smith Use MatMFFDSetHHistory() to create the original history counter. 1060e884886fSBarry Smith 1061e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1062e884886fSBarry Smith 1063e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 10641d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDSetFunctionError() 1065e884886fSBarry Smith 1066e884886fSBarry Smith @*/ 10677087cfbeSBarry Smith PetscErrorCode MatMFFDResetHHistory(Mat J) 1068e884886fSBarry Smith { 1069bc0f08ceSBarry Smith PetscErrorCode ierr; 1070e884886fSBarry Smith 1071e884886fSBarry Smith PetscFunctionBegin; 107288b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1073bc0f08ceSBarry Smith ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr); 1074e884886fSBarry Smith PetscFunctionReturn(0); 1075e884886fSBarry Smith } 1076e884886fSBarry Smith 1077e884886fSBarry Smith /*@ 1078e884886fSBarry Smith MatMFFDSetBase - Sets the vector U at which matrix vector products of the 1079e884886fSBarry Smith Jacobian are computed 1080e884886fSBarry Smith 10813f9fe445SBarry Smith Logically Collective on Mat 1082e884886fSBarry Smith 1083e884886fSBarry Smith Input Parameters: 1084e884886fSBarry Smith + J - the MatMFFD matrix 10853ec795f1SBarry Smith . U - the vector 1086bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed 1087e884886fSBarry Smith 1088e884886fSBarry Smith Notes: This is rarely used directly 1089e884886fSBarry Smith 10908af5ae88SBarry Smith If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 10918af5ae88SBarry Smith point during the first MatMult() after each call to MatMFFDSetBase(). 1092dff2f722SBarry Smith 1093e884886fSBarry Smith Level: advanced 1094e884886fSBarry Smith 1095e884886fSBarry Smith @*/ 10967087cfbeSBarry Smith PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 1097e884886fSBarry Smith { 10984ac538c5SBarry Smith PetscErrorCode ierr; 1099e884886fSBarry Smith 1100e884886fSBarry Smith PetscFunctionBegin; 11010700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 11020700a824SBarry Smith PetscValidHeaderSpecific(U,VEC_CLASSID,2); 11030700a824SBarry Smith if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 11044ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr); 1105e884886fSBarry Smith PetscFunctionReturn(0); 1106e884886fSBarry Smith } 1107e884886fSBarry Smith 1108e884886fSBarry Smith /*@C 1109e884886fSBarry Smith MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 1110e884886fSBarry Smith it to satisfy some criteria 1111e884886fSBarry Smith 11123f9fe445SBarry Smith Logically Collective on Mat 1113e884886fSBarry Smith 1114e884886fSBarry Smith Input Parameters: 1115e884886fSBarry Smith + J - the MatMFFD matrix 1116e884886fSBarry Smith . fun - the function that checks h 1117e884886fSBarry Smith - ctx - any context needed by the function 1118e884886fSBarry Smith 1119e884886fSBarry Smith Options Database Keys: 1120e884886fSBarry Smith . -mat_mffd_check_positivity 1121e884886fSBarry Smith 1122e884886fSBarry Smith Level: advanced 1123e884886fSBarry Smith 1124755b7f64SBarry Smith Notes: For example, MatMFFDCheckPositivity() insures that all entries 1125e884886fSBarry Smith of U + h*a are non-negative 1126e884886fSBarry Smith 1127755b7f64SBarry Smith The function you provide is called after the default h has been computed and allows you to 1128755b7f64SBarry Smith modify it. 1129755b7f64SBarry Smith 1130755b7f64SBarry Smith .seealso: MatMFFDCheckPositivity() 1131e884886fSBarry Smith @*/ 11327087cfbeSBarry Smith PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx) 1133e884886fSBarry Smith { 11344ac538c5SBarry Smith PetscErrorCode ierr; 1135e884886fSBarry Smith 1136e884886fSBarry Smith PetscFunctionBegin; 11370700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 11384ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr); 1139e884886fSBarry Smith PetscFunctionReturn(0); 1140e884886fSBarry Smith } 1141e884886fSBarry Smith 1142e884886fSBarry Smith /*@ 1143e884886fSBarry Smith MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1144e884886fSBarry Smith zero, decreases h until this is satisfied. 1145e884886fSBarry Smith 11463f9fe445SBarry Smith Logically Collective on Vec 1147e884886fSBarry Smith 1148e884886fSBarry Smith Input Parameters: 1149e884886fSBarry Smith + U - base vector that is added to 1150e884886fSBarry Smith . a - vector that is added 1151e884886fSBarry Smith . h - scaling factor on a 1152e884886fSBarry Smith - dummy - context variable (unused) 1153e884886fSBarry Smith 1154e884886fSBarry Smith Options Database Keys: 1155e884886fSBarry Smith . -mat_mffd_check_positivity 1156e884886fSBarry Smith 1157e884886fSBarry Smith Level: advanced 1158e884886fSBarry Smith 1159e884886fSBarry Smith Notes: This is rarely used directly, rather it is passed as an argument to 1160e884886fSBarry Smith MatMFFDSetCheckh() 1161e884886fSBarry Smith 1162e884886fSBarry Smith .seealso: MatMFFDSetCheckh() 1163e884886fSBarry Smith @*/ 11647087cfbeSBarry Smith PetscErrorCode MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h) 1165e884886fSBarry Smith { 1166e884886fSBarry Smith PetscReal val, minval; 1167e884886fSBarry Smith PetscScalar *u_vec, *a_vec; 1168e884886fSBarry Smith PetscErrorCode ierr; 1169e884886fSBarry Smith PetscInt i,n; 1170e884886fSBarry Smith MPI_Comm comm; 1171e884886fSBarry Smith 1172e884886fSBarry Smith PetscFunctionBegin; 117388b4c220SStefano Zampini PetscValidHeaderSpecific(U,VEC_CLASSID,2); 117488b4c220SStefano Zampini PetscValidHeaderSpecific(a,VEC_CLASSID,3); 117588b4c220SStefano Zampini PetscValidPointer(h,4); 1176e884886fSBarry Smith ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); 1177e884886fSBarry Smith ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); 1178e884886fSBarry Smith ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); 1179e884886fSBarry Smith ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); 1180e884886fSBarry Smith minval = PetscAbsScalar(*h*1.01); 1181e884886fSBarry Smith for (i=0; i<n; i++) { 1182e884886fSBarry Smith if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1183e884886fSBarry Smith val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1184e884886fSBarry Smith if (val < minval) minval = val; 1185e884886fSBarry Smith } 1186e884886fSBarry Smith } 1187e884886fSBarry Smith ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); 1188e884886fSBarry Smith ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); 1189b2566f29SBarry Smith ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr); 1190e884886fSBarry Smith if (val <= PetscAbsScalar(*h)) { 11916712e2f1SBarry Smith ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr); 1192e884886fSBarry Smith if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1193e884886fSBarry Smith else *h = -0.99*val; 1194e884886fSBarry Smith } 1195e884886fSBarry Smith PetscFunctionReturn(0); 1196e884886fSBarry Smith } 1197