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 1295d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat,Vec); 1305d21e1f8SStefano Zampini 131e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/ 13240244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func) 133e884886fSBarry Smith { 134e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 135e884886fSBarry Smith 136e884886fSBarry Smith PetscFunctionBegin; 137e884886fSBarry Smith ctx->funcisetbase = func; 138*486afcceSStefano Zampini /* allow users to compose their own getdiagonal and allow MatHasOperation 139*486afcceSStefano Zampini to return false if the two functions pointers are not set */ 140*486afcceSStefano Zampini if (!mat->ops->getdiagonal && func) { 141*486afcceSStefano Zampini mat->ops->getdiagonal == MatGetDiagonal_MFFD; 1425d21e1f8SStefano Zampini } 143e884886fSBarry Smith PetscFunctionReturn(0); 144e884886fSBarry Smith } 145e884886fSBarry Smith 146e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 14740244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci) 148e884886fSBarry Smith { 149e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 150e884886fSBarry Smith 151e884886fSBarry Smith PetscFunctionBegin; 152e884886fSBarry Smith ctx->funci = funci; 153*486afcceSStefano Zampini /* allow users to compose their own getdiagonal and allow MatHasOperation 154*486afcceSStefano Zampini to return false if the two functions pointers are not set */ 155*486afcceSStefano Zampini if (!mat->ops->getdiagonal && funci) { 156*486afcceSStefano Zampini mat->ops->getdiagonal == MatGetDiagonal_MFFD; 1575d21e1f8SStefano Zampini } 158e884886fSBarry Smith PetscFunctionReturn(0); 159e884886fSBarry Smith } 160e884886fSBarry Smith 16140244768SBarry Smith static 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 1701c84c290SBarry Smith /*@C 1711c84c290SBarry Smith MatMFFDRegister - Adds a method to the MatMFFD registry. 1721c84c290SBarry Smith 1731c84c290SBarry Smith Not Collective 1741c84c290SBarry Smith 1751c84c290SBarry Smith Input Parameters: 1761c84c290SBarry Smith + name_solver - name of a new user-defined compute-h module 1771c84c290SBarry Smith - routine_create - routine to create method context 1781c84c290SBarry Smith 1791c84c290SBarry Smith Level: developer 1801c84c290SBarry Smith 1811c84c290SBarry Smith Notes: 1821c84c290SBarry Smith MatMFFDRegister() may be called multiple times to add several user-defined solvers. 1831c84c290SBarry Smith 1841c84c290SBarry Smith Sample usage: 1851c84c290SBarry Smith .vb 186bdf89e91SBarry Smith MatMFFDRegister("my_h",MyHCreate); 1871c84c290SBarry Smith .ve 1881c84c290SBarry Smith 1891c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 1901c84c290SBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1911c84c290SBarry Smith or at runtime via the option 192be7a6d03SBarry Smith $ -mat_mffd_type my_h 1931c84c290SBarry Smith 1941c84c290SBarry Smith .keywords: MatMFFD, register 1951c84c290SBarry Smith 1961c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1971c84c290SBarry Smith @*/ 198bdf89e91SBarry Smith PetscErrorCode MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD)) 199e884886fSBarry Smith { 200e884886fSBarry Smith PetscErrorCode ierr; 201e884886fSBarry Smith 202e884886fSBarry Smith PetscFunctionBegin; 203a240a19fSJed Brown ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr); 204e884886fSBarry Smith PetscFunctionReturn(0); 205e884886fSBarry Smith } 206e884886fSBarry Smith 207e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/ 20840244768SBarry Smith static PetscErrorCode MatDestroy_MFFD(Mat mat) 209e884886fSBarry Smith { 210e884886fSBarry Smith PetscErrorCode ierr; 211e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 212e884886fSBarry Smith 213e884886fSBarry Smith PetscFunctionBegin; 2146bf464f9SBarry Smith ierr = VecDestroy(&ctx->w);CHKERRQ(ierr); 2156bf464f9SBarry Smith ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr); 2166bf464f9SBarry Smith ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr); 2176bf464f9SBarry Smith ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr); 218c51ad4d4SStefano Zampini ierr = VecDestroy(&ctx->dshiftw);CHKERRQ(ierr); 219c51ad4d4SStefano Zampini ierr = VecDestroy(&ctx->current_u);CHKERRQ(ierr); 220cfe22f5eSBarry Smith if (ctx->current_f_allocated) { 2216bf464f9SBarry Smith ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr); 222cfe22f5eSBarry Smith } 223e884886fSBarry Smith if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);} 2246bf464f9SBarry Smith ierr = PetscHeaderDestroy(&ctx);CHKERRQ(ierr); 225bf0cc555SLisandro Dalcin mat->data = 0; 226e884886fSBarry Smith 227bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr); 228bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr); 229bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr); 230bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr); 231bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr); 232bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr); 233bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr); 234bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr); 235e884886fSBarry Smith PetscFunctionReturn(0); 236e884886fSBarry Smith } 237e884886fSBarry Smith 238e884886fSBarry Smith /* 239e884886fSBarry Smith MatMFFDView_MFFD - Views matrix-free parameters. 240e884886fSBarry Smith 241e884886fSBarry Smith */ 24240244768SBarry Smith static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer) 243e884886fSBarry Smith { 244e884886fSBarry Smith PetscErrorCode ierr; 245e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 246a283635eSDmitry Karpeev PetscBool iascii, viewbase, viewfunction; 247a283635eSDmitry Karpeev const char *prefix; 248e884886fSBarry Smith 249e884886fSBarry Smith PetscFunctionBegin; 250251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 251e884886fSBarry Smith if (iascii) { 252a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr); 253a283635eSDmitry Karpeev ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 25457622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr); 2557adad957SLisandro Dalcin if (!((PetscObject)ctx)->type_name) { 256e884886fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr); 257e884886fSBarry Smith } else { 2587adad957SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr); 259e884886fSBarry Smith } 260e884886fSBarry Smith if (ctx->ops->view) { 261e884886fSBarry Smith ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr); 262e884886fSBarry Smith } 263a283635eSDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr); 264a283635eSDmitry Karpeev 265c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr); 266a283635eSDmitry Karpeev if (viewbase) { 267a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr); 268a283635eSDmitry Karpeev ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr); 269a283635eSDmitry Karpeev } 270c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr); 271a283635eSDmitry Karpeev if (viewfunction) { 272a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr); 273a283635eSDmitry Karpeev ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr); 274a283635eSDmitry Karpeev } 275a283635eSDmitry Karpeev ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 27611aeaf0aSBarry Smith } 277e884886fSBarry Smith PetscFunctionReturn(0); 278e884886fSBarry Smith } 279e884886fSBarry Smith 280e884886fSBarry Smith /* 281e884886fSBarry Smith MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 282e884886fSBarry Smith allows the user to indicate the beginning of a new linear solve by calling 283e884886fSBarry Smith MatAssemblyXXX() on the matrix free matrix. This then allows the 2841d0fab5eSBarry Smith MatCreateMFFD_WP() to properly compute ||U|| only the first time 285e884886fSBarry Smith in the linear solver rather than every time. 2865a576424SJed Brown 287cc2e6a90SBarry Smith This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence 288cc2e6a90SBarry Smith it must be labeled as PETSC_EXTERN 289e884886fSBarry Smith */ 2905a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt) 291e884886fSBarry Smith { 292e884886fSBarry Smith PetscErrorCode ierr; 293e884886fSBarry Smith MatMFFD j = (MatMFFD)J->data; 294e884886fSBarry Smith 295e884886fSBarry Smith PetscFunctionBegin; 296e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 297e884886fSBarry Smith j->vshift = 0.0; 298e884886fSBarry Smith j->vscale = 1.0; 299e884886fSBarry Smith PetscFunctionReturn(0); 300e884886fSBarry Smith } 301e884886fSBarry Smith 302e884886fSBarry Smith /* 303e884886fSBarry Smith MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 304e884886fSBarry Smith 305e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 306e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 307e884886fSBarry Smith u = current iterate 308e884886fSBarry Smith h = difference interval 309e884886fSBarry Smith */ 31040244768SBarry Smith static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y) 311e884886fSBarry Smith { 312e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 313e884886fSBarry Smith PetscScalar h; 314e884886fSBarry Smith Vec w,U,F; 315e884886fSBarry Smith PetscErrorCode ierr; 316ace3abfcSBarry Smith PetscBool zeroa; 317e884886fSBarry Smith 318e884886fSBarry Smith PetscFunctionBegin; 319ce94432eSBarry 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"); 320e884886fSBarry Smith /* We log matrix-free matrix-vector products separately, so that we can 321e884886fSBarry Smith separate the performance monitoring from the cases that use conventional 322e884886fSBarry Smith storage. We may eventually modify event logging to associate events 323e884886fSBarry Smith with particular objects, hence alleviating the more general problem. */ 324e884886fSBarry Smith ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 325e884886fSBarry Smith 326e884886fSBarry Smith w = ctx->w; 327e884886fSBarry Smith U = ctx->current_u; 3283ec795f1SBarry Smith F = ctx->current_f; 329e884886fSBarry Smith /* 330e884886fSBarry Smith Compute differencing parameter 331e884886fSBarry Smith */ 3324a2f8832SBarry Smith if (!((PetscObject)ctx)->type_name) { 333e884886fSBarry Smith ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr); 3349c6ac3b3SBarry Smith ierr = MatSetFromOptions(mat);CHKERRQ(ierr); 335e884886fSBarry Smith } 336e884886fSBarry Smith ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr); 337e884886fSBarry Smith if (zeroa) { 338e884886fSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 339e884886fSBarry Smith PetscFunctionReturn(0); 340e884886fSBarry Smith } 341e884886fSBarry Smith 34284d44b13SHong Zhang if (mat->erroriffailure && PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h"); 343e884886fSBarry Smith if (ctx->checkh) { 344e884886fSBarry Smith ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr); 345e884886fSBarry Smith } 346e884886fSBarry Smith 347e884886fSBarry Smith /* keep a record of the current differencing parameter h */ 348e884886fSBarry Smith ctx->currenth = h; 349e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX) 35057622a8eSBarry Smith ierr = PetscInfo2(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h));CHKERRQ(ierr); 351e884886fSBarry Smith #else 352e884886fSBarry Smith ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr); 353e884886fSBarry Smith #endif 354e884886fSBarry Smith if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 355e884886fSBarry Smith ctx->historyh[ctx->ncurrenth] = h; 356e884886fSBarry Smith } 357e884886fSBarry Smith ctx->ncurrenth++; 358e884886fSBarry Smith 359e884886fSBarry Smith /* w = u + ha */ 360c3bb7e23SBarry Smith if (ctx->drscale) { 361c3bb7e23SBarry Smith ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr); 362c3bb7e23SBarry Smith ierr = VecAYPX(U,h,w);CHKERRQ(ierr); 363c3bb7e23SBarry Smith } else { 364e884886fSBarry Smith ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr); 365c3bb7e23SBarry Smith } 366e884886fSBarry Smith 3671b797266SDmitry Karpeev /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 3681b797266SDmitry Karpeev if (ctx->ncurrenth == 1 && ctx->current_f_allocated) { 369e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr); 37052121784SBarry Smith } 371e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr); 372e884886fSBarry Smith 373e884886fSBarry Smith ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 374e884886fSBarry Smith ierr = VecScale(y,1.0/h);CHKERRQ(ierr); 375e884886fSBarry Smith 376c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 377e884886fSBarry Smith ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr); 378c3bb7e23SBarry Smith } 379c3bb7e23SBarry Smith if (ctx->dlscale) { 380c3bb7e23SBarry Smith ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr); 381c3bb7e23SBarry Smith } 382c3bb7e23SBarry Smith if (ctx->dshift) { 383c51ad4d4SStefano Zampini if (!ctx->dshiftw) { 384c51ad4d4SStefano Zampini ierr = VecDuplicate(y,&ctx->dshiftw);CHKERRQ(ierr); 385c51ad4d4SStefano Zampini } 386c51ad4d4SStefano Zampini ierr = VecPointwiseMult(ctx->dshift,a,ctx->dshiftw);CHKERRQ(ierr); 387c51ad4d4SStefano Zampini ierr = VecAXPY(y,1.0,ctx->dshiftw);CHKERRQ(ierr); 388c3bb7e23SBarry Smith } 389e884886fSBarry Smith 39039601f4eSBarry Smith if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);} 391e884886fSBarry Smith 392e884886fSBarry Smith ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 393e884886fSBarry Smith PetscFunctionReturn(0); 394e884886fSBarry Smith } 395e884886fSBarry Smith 396e884886fSBarry Smith /* 397e884886fSBarry Smith MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 398e884886fSBarry Smith 399e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 400e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 401e884886fSBarry Smith u = current iterate 402e884886fSBarry Smith h = difference interval 403e884886fSBarry Smith */ 4045d21e1f8SStefano Zampini PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 405e884886fSBarry Smith { 406e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 407e884886fSBarry Smith PetscScalar h,*aa,*ww,v; 408e884886fSBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 409e884886fSBarry Smith Vec w,U; 410e884886fSBarry Smith PetscErrorCode ierr; 411e884886fSBarry Smith PetscInt i,rstart,rend; 412e884886fSBarry Smith 413e884886fSBarry Smith PetscFunctionBegin; 414*486afcceSStefano Zampini if (!ctx->func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunction() first"); 415*486afcceSStefano Zampini if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioni() first"); 416*486afcceSStefano Zampini if (!ctx->funcisetbase) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioniBase() first"); 417e884886fSBarry Smith w = ctx->w; 418e884886fSBarry Smith U = ctx->current_u; 419e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr); 420e884886fSBarry Smith ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr); 421e884886fSBarry Smith ierr = VecCopy(U,w);CHKERRQ(ierr); 422e884886fSBarry Smith 423e884886fSBarry Smith ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr); 424e884886fSBarry Smith ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 425e884886fSBarry Smith for (i=rstart; i<rend; i++) { 426e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 427e884886fSBarry Smith h = ww[i-rstart]; 428e884886fSBarry Smith if (h == 0.0) h = 1.0; 429e884886fSBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 430e884886fSBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 431e884886fSBarry Smith h *= epsilon; 432e884886fSBarry Smith 433e884886fSBarry Smith ww[i-rstart] += h; 434e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 435e884886fSBarry Smith ierr = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr); 436e884886fSBarry Smith aa[i-rstart] = (v - aa[i-rstart])/h; 437e884886fSBarry Smith 438e884886fSBarry Smith /* possibly shift and scale result */ 439c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 440e884886fSBarry Smith aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart]; 441c3bb7e23SBarry Smith } 442e884886fSBarry Smith 443e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 444e884886fSBarry Smith ww[i-rstart] -= h; 445e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 446e884886fSBarry Smith } 447e884886fSBarry Smith ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 448e884886fSBarry Smith PetscFunctionReturn(0); 449e884886fSBarry Smith } 450e884886fSBarry Smith 45140244768SBarry Smith static PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr) 452c3bb7e23SBarry Smith { 453c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 454c3bb7e23SBarry Smith PetscErrorCode ierr; 455c3bb7e23SBarry Smith 456c3bb7e23SBarry Smith PetscFunctionBegin; 457c3bb7e23SBarry Smith if (ll && !aij->dlscale) { 458c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr); 459c3bb7e23SBarry Smith } 460c3bb7e23SBarry Smith if (rr && !aij->drscale) { 461c3bb7e23SBarry Smith ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr); 462c3bb7e23SBarry Smith } 463c3bb7e23SBarry Smith if (ll) { 464c3bb7e23SBarry Smith ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr); 465c3bb7e23SBarry Smith } 466c3bb7e23SBarry Smith if (rr) { 467c3bb7e23SBarry Smith ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr); 468c3bb7e23SBarry Smith } 469c3bb7e23SBarry Smith PetscFunctionReturn(0); 470c3bb7e23SBarry Smith } 471c3bb7e23SBarry Smith 47240244768SBarry Smith static PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode) 473c3bb7e23SBarry Smith { 474c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 475c3bb7e23SBarry Smith PetscErrorCode ierr; 476c3bb7e23SBarry Smith 477c3bb7e23SBarry Smith PetscFunctionBegin; 478ce94432eSBarry Smith if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES"); 479c3bb7e23SBarry Smith if (!aij->dshift) { 480c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr); 481c3bb7e23SBarry Smith } 482c3bb7e23SBarry Smith ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr); 483c3bb7e23SBarry Smith PetscFunctionReturn(0); 484c3bb7e23SBarry Smith } 485c3bb7e23SBarry Smith 48640244768SBarry Smith static PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a) 487e884886fSBarry Smith { 488e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 4895fd66863SKarl Rupp 490e884886fSBarry Smith PetscFunctionBegin; 491e884886fSBarry Smith shell->vshift += a; 492e884886fSBarry Smith PetscFunctionReturn(0); 493e884886fSBarry Smith } 494e884886fSBarry Smith 49540244768SBarry Smith static PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a) 496e884886fSBarry Smith { 497e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 4985fd66863SKarl Rupp 499e884886fSBarry Smith PetscFunctionBegin; 500e884886fSBarry Smith shell->vscale *= a; 501e884886fSBarry Smith PetscFunctionReturn(0); 502e884886fSBarry Smith } 503e884886fSBarry Smith 504d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F) 505e884886fSBarry Smith { 506e884886fSBarry Smith PetscErrorCode ierr; 507e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 508e884886fSBarry Smith 509e884886fSBarry Smith PetscFunctionBegin; 510e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 511c51ad4d4SStefano Zampini if (!ctx->current_u) { 512c51ad4d4SStefano Zampini ierr = VecDuplicate(U,&ctx->current_u);CHKERRQ(ierr); 513c51ad4d4SStefano Zampini ierr = VecLockPush(ctx->current_u);CHKERRQ(ierr); 514c51ad4d4SStefano Zampini } 515c51ad4d4SStefano Zampini ierr = VecLockPop(ctx->current_u);CHKERRQ(ierr); 516c51ad4d4SStefano Zampini ierr = VecCopy(U,ctx->current_u);CHKERRQ(ierr); 517c51ad4d4SStefano Zampini ierr = VecLockPush(ctx->current_u);CHKERRQ(ierr); 51852121784SBarry Smith if (F) { 5196bf464f9SBarry Smith if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);} 5203ec795f1SBarry Smith ctx->current_f = F; 521cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_FALSE; 522cfe22f5eSBarry Smith } else if (!ctx->current_f_allocated) { 52306c4b6baSBarry Smith ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr); 5242205254eSKarl Rupp 525cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_TRUE; 52652121784SBarry Smith } 527e884886fSBarry Smith if (!ctx->w) { 528e884886fSBarry Smith ierr = VecDuplicate(ctx->current_u,&ctx->w);CHKERRQ(ierr); 529e884886fSBarry Smith } 530e884886fSBarry Smith J->assembled = PETSC_TRUE; 531e884886fSBarry Smith PetscFunctionReturn(0); 532e884886fSBarry Smith } 5335a576424SJed Brown 534e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 535b2573a8aSBarry Smith 53640244768SBarry Smith static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx) 537e884886fSBarry Smith { 538e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 539e884886fSBarry Smith 540e884886fSBarry Smith PetscFunctionBegin; 541e884886fSBarry Smith ctx->checkh = fun; 542e884886fSBarry Smith ctx->checkhctx = ectx; 543e884886fSBarry Smith PetscFunctionReturn(0); 544e884886fSBarry Smith } 545e884886fSBarry Smith 5466aa9148fSLisandro Dalcin /*@C 5476aa9148fSLisandro Dalcin MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 5486aa9148fSLisandro Dalcin MatMFFD options in the database. 5496aa9148fSLisandro Dalcin 5506aa9148fSLisandro Dalcin Collective on Mat 5516aa9148fSLisandro Dalcin 5526aa9148fSLisandro Dalcin Input Parameter: 5536aa9148fSLisandro Dalcin + A - the Mat context 5546aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names 5556aa9148fSLisandro Dalcin 5566aa9148fSLisandro Dalcin Notes: 5576aa9148fSLisandro Dalcin A hyphen (-) must NOT be given at the beginning of the prefix name. 5586aa9148fSLisandro Dalcin The first character of all runtime options is AUTOMATICALLY the hyphen. 5596aa9148fSLisandro Dalcin 5606aa9148fSLisandro Dalcin Level: advanced 5616aa9148fSLisandro Dalcin 5626aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters 5636aa9148fSLisandro Dalcin 564755b7f64SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD() 5656aa9148fSLisandro Dalcin @*/ 5667087cfbeSBarry Smith PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[]) 5676aa9148fSLisandro Dalcin 5686aa9148fSLisandro Dalcin { 5690298fd71SBarry Smith MatMFFD mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL; 5706aa9148fSLisandro Dalcin PetscErrorCode ierr; 5715fd66863SKarl Rupp 5726aa9148fSLisandro Dalcin PetscFunctionBegin; 5730700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5740700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 5756aa9148fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr); 5766aa9148fSLisandro Dalcin PetscFunctionReturn(0); 5776aa9148fSLisandro Dalcin } 5786aa9148fSLisandro Dalcin 57940244768SBarry Smith static PetscErrorCode MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat) 580e884886fSBarry Smith { 58171f08665SBarry Smith MatMFFD mfctx = (MatMFFD)mat->data; 582e884886fSBarry Smith PetscErrorCode ierr; 583ace3abfcSBarry Smith PetscBool flg; 584e884886fSBarry Smith char ftype[256]; 585e884886fSBarry Smith 586e884886fSBarry Smith PetscFunctionBegin; 5870700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5880700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 5893194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr); 590a264d7a6SBarry Smith ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr); 591e884886fSBarry Smith if (flg) { 592e884886fSBarry Smith ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr); 593e884886fSBarry Smith } 594e884886fSBarry Smith 595e884886fSBarry Smith ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 596e884886fSBarry Smith ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 597e884886fSBarry Smith 59890d69ab7SBarry Smith flg = PETSC_FALSE; 5990298fd71SBarry Smith ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr); 600e884886fSBarry Smith if (flg) { 601e884886fSBarry Smith ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr); 602e884886fSBarry Smith } 603e884886fSBarry Smith if (mfctx->ops->setfromoptions) { 604e55864a3SBarry Smith ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr); 605e884886fSBarry Smith } 606e884886fSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 607e884886fSBarry Smith PetscFunctionReturn(0); 608e884886fSBarry Smith } 609e884886fSBarry Smith 61040244768SBarry Smith static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period) 611bc0f08ceSBarry Smith { 612bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 613bc0f08ceSBarry Smith 614bc0f08ceSBarry Smith PetscFunctionBegin; 615bc0f08ceSBarry Smith ctx->recomputeperiod = period; 616bc0f08ceSBarry Smith PetscFunctionReturn(0); 617bc0f08ceSBarry Smith } 618bc0f08ceSBarry Smith 61940244768SBarry Smith static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 620bc0f08ceSBarry Smith { 621bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 622bc0f08ceSBarry Smith 623bc0f08ceSBarry Smith PetscFunctionBegin; 624bc0f08ceSBarry Smith ctx->func = func; 625bc0f08ceSBarry Smith ctx->funcctx = funcctx; 626bc0f08ceSBarry Smith PetscFunctionReturn(0); 627bc0f08ceSBarry Smith } 628bc0f08ceSBarry Smith 62940244768SBarry Smith static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error) 630bc0f08ceSBarry Smith { 631bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 632bc0f08ceSBarry Smith 633bc0f08ceSBarry Smith PetscFunctionBegin; 634bc0f08ceSBarry Smith if (error != PETSC_DEFAULT) ctx->error_rel = error; 635bc0f08ceSBarry Smith PetscFunctionReturn(0); 636bc0f08ceSBarry Smith } 637bc0f08ceSBarry Smith 6383b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MFFD(Mat A,PetscBool *missing,PetscInt *d) 6393b49f96aSBarry Smith { 6403b49f96aSBarry Smith PetscFunctionBegin; 6413b49f96aSBarry Smith *missing = PETSC_FALSE; 6423b49f96aSBarry Smith PetscFunctionReturn(0); 6433b49f96aSBarry Smith } 6443b49f96aSBarry Smith 645e884886fSBarry Smith /*MC 646e884886fSBarry Smith MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 647e884886fSBarry Smith 648e884886fSBarry Smith Level: advanced 649e884886fSBarry Smith 650755b7f64SBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(), 651755b7f64SBarry Smith MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 652755b7f64SBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 653755b7f64SBarry Smith MatMFFDGetH(), 654e884886fSBarry Smith M*/ 6558cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 656e884886fSBarry Smith { 657e884886fSBarry Smith MatMFFD mfctx; 658e884886fSBarry Smith PetscErrorCode ierr; 659e884886fSBarry Smith 660e884886fSBarry Smith PetscFunctionBegin; 661607a6623SBarry Smith ierr = MatMFFDInitializePackage();CHKERRQ(ierr); 662e884886fSBarry Smith 66373107ff1SLisandro Dalcin ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr); 6642205254eSKarl Rupp 665e884886fSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 666e884886fSBarry Smith mfctx->recomputeperiod = 1; 667e884886fSBarry Smith mfctx->count = 0; 668e884886fSBarry Smith mfctx->currenth = 0.0; 6690298fd71SBarry Smith mfctx->historyh = NULL; 670e884886fSBarry Smith mfctx->ncurrenth = 0; 671e884886fSBarry Smith mfctx->maxcurrenth = 0; 6727adad957SLisandro Dalcin ((PetscObject)mfctx)->type_name = 0; 673e884886fSBarry Smith 674e884886fSBarry Smith mfctx->vshift = 0.0; 675e884886fSBarry Smith mfctx->vscale = 1.0; 676e884886fSBarry Smith 677e884886fSBarry Smith /* 678e884886fSBarry Smith Create the empty data structure to contain compute-h routines. 679e884886fSBarry Smith These will be filled in below from the command line options or 680e884886fSBarry Smith a later call with MatMFFDSetType() or if that is not called 681e884886fSBarry Smith then it will default in the first use of MatMult_MFFD() 682e884886fSBarry Smith */ 683e884886fSBarry Smith mfctx->ops->compute = 0; 684e884886fSBarry Smith mfctx->ops->destroy = 0; 685e884886fSBarry Smith mfctx->ops->view = 0; 686e884886fSBarry Smith mfctx->ops->setfromoptions = 0; 687e884886fSBarry Smith mfctx->hctx = 0; 688e884886fSBarry Smith 689e884886fSBarry Smith mfctx->func = 0; 690e884886fSBarry Smith mfctx->funcctx = 0; 6910298fd71SBarry Smith mfctx->w = NULL; 692e884886fSBarry Smith 693e884886fSBarry Smith A->data = mfctx; 694e884886fSBarry Smith 695e884886fSBarry Smith A->ops->mult = MatMult_MFFD; 696e884886fSBarry Smith A->ops->destroy = MatDestroy_MFFD; 697e884886fSBarry Smith A->ops->view = MatView_MFFD; 698e884886fSBarry Smith A->ops->assemblyend = MatAssemblyEnd_MFFD; 699e884886fSBarry Smith A->ops->scale = MatScale_MFFD; 700e884886fSBarry Smith A->ops->shift = MatShift_MFFD; 701c3bb7e23SBarry Smith A->ops->diagonalscale = MatDiagonalScale_MFFD; 702c3bb7e23SBarry Smith A->ops->diagonalset = MatDiagonalSet_MFFD; 7039c6ac3b3SBarry Smith A->ops->setfromoptions = MatSetFromOptions_MFFD; 7043b49f96aSBarry Smith A->ops->missingdiagonal = MatMissingDiagonal_MFFD; 705e884886fSBarry Smith A->assembled = PETSC_TRUE; 706e884886fSBarry Smith 707bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr); 708bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr); 709bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr); 710bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr); 711bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr); 712bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr); 713bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr); 714bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr); 7152205254eSKarl Rupp 716e884886fSBarry Smith mfctx->mat = A; 7172205254eSKarl Rupp 718e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr); 719e884886fSBarry Smith PetscFunctionReturn(0); 720e884886fSBarry Smith } 721e884886fSBarry Smith 722e884886fSBarry Smith /*@ 723e884886fSBarry Smith MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 724e884886fSBarry Smith 725e884886fSBarry Smith Collective on Vec 726e884886fSBarry Smith 727e884886fSBarry Smith Input Parameters: 728fef1beadSBarry Smith + comm - MPI communicator 729fef1beadSBarry Smith . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 730fef1beadSBarry Smith This value should be the same as the local size used in creating the 731fef1beadSBarry Smith y vector for the matrix-vector product y = Ax. 732fef1beadSBarry Smith . n - This value should be the same as the local size used in creating the 733fef1beadSBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 734fef1beadSBarry Smith calculated if N is given) For square matrices n is almost always m. 735fef1beadSBarry Smith . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 736fef1beadSBarry Smith - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 737fef1beadSBarry Smith 738e884886fSBarry Smith 739e884886fSBarry Smith Output Parameter: 740e884886fSBarry Smith . J - the matrix-free matrix 741e884886fSBarry Smith 7429c6ac3b3SBarry Smith Options Database Keys: call MatSetFromOptions() to trigger these 7439c6ac3b3SBarry Smith + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 7449c6ac3b3SBarry Smith - -mat_mffd_err - square root of estimated relative error in function evaluation 7459c6ac3b3SBarry Smith - -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 7469c6ac3b3SBarry Smith 7479c6ac3b3SBarry Smith 748e884886fSBarry Smith Level: advanced 749e884886fSBarry Smith 750e884886fSBarry Smith Notes: 751e884886fSBarry Smith The matrix-free matrix context merely contains the function pointers 752e884886fSBarry Smith and work space for performing finite difference approximations of 753e884886fSBarry Smith Jacobian-vector products, F'(u)*a, 754e884886fSBarry Smith 755e884886fSBarry Smith The default code uses the following approach to compute h 756e884886fSBarry Smith 757e884886fSBarry Smith .vb 758e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 759e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 760e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 761e884886fSBarry Smith where 762e884886fSBarry Smith error_rel = square root of relative error in function evaluation 763e884886fSBarry Smith umin = minimum iterate parameter 764e884886fSBarry Smith .ve 765e884886fSBarry Smith 766ca93e954SBarry Smith You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different 767ca93e954SBarry Smith preconditioner matrix 768ca93e954SBarry Smith 769e884886fSBarry Smith The user can set the error_rel via MatMFFDSetFunctionError() and 770a7f22e61SSatish Balay umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details. 771e884886fSBarry Smith 772e884886fSBarry Smith The user should call MatDestroy() when finished with the matrix-free 773e884886fSBarry Smith matrix context. 774e884886fSBarry Smith 775e884886fSBarry Smith Options Database Keys: 776e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 777e884886fSBarry Smith . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only) 778e884886fSBarry Smith - -mat_mffd_check_positivity 779e884886fSBarry Smith 780e884886fSBarry Smith .keywords: default, matrix-free, create, matrix 781e884886fSBarry Smith 7821d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 783e884886fSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 78481242352SJed Brown MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian() 785e884886fSBarry Smith 786e884886fSBarry Smith @*/ 7877087cfbeSBarry Smith PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 788e884886fSBarry Smith { 789e884886fSBarry Smith PetscErrorCode ierr; 790e884886fSBarry Smith 791e884886fSBarry Smith PetscFunctionBegin; 792e884886fSBarry Smith ierr = MatCreate(comm,J);CHKERRQ(ierr); 793fef1beadSBarry Smith ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr); 794e884886fSBarry Smith ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 795be7c243fSJed Brown ierr = MatSetUp(*J);CHKERRQ(ierr); 796e884886fSBarry Smith PetscFunctionReturn(0); 797e884886fSBarry Smith } 798e884886fSBarry Smith 799e884886fSBarry Smith /*@ 800e884886fSBarry Smith MatMFFDGetH - Gets the last value that was used as the differencing 801e884886fSBarry Smith parameter. 802e884886fSBarry Smith 803e884886fSBarry Smith Not Collective 804e884886fSBarry Smith 805e884886fSBarry Smith Input Parameters: 806e884886fSBarry Smith . mat - the matrix obtained with MatCreateSNESMF() 807e884886fSBarry Smith 808e884886fSBarry Smith Output Paramter: 809e884886fSBarry Smith . h - the differencing step size 810e884886fSBarry Smith 811e884886fSBarry Smith Level: advanced 812e884886fSBarry Smith 813e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 814e884886fSBarry Smith 8151d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory() 816e884886fSBarry Smith @*/ 8177087cfbeSBarry Smith PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 818e884886fSBarry Smith { 819e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 820bc0f08ceSBarry Smith PetscErrorCode ierr; 821bc0f08ceSBarry Smith PetscBool match; 822e884886fSBarry Smith 823e884886fSBarry Smith PetscFunctionBegin; 82488b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 82588b4c220SStefano Zampini PetscValidPointer(h,2); 826251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 827ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 828bc0f08ceSBarry Smith 829e884886fSBarry Smith *h = ctx->currenth; 830e884886fSBarry Smith PetscFunctionReturn(0); 831e884886fSBarry Smith } 832e884886fSBarry Smith 833e884886fSBarry Smith /*@C 834e884886fSBarry Smith MatMFFDSetFunction - Sets the function used in applying the matrix free. 835e884886fSBarry Smith 8363f9fe445SBarry Smith Logically Collective on Mat 837e884886fSBarry Smith 838e884886fSBarry Smith Input Parameters: 83914f633e2SBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD() 840e884886fSBarry Smith . func - the function to use 841e884886fSBarry Smith - funcctx - optional function context passed to function 842e884886fSBarry Smith 84314f633e2SBarry Smith Calling Sequence of func: 84414f633e2SBarry Smith $ func (void *funcctx, Vec x, Vec f) 84514f633e2SBarry Smith 84614f633e2SBarry Smith + funcctx - user provided context 84714f633e2SBarry Smith . x - input vector 84814f633e2SBarry Smith - f - computed output function 84914f633e2SBarry Smith 850e884886fSBarry Smith Level: advanced 851e884886fSBarry Smith 852e884886fSBarry Smith Notes: 853e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 854e884886fSBarry Smith matrix inside your compute Jacobian routine 855e884886fSBarry Smith 8563ec795f1SBarry Smith If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 857e884886fSBarry Smith 858e884886fSBarry Smith .keywords: SNES, matrix-free, function 859e884886fSBarry Smith 8601d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD, 8611d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 862e884886fSBarry Smith @*/ 8637087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 864e884886fSBarry Smith { 865bc0f08ceSBarry Smith PetscErrorCode ierr; 866e884886fSBarry Smith 867e884886fSBarry Smith PetscFunctionBegin; 86888b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 869bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr); 870e884886fSBarry Smith PetscFunctionReturn(0); 871e884886fSBarry Smith } 872e884886fSBarry Smith 873e884886fSBarry Smith /*@C 874e884886fSBarry Smith MatMFFDSetFunctioni - Sets the function for a single component 875e884886fSBarry Smith 8763f9fe445SBarry Smith Logically Collective on Mat 877e884886fSBarry Smith 878e884886fSBarry Smith Input Parameters: 879e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 880e884886fSBarry Smith - funci - the function to use 881e884886fSBarry Smith 882e884886fSBarry Smith Level: advanced 883e884886fSBarry Smith 884e884886fSBarry Smith Notes: 885e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 88694eb4328SStefano Zampini matrix inside your compute Jacobian routine. 88794eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 888*486afcceSStefano Zampini funci must not contain any MPI call as it is called inside a loop on the local portion of the vector. 889e884886fSBarry Smith 890e884886fSBarry Smith .keywords: SNES, matrix-free, function 891e884886fSBarry Smith 89294eb4328SStefano Zampini .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 8931d0fab5eSBarry Smith 894e884886fSBarry Smith @*/ 8957087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 896e884886fSBarry Smith { 8974ac538c5SBarry Smith PetscErrorCode ierr; 898e884886fSBarry Smith 899e884886fSBarry Smith PetscFunctionBegin; 9000700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 9014ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr); 902e884886fSBarry Smith PetscFunctionReturn(0); 903e884886fSBarry Smith } 904e884886fSBarry Smith 905e884886fSBarry Smith /*@C 906e884886fSBarry Smith MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 907e884886fSBarry Smith 9083f9fe445SBarry Smith Logically Collective on Mat 909e884886fSBarry Smith 910e884886fSBarry Smith Input Parameters: 911e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 912e884886fSBarry Smith - func - the function to use 913e884886fSBarry Smith 914e884886fSBarry Smith Level: advanced 915e884886fSBarry Smith 916e884886fSBarry Smith Notes: 917e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 91894eb4328SStefano Zampini matrix inside your compute Jacobian routine. 91994eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 920e884886fSBarry Smith 921e884886fSBarry Smith 922e884886fSBarry Smith .keywords: SNES, matrix-free, function 923e884886fSBarry Smith 924e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 92594eb4328SStefano Zampini MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 926e884886fSBarry Smith @*/ 9277087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 928e884886fSBarry Smith { 9294ac538c5SBarry Smith PetscErrorCode ierr; 930e884886fSBarry Smith 931e884886fSBarry Smith PetscFunctionBegin; 9320700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 9334ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr); 934e884886fSBarry Smith PetscFunctionReturn(0); 935e884886fSBarry Smith } 936e884886fSBarry Smith 937e884886fSBarry Smith /*@ 938e884886fSBarry Smith MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 939e884886fSBarry Smith 9403f9fe445SBarry Smith Logically Collective on Mat 941e884886fSBarry Smith 942e884886fSBarry Smith Input Parameters: 943e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 944e884886fSBarry Smith - period - 1 for everytime, 2 for every second etc 945e884886fSBarry Smith 946e884886fSBarry Smith Options Database Keys: 947e884886fSBarry Smith + -mat_mffd_period <period> 948e884886fSBarry Smith 949e884886fSBarry Smith Level: advanced 950e884886fSBarry Smith 951e884886fSBarry Smith 952e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 953e884886fSBarry Smith 954e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), 9551d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 956e884886fSBarry Smith @*/ 9577087cfbeSBarry Smith PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 958e884886fSBarry Smith { 959bc0f08ceSBarry Smith PetscErrorCode ierr; 960e884886fSBarry Smith 961e884886fSBarry Smith PetscFunctionBegin; 96288b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 96388b4c220SStefano Zampini PetscValidLogicalCollectiveInt(mat,period,2); 964bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr); 965e884886fSBarry Smith PetscFunctionReturn(0); 966e884886fSBarry Smith } 967e884886fSBarry Smith 968e884886fSBarry Smith /*@ 969e884886fSBarry Smith MatMFFDSetFunctionError - Sets the error_rel for the approximation of 970e884886fSBarry Smith matrix-vector products using finite differences. 971e884886fSBarry Smith 9723f9fe445SBarry Smith Logically Collective on Mat 973e884886fSBarry Smith 974e884886fSBarry Smith Input Parameters: 975e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 976e884886fSBarry Smith - error_rel - relative error (should be set to the square root of 977e884886fSBarry Smith the relative error in the function evaluations) 978e884886fSBarry Smith 979e884886fSBarry Smith Options Database Keys: 980e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 981e884886fSBarry Smith 982e884886fSBarry Smith Level: advanced 983e884886fSBarry Smith 984e884886fSBarry Smith Notes: 985e884886fSBarry Smith The default matrix-free matrix-vector product routine computes 986e884886fSBarry Smith .vb 987e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 988e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 989e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 990e884886fSBarry Smith .ve 991e884886fSBarry Smith 992e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 993e884886fSBarry Smith 994e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 9951d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 996e884886fSBarry Smith @*/ 9977087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 998e884886fSBarry Smith { 999bc0f08ceSBarry Smith PetscErrorCode ierr; 1000e884886fSBarry Smith 1001e884886fSBarry Smith PetscFunctionBegin; 100288b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 100388b4c220SStefano Zampini PetscValidLogicalCollectiveReal(mat,error,2); 1004bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr); 1005e884886fSBarry Smith PetscFunctionReturn(0); 1006e884886fSBarry Smith } 1007e884886fSBarry Smith 1008e884886fSBarry Smith /*@ 1009e884886fSBarry Smith MatMFFDSetHHistory - Sets an array to collect a history of the 1010e884886fSBarry Smith differencing values (h) computed for the matrix-free product. 1011e884886fSBarry Smith 10123f9fe445SBarry Smith Logically Collective on Mat 1013e884886fSBarry Smith 1014e884886fSBarry Smith Input Parameters: 1015e884886fSBarry Smith + J - the matrix-free matrix context 1016e884886fSBarry Smith . histroy - space to hold the history 1017e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than 1018e884886fSBarry Smith nhistory, then the later ones are discarded 1019e884886fSBarry Smith 1020e884886fSBarry Smith Level: advanced 1021e884886fSBarry Smith 1022e884886fSBarry Smith Notes: 1023e884886fSBarry Smith Use MatMFFDResetHHistory() to reset the history counter and collect 1024e884886fSBarry Smith a new batch of differencing parameters, h. 1025e884886fSBarry Smith 1026e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1027e884886fSBarry Smith 1028e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 10291d0fab5eSBarry Smith MatMFFDResetHHistory(), MatMFFDSetFunctionError() 1030e884886fSBarry Smith 1031e884886fSBarry Smith @*/ 10327087cfbeSBarry Smith PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 1033e884886fSBarry Smith { 1034e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 1035bc0f08ceSBarry Smith PetscErrorCode ierr; 1036bc0f08ceSBarry Smith PetscBool match; 1037e884886fSBarry Smith 1038e884886fSBarry Smith PetscFunctionBegin; 103988b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 104088b4c220SStefano Zampini if (history) PetscValidPointer(history,2); 104188b4c220SStefano Zampini PetscValidLogicalCollectiveInt(J,nhistory,3); 1042251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr); 1043ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 1044e884886fSBarry Smith ctx->historyh = history; 1045e884886fSBarry Smith ctx->maxcurrenth = nhistory; 104675567043SBarry Smith ctx->currenth = 0.; 1047e884886fSBarry Smith PetscFunctionReturn(0); 1048e884886fSBarry Smith } 1049e884886fSBarry Smith 1050e884886fSBarry Smith /*@ 1051e884886fSBarry Smith MatMFFDResetHHistory - Resets the counter to zero to begin 1052e884886fSBarry Smith collecting a new set of differencing histories. 1053e884886fSBarry Smith 10543f9fe445SBarry Smith Logically Collective on Mat 1055e884886fSBarry Smith 1056e884886fSBarry Smith Input Parameters: 1057e884886fSBarry Smith . J - the matrix-free matrix context 1058e884886fSBarry Smith 1059e884886fSBarry Smith Level: advanced 1060e884886fSBarry Smith 1061e884886fSBarry Smith Notes: 1062e884886fSBarry Smith Use MatMFFDSetHHistory() to create the original history counter. 1063e884886fSBarry Smith 1064e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1065e884886fSBarry Smith 1066e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 10671d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDSetFunctionError() 1068e884886fSBarry Smith 1069e884886fSBarry Smith @*/ 10707087cfbeSBarry Smith PetscErrorCode MatMFFDResetHHistory(Mat J) 1071e884886fSBarry Smith { 1072bc0f08ceSBarry Smith PetscErrorCode ierr; 1073e884886fSBarry Smith 1074e884886fSBarry Smith PetscFunctionBegin; 107588b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1076bc0f08ceSBarry Smith ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr); 1077e884886fSBarry Smith PetscFunctionReturn(0); 1078e884886fSBarry Smith } 1079e884886fSBarry Smith 1080e884886fSBarry Smith /*@ 1081e884886fSBarry Smith MatMFFDSetBase - Sets the vector U at which matrix vector products of the 1082e884886fSBarry Smith Jacobian are computed 1083e884886fSBarry Smith 10843f9fe445SBarry Smith Logically Collective on Mat 1085e884886fSBarry Smith 1086e884886fSBarry Smith Input Parameters: 1087e884886fSBarry Smith + J - the MatMFFD matrix 10883ec795f1SBarry Smith . U - the vector 1089bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed 1090e884886fSBarry Smith 1091e884886fSBarry Smith Notes: This is rarely used directly 1092e884886fSBarry Smith 10938af5ae88SBarry Smith If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 10948af5ae88SBarry Smith point during the first MatMult() after each call to MatMFFDSetBase(). 1095dff2f722SBarry Smith 1096e884886fSBarry Smith Level: advanced 1097e884886fSBarry Smith 1098e884886fSBarry Smith @*/ 10997087cfbeSBarry Smith PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 1100e884886fSBarry Smith { 11014ac538c5SBarry Smith PetscErrorCode ierr; 1102e884886fSBarry Smith 1103e884886fSBarry Smith PetscFunctionBegin; 11040700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 11050700a824SBarry Smith PetscValidHeaderSpecific(U,VEC_CLASSID,2); 11060700a824SBarry Smith if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 11074ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr); 1108e884886fSBarry Smith PetscFunctionReturn(0); 1109e884886fSBarry Smith } 1110e884886fSBarry Smith 1111e884886fSBarry Smith /*@C 1112e884886fSBarry Smith MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 1113e884886fSBarry Smith it to satisfy some criteria 1114e884886fSBarry Smith 11153f9fe445SBarry Smith Logically Collective on Mat 1116e884886fSBarry Smith 1117e884886fSBarry Smith Input Parameters: 1118e884886fSBarry Smith + J - the MatMFFD matrix 1119e884886fSBarry Smith . fun - the function that checks h 1120e884886fSBarry Smith - ctx - any context needed by the function 1121e884886fSBarry Smith 1122e884886fSBarry Smith Options Database Keys: 1123e884886fSBarry Smith . -mat_mffd_check_positivity 1124e884886fSBarry Smith 1125e884886fSBarry Smith Level: advanced 1126e884886fSBarry Smith 1127755b7f64SBarry Smith Notes: For example, MatMFFDCheckPositivity() insures that all entries 1128e884886fSBarry Smith of U + h*a are non-negative 1129e884886fSBarry Smith 1130755b7f64SBarry Smith The function you provide is called after the default h has been computed and allows you to 1131755b7f64SBarry Smith modify it. 1132755b7f64SBarry Smith 1133755b7f64SBarry Smith .seealso: MatMFFDCheckPositivity() 1134e884886fSBarry Smith @*/ 11357087cfbeSBarry Smith PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx) 1136e884886fSBarry Smith { 11374ac538c5SBarry Smith PetscErrorCode ierr; 1138e884886fSBarry Smith 1139e884886fSBarry Smith PetscFunctionBegin; 11400700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 11414ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr); 1142e884886fSBarry Smith PetscFunctionReturn(0); 1143e884886fSBarry Smith } 1144e884886fSBarry Smith 1145e884886fSBarry Smith /*@ 1146e884886fSBarry Smith MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1147e884886fSBarry Smith zero, decreases h until this is satisfied. 1148e884886fSBarry Smith 11493f9fe445SBarry Smith Logically Collective on Vec 1150e884886fSBarry Smith 1151e884886fSBarry Smith Input Parameters: 1152e884886fSBarry Smith + U - base vector that is added to 1153e884886fSBarry Smith . a - vector that is added 1154e884886fSBarry Smith . h - scaling factor on a 1155e884886fSBarry Smith - dummy - context variable (unused) 1156e884886fSBarry Smith 1157e884886fSBarry Smith Options Database Keys: 1158e884886fSBarry Smith . -mat_mffd_check_positivity 1159e884886fSBarry Smith 1160e884886fSBarry Smith Level: advanced 1161e884886fSBarry Smith 1162e884886fSBarry Smith Notes: This is rarely used directly, rather it is passed as an argument to 1163e884886fSBarry Smith MatMFFDSetCheckh() 1164e884886fSBarry Smith 1165e884886fSBarry Smith .seealso: MatMFFDSetCheckh() 1166e884886fSBarry Smith @*/ 11677087cfbeSBarry Smith PetscErrorCode MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h) 1168e884886fSBarry Smith { 1169e884886fSBarry Smith PetscReal val, minval; 1170e884886fSBarry Smith PetscScalar *u_vec, *a_vec; 1171e884886fSBarry Smith PetscErrorCode ierr; 1172e884886fSBarry Smith PetscInt i,n; 1173e884886fSBarry Smith MPI_Comm comm; 1174e884886fSBarry Smith 1175e884886fSBarry Smith PetscFunctionBegin; 117688b4c220SStefano Zampini PetscValidHeaderSpecific(U,VEC_CLASSID,2); 117788b4c220SStefano Zampini PetscValidHeaderSpecific(a,VEC_CLASSID,3); 117888b4c220SStefano Zampini PetscValidPointer(h,4); 1179e884886fSBarry Smith ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); 1180e884886fSBarry Smith ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); 1181e884886fSBarry Smith ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); 1182e884886fSBarry Smith ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); 1183e884886fSBarry Smith minval = PetscAbsScalar(*h*1.01); 1184e884886fSBarry Smith for (i=0; i<n; i++) { 1185e884886fSBarry Smith if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1186e884886fSBarry Smith val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1187e884886fSBarry Smith if (val < minval) minval = val; 1188e884886fSBarry Smith } 1189e884886fSBarry Smith } 1190e884886fSBarry Smith ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); 1191e884886fSBarry Smith ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); 1192b2566f29SBarry Smith ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr); 1193e884886fSBarry Smith if (val <= PetscAbsScalar(*h)) { 11946712e2f1SBarry Smith ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr); 1195e884886fSBarry Smith if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1196e884886fSBarry Smith else *h = -0.99*val; 1197e884886fSBarry Smith } 1198e884886fSBarry Smith PetscFunctionReturn(0); 1199e884886fSBarry Smith } 1200