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 348a690491SBarry Smith from MatInitializePackage(). 353ec795f1SBarry Smith 363ec795f1SBarry Smith Level: developer 373ec795f1SBarry Smith 383ec795f1SBarry Smith .keywords: Vec, initialize, package 393ec795f1SBarry Smith .seealso: PetscInitialize() 403ec795f1SBarry Smith @*/ 41607a6623SBarry Smith PetscErrorCode MatMFFDInitializePackage(void) 423ec795f1SBarry Smith { 433ec795f1SBarry Smith char logList[256]; 448e81d068SLisandro Dalcin PetscBool opt,pkg; 453ec795f1SBarry Smith PetscErrorCode ierr; 463ec795f1SBarry Smith 473ec795f1SBarry Smith PetscFunctionBegin; 48b022a5c1SBarry Smith if (MatMFFDPackageInitialized) PetscFunctionReturn(0); 49b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_TRUE; 503ec795f1SBarry Smith /* Register Classes */ 510700a824SBarry Smith ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr); 523ec795f1SBarry Smith /* Register Constructors */ 53607a6623SBarry Smith ierr = MatMFFDRegisterAll();CHKERRQ(ierr); 543ec795f1SBarry Smith /* Register Events */ 550700a824SBarry Smith ierr = PetscLogEventRegister("MatMult MF",MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr); 563ec795f1SBarry Smith /* Process info exclusions */ 578e81d068SLisandro Dalcin ierr = PetscOptionsGetString(NULL,NULL,"-info_exclude",logList,sizeof(logList),&opt);CHKERRQ(ierr); 583ec795f1SBarry Smith if (opt) { 598e81d068SLisandro Dalcin ierr = PetscStrInList("matmffd",logList,',',&pkg);CHKERRQ(ierr); 608e81d068SLisandro Dalcin if (pkg) {ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);} 613ec795f1SBarry Smith } 623ec795f1SBarry Smith /* Process summary exclusions */ 638e81d068SLisandro Dalcin ierr = PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt);CHKERRQ(ierr); 643ec795f1SBarry Smith if (opt) { 658e81d068SLisandro Dalcin ierr = PetscStrInList("matmffd",logList,',',&pkg);CHKERRQ(ierr); 66fa2bb9feSLisandro Dalcin if (pkg) {ierr = PetscLogEventExcludeClass(MATMFFD_CLASSID);CHKERRQ(ierr);} 673ec795f1SBarry Smith } 688e81d068SLisandro Dalcin /* Register package finalizer */ 69b022a5c1SBarry Smith ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr); 703ec795f1SBarry Smith PetscFunctionReturn(0); 713ec795f1SBarry Smith } 723ec795f1SBarry Smith 73e884886fSBarry Smith /*@C 74e884886fSBarry Smith MatMFFDSetType - Sets the method that is used to compute the 75e884886fSBarry Smith differencing parameter for finite differene matrix-free formulations. 76e884886fSBarry Smith 77e884886fSBarry Smith Input Parameters: 78e884886fSBarry Smith + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD() 79e884886fSBarry Smith or MatSetType(mat,MATMFFD); 80e884886fSBarry Smith - ftype - the type requested, either MATMFFD_WP or MATMFFD_DS 81e884886fSBarry Smith 82e884886fSBarry Smith Level: advanced 83e884886fSBarry Smith 84e884886fSBarry Smith Notes: 85e884886fSBarry Smith For example, such routines can compute h for use in 86e884886fSBarry Smith Jacobian-vector products of the form 87e884886fSBarry Smith 88e884886fSBarry Smith F(x+ha) - F(x) 89e884886fSBarry Smith F'(u)a ~= ---------------- 90e884886fSBarry Smith h 91e884886fSBarry Smith 92755b7f64SBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction(), MatCreateMFFD() 93e884886fSBarry Smith @*/ 9419fd82e9SBarry Smith PetscErrorCode MatMFFDSetType(Mat mat,MatMFFDType ftype) 95e884886fSBarry Smith { 96e884886fSBarry Smith PetscErrorCode ierr,(*r)(MatMFFD); 97e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 98ace3abfcSBarry Smith PetscBool match; 99e884886fSBarry Smith 100e884886fSBarry Smith PetscFunctionBegin; 1010700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 102e884886fSBarry Smith PetscValidCharPointer(ftype,2); 103e884886fSBarry Smith 104251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 1059a13eae7SJed Brown if (!match) PetscFunctionReturn(0); 1069a13eae7SJed Brown 107e884886fSBarry Smith /* already set, so just return */ 108251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr); 109e884886fSBarry Smith if (match) PetscFunctionReturn(0); 110e884886fSBarry Smith 111e884886fSBarry Smith /* destroy the old one if it exists */ 112e884886fSBarry Smith if (ctx->ops->destroy) { 113e884886fSBarry Smith ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr); 114e884886fSBarry Smith } 115e884886fSBarry Smith 1161c9cd337SJed Brown ierr = PetscFunctionListFind(MatMFFDList,ftype,&r);CHKERRQ(ierr); 117e32f2f54SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype); 118e884886fSBarry Smith ierr = (*r)(ctx);CHKERRQ(ierr); 119e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr); 120e884886fSBarry Smith PetscFunctionReturn(0); 121e884886fSBarry Smith } 122e884886fSBarry Smith 1235d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat,Vec); 1245d21e1f8SStefano Zampini 125e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/ 12640244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func) 127e884886fSBarry Smith { 128e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 129e884886fSBarry Smith 130e884886fSBarry Smith PetscFunctionBegin; 131e884886fSBarry Smith ctx->funcisetbase = func; 132486afcceSStefano Zampini /* allow users to compose their own getdiagonal and allow MatHasOperation 133486afcceSStefano Zampini to return false if the two functions pointers are not set */ 134486afcceSStefano Zampini if (!mat->ops->getdiagonal && func) { 135dbfa1270SStefano Zampini mat->ops->getdiagonal = MatGetDiagonal_MFFD; 1365d21e1f8SStefano Zampini } 137e884886fSBarry Smith PetscFunctionReturn(0); 138e884886fSBarry Smith } 139e884886fSBarry Smith 140e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 14140244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci) 142e884886fSBarry Smith { 143e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 144e884886fSBarry Smith 145e884886fSBarry Smith PetscFunctionBegin; 146e884886fSBarry Smith ctx->funci = funci; 147486afcceSStefano Zampini /* allow users to compose their own getdiagonal and allow MatHasOperation 148486afcceSStefano Zampini to return false if the two functions pointers are not set */ 149486afcceSStefano Zampini if (!mat->ops->getdiagonal && funci) { 150dbfa1270SStefano Zampini mat->ops->getdiagonal = MatGetDiagonal_MFFD; 1515d21e1f8SStefano Zampini } 152e884886fSBarry Smith PetscFunctionReturn(0); 153e884886fSBarry Smith } 154e884886fSBarry Smith 15540244768SBarry Smith static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) 156bc0f08ceSBarry Smith { 157bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 158bc0f08ceSBarry Smith 159bc0f08ceSBarry Smith PetscFunctionBegin; 160bc0f08ceSBarry Smith ctx->ncurrenth = 0; 161bc0f08ceSBarry Smith PetscFunctionReturn(0); 162bc0f08ceSBarry Smith } 163e884886fSBarry Smith 1641c84c290SBarry Smith /*@C 1651c84c290SBarry Smith MatMFFDRegister - Adds a method to the MatMFFD registry. 1661c84c290SBarry Smith 1671c84c290SBarry Smith Not Collective 1681c84c290SBarry Smith 1691c84c290SBarry Smith Input Parameters: 1701c84c290SBarry Smith + name_solver - name of a new user-defined compute-h module 1711c84c290SBarry Smith - routine_create - routine to create method context 1721c84c290SBarry Smith 1731c84c290SBarry Smith Level: developer 1741c84c290SBarry Smith 1751c84c290SBarry Smith Notes: 1761c84c290SBarry Smith MatMFFDRegister() may be called multiple times to add several user-defined solvers. 1771c84c290SBarry Smith 1781c84c290SBarry Smith Sample usage: 1791c84c290SBarry Smith .vb 180bdf89e91SBarry Smith MatMFFDRegister("my_h",MyHCreate); 1811c84c290SBarry Smith .ve 1821c84c290SBarry Smith 1831c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 1841c84c290SBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1851c84c290SBarry Smith or at runtime via the option 186be7a6d03SBarry Smith $ -mat_mffd_type my_h 1871c84c290SBarry Smith 1881c84c290SBarry Smith .keywords: MatMFFD, register 1891c84c290SBarry Smith 1901c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1911c84c290SBarry Smith @*/ 192bdf89e91SBarry Smith PetscErrorCode MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD)) 193e884886fSBarry Smith { 194e884886fSBarry Smith PetscErrorCode ierr; 195e884886fSBarry Smith 196e884886fSBarry Smith PetscFunctionBegin; 1979cc31a68SJed Brown ierr = MatInitializePackage();CHKERRQ(ierr); 198a240a19fSJed Brown ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr); 199e884886fSBarry Smith PetscFunctionReturn(0); 200e884886fSBarry Smith } 201e884886fSBarry Smith 202e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/ 20340244768SBarry Smith static PetscErrorCode MatDestroy_MFFD(Mat mat) 204e884886fSBarry Smith { 205e884886fSBarry Smith PetscErrorCode ierr; 206e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 207e884886fSBarry Smith 208e884886fSBarry Smith PetscFunctionBegin; 2096bf464f9SBarry Smith ierr = VecDestroy(&ctx->w);CHKERRQ(ierr); 2106bf464f9SBarry Smith ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr); 2116bf464f9SBarry Smith ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr); 2126bf464f9SBarry Smith ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr); 213c51ad4d4SStefano Zampini ierr = VecDestroy(&ctx->dshiftw);CHKERRQ(ierr); 214c51ad4d4SStefano Zampini ierr = VecDestroy(&ctx->current_u);CHKERRQ(ierr); 215cfe22f5eSBarry Smith if (ctx->current_f_allocated) { 2166bf464f9SBarry Smith ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr); 217cfe22f5eSBarry Smith } 218e884886fSBarry Smith if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);} 2196bf464f9SBarry Smith ierr = PetscHeaderDestroy(&ctx);CHKERRQ(ierr); 220bf0cc555SLisandro Dalcin mat->data = 0; 221e884886fSBarry Smith 222bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr); 223bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr); 224bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr); 225bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr); 226bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr); 227bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr); 228bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr); 229bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr); 230e884886fSBarry Smith PetscFunctionReturn(0); 231e884886fSBarry Smith } 232e884886fSBarry Smith 233e884886fSBarry Smith /* 234e884886fSBarry Smith MatMFFDView_MFFD - Views matrix-free parameters. 235e884886fSBarry Smith 236e884886fSBarry Smith */ 23740244768SBarry Smith static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer) 238e884886fSBarry Smith { 239e884886fSBarry Smith PetscErrorCode ierr; 240e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 241a283635eSDmitry Karpeev PetscBool iascii, viewbase, viewfunction; 242a283635eSDmitry Karpeev const char *prefix; 243e884886fSBarry Smith 244e884886fSBarry Smith PetscFunctionBegin; 245251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 246e884886fSBarry Smith if (iascii) { 247a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr); 248a283635eSDmitry Karpeev ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 24957622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr); 2507adad957SLisandro Dalcin if (!((PetscObject)ctx)->type_name) { 251e884886fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr); 252e884886fSBarry Smith } else { 2537adad957SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr); 254e884886fSBarry Smith } 255c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 256c92e3469SBarry Smith if (ctx->usecomplex) { 257c92e3469SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Using Lyness complex number trick to compute the matrix-vector product\n");CHKERRQ(ierr); 258c92e3469SBarry Smith } 259c92e3469SBarry Smith #endif 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 359c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 360c92e3469SBarry Smith if (ctx->usecomplex) h = PETSC_i*h; 361c92e3469SBarry Smith #endif 362c92e3469SBarry Smith 363e884886fSBarry Smith /* w = u + ha */ 364c3bb7e23SBarry Smith if (ctx->drscale) { 365c3bb7e23SBarry Smith ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr); 366c3bb7e23SBarry Smith ierr = VecAYPX(U,h,w);CHKERRQ(ierr); 367c3bb7e23SBarry Smith } else { 368e884886fSBarry Smith ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr); 369c3bb7e23SBarry Smith } 370e884886fSBarry Smith 3711b797266SDmitry Karpeev /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 3721b797266SDmitry Karpeev if (ctx->ncurrenth == 1 && ctx->current_f_allocated) { 373e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr); 37452121784SBarry Smith } 375e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr); 376e884886fSBarry Smith 377c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 378c92e3469SBarry Smith if (ctx->usecomplex) { 379c92e3469SBarry Smith ierr = VecImaginaryPart(y);CHKERRQ(ierr); 380c92e3469SBarry Smith h = PetscImaginaryPart(h); 381c92e3469SBarry Smith } else { 382e884886fSBarry Smith ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 383c92e3469SBarry Smith } 384c92e3469SBarry Smith #else 385c92e3469SBarry Smith ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 386c92e3469SBarry Smith #endif 387e884886fSBarry Smith ierr = VecScale(y,1.0/h);CHKERRQ(ierr); 388e884886fSBarry Smith 389c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 390e884886fSBarry Smith ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr); 391c3bb7e23SBarry Smith } 392c3bb7e23SBarry Smith if (ctx->dlscale) { 393c3bb7e23SBarry Smith ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr); 394c3bb7e23SBarry Smith } 395c3bb7e23SBarry Smith if (ctx->dshift) { 396c51ad4d4SStefano Zampini if (!ctx->dshiftw) { 397c51ad4d4SStefano Zampini ierr = VecDuplicate(y,&ctx->dshiftw);CHKERRQ(ierr); 398c51ad4d4SStefano Zampini } 399c51ad4d4SStefano Zampini ierr = VecPointwiseMult(ctx->dshift,a,ctx->dshiftw);CHKERRQ(ierr); 400c51ad4d4SStefano Zampini ierr = VecAXPY(y,1.0,ctx->dshiftw);CHKERRQ(ierr); 401c3bb7e23SBarry Smith } 402e884886fSBarry Smith 40339601f4eSBarry Smith if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);} 404e884886fSBarry Smith 405e884886fSBarry Smith ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 406e884886fSBarry Smith PetscFunctionReturn(0); 407e884886fSBarry Smith } 408e884886fSBarry Smith 409e884886fSBarry Smith /* 410e884886fSBarry Smith MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 411e884886fSBarry Smith 412e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 413e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 414e884886fSBarry Smith u = current iterate 415e884886fSBarry Smith h = difference interval 416e884886fSBarry Smith */ 4175d21e1f8SStefano Zampini PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 418e884886fSBarry Smith { 419e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 420e884886fSBarry Smith PetscScalar h,*aa,*ww,v; 421e884886fSBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 422e884886fSBarry Smith Vec w,U; 423e884886fSBarry Smith PetscErrorCode ierr; 424e884886fSBarry Smith PetscInt i,rstart,rend; 425e884886fSBarry Smith 426e884886fSBarry Smith PetscFunctionBegin; 427486afcceSStefano Zampini if (!ctx->func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunction() first"); 428486afcceSStefano Zampini if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioni() first"); 429486afcceSStefano Zampini if (!ctx->funcisetbase) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioniBase() first"); 430e884886fSBarry Smith w = ctx->w; 431e884886fSBarry Smith U = ctx->current_u; 432e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr); 433e884886fSBarry Smith ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr); 434e884886fSBarry Smith ierr = VecCopy(U,w);CHKERRQ(ierr); 435e884886fSBarry Smith 436e884886fSBarry Smith ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr); 437e884886fSBarry Smith ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 438e884886fSBarry Smith for (i=rstart; i<rend; i++) { 439e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 440e884886fSBarry Smith h = ww[i-rstart]; 441e884886fSBarry Smith if (h == 0.0) h = 1.0; 442e884886fSBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 443e884886fSBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 444e884886fSBarry Smith h *= epsilon; 445e884886fSBarry Smith 446e884886fSBarry Smith ww[i-rstart] += h; 447e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 448e884886fSBarry Smith ierr = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr); 449e884886fSBarry Smith aa[i-rstart] = (v - aa[i-rstart])/h; 450e884886fSBarry Smith 451e884886fSBarry Smith /* possibly shift and scale result */ 452c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 453e884886fSBarry Smith aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart]; 454c3bb7e23SBarry Smith } 455e884886fSBarry Smith 456e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 457e884886fSBarry Smith ww[i-rstart] -= h; 458e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 459e884886fSBarry Smith } 460e884886fSBarry Smith ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 461e884886fSBarry Smith PetscFunctionReturn(0); 462e884886fSBarry Smith } 463e884886fSBarry Smith 46440244768SBarry Smith static PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr) 465c3bb7e23SBarry Smith { 466c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 467c3bb7e23SBarry Smith PetscErrorCode ierr; 468c3bb7e23SBarry Smith 469c3bb7e23SBarry Smith PetscFunctionBegin; 470c3bb7e23SBarry Smith if (ll && !aij->dlscale) { 471c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr); 472c3bb7e23SBarry Smith } 473c3bb7e23SBarry Smith if (rr && !aij->drscale) { 474c3bb7e23SBarry Smith ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr); 475c3bb7e23SBarry Smith } 476c3bb7e23SBarry Smith if (ll) { 477c3bb7e23SBarry Smith ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr); 478c3bb7e23SBarry Smith } 479c3bb7e23SBarry Smith if (rr) { 480c3bb7e23SBarry Smith ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr); 481c3bb7e23SBarry Smith } 482c3bb7e23SBarry Smith PetscFunctionReturn(0); 483c3bb7e23SBarry Smith } 484c3bb7e23SBarry Smith 48540244768SBarry Smith static PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode) 486c3bb7e23SBarry Smith { 487c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 488c3bb7e23SBarry Smith PetscErrorCode ierr; 489c3bb7e23SBarry Smith 490c3bb7e23SBarry Smith PetscFunctionBegin; 491ce94432eSBarry Smith if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES"); 492c3bb7e23SBarry Smith if (!aij->dshift) { 493c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr); 494c3bb7e23SBarry Smith } 495c3bb7e23SBarry Smith ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr); 496c3bb7e23SBarry Smith PetscFunctionReturn(0); 497c3bb7e23SBarry Smith } 498c3bb7e23SBarry Smith 49940244768SBarry Smith static PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a) 500e884886fSBarry Smith { 501e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 5025fd66863SKarl Rupp 503e884886fSBarry Smith PetscFunctionBegin; 504e884886fSBarry Smith shell->vshift += a; 505e884886fSBarry Smith PetscFunctionReturn(0); 506e884886fSBarry Smith } 507e884886fSBarry Smith 50840244768SBarry Smith static PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a) 509e884886fSBarry Smith { 510e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 5115fd66863SKarl Rupp 512e884886fSBarry Smith PetscFunctionBegin; 513e884886fSBarry Smith shell->vscale *= a; 514e884886fSBarry Smith PetscFunctionReturn(0); 515e884886fSBarry Smith } 516e884886fSBarry Smith 517d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F) 518e884886fSBarry Smith { 519e884886fSBarry Smith PetscErrorCode ierr; 520e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 521e884886fSBarry Smith 522e884886fSBarry Smith PetscFunctionBegin; 523e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 524c51ad4d4SStefano Zampini if (!ctx->current_u) { 525c51ad4d4SStefano Zampini ierr = VecDuplicate(U,&ctx->current_u);CHKERRQ(ierr); 5268860a134SJunchao Zhang ierr = VecLockReadPush(ctx->current_u);CHKERRQ(ierr); 527c51ad4d4SStefano Zampini } 5288860a134SJunchao Zhang ierr = VecLockReadPop(ctx->current_u);CHKERRQ(ierr); 529c51ad4d4SStefano Zampini ierr = VecCopy(U,ctx->current_u);CHKERRQ(ierr); 5308860a134SJunchao Zhang ierr = VecLockReadPush(ctx->current_u);CHKERRQ(ierr); 53152121784SBarry Smith if (F) { 5326bf464f9SBarry Smith if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);} 5333ec795f1SBarry Smith ctx->current_f = F; 534cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_FALSE; 535cfe22f5eSBarry Smith } else if (!ctx->current_f_allocated) { 53606c4b6baSBarry Smith ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr); 5372205254eSKarl Rupp 538cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_TRUE; 53952121784SBarry Smith } 540e884886fSBarry Smith if (!ctx->w) { 541e884886fSBarry Smith ierr = VecDuplicate(ctx->current_u,&ctx->w);CHKERRQ(ierr); 542e884886fSBarry Smith } 543e884886fSBarry Smith J->assembled = PETSC_TRUE; 544e884886fSBarry Smith PetscFunctionReturn(0); 545e884886fSBarry Smith } 5465a576424SJed Brown 547e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 548b2573a8aSBarry Smith 54940244768SBarry Smith static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx) 550e884886fSBarry Smith { 551e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 552e884886fSBarry Smith 553e884886fSBarry Smith PetscFunctionBegin; 554e884886fSBarry Smith ctx->checkh = fun; 555e884886fSBarry Smith ctx->checkhctx = ectx; 556e884886fSBarry Smith PetscFunctionReturn(0); 557e884886fSBarry Smith } 558e884886fSBarry Smith 5596aa9148fSLisandro Dalcin /*@C 5606aa9148fSLisandro Dalcin MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 5616aa9148fSLisandro Dalcin MatMFFD options in the database. 5626aa9148fSLisandro Dalcin 5636aa9148fSLisandro Dalcin Collective on Mat 5646aa9148fSLisandro Dalcin 5656aa9148fSLisandro Dalcin Input Parameter: 5666aa9148fSLisandro Dalcin + A - the Mat context 5676aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names 5686aa9148fSLisandro Dalcin 5696aa9148fSLisandro Dalcin Notes: 5706aa9148fSLisandro Dalcin A hyphen (-) must NOT be given at the beginning of the prefix name. 5716aa9148fSLisandro Dalcin The first character of all runtime options is AUTOMATICALLY the hyphen. 5726aa9148fSLisandro Dalcin 5736aa9148fSLisandro Dalcin Level: advanced 5746aa9148fSLisandro Dalcin 5756aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters 5766aa9148fSLisandro Dalcin 577755b7f64SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD() 5786aa9148fSLisandro Dalcin @*/ 5797087cfbeSBarry Smith PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[]) 5806aa9148fSLisandro Dalcin 5816aa9148fSLisandro Dalcin { 5820298fd71SBarry Smith MatMFFD mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL; 5836aa9148fSLisandro Dalcin PetscErrorCode ierr; 5845fd66863SKarl Rupp 5856aa9148fSLisandro Dalcin PetscFunctionBegin; 5860700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5870700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 5886aa9148fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr); 5896aa9148fSLisandro Dalcin PetscFunctionReturn(0); 5906aa9148fSLisandro Dalcin } 5916aa9148fSLisandro Dalcin 59240244768SBarry Smith static PetscErrorCode MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat) 593e884886fSBarry Smith { 59471f08665SBarry Smith MatMFFD mfctx = (MatMFFD)mat->data; 595e884886fSBarry Smith PetscErrorCode ierr; 596ace3abfcSBarry Smith PetscBool flg; 597e884886fSBarry Smith char ftype[256]; 598e884886fSBarry Smith 599e884886fSBarry Smith PetscFunctionBegin; 6000700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 6010700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 6023194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr); 603a264d7a6SBarry Smith ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr); 604e884886fSBarry Smith if (flg) { 605e884886fSBarry Smith ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr); 606e884886fSBarry Smith } 607e884886fSBarry Smith 608e884886fSBarry Smith ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 609e884886fSBarry Smith ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 610e884886fSBarry Smith 61190d69ab7SBarry Smith flg = PETSC_FALSE; 6120298fd71SBarry Smith ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr); 613e884886fSBarry Smith if (flg) { 614e884886fSBarry Smith ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr); 615e884886fSBarry Smith } 616c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 617c92e3469SBarry Smith ierr = PetscOptionsBool("-mat_mffd_complex","Use Lyness complex number trick to compute the matrix-vector product","None",mfctx->usecomplex,&mfctx->usecomplex,NULL);CHKERRQ(ierr); 618c92e3469SBarry Smith #endif 619e884886fSBarry Smith if (mfctx->ops->setfromoptions) { 620e55864a3SBarry Smith ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr); 621e884886fSBarry Smith } 622e884886fSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 623e884886fSBarry Smith PetscFunctionReturn(0); 624e884886fSBarry Smith } 625e884886fSBarry Smith 62640244768SBarry Smith static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period) 627bc0f08ceSBarry Smith { 628bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 629bc0f08ceSBarry Smith 630bc0f08ceSBarry Smith PetscFunctionBegin; 631bc0f08ceSBarry Smith ctx->recomputeperiod = period; 632bc0f08ceSBarry Smith PetscFunctionReturn(0); 633bc0f08ceSBarry Smith } 634bc0f08ceSBarry Smith 63540244768SBarry Smith static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 636bc0f08ceSBarry Smith { 637bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 638bc0f08ceSBarry Smith 639bc0f08ceSBarry Smith PetscFunctionBegin; 640bc0f08ceSBarry Smith ctx->func = func; 641bc0f08ceSBarry Smith ctx->funcctx = funcctx; 642bc0f08ceSBarry Smith PetscFunctionReturn(0); 643bc0f08ceSBarry Smith } 644bc0f08ceSBarry Smith 64540244768SBarry Smith static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error) 646bc0f08ceSBarry Smith { 647bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 648bc0f08ceSBarry Smith 649bc0f08ceSBarry Smith PetscFunctionBegin; 650bc0f08ceSBarry Smith if (error != PETSC_DEFAULT) ctx->error_rel = error; 651bc0f08ceSBarry Smith PetscFunctionReturn(0); 652bc0f08ceSBarry Smith } 653bc0f08ceSBarry Smith 6543b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MFFD(Mat A,PetscBool *missing,PetscInt *d) 6553b49f96aSBarry Smith { 6563b49f96aSBarry Smith PetscFunctionBegin; 6573b49f96aSBarry Smith *missing = PETSC_FALSE; 6583b49f96aSBarry Smith PetscFunctionReturn(0); 6593b49f96aSBarry Smith } 6603b49f96aSBarry Smith 661e884886fSBarry Smith /*MC 662e884886fSBarry Smith MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 663e884886fSBarry Smith 664e884886fSBarry Smith Level: advanced 665e884886fSBarry Smith 666755b7f64SBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(), 667755b7f64SBarry Smith MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 668755b7f64SBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 669755b7f64SBarry Smith MatMFFDGetH(), 670e884886fSBarry Smith M*/ 6718cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 672e884886fSBarry Smith { 673e884886fSBarry Smith MatMFFD mfctx; 674e884886fSBarry Smith PetscErrorCode ierr; 675e884886fSBarry Smith 676e884886fSBarry Smith PetscFunctionBegin; 677607a6623SBarry Smith ierr = MatMFFDInitializePackage();CHKERRQ(ierr); 678e884886fSBarry Smith 67973107ff1SLisandro Dalcin ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr); 6802205254eSKarl Rupp 681e884886fSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 682e884886fSBarry Smith mfctx->recomputeperiod = 1; 683e884886fSBarry Smith mfctx->count = 0; 684e884886fSBarry Smith mfctx->currenth = 0.0; 6850298fd71SBarry Smith mfctx->historyh = NULL; 686e884886fSBarry Smith mfctx->ncurrenth = 0; 687e884886fSBarry Smith mfctx->maxcurrenth = 0; 6887adad957SLisandro Dalcin ((PetscObject)mfctx)->type_name = 0; 689e884886fSBarry Smith 690e884886fSBarry Smith mfctx->vshift = 0.0; 691e884886fSBarry Smith mfctx->vscale = 1.0; 692e884886fSBarry Smith 693e884886fSBarry Smith /* 694e884886fSBarry Smith Create the empty data structure to contain compute-h routines. 695e884886fSBarry Smith These will be filled in below from the command line options or 696e884886fSBarry Smith a later call with MatMFFDSetType() or if that is not called 697e884886fSBarry Smith then it will default in the first use of MatMult_MFFD() 698e884886fSBarry Smith */ 699e884886fSBarry Smith mfctx->ops->compute = 0; 700e884886fSBarry Smith mfctx->ops->destroy = 0; 701e884886fSBarry Smith mfctx->ops->view = 0; 702e884886fSBarry Smith mfctx->ops->setfromoptions = 0; 703e884886fSBarry Smith mfctx->hctx = 0; 704e884886fSBarry Smith 705e884886fSBarry Smith mfctx->func = 0; 706e884886fSBarry Smith mfctx->funcctx = 0; 7070298fd71SBarry Smith mfctx->w = NULL; 708e884886fSBarry Smith 709e884886fSBarry Smith A->data = mfctx; 710e884886fSBarry Smith 711e884886fSBarry Smith A->ops->mult = MatMult_MFFD; 712e884886fSBarry Smith A->ops->destroy = MatDestroy_MFFD; 713e884886fSBarry Smith A->ops->view = MatView_MFFD; 714e884886fSBarry Smith A->ops->assemblyend = MatAssemblyEnd_MFFD; 715e884886fSBarry Smith A->ops->scale = MatScale_MFFD; 716e884886fSBarry Smith A->ops->shift = MatShift_MFFD; 717c3bb7e23SBarry Smith A->ops->diagonalscale = MatDiagonalScale_MFFD; 718c3bb7e23SBarry Smith A->ops->diagonalset = MatDiagonalSet_MFFD; 7199c6ac3b3SBarry Smith A->ops->setfromoptions = MatSetFromOptions_MFFD; 7203b49f96aSBarry Smith A->ops->missingdiagonal = MatMissingDiagonal_MFFD; 721e884886fSBarry Smith A->assembled = PETSC_TRUE; 722e884886fSBarry Smith 723bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr); 724bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr); 725bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr); 726bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr); 727bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr); 728bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr); 729bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr); 730bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr); 7312205254eSKarl Rupp 732e884886fSBarry Smith mfctx->mat = A; 7332205254eSKarl Rupp 734e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr); 735e884886fSBarry Smith PetscFunctionReturn(0); 736e884886fSBarry Smith } 737e884886fSBarry Smith 738e884886fSBarry Smith /*@ 739e884886fSBarry Smith MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 740e884886fSBarry Smith 741e884886fSBarry Smith Collective on Vec 742e884886fSBarry Smith 743e884886fSBarry Smith Input Parameters: 744fef1beadSBarry Smith + comm - MPI communicator 745fef1beadSBarry Smith . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 746fef1beadSBarry Smith This value should be the same as the local size used in creating the 747fef1beadSBarry Smith y vector for the matrix-vector product y = Ax. 748fef1beadSBarry Smith . n - This value should be the same as the local size used in creating the 749fef1beadSBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 750fef1beadSBarry Smith calculated if N is given) For square matrices n is almost always m. 751fef1beadSBarry Smith . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 752fef1beadSBarry Smith - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 753fef1beadSBarry Smith 754e884886fSBarry Smith 755e884886fSBarry Smith Output Parameter: 756e884886fSBarry Smith . J - the matrix-free matrix 757e884886fSBarry Smith 7589c6ac3b3SBarry Smith Options Database Keys: call MatSetFromOptions() to trigger these 7599c6ac3b3SBarry Smith + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 760c92e3469SBarry Smith . -mat_mffd_err - square root of estimated relative error in function evaluation 761c92e3469SBarry Smith . -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 762c92e3469SBarry Smith . -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values 763c92e3469SBarry Smith - -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing 764c92e3469SBarry Smith (requires real valued functions but that PETSc be configured for complex numbers) 7659c6ac3b3SBarry Smith 7669c6ac3b3SBarry Smith 767e884886fSBarry Smith Level: advanced 768e884886fSBarry Smith 769e884886fSBarry Smith Notes: 770e884886fSBarry Smith The matrix-free matrix context merely contains the function pointers 771e884886fSBarry Smith and work space for performing finite difference approximations of 772e884886fSBarry Smith Jacobian-vector products, F'(u)*a, 773e884886fSBarry Smith 774e884886fSBarry Smith The default code uses the following approach to compute h 775e884886fSBarry Smith 776e884886fSBarry Smith .vb 777e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 778e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 779e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 780e884886fSBarry Smith where 781e884886fSBarry Smith error_rel = square root of relative error in function evaluation 782e884886fSBarry Smith umin = minimum iterate parameter 783e884886fSBarry Smith .ve 784e884886fSBarry Smith 785ca93e954SBarry Smith You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different 786ca93e954SBarry Smith preconditioner matrix 787ca93e954SBarry Smith 788e884886fSBarry Smith The user can set the error_rel via MatMFFDSetFunctionError() and 789a7f22e61SSatish Balay umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details. 790e884886fSBarry Smith 791e884886fSBarry Smith The user should call MatDestroy() when finished with the matrix-free 792e884886fSBarry Smith matrix context. 793e884886fSBarry Smith 794e884886fSBarry Smith Options Database Keys: 795e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 796e884886fSBarry Smith . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only) 797e884886fSBarry Smith - -mat_mffd_check_positivity 798e884886fSBarry Smith 799e884886fSBarry Smith .keywords: default, matrix-free, create, matrix 800e884886fSBarry Smith 8011d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 802e884886fSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 80381242352SJed Brown MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian() 804e884886fSBarry Smith 805e884886fSBarry Smith @*/ 8067087cfbeSBarry Smith PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 807e884886fSBarry Smith { 808e884886fSBarry Smith PetscErrorCode ierr; 809e884886fSBarry Smith 810e884886fSBarry Smith PetscFunctionBegin; 811e884886fSBarry Smith ierr = MatCreate(comm,J);CHKERRQ(ierr); 812fef1beadSBarry Smith ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr); 813e884886fSBarry Smith ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 814be7c243fSJed Brown ierr = MatSetUp(*J);CHKERRQ(ierr); 815e884886fSBarry Smith PetscFunctionReturn(0); 816e884886fSBarry Smith } 817e884886fSBarry Smith 818e884886fSBarry Smith /*@ 819e884886fSBarry Smith MatMFFDGetH - Gets the last value that was used as the differencing 820e884886fSBarry Smith parameter. 821e884886fSBarry Smith 822e884886fSBarry Smith Not Collective 823e884886fSBarry Smith 824e884886fSBarry Smith Input Parameters: 825e884886fSBarry Smith . mat - the matrix obtained with MatCreateSNESMF() 826e884886fSBarry Smith 827e884886fSBarry Smith Output Paramter: 828e884886fSBarry Smith . h - the differencing step size 829e884886fSBarry Smith 830e884886fSBarry Smith Level: advanced 831e884886fSBarry Smith 832e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 833e884886fSBarry Smith 8341d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory() 835e884886fSBarry Smith @*/ 8367087cfbeSBarry Smith PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 837e884886fSBarry Smith { 838e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 839bc0f08ceSBarry Smith PetscErrorCode ierr; 840bc0f08ceSBarry Smith PetscBool match; 841e884886fSBarry Smith 842e884886fSBarry Smith PetscFunctionBegin; 84388b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 84488b4c220SStefano Zampini PetscValidPointer(h,2); 845251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 846ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 847bc0f08ceSBarry Smith 848e884886fSBarry Smith *h = ctx->currenth; 849e884886fSBarry Smith PetscFunctionReturn(0); 850e884886fSBarry Smith } 851e884886fSBarry Smith 852e884886fSBarry Smith /*@C 853e884886fSBarry Smith MatMFFDSetFunction - Sets the function used in applying the matrix free. 854e884886fSBarry Smith 8553f9fe445SBarry Smith Logically Collective on Mat 856e884886fSBarry Smith 857e884886fSBarry Smith Input Parameters: 85814f633e2SBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD() 859e884886fSBarry Smith . func - the function to use 860e884886fSBarry Smith - funcctx - optional function context passed to function 861e884886fSBarry Smith 86214f633e2SBarry Smith Calling Sequence of func: 86314f633e2SBarry Smith $ func (void *funcctx, Vec x, Vec f) 86414f633e2SBarry Smith 86514f633e2SBarry Smith + funcctx - user provided context 86614f633e2SBarry Smith . x - input vector 86714f633e2SBarry Smith - f - computed output function 86814f633e2SBarry Smith 869e884886fSBarry Smith Level: advanced 870e884886fSBarry Smith 871e884886fSBarry Smith Notes: 872e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 873e884886fSBarry Smith matrix inside your compute Jacobian routine 874e884886fSBarry Smith 8753ec795f1SBarry Smith If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 876e884886fSBarry Smith 877e884886fSBarry Smith .keywords: SNES, matrix-free, function 878e884886fSBarry Smith 8791d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD, 8801d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 881e884886fSBarry Smith @*/ 8827087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 883e884886fSBarry Smith { 884bc0f08ceSBarry Smith PetscErrorCode ierr; 885e884886fSBarry Smith 886e884886fSBarry Smith PetscFunctionBegin; 88788b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 888bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr); 889e884886fSBarry Smith PetscFunctionReturn(0); 890e884886fSBarry Smith } 891e884886fSBarry Smith 892e884886fSBarry Smith /*@C 893e884886fSBarry Smith MatMFFDSetFunctioni - Sets the function for a single component 894e884886fSBarry Smith 8953f9fe445SBarry Smith Logically Collective on Mat 896e884886fSBarry Smith 897e884886fSBarry Smith Input Parameters: 898e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 899e884886fSBarry Smith - funci - the function to use 900e884886fSBarry Smith 901e884886fSBarry Smith Level: advanced 902e884886fSBarry Smith 903e884886fSBarry Smith Notes: 904e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 90594eb4328SStefano Zampini matrix inside your compute Jacobian routine. 90694eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 907486afcceSStefano Zampini funci must not contain any MPI call as it is called inside a loop on the local portion of the vector. 908e884886fSBarry Smith 909e884886fSBarry Smith .keywords: SNES, matrix-free, function 910e884886fSBarry Smith 91194eb4328SStefano Zampini .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 9121d0fab5eSBarry Smith 913e884886fSBarry Smith @*/ 9147087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 915e884886fSBarry Smith { 9164ac538c5SBarry Smith PetscErrorCode ierr; 917e884886fSBarry Smith 918e884886fSBarry Smith PetscFunctionBegin; 9190700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 9204ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr); 921e884886fSBarry Smith PetscFunctionReturn(0); 922e884886fSBarry Smith } 923e884886fSBarry Smith 924e884886fSBarry Smith /*@C 925e884886fSBarry Smith MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 926e884886fSBarry Smith 9273f9fe445SBarry Smith Logically Collective on Mat 928e884886fSBarry Smith 929e884886fSBarry Smith Input Parameters: 930e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 931e884886fSBarry Smith - func - the function to use 932e884886fSBarry Smith 933e884886fSBarry Smith Level: advanced 934e884886fSBarry Smith 935e884886fSBarry Smith Notes: 936e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 93794eb4328SStefano Zampini matrix inside your compute Jacobian routine. 93894eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 939e884886fSBarry Smith 940e884886fSBarry Smith 941e884886fSBarry Smith .keywords: SNES, matrix-free, function 942e884886fSBarry Smith 943e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 94494eb4328SStefano Zampini MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 945e884886fSBarry Smith @*/ 9467087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 947e884886fSBarry Smith { 9484ac538c5SBarry Smith PetscErrorCode ierr; 949e884886fSBarry Smith 950e884886fSBarry Smith PetscFunctionBegin; 9510700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 9524ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr); 953e884886fSBarry Smith PetscFunctionReturn(0); 954e884886fSBarry Smith } 955e884886fSBarry Smith 956e884886fSBarry Smith /*@ 957e884886fSBarry Smith MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 958e884886fSBarry Smith 9593f9fe445SBarry Smith Logically Collective on Mat 960e884886fSBarry Smith 961e884886fSBarry Smith Input Parameters: 962e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 963e884886fSBarry Smith - period - 1 for everytime, 2 for every second etc 964e884886fSBarry Smith 965e884886fSBarry Smith Options Database Keys: 966*a2b725a8SWilliam Gropp . -mat_mffd_period <period> 967e884886fSBarry Smith 968e884886fSBarry Smith Level: advanced 969e884886fSBarry Smith 970e884886fSBarry Smith 971e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 972e884886fSBarry Smith 973e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), 9741d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 975e884886fSBarry Smith @*/ 9767087cfbeSBarry Smith PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 977e884886fSBarry Smith { 978bc0f08ceSBarry Smith PetscErrorCode ierr; 979e884886fSBarry Smith 980e884886fSBarry Smith PetscFunctionBegin; 98188b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 98288b4c220SStefano Zampini PetscValidLogicalCollectiveInt(mat,period,2); 983bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr); 984e884886fSBarry Smith PetscFunctionReturn(0); 985e884886fSBarry Smith } 986e884886fSBarry Smith 987e884886fSBarry Smith /*@ 988e884886fSBarry Smith MatMFFDSetFunctionError - Sets the error_rel for the approximation of 989e884886fSBarry Smith matrix-vector products using finite differences. 990e884886fSBarry Smith 9913f9fe445SBarry Smith Logically Collective on Mat 992e884886fSBarry Smith 993e884886fSBarry Smith Input Parameters: 994e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 995e884886fSBarry Smith - error_rel - relative error (should be set to the square root of 996e884886fSBarry Smith the relative error in the function evaluations) 997e884886fSBarry Smith 998e884886fSBarry Smith Options Database Keys: 999*a2b725a8SWilliam Gropp . -mat_mffd_err <error_rel> - Sets error_rel 1000e884886fSBarry Smith 1001e884886fSBarry Smith Level: advanced 1002e884886fSBarry Smith 1003e884886fSBarry Smith Notes: 1004e884886fSBarry Smith The default matrix-free matrix-vector product routine computes 1005e884886fSBarry Smith .vb 1006e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 1007e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 1008e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 1009e884886fSBarry Smith .ve 1010e884886fSBarry Smith 1011e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 1012e884886fSBarry Smith 1013e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 10141d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 1015e884886fSBarry Smith @*/ 10167087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 1017e884886fSBarry Smith { 1018bc0f08ceSBarry Smith PetscErrorCode ierr; 1019e884886fSBarry Smith 1020e884886fSBarry Smith PetscFunctionBegin; 102188b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 102288b4c220SStefano Zampini PetscValidLogicalCollectiveReal(mat,error,2); 1023bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr); 1024e884886fSBarry Smith PetscFunctionReturn(0); 1025e884886fSBarry Smith } 1026e884886fSBarry Smith 1027e884886fSBarry Smith /*@ 1028e884886fSBarry Smith MatMFFDSetHHistory - Sets an array to collect a history of the 1029e884886fSBarry Smith differencing values (h) computed for the matrix-free product. 1030e884886fSBarry Smith 10313f9fe445SBarry Smith Logically Collective on Mat 1032e884886fSBarry Smith 1033e884886fSBarry Smith Input Parameters: 1034e884886fSBarry Smith + J - the matrix-free matrix context 1035e884886fSBarry Smith . histroy - space to hold the history 1036e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than 1037e884886fSBarry Smith nhistory, then the later ones are discarded 1038e884886fSBarry Smith 1039e884886fSBarry Smith Level: advanced 1040e884886fSBarry Smith 1041e884886fSBarry Smith Notes: 1042e884886fSBarry Smith Use MatMFFDResetHHistory() to reset the history counter and collect 1043e884886fSBarry Smith a new batch of differencing parameters, h. 1044e884886fSBarry Smith 1045e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1046e884886fSBarry Smith 1047e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 10481d0fab5eSBarry Smith MatMFFDResetHHistory(), MatMFFDSetFunctionError() 1049e884886fSBarry Smith 1050e884886fSBarry Smith @*/ 10517087cfbeSBarry Smith PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 1052e884886fSBarry Smith { 1053e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 1054bc0f08ceSBarry Smith PetscErrorCode ierr; 1055bc0f08ceSBarry Smith PetscBool match; 1056e884886fSBarry Smith 1057e884886fSBarry Smith PetscFunctionBegin; 105888b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 105988b4c220SStefano Zampini if (history) PetscValidPointer(history,2); 106088b4c220SStefano Zampini PetscValidLogicalCollectiveInt(J,nhistory,3); 1061251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr); 1062ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 1063e884886fSBarry Smith ctx->historyh = history; 1064e884886fSBarry Smith ctx->maxcurrenth = nhistory; 106575567043SBarry Smith ctx->currenth = 0.; 1066e884886fSBarry Smith PetscFunctionReturn(0); 1067e884886fSBarry Smith } 1068e884886fSBarry Smith 1069e884886fSBarry Smith /*@ 1070e884886fSBarry Smith MatMFFDResetHHistory - Resets the counter to zero to begin 1071e884886fSBarry Smith collecting a new set of differencing histories. 1072e884886fSBarry Smith 10733f9fe445SBarry Smith Logically Collective on Mat 1074e884886fSBarry Smith 1075e884886fSBarry Smith Input Parameters: 1076e884886fSBarry Smith . J - the matrix-free matrix context 1077e884886fSBarry Smith 1078e884886fSBarry Smith Level: advanced 1079e884886fSBarry Smith 1080e884886fSBarry Smith Notes: 1081e884886fSBarry Smith Use MatMFFDSetHHistory() to create the original history counter. 1082e884886fSBarry Smith 1083e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1084e884886fSBarry Smith 1085e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 10861d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDSetFunctionError() 1087e884886fSBarry Smith 1088e884886fSBarry Smith @*/ 10897087cfbeSBarry Smith PetscErrorCode MatMFFDResetHHistory(Mat J) 1090e884886fSBarry Smith { 1091bc0f08ceSBarry Smith PetscErrorCode ierr; 1092e884886fSBarry Smith 1093e884886fSBarry Smith PetscFunctionBegin; 109488b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1095bc0f08ceSBarry Smith ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr); 1096e884886fSBarry Smith PetscFunctionReturn(0); 1097e884886fSBarry Smith } 1098e884886fSBarry Smith 1099487a658cSBarry Smith /*@ 1100e884886fSBarry Smith MatMFFDSetBase - Sets the vector U at which matrix vector products of the 1101e884886fSBarry Smith Jacobian are computed 1102e884886fSBarry Smith 11033f9fe445SBarry Smith Logically Collective on Mat 1104e884886fSBarry Smith 1105e884886fSBarry Smith Input Parameters: 1106e884886fSBarry Smith + J - the MatMFFD matrix 11073ec795f1SBarry Smith . U - the vector 1108bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed 1109e884886fSBarry Smith 111095452b02SPatrick Sanan Notes: 111195452b02SPatrick Sanan This is rarely used directly 1112e884886fSBarry Smith 11138af5ae88SBarry Smith If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 11148af5ae88SBarry Smith point during the first MatMult() after each call to MatMFFDSetBase(). 1115dff2f722SBarry Smith 1116e884886fSBarry Smith Level: advanced 1117e884886fSBarry Smith 1118e884886fSBarry Smith @*/ 11197087cfbeSBarry Smith PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 1120e884886fSBarry Smith { 11214ac538c5SBarry Smith PetscErrorCode ierr; 1122e884886fSBarry Smith 1123e884886fSBarry Smith PetscFunctionBegin; 11240700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 11250700a824SBarry Smith PetscValidHeaderSpecific(U,VEC_CLASSID,2); 11260700a824SBarry Smith if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 11274ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr); 1128e884886fSBarry Smith PetscFunctionReturn(0); 1129e884886fSBarry Smith } 1130e884886fSBarry Smith 1131e884886fSBarry Smith /*@C 1132e884886fSBarry Smith MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 1133e884886fSBarry Smith it to satisfy some criteria 1134e884886fSBarry Smith 11353f9fe445SBarry Smith Logically Collective on Mat 1136e884886fSBarry Smith 1137e884886fSBarry Smith Input Parameters: 1138e884886fSBarry Smith + J - the MatMFFD matrix 1139e884886fSBarry Smith . fun - the function that checks h 1140e884886fSBarry Smith - ctx - any context needed by the function 1141e884886fSBarry Smith 1142e884886fSBarry Smith Options Database Keys: 1143e884886fSBarry Smith . -mat_mffd_check_positivity 1144e884886fSBarry Smith 1145e884886fSBarry Smith Level: advanced 1146e884886fSBarry Smith 114795452b02SPatrick Sanan Notes: 114895452b02SPatrick Sanan For example, MatMFFDCheckPositivity() insures that all entries 1149e884886fSBarry Smith of U + h*a are non-negative 1150e884886fSBarry Smith 1151755b7f64SBarry Smith The function you provide is called after the default h has been computed and allows you to 1152755b7f64SBarry Smith modify it. 1153755b7f64SBarry Smith 1154755b7f64SBarry Smith .seealso: MatMFFDCheckPositivity() 1155e884886fSBarry Smith @*/ 11567087cfbeSBarry Smith PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx) 1157e884886fSBarry Smith { 11584ac538c5SBarry Smith PetscErrorCode ierr; 1159e884886fSBarry Smith 1160e884886fSBarry Smith PetscFunctionBegin; 11610700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 11624ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr); 1163e884886fSBarry Smith PetscFunctionReturn(0); 1164e884886fSBarry Smith } 1165e884886fSBarry Smith 1166e884886fSBarry Smith /*@ 1167e884886fSBarry Smith MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1168e884886fSBarry Smith zero, decreases h until this is satisfied. 1169e884886fSBarry Smith 11703f9fe445SBarry Smith Logically Collective on Vec 1171e884886fSBarry Smith 1172e884886fSBarry Smith Input Parameters: 1173e884886fSBarry Smith + U - base vector that is added to 1174e884886fSBarry Smith . a - vector that is added 1175e884886fSBarry Smith . h - scaling factor on a 1176e884886fSBarry Smith - dummy - context variable (unused) 1177e884886fSBarry Smith 1178e884886fSBarry Smith Options Database Keys: 1179e884886fSBarry Smith . -mat_mffd_check_positivity 1180e884886fSBarry Smith 1181e884886fSBarry Smith Level: advanced 1182e884886fSBarry Smith 118395452b02SPatrick Sanan Notes: 118495452b02SPatrick Sanan This is rarely used directly, rather it is passed as an argument to 1185e884886fSBarry Smith MatMFFDSetCheckh() 1186e884886fSBarry Smith 1187e884886fSBarry Smith .seealso: MatMFFDSetCheckh() 1188e884886fSBarry Smith @*/ 11897087cfbeSBarry Smith PetscErrorCode MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h) 1190e884886fSBarry Smith { 1191e884886fSBarry Smith PetscReal val, minval; 1192e884886fSBarry Smith PetscScalar *u_vec, *a_vec; 1193e884886fSBarry Smith PetscErrorCode ierr; 1194e884886fSBarry Smith PetscInt i,n; 1195e884886fSBarry Smith MPI_Comm comm; 1196e884886fSBarry Smith 1197e884886fSBarry Smith PetscFunctionBegin; 119888b4c220SStefano Zampini PetscValidHeaderSpecific(U,VEC_CLASSID,2); 119988b4c220SStefano Zampini PetscValidHeaderSpecific(a,VEC_CLASSID,3); 120088b4c220SStefano Zampini PetscValidPointer(h,4); 1201e884886fSBarry Smith ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); 1202e884886fSBarry Smith ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); 1203e884886fSBarry Smith ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); 1204e884886fSBarry Smith ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); 1205c068d9bbSLisandro Dalcin minval = PetscAbsScalar(*h)*PetscRealConstant(1.01); 1206e884886fSBarry Smith for (i=0; i<n; i++) { 1207e884886fSBarry Smith if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1208e884886fSBarry Smith val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1209e884886fSBarry Smith if (val < minval) minval = val; 1210e884886fSBarry Smith } 1211e884886fSBarry Smith } 1212e884886fSBarry Smith ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); 1213e884886fSBarry Smith ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); 1214b2566f29SBarry Smith ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr); 1215e884886fSBarry Smith if (val <= PetscAbsScalar(*h)) { 12166712e2f1SBarry Smith ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr); 1217e884886fSBarry Smith if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1218e884886fSBarry Smith else *h = -0.99*val; 1219e884886fSBarry Smith } 1220e884886fSBarry Smith PetscFunctionReturn(0); 1221e884886fSBarry Smith } 1222