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]; 458e81d068SLisandro Dalcin PetscBool opt,pkg; 463ec795f1SBarry Smith PetscErrorCode ierr; 473ec795f1SBarry Smith 483ec795f1SBarry Smith PetscFunctionBegin; 49b022a5c1SBarry Smith if (MatMFFDPackageInitialized) PetscFunctionReturn(0); 50b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_TRUE; 513ec795f1SBarry Smith /* Register Classes */ 520700a824SBarry Smith ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr); 533ec795f1SBarry Smith /* Register Constructors */ 54607a6623SBarry Smith ierr = MatMFFDRegisterAll();CHKERRQ(ierr); 553ec795f1SBarry Smith /* Register Events */ 560700a824SBarry Smith ierr = PetscLogEventRegister("MatMult MF",MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr); 573ec795f1SBarry Smith /* Process info exclusions */ 588e81d068SLisandro Dalcin ierr = PetscOptionsGetString(NULL,NULL,"-info_exclude",logList,sizeof(logList),&opt);CHKERRQ(ierr); 593ec795f1SBarry Smith if (opt) { 608e81d068SLisandro Dalcin ierr = PetscStrInList("matmffd",logList,',',&pkg);CHKERRQ(ierr); 618e81d068SLisandro Dalcin if (pkg) {ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);} 623ec795f1SBarry Smith } 633ec795f1SBarry Smith /* Process summary exclusions */ 648e81d068SLisandro Dalcin ierr = PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt);CHKERRQ(ierr); 653ec795f1SBarry Smith if (opt) { 668e81d068SLisandro Dalcin ierr = PetscStrInList("matmffd",logList,',',&pkg);CHKERRQ(ierr); 67fa2bb9feSLisandro Dalcin if (pkg) {ierr = PetscLogEventExcludeClass(MATMFFD_CLASSID);CHKERRQ(ierr);} 683ec795f1SBarry Smith } 698e81d068SLisandro Dalcin /* Register package finalizer */ 70b022a5c1SBarry Smith ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr); 713ec795f1SBarry Smith PetscFunctionReturn(0); 723ec795f1SBarry Smith } 733ec795f1SBarry Smith 74e884886fSBarry Smith /*@C 75e884886fSBarry Smith MatMFFDSetType - Sets the method that is used to compute the 76e884886fSBarry Smith differencing parameter for finite differene matrix-free formulations. 77e884886fSBarry Smith 78e884886fSBarry Smith Input Parameters: 79e884886fSBarry Smith + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD() 80e884886fSBarry Smith or MatSetType(mat,MATMFFD); 81e884886fSBarry Smith - ftype - the type requested, either MATMFFD_WP or MATMFFD_DS 82e884886fSBarry Smith 83e884886fSBarry Smith Level: advanced 84e884886fSBarry Smith 85e884886fSBarry Smith Notes: 86e884886fSBarry Smith For example, such routines can compute h for use in 87e884886fSBarry Smith Jacobian-vector products of the form 88e884886fSBarry Smith 89e884886fSBarry Smith F(x+ha) - F(x) 90e884886fSBarry Smith F'(u)a ~= ---------------- 91e884886fSBarry Smith h 92e884886fSBarry Smith 93755b7f64SBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction(), MatCreateMFFD() 94e884886fSBarry Smith @*/ 9519fd82e9SBarry Smith PetscErrorCode MatMFFDSetType(Mat mat,MatMFFDType ftype) 96e884886fSBarry Smith { 97e884886fSBarry Smith PetscErrorCode ierr,(*r)(MatMFFD); 98e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 99ace3abfcSBarry Smith PetscBool match; 100e884886fSBarry Smith 101e884886fSBarry Smith PetscFunctionBegin; 1020700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 103e884886fSBarry Smith PetscValidCharPointer(ftype,2); 104e884886fSBarry Smith 105251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 1069a13eae7SJed Brown if (!match) PetscFunctionReturn(0); 1079a13eae7SJed Brown 108e884886fSBarry Smith /* already set, so just return */ 109251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr); 110e884886fSBarry Smith if (match) PetscFunctionReturn(0); 111e884886fSBarry Smith 112e884886fSBarry Smith /* destroy the old one if it exists */ 113e884886fSBarry Smith if (ctx->ops->destroy) { 114e884886fSBarry Smith ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr); 115e884886fSBarry Smith } 116e884886fSBarry Smith 1171c9cd337SJed Brown ierr = PetscFunctionListFind(MatMFFDList,ftype,&r);CHKERRQ(ierr); 118e32f2f54SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype); 119e884886fSBarry Smith ierr = (*r)(ctx);CHKERRQ(ierr); 120e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr); 121e884886fSBarry Smith PetscFunctionReturn(0); 122e884886fSBarry Smith } 123e884886fSBarry Smith 1245d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat,Vec); 1255d21e1f8SStefano Zampini 126e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/ 12740244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func) 128e884886fSBarry Smith { 129e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 130e884886fSBarry Smith 131e884886fSBarry Smith PetscFunctionBegin; 132e884886fSBarry Smith ctx->funcisetbase = func; 133486afcceSStefano Zampini /* allow users to compose their own getdiagonal and allow MatHasOperation 134486afcceSStefano Zampini to return false if the two functions pointers are not set */ 135486afcceSStefano Zampini if (!mat->ops->getdiagonal && func) { 136dbfa1270SStefano Zampini mat->ops->getdiagonal = MatGetDiagonal_MFFD; 1375d21e1f8SStefano Zampini } 138e884886fSBarry Smith PetscFunctionReturn(0); 139e884886fSBarry Smith } 140e884886fSBarry Smith 141e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 14240244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci) 143e884886fSBarry Smith { 144e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 145e884886fSBarry Smith 146e884886fSBarry Smith PetscFunctionBegin; 147e884886fSBarry Smith ctx->funci = funci; 148486afcceSStefano Zampini /* allow users to compose their own getdiagonal and allow MatHasOperation 149486afcceSStefano Zampini to return false if the two functions pointers are not set */ 150486afcceSStefano Zampini if (!mat->ops->getdiagonal && funci) { 151dbfa1270SStefano Zampini mat->ops->getdiagonal = MatGetDiagonal_MFFD; 1525d21e1f8SStefano Zampini } 153e884886fSBarry Smith PetscFunctionReturn(0); 154e884886fSBarry Smith } 155e884886fSBarry Smith 15640244768SBarry Smith static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) 157bc0f08ceSBarry Smith { 158bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 159bc0f08ceSBarry Smith 160bc0f08ceSBarry Smith PetscFunctionBegin; 161bc0f08ceSBarry Smith ctx->ncurrenth = 0; 162bc0f08ceSBarry Smith PetscFunctionReturn(0); 163bc0f08ceSBarry Smith } 164e884886fSBarry Smith 1651c84c290SBarry Smith /*@C 1661c84c290SBarry Smith MatMFFDRegister - Adds a method to the MatMFFD registry. 1671c84c290SBarry Smith 1681c84c290SBarry Smith Not Collective 1691c84c290SBarry Smith 1701c84c290SBarry Smith Input Parameters: 1711c84c290SBarry Smith + name_solver - name of a new user-defined compute-h module 1721c84c290SBarry Smith - routine_create - routine to create method context 1731c84c290SBarry Smith 1741c84c290SBarry Smith Level: developer 1751c84c290SBarry Smith 1761c84c290SBarry Smith Notes: 1771c84c290SBarry Smith MatMFFDRegister() may be called multiple times to add several user-defined solvers. 1781c84c290SBarry Smith 1791c84c290SBarry Smith Sample usage: 1801c84c290SBarry Smith .vb 181bdf89e91SBarry Smith MatMFFDRegister("my_h",MyHCreate); 1821c84c290SBarry Smith .ve 1831c84c290SBarry Smith 1841c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 1851c84c290SBarry Smith $ MatMFFDSetType(mfctx,"my_h") 1861c84c290SBarry Smith or at runtime via the option 187be7a6d03SBarry Smith $ -mat_mffd_type my_h 1881c84c290SBarry Smith 1891c84c290SBarry Smith .keywords: MatMFFD, register 1901c84c290SBarry Smith 1911c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 1921c84c290SBarry Smith @*/ 193bdf89e91SBarry Smith PetscErrorCode MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD)) 194e884886fSBarry Smith { 195e884886fSBarry Smith PetscErrorCode ierr; 196e884886fSBarry Smith 197e884886fSBarry Smith PetscFunctionBegin; 1989cc31a68SJed Brown ierr = MatInitializePackage();CHKERRQ(ierr); 199a240a19fSJed Brown ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr); 200e884886fSBarry Smith PetscFunctionReturn(0); 201e884886fSBarry Smith } 202e884886fSBarry Smith 203e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/ 20440244768SBarry Smith static PetscErrorCode MatDestroy_MFFD(Mat mat) 205e884886fSBarry Smith { 206e884886fSBarry Smith PetscErrorCode ierr; 207e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 208e884886fSBarry Smith 209e884886fSBarry Smith PetscFunctionBegin; 2106bf464f9SBarry Smith ierr = VecDestroy(&ctx->w);CHKERRQ(ierr); 2116bf464f9SBarry Smith ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr); 2126bf464f9SBarry Smith ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr); 2136bf464f9SBarry Smith ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr); 214c51ad4d4SStefano Zampini ierr = VecDestroy(&ctx->dshiftw);CHKERRQ(ierr); 215c51ad4d4SStefano Zampini ierr = VecDestroy(&ctx->current_u);CHKERRQ(ierr); 216cfe22f5eSBarry Smith if (ctx->current_f_allocated) { 2176bf464f9SBarry Smith ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr); 218cfe22f5eSBarry Smith } 219e884886fSBarry Smith if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);} 2206bf464f9SBarry Smith ierr = PetscHeaderDestroy(&ctx);CHKERRQ(ierr); 221bf0cc555SLisandro Dalcin mat->data = 0; 222e884886fSBarry Smith 223bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr); 224bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr); 225bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr); 226bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr); 227bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr); 228bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr); 229bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr); 230bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr); 231e884886fSBarry Smith PetscFunctionReturn(0); 232e884886fSBarry Smith } 233e884886fSBarry Smith 234e884886fSBarry Smith /* 235e884886fSBarry Smith MatMFFDView_MFFD - Views matrix-free parameters. 236e884886fSBarry Smith 237e884886fSBarry Smith */ 23840244768SBarry Smith static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer) 239e884886fSBarry Smith { 240e884886fSBarry Smith PetscErrorCode ierr; 241e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 242a283635eSDmitry Karpeev PetscBool iascii, viewbase, viewfunction; 243a283635eSDmitry Karpeev const char *prefix; 244e884886fSBarry Smith 245e884886fSBarry Smith PetscFunctionBegin; 246251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 247e884886fSBarry Smith if (iascii) { 248a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr); 249a283635eSDmitry Karpeev ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 25057622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr); 2517adad957SLisandro Dalcin if (!((PetscObject)ctx)->type_name) { 252e884886fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr); 253e884886fSBarry Smith } else { 2547adad957SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr); 255e884886fSBarry Smith } 256*c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 257*c92e3469SBarry Smith if (ctx->usecomplex) { 258*c92e3469SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Using Lyness complex number trick to compute the matrix-vector product\n");CHKERRQ(ierr); 259*c92e3469SBarry Smith } 260*c92e3469SBarry Smith #endif 261e884886fSBarry Smith if (ctx->ops->view) { 262e884886fSBarry Smith ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr); 263e884886fSBarry Smith } 264a283635eSDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr); 265a283635eSDmitry Karpeev 266c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr); 267a283635eSDmitry Karpeev if (viewbase) { 268a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr); 269a283635eSDmitry Karpeev ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr); 270a283635eSDmitry Karpeev } 271c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr); 272a283635eSDmitry Karpeev if (viewfunction) { 273a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr); 274a283635eSDmitry Karpeev ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr); 275a283635eSDmitry Karpeev } 276a283635eSDmitry Karpeev ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 27711aeaf0aSBarry Smith } 278e884886fSBarry Smith PetscFunctionReturn(0); 279e884886fSBarry Smith } 280e884886fSBarry Smith 281e884886fSBarry Smith /* 282e884886fSBarry Smith MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 283e884886fSBarry Smith allows the user to indicate the beginning of a new linear solve by calling 284e884886fSBarry Smith MatAssemblyXXX() on the matrix free matrix. This then allows the 2851d0fab5eSBarry Smith MatCreateMFFD_WP() to properly compute ||U|| only the first time 286e884886fSBarry Smith in the linear solver rather than every time. 2875a576424SJed Brown 288cc2e6a90SBarry Smith This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence 289cc2e6a90SBarry Smith it must be labeled as PETSC_EXTERN 290e884886fSBarry Smith */ 2915a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt) 292e884886fSBarry Smith { 293e884886fSBarry Smith PetscErrorCode ierr; 294e884886fSBarry Smith MatMFFD j = (MatMFFD)J->data; 295e884886fSBarry Smith 296e884886fSBarry Smith PetscFunctionBegin; 297e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 298e884886fSBarry Smith j->vshift = 0.0; 299e884886fSBarry Smith j->vscale = 1.0; 300e884886fSBarry Smith PetscFunctionReturn(0); 301e884886fSBarry Smith } 302e884886fSBarry Smith 303e884886fSBarry Smith /* 304e884886fSBarry Smith MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 305e884886fSBarry Smith 306e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 307e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 308e884886fSBarry Smith u = current iterate 309e884886fSBarry Smith h = difference interval 310e884886fSBarry Smith */ 31140244768SBarry Smith static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y) 312e884886fSBarry Smith { 313e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 314e884886fSBarry Smith PetscScalar h; 315e884886fSBarry Smith Vec w,U,F; 316e884886fSBarry Smith PetscErrorCode ierr; 317ace3abfcSBarry Smith PetscBool zeroa; 318e884886fSBarry Smith 319e884886fSBarry Smith PetscFunctionBegin; 320ce94432eSBarry 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"); 321e884886fSBarry Smith /* We log matrix-free matrix-vector products separately, so that we can 322e884886fSBarry Smith separate the performance monitoring from the cases that use conventional 323e884886fSBarry Smith storage. We may eventually modify event logging to associate events 324e884886fSBarry Smith with particular objects, hence alleviating the more general problem. */ 325e884886fSBarry Smith ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 326e884886fSBarry Smith 327e884886fSBarry Smith w = ctx->w; 328e884886fSBarry Smith U = ctx->current_u; 3293ec795f1SBarry Smith F = ctx->current_f; 330e884886fSBarry Smith /* 331e884886fSBarry Smith Compute differencing parameter 332e884886fSBarry Smith */ 3334a2f8832SBarry Smith if (!((PetscObject)ctx)->type_name) { 334e884886fSBarry Smith ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr); 3359c6ac3b3SBarry Smith ierr = MatSetFromOptions(mat);CHKERRQ(ierr); 336e884886fSBarry Smith } 337e884886fSBarry Smith ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr); 338e884886fSBarry Smith if (zeroa) { 339e884886fSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 340e884886fSBarry Smith PetscFunctionReturn(0); 341e884886fSBarry Smith } 342e884886fSBarry Smith 34384d44b13SHong Zhang if (mat->erroriffailure && PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h"); 344e884886fSBarry Smith if (ctx->checkh) { 345e884886fSBarry Smith ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr); 346e884886fSBarry Smith } 347e884886fSBarry Smith 348e884886fSBarry Smith /* keep a record of the current differencing parameter h */ 349e884886fSBarry Smith ctx->currenth = h; 350e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX) 35157622a8eSBarry Smith ierr = PetscInfo2(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h));CHKERRQ(ierr); 352e884886fSBarry Smith #else 353e884886fSBarry Smith ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr); 354e884886fSBarry Smith #endif 355e884886fSBarry Smith if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 356e884886fSBarry Smith ctx->historyh[ctx->ncurrenth] = h; 357e884886fSBarry Smith } 358e884886fSBarry Smith ctx->ncurrenth++; 359e884886fSBarry Smith 360*c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 361*c92e3469SBarry Smith if (ctx->usecomplex) h = PETSC_i*h; 362*c92e3469SBarry Smith #endif 363*c92e3469SBarry Smith 364e884886fSBarry Smith /* w = u + ha */ 365c3bb7e23SBarry Smith if (ctx->drscale) { 366c3bb7e23SBarry Smith ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr); 367c3bb7e23SBarry Smith ierr = VecAYPX(U,h,w);CHKERRQ(ierr); 368c3bb7e23SBarry Smith } else { 369e884886fSBarry Smith ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr); 370c3bb7e23SBarry Smith } 371e884886fSBarry Smith 3721b797266SDmitry Karpeev /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 3731b797266SDmitry Karpeev if (ctx->ncurrenth == 1 && ctx->current_f_allocated) { 374e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr); 37552121784SBarry Smith } 376e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr); 377e884886fSBarry Smith 378*c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 379*c92e3469SBarry Smith if (ctx->usecomplex) { 380*c92e3469SBarry Smith ierr = VecImaginaryPart(y);CHKERRQ(ierr); 381*c92e3469SBarry Smith h = PetscImaginaryPart(h); 382*c92e3469SBarry Smith } else { 383e884886fSBarry Smith ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 384*c92e3469SBarry Smith } 385*c92e3469SBarry Smith #else 386*c92e3469SBarry Smith ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 387*c92e3469SBarry Smith #endif 388e884886fSBarry Smith ierr = VecScale(y,1.0/h);CHKERRQ(ierr); 389e884886fSBarry Smith 390c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 391e884886fSBarry Smith ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr); 392c3bb7e23SBarry Smith } 393c3bb7e23SBarry Smith if (ctx->dlscale) { 394c3bb7e23SBarry Smith ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr); 395c3bb7e23SBarry Smith } 396c3bb7e23SBarry Smith if (ctx->dshift) { 397c51ad4d4SStefano Zampini if (!ctx->dshiftw) { 398c51ad4d4SStefano Zampini ierr = VecDuplicate(y,&ctx->dshiftw);CHKERRQ(ierr); 399c51ad4d4SStefano Zampini } 400c51ad4d4SStefano Zampini ierr = VecPointwiseMult(ctx->dshift,a,ctx->dshiftw);CHKERRQ(ierr); 401c51ad4d4SStefano Zampini ierr = VecAXPY(y,1.0,ctx->dshiftw);CHKERRQ(ierr); 402c3bb7e23SBarry Smith } 403e884886fSBarry Smith 40439601f4eSBarry Smith if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);} 405e884886fSBarry Smith 406e884886fSBarry Smith ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 407e884886fSBarry Smith PetscFunctionReturn(0); 408e884886fSBarry Smith } 409e884886fSBarry Smith 410e884886fSBarry Smith /* 411e884886fSBarry Smith MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 412e884886fSBarry Smith 413e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 414e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 415e884886fSBarry Smith u = current iterate 416e884886fSBarry Smith h = difference interval 417e884886fSBarry Smith */ 4185d21e1f8SStefano Zampini PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 419e884886fSBarry Smith { 420e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 421e884886fSBarry Smith PetscScalar h,*aa,*ww,v; 422e884886fSBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 423e884886fSBarry Smith Vec w,U; 424e884886fSBarry Smith PetscErrorCode ierr; 425e884886fSBarry Smith PetscInt i,rstart,rend; 426e884886fSBarry Smith 427e884886fSBarry Smith PetscFunctionBegin; 428486afcceSStefano Zampini if (!ctx->func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunction() first"); 429486afcceSStefano Zampini if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioni() first"); 430486afcceSStefano Zampini if (!ctx->funcisetbase) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioniBase() first"); 431e884886fSBarry Smith w = ctx->w; 432e884886fSBarry Smith U = ctx->current_u; 433e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr); 434e884886fSBarry Smith ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr); 435e884886fSBarry Smith ierr = VecCopy(U,w);CHKERRQ(ierr); 436e884886fSBarry Smith 437e884886fSBarry Smith ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr); 438e884886fSBarry Smith ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 439e884886fSBarry Smith for (i=rstart; i<rend; i++) { 440e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 441e884886fSBarry Smith h = ww[i-rstart]; 442e884886fSBarry Smith if (h == 0.0) h = 1.0; 443e884886fSBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 444e884886fSBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 445e884886fSBarry Smith h *= epsilon; 446e884886fSBarry Smith 447e884886fSBarry Smith ww[i-rstart] += h; 448e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 449e884886fSBarry Smith ierr = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr); 450e884886fSBarry Smith aa[i-rstart] = (v - aa[i-rstart])/h; 451e884886fSBarry Smith 452e884886fSBarry Smith /* possibly shift and scale result */ 453c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 454e884886fSBarry Smith aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart]; 455c3bb7e23SBarry Smith } 456e884886fSBarry Smith 457e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 458e884886fSBarry Smith ww[i-rstart] -= h; 459e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 460e884886fSBarry Smith } 461e884886fSBarry Smith ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 462e884886fSBarry Smith PetscFunctionReturn(0); 463e884886fSBarry Smith } 464e884886fSBarry Smith 46540244768SBarry Smith static PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr) 466c3bb7e23SBarry Smith { 467c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 468c3bb7e23SBarry Smith PetscErrorCode ierr; 469c3bb7e23SBarry Smith 470c3bb7e23SBarry Smith PetscFunctionBegin; 471c3bb7e23SBarry Smith if (ll && !aij->dlscale) { 472c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr); 473c3bb7e23SBarry Smith } 474c3bb7e23SBarry Smith if (rr && !aij->drscale) { 475c3bb7e23SBarry Smith ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr); 476c3bb7e23SBarry Smith } 477c3bb7e23SBarry Smith if (ll) { 478c3bb7e23SBarry Smith ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr); 479c3bb7e23SBarry Smith } 480c3bb7e23SBarry Smith if (rr) { 481c3bb7e23SBarry Smith ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr); 482c3bb7e23SBarry Smith } 483c3bb7e23SBarry Smith PetscFunctionReturn(0); 484c3bb7e23SBarry Smith } 485c3bb7e23SBarry Smith 48640244768SBarry Smith static PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode) 487c3bb7e23SBarry Smith { 488c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 489c3bb7e23SBarry Smith PetscErrorCode ierr; 490c3bb7e23SBarry Smith 491c3bb7e23SBarry Smith PetscFunctionBegin; 492ce94432eSBarry Smith if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES"); 493c3bb7e23SBarry Smith if (!aij->dshift) { 494c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr); 495c3bb7e23SBarry Smith } 496c3bb7e23SBarry Smith ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr); 497c3bb7e23SBarry Smith PetscFunctionReturn(0); 498c3bb7e23SBarry Smith } 499c3bb7e23SBarry Smith 50040244768SBarry Smith static PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a) 501e884886fSBarry Smith { 502e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 5035fd66863SKarl Rupp 504e884886fSBarry Smith PetscFunctionBegin; 505e884886fSBarry Smith shell->vshift += a; 506e884886fSBarry Smith PetscFunctionReturn(0); 507e884886fSBarry Smith } 508e884886fSBarry Smith 50940244768SBarry Smith static PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a) 510e884886fSBarry Smith { 511e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 5125fd66863SKarl Rupp 513e884886fSBarry Smith PetscFunctionBegin; 514e884886fSBarry Smith shell->vscale *= a; 515e884886fSBarry Smith PetscFunctionReturn(0); 516e884886fSBarry Smith } 517e884886fSBarry Smith 518d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F) 519e884886fSBarry Smith { 520e884886fSBarry Smith PetscErrorCode ierr; 521e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 522e884886fSBarry Smith 523e884886fSBarry Smith PetscFunctionBegin; 524e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 525c51ad4d4SStefano Zampini if (!ctx->current_u) { 526c51ad4d4SStefano Zampini ierr = VecDuplicate(U,&ctx->current_u);CHKERRQ(ierr); 527c51ad4d4SStefano Zampini ierr = VecLockPush(ctx->current_u);CHKERRQ(ierr); 528c51ad4d4SStefano Zampini } 529c51ad4d4SStefano Zampini ierr = VecLockPop(ctx->current_u);CHKERRQ(ierr); 530c51ad4d4SStefano Zampini ierr = VecCopy(U,ctx->current_u);CHKERRQ(ierr); 531c51ad4d4SStefano Zampini ierr = VecLockPush(ctx->current_u);CHKERRQ(ierr); 53252121784SBarry Smith if (F) { 5336bf464f9SBarry Smith if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);} 5343ec795f1SBarry Smith ctx->current_f = F; 535cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_FALSE; 536cfe22f5eSBarry Smith } else if (!ctx->current_f_allocated) { 53706c4b6baSBarry Smith ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr); 5382205254eSKarl Rupp 539cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_TRUE; 54052121784SBarry Smith } 541e884886fSBarry Smith if (!ctx->w) { 542e884886fSBarry Smith ierr = VecDuplicate(ctx->current_u,&ctx->w);CHKERRQ(ierr); 543e884886fSBarry Smith } 544e884886fSBarry Smith J->assembled = PETSC_TRUE; 545e884886fSBarry Smith PetscFunctionReturn(0); 546e884886fSBarry Smith } 5475a576424SJed Brown 548e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 549b2573a8aSBarry Smith 55040244768SBarry Smith static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx) 551e884886fSBarry Smith { 552e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 553e884886fSBarry Smith 554e884886fSBarry Smith PetscFunctionBegin; 555e884886fSBarry Smith ctx->checkh = fun; 556e884886fSBarry Smith ctx->checkhctx = ectx; 557e884886fSBarry Smith PetscFunctionReturn(0); 558e884886fSBarry Smith } 559e884886fSBarry Smith 5606aa9148fSLisandro Dalcin /*@C 5616aa9148fSLisandro Dalcin MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 5626aa9148fSLisandro Dalcin MatMFFD options in the database. 5636aa9148fSLisandro Dalcin 5646aa9148fSLisandro Dalcin Collective on Mat 5656aa9148fSLisandro Dalcin 5666aa9148fSLisandro Dalcin Input Parameter: 5676aa9148fSLisandro Dalcin + A - the Mat context 5686aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names 5696aa9148fSLisandro Dalcin 5706aa9148fSLisandro Dalcin Notes: 5716aa9148fSLisandro Dalcin A hyphen (-) must NOT be given at the beginning of the prefix name. 5726aa9148fSLisandro Dalcin The first character of all runtime options is AUTOMATICALLY the hyphen. 5736aa9148fSLisandro Dalcin 5746aa9148fSLisandro Dalcin Level: advanced 5756aa9148fSLisandro Dalcin 5766aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters 5776aa9148fSLisandro Dalcin 578755b7f64SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD() 5796aa9148fSLisandro Dalcin @*/ 5807087cfbeSBarry Smith PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[]) 5816aa9148fSLisandro Dalcin 5826aa9148fSLisandro Dalcin { 5830298fd71SBarry Smith MatMFFD mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL; 5846aa9148fSLisandro Dalcin PetscErrorCode ierr; 5855fd66863SKarl Rupp 5866aa9148fSLisandro Dalcin PetscFunctionBegin; 5870700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 5880700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 5896aa9148fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr); 5906aa9148fSLisandro Dalcin PetscFunctionReturn(0); 5916aa9148fSLisandro Dalcin } 5926aa9148fSLisandro Dalcin 59340244768SBarry Smith static PetscErrorCode MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat) 594e884886fSBarry Smith { 59571f08665SBarry Smith MatMFFD mfctx = (MatMFFD)mat->data; 596e884886fSBarry Smith PetscErrorCode ierr; 597ace3abfcSBarry Smith PetscBool flg; 598e884886fSBarry Smith char ftype[256]; 599e884886fSBarry Smith 600e884886fSBarry Smith PetscFunctionBegin; 6010700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 6020700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 6033194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr); 604a264d7a6SBarry Smith ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr); 605e884886fSBarry Smith if (flg) { 606e884886fSBarry Smith ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr); 607e884886fSBarry Smith } 608e884886fSBarry Smith 609e884886fSBarry Smith ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 610e884886fSBarry Smith ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 611e884886fSBarry Smith 61290d69ab7SBarry Smith flg = PETSC_FALSE; 6130298fd71SBarry Smith ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr); 614e884886fSBarry Smith if (flg) { 615e884886fSBarry Smith ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr); 616e884886fSBarry Smith } 617*c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 618*c92e3469SBarry 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); 619*c92e3469SBarry Smith #endif 620e884886fSBarry Smith if (mfctx->ops->setfromoptions) { 621e55864a3SBarry Smith ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr); 622e884886fSBarry Smith } 623e884886fSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 624e884886fSBarry Smith PetscFunctionReturn(0); 625e884886fSBarry Smith } 626e884886fSBarry Smith 62740244768SBarry Smith static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period) 628bc0f08ceSBarry Smith { 629bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 630bc0f08ceSBarry Smith 631bc0f08ceSBarry Smith PetscFunctionBegin; 632bc0f08ceSBarry Smith ctx->recomputeperiod = period; 633bc0f08ceSBarry Smith PetscFunctionReturn(0); 634bc0f08ceSBarry Smith } 635bc0f08ceSBarry Smith 63640244768SBarry Smith static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 637bc0f08ceSBarry Smith { 638bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 639bc0f08ceSBarry Smith 640bc0f08ceSBarry Smith PetscFunctionBegin; 641bc0f08ceSBarry Smith ctx->func = func; 642bc0f08ceSBarry Smith ctx->funcctx = funcctx; 643bc0f08ceSBarry Smith PetscFunctionReturn(0); 644bc0f08ceSBarry Smith } 645bc0f08ceSBarry Smith 64640244768SBarry Smith static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error) 647bc0f08ceSBarry Smith { 648bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 649bc0f08ceSBarry Smith 650bc0f08ceSBarry Smith PetscFunctionBegin; 651bc0f08ceSBarry Smith if (error != PETSC_DEFAULT) ctx->error_rel = error; 652bc0f08ceSBarry Smith PetscFunctionReturn(0); 653bc0f08ceSBarry Smith } 654bc0f08ceSBarry Smith 6553b49f96aSBarry Smith static PetscErrorCode MatMissingDiagonal_MFFD(Mat A,PetscBool *missing,PetscInt *d) 6563b49f96aSBarry Smith { 6573b49f96aSBarry Smith PetscFunctionBegin; 6583b49f96aSBarry Smith *missing = PETSC_FALSE; 6593b49f96aSBarry Smith PetscFunctionReturn(0); 6603b49f96aSBarry Smith } 6613b49f96aSBarry Smith 662e884886fSBarry Smith /*MC 663e884886fSBarry Smith MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 664e884886fSBarry Smith 665e884886fSBarry Smith Level: advanced 666e884886fSBarry Smith 667755b7f64SBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(), 668755b7f64SBarry Smith MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 669755b7f64SBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 670755b7f64SBarry Smith MatMFFDGetH(), 671e884886fSBarry Smith M*/ 6728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 673e884886fSBarry Smith { 674e884886fSBarry Smith MatMFFD mfctx; 675e884886fSBarry Smith PetscErrorCode ierr; 676e884886fSBarry Smith 677e884886fSBarry Smith PetscFunctionBegin; 678607a6623SBarry Smith ierr = MatMFFDInitializePackage();CHKERRQ(ierr); 679e884886fSBarry Smith 68073107ff1SLisandro Dalcin ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr); 6812205254eSKarl Rupp 682e884886fSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 683e884886fSBarry Smith mfctx->recomputeperiod = 1; 684e884886fSBarry Smith mfctx->count = 0; 685e884886fSBarry Smith mfctx->currenth = 0.0; 6860298fd71SBarry Smith mfctx->historyh = NULL; 687e884886fSBarry Smith mfctx->ncurrenth = 0; 688e884886fSBarry Smith mfctx->maxcurrenth = 0; 6897adad957SLisandro Dalcin ((PetscObject)mfctx)->type_name = 0; 690e884886fSBarry Smith 691e884886fSBarry Smith mfctx->vshift = 0.0; 692e884886fSBarry Smith mfctx->vscale = 1.0; 693e884886fSBarry Smith 694e884886fSBarry Smith /* 695e884886fSBarry Smith Create the empty data structure to contain compute-h routines. 696e884886fSBarry Smith These will be filled in below from the command line options or 697e884886fSBarry Smith a later call with MatMFFDSetType() or if that is not called 698e884886fSBarry Smith then it will default in the first use of MatMult_MFFD() 699e884886fSBarry Smith */ 700e884886fSBarry Smith mfctx->ops->compute = 0; 701e884886fSBarry Smith mfctx->ops->destroy = 0; 702e884886fSBarry Smith mfctx->ops->view = 0; 703e884886fSBarry Smith mfctx->ops->setfromoptions = 0; 704e884886fSBarry Smith mfctx->hctx = 0; 705e884886fSBarry Smith 706e884886fSBarry Smith mfctx->func = 0; 707e884886fSBarry Smith mfctx->funcctx = 0; 7080298fd71SBarry Smith mfctx->w = NULL; 709e884886fSBarry Smith 710e884886fSBarry Smith A->data = mfctx; 711e884886fSBarry Smith 712e884886fSBarry Smith A->ops->mult = MatMult_MFFD; 713e884886fSBarry Smith A->ops->destroy = MatDestroy_MFFD; 714e884886fSBarry Smith A->ops->view = MatView_MFFD; 715e884886fSBarry Smith A->ops->assemblyend = MatAssemblyEnd_MFFD; 716e884886fSBarry Smith A->ops->scale = MatScale_MFFD; 717e884886fSBarry Smith A->ops->shift = MatShift_MFFD; 718c3bb7e23SBarry Smith A->ops->diagonalscale = MatDiagonalScale_MFFD; 719c3bb7e23SBarry Smith A->ops->diagonalset = MatDiagonalSet_MFFD; 7209c6ac3b3SBarry Smith A->ops->setfromoptions = MatSetFromOptions_MFFD; 7213b49f96aSBarry Smith A->ops->missingdiagonal = MatMissingDiagonal_MFFD; 722e884886fSBarry Smith A->assembled = PETSC_TRUE; 723e884886fSBarry Smith 724bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr); 725bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr); 726bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr); 727bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr); 728bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr); 729bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr); 730bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr); 731bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr); 7322205254eSKarl Rupp 733e884886fSBarry Smith mfctx->mat = A; 7342205254eSKarl Rupp 735e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr); 736e884886fSBarry Smith PetscFunctionReturn(0); 737e884886fSBarry Smith } 738e884886fSBarry Smith 739e884886fSBarry Smith /*@ 740e884886fSBarry Smith MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 741e884886fSBarry Smith 742e884886fSBarry Smith Collective on Vec 743e884886fSBarry Smith 744e884886fSBarry Smith Input Parameters: 745fef1beadSBarry Smith + comm - MPI communicator 746fef1beadSBarry Smith . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 747fef1beadSBarry Smith This value should be the same as the local size used in creating the 748fef1beadSBarry Smith y vector for the matrix-vector product y = Ax. 749fef1beadSBarry Smith . n - This value should be the same as the local size used in creating the 750fef1beadSBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 751fef1beadSBarry Smith calculated if N is given) For square matrices n is almost always m. 752fef1beadSBarry Smith . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 753fef1beadSBarry Smith - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 754fef1beadSBarry Smith 755e884886fSBarry Smith 756e884886fSBarry Smith Output Parameter: 757e884886fSBarry Smith . J - the matrix-free matrix 758e884886fSBarry Smith 7599c6ac3b3SBarry Smith Options Database Keys: call MatSetFromOptions() to trigger these 7609c6ac3b3SBarry Smith + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 761*c92e3469SBarry Smith . -mat_mffd_err - square root of estimated relative error in function evaluation 762*c92e3469SBarry Smith . -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 763*c92e3469SBarry Smith . -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values 764*c92e3469SBarry Smith - -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing 765*c92e3469SBarry Smith (requires real valued functions but that PETSc be configured for complex numbers) 7669c6ac3b3SBarry Smith 7679c6ac3b3SBarry Smith 768e884886fSBarry Smith Level: advanced 769e884886fSBarry Smith 770e884886fSBarry Smith Notes: 771e884886fSBarry Smith The matrix-free matrix context merely contains the function pointers 772e884886fSBarry Smith and work space for performing finite difference approximations of 773e884886fSBarry Smith Jacobian-vector products, F'(u)*a, 774e884886fSBarry Smith 775e884886fSBarry Smith The default code uses the following approach to compute h 776e884886fSBarry Smith 777e884886fSBarry Smith .vb 778e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 779e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 780e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 781e884886fSBarry Smith where 782e884886fSBarry Smith error_rel = square root of relative error in function evaluation 783e884886fSBarry Smith umin = minimum iterate parameter 784e884886fSBarry Smith .ve 785e884886fSBarry Smith 786ca93e954SBarry Smith You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different 787ca93e954SBarry Smith preconditioner matrix 788ca93e954SBarry Smith 789e884886fSBarry Smith The user can set the error_rel via MatMFFDSetFunctionError() and 790a7f22e61SSatish Balay umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details. 791e884886fSBarry Smith 792e884886fSBarry Smith The user should call MatDestroy() when finished with the matrix-free 793e884886fSBarry Smith matrix context. 794e884886fSBarry Smith 795e884886fSBarry Smith Options Database Keys: 796e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 797e884886fSBarry Smith . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only) 798e884886fSBarry Smith - -mat_mffd_check_positivity 799e884886fSBarry Smith 800e884886fSBarry Smith .keywords: default, matrix-free, create, matrix 801e884886fSBarry Smith 8021d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 803e884886fSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 80481242352SJed Brown MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian() 805e884886fSBarry Smith 806e884886fSBarry Smith @*/ 8077087cfbeSBarry Smith PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 808e884886fSBarry Smith { 809e884886fSBarry Smith PetscErrorCode ierr; 810e884886fSBarry Smith 811e884886fSBarry Smith PetscFunctionBegin; 812e884886fSBarry Smith ierr = MatCreate(comm,J);CHKERRQ(ierr); 813fef1beadSBarry Smith ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr); 814e884886fSBarry Smith ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 815be7c243fSJed Brown ierr = MatSetUp(*J);CHKERRQ(ierr); 816e884886fSBarry Smith PetscFunctionReturn(0); 817e884886fSBarry Smith } 818e884886fSBarry Smith 819e884886fSBarry Smith /*@ 820e884886fSBarry Smith MatMFFDGetH - Gets the last value that was used as the differencing 821e884886fSBarry Smith parameter. 822e884886fSBarry Smith 823e884886fSBarry Smith Not Collective 824e884886fSBarry Smith 825e884886fSBarry Smith Input Parameters: 826e884886fSBarry Smith . mat - the matrix obtained with MatCreateSNESMF() 827e884886fSBarry Smith 828e884886fSBarry Smith Output Paramter: 829e884886fSBarry Smith . h - the differencing step size 830e884886fSBarry Smith 831e884886fSBarry Smith Level: advanced 832e884886fSBarry Smith 833e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 834e884886fSBarry Smith 8351d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory() 836e884886fSBarry Smith @*/ 8377087cfbeSBarry Smith PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 838e884886fSBarry Smith { 839e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 840bc0f08ceSBarry Smith PetscErrorCode ierr; 841bc0f08ceSBarry Smith PetscBool match; 842e884886fSBarry Smith 843e884886fSBarry Smith PetscFunctionBegin; 84488b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 84588b4c220SStefano Zampini PetscValidPointer(h,2); 846251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 847ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 848bc0f08ceSBarry Smith 849e884886fSBarry Smith *h = ctx->currenth; 850e884886fSBarry Smith PetscFunctionReturn(0); 851e884886fSBarry Smith } 852e884886fSBarry Smith 853e884886fSBarry Smith /*@C 854e884886fSBarry Smith MatMFFDSetFunction - Sets the function used in applying the matrix free. 855e884886fSBarry Smith 8563f9fe445SBarry Smith Logically Collective on Mat 857e884886fSBarry Smith 858e884886fSBarry Smith Input Parameters: 85914f633e2SBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD() 860e884886fSBarry Smith . func - the function to use 861e884886fSBarry Smith - funcctx - optional function context passed to function 862e884886fSBarry Smith 86314f633e2SBarry Smith Calling Sequence of func: 86414f633e2SBarry Smith $ func (void *funcctx, Vec x, Vec f) 86514f633e2SBarry Smith 86614f633e2SBarry Smith + funcctx - user provided context 86714f633e2SBarry Smith . x - input vector 86814f633e2SBarry Smith - f - computed output function 86914f633e2SBarry Smith 870e884886fSBarry Smith Level: advanced 871e884886fSBarry Smith 872e884886fSBarry Smith Notes: 873e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 874e884886fSBarry Smith matrix inside your compute Jacobian routine 875e884886fSBarry Smith 8763ec795f1SBarry Smith If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 877e884886fSBarry Smith 878e884886fSBarry Smith .keywords: SNES, matrix-free, function 879e884886fSBarry Smith 8801d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD, 8811d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 882e884886fSBarry Smith @*/ 8837087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 884e884886fSBarry Smith { 885bc0f08ceSBarry Smith PetscErrorCode ierr; 886e884886fSBarry Smith 887e884886fSBarry Smith PetscFunctionBegin; 88888b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 889bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr); 890e884886fSBarry Smith PetscFunctionReturn(0); 891e884886fSBarry Smith } 892e884886fSBarry Smith 893e884886fSBarry Smith /*@C 894e884886fSBarry Smith MatMFFDSetFunctioni - Sets the function for a single component 895e884886fSBarry Smith 8963f9fe445SBarry Smith Logically Collective on Mat 897e884886fSBarry Smith 898e884886fSBarry Smith Input Parameters: 899e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 900e884886fSBarry Smith - funci - the function to use 901e884886fSBarry Smith 902e884886fSBarry Smith Level: advanced 903e884886fSBarry Smith 904e884886fSBarry Smith Notes: 905e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 90694eb4328SStefano Zampini matrix inside your compute Jacobian routine. 90794eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 908486afcceSStefano Zampini funci must not contain any MPI call as it is called inside a loop on the local portion of the vector. 909e884886fSBarry Smith 910e884886fSBarry Smith .keywords: SNES, matrix-free, function 911e884886fSBarry Smith 91294eb4328SStefano Zampini .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 9131d0fab5eSBarry Smith 914e884886fSBarry Smith @*/ 9157087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 916e884886fSBarry Smith { 9174ac538c5SBarry Smith PetscErrorCode ierr; 918e884886fSBarry Smith 919e884886fSBarry Smith PetscFunctionBegin; 9200700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 9214ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr); 922e884886fSBarry Smith PetscFunctionReturn(0); 923e884886fSBarry Smith } 924e884886fSBarry Smith 925e884886fSBarry Smith /*@C 926e884886fSBarry Smith MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 927e884886fSBarry Smith 9283f9fe445SBarry Smith Logically Collective on Mat 929e884886fSBarry Smith 930e884886fSBarry Smith Input Parameters: 931e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 932e884886fSBarry Smith - func - the function to use 933e884886fSBarry Smith 934e884886fSBarry Smith Level: advanced 935e884886fSBarry Smith 936e884886fSBarry Smith Notes: 937e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 93894eb4328SStefano Zampini matrix inside your compute Jacobian routine. 93994eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 940e884886fSBarry Smith 941e884886fSBarry Smith 942e884886fSBarry Smith .keywords: SNES, matrix-free, function 943e884886fSBarry Smith 944e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 94594eb4328SStefano Zampini MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 946e884886fSBarry Smith @*/ 9477087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 948e884886fSBarry Smith { 9494ac538c5SBarry Smith PetscErrorCode ierr; 950e884886fSBarry Smith 951e884886fSBarry Smith PetscFunctionBegin; 9520700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 9534ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr); 954e884886fSBarry Smith PetscFunctionReturn(0); 955e884886fSBarry Smith } 956e884886fSBarry Smith 957e884886fSBarry Smith /*@ 958e884886fSBarry Smith MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 959e884886fSBarry Smith 9603f9fe445SBarry Smith Logically Collective on Mat 961e884886fSBarry Smith 962e884886fSBarry Smith Input Parameters: 963e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 964e884886fSBarry Smith - period - 1 for everytime, 2 for every second etc 965e884886fSBarry Smith 966e884886fSBarry Smith Options Database Keys: 967e884886fSBarry Smith + -mat_mffd_period <period> 968e884886fSBarry Smith 969e884886fSBarry Smith Level: advanced 970e884886fSBarry Smith 971e884886fSBarry Smith 972e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 973e884886fSBarry Smith 974e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), 9751d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 976e884886fSBarry Smith @*/ 9777087cfbeSBarry Smith PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 978e884886fSBarry Smith { 979bc0f08ceSBarry Smith PetscErrorCode ierr; 980e884886fSBarry Smith 981e884886fSBarry Smith PetscFunctionBegin; 98288b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 98388b4c220SStefano Zampini PetscValidLogicalCollectiveInt(mat,period,2); 984bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr); 985e884886fSBarry Smith PetscFunctionReturn(0); 986e884886fSBarry Smith } 987e884886fSBarry Smith 988e884886fSBarry Smith /*@ 989e884886fSBarry Smith MatMFFDSetFunctionError - Sets the error_rel for the approximation of 990e884886fSBarry Smith matrix-vector products using finite differences. 991e884886fSBarry Smith 9923f9fe445SBarry Smith Logically Collective on Mat 993e884886fSBarry Smith 994e884886fSBarry Smith Input Parameters: 995e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 996e884886fSBarry Smith - error_rel - relative error (should be set to the square root of 997e884886fSBarry Smith the relative error in the function evaluations) 998e884886fSBarry Smith 999e884886fSBarry Smith Options Database Keys: 1000e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 1001e884886fSBarry Smith 1002e884886fSBarry Smith Level: advanced 1003e884886fSBarry Smith 1004e884886fSBarry Smith Notes: 1005e884886fSBarry Smith The default matrix-free matrix-vector product routine computes 1006e884886fSBarry Smith .vb 1007e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 1008e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 1009e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 1010e884886fSBarry Smith .ve 1011e884886fSBarry Smith 1012e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 1013e884886fSBarry Smith 1014e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 10151d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 1016e884886fSBarry Smith @*/ 10177087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 1018e884886fSBarry Smith { 1019bc0f08ceSBarry Smith PetscErrorCode ierr; 1020e884886fSBarry Smith 1021e884886fSBarry Smith PetscFunctionBegin; 102288b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 102388b4c220SStefano Zampini PetscValidLogicalCollectiveReal(mat,error,2); 1024bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr); 1025e884886fSBarry Smith PetscFunctionReturn(0); 1026e884886fSBarry Smith } 1027e884886fSBarry Smith 1028e884886fSBarry Smith /*@ 1029e884886fSBarry Smith MatMFFDSetHHistory - Sets an array to collect a history of the 1030e884886fSBarry Smith differencing values (h) computed for the matrix-free product. 1031e884886fSBarry Smith 10323f9fe445SBarry Smith Logically Collective on Mat 1033e884886fSBarry Smith 1034e884886fSBarry Smith Input Parameters: 1035e884886fSBarry Smith + J - the matrix-free matrix context 1036e884886fSBarry Smith . histroy - space to hold the history 1037e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than 1038e884886fSBarry Smith nhistory, then the later ones are discarded 1039e884886fSBarry Smith 1040e884886fSBarry Smith Level: advanced 1041e884886fSBarry Smith 1042e884886fSBarry Smith Notes: 1043e884886fSBarry Smith Use MatMFFDResetHHistory() to reset the history counter and collect 1044e884886fSBarry Smith a new batch of differencing parameters, h. 1045e884886fSBarry Smith 1046e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1047e884886fSBarry Smith 1048e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 10491d0fab5eSBarry Smith MatMFFDResetHHistory(), MatMFFDSetFunctionError() 1050e884886fSBarry Smith 1051e884886fSBarry Smith @*/ 10527087cfbeSBarry Smith PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 1053e884886fSBarry Smith { 1054e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 1055bc0f08ceSBarry Smith PetscErrorCode ierr; 1056bc0f08ceSBarry Smith PetscBool match; 1057e884886fSBarry Smith 1058e884886fSBarry Smith PetscFunctionBegin; 105988b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 106088b4c220SStefano Zampini if (history) PetscValidPointer(history,2); 106188b4c220SStefano Zampini PetscValidLogicalCollectiveInt(J,nhistory,3); 1062251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr); 1063ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 1064e884886fSBarry Smith ctx->historyh = history; 1065e884886fSBarry Smith ctx->maxcurrenth = nhistory; 106675567043SBarry Smith ctx->currenth = 0.; 1067e884886fSBarry Smith PetscFunctionReturn(0); 1068e884886fSBarry Smith } 1069e884886fSBarry Smith 1070e884886fSBarry Smith /*@ 1071e884886fSBarry Smith MatMFFDResetHHistory - Resets the counter to zero to begin 1072e884886fSBarry Smith collecting a new set of differencing histories. 1073e884886fSBarry Smith 10743f9fe445SBarry Smith Logically Collective on Mat 1075e884886fSBarry Smith 1076e884886fSBarry Smith Input Parameters: 1077e884886fSBarry Smith . J - the matrix-free matrix context 1078e884886fSBarry Smith 1079e884886fSBarry Smith Level: advanced 1080e884886fSBarry Smith 1081e884886fSBarry Smith Notes: 1082e884886fSBarry Smith Use MatMFFDSetHHistory() to create the original history counter. 1083e884886fSBarry Smith 1084e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1085e884886fSBarry Smith 1086e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 10871d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDSetFunctionError() 1088e884886fSBarry Smith 1089e884886fSBarry Smith @*/ 10907087cfbeSBarry Smith PetscErrorCode MatMFFDResetHHistory(Mat J) 1091e884886fSBarry Smith { 1092bc0f08ceSBarry Smith PetscErrorCode ierr; 1093e884886fSBarry Smith 1094e884886fSBarry Smith PetscFunctionBegin; 109588b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1096bc0f08ceSBarry Smith ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr); 1097e884886fSBarry Smith PetscFunctionReturn(0); 1098e884886fSBarry Smith } 1099e884886fSBarry Smith 1100487a658cSBarry Smith /*@ 1101e884886fSBarry Smith MatMFFDSetBase - Sets the vector U at which matrix vector products of the 1102e884886fSBarry Smith Jacobian are computed 1103e884886fSBarry Smith 11043f9fe445SBarry Smith Logically Collective on Mat 1105e884886fSBarry Smith 1106e884886fSBarry Smith Input Parameters: 1107e884886fSBarry Smith + J - the MatMFFD matrix 11083ec795f1SBarry Smith . U - the vector 1109bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed 1110e884886fSBarry Smith 111195452b02SPatrick Sanan Notes: 111295452b02SPatrick Sanan This is rarely used directly 1113e884886fSBarry Smith 11148af5ae88SBarry Smith If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 11158af5ae88SBarry Smith point during the first MatMult() after each call to MatMFFDSetBase(). 1116dff2f722SBarry Smith 1117e884886fSBarry Smith Level: advanced 1118e884886fSBarry Smith 1119e884886fSBarry Smith @*/ 11207087cfbeSBarry Smith PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 1121e884886fSBarry Smith { 11224ac538c5SBarry Smith PetscErrorCode ierr; 1123e884886fSBarry Smith 1124e884886fSBarry Smith PetscFunctionBegin; 11250700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 11260700a824SBarry Smith PetscValidHeaderSpecific(U,VEC_CLASSID,2); 11270700a824SBarry Smith if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 11284ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr); 1129e884886fSBarry Smith PetscFunctionReturn(0); 1130e884886fSBarry Smith } 1131e884886fSBarry Smith 1132e884886fSBarry Smith /*@C 1133e884886fSBarry Smith MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 1134e884886fSBarry Smith it to satisfy some criteria 1135e884886fSBarry Smith 11363f9fe445SBarry Smith Logically Collective on Mat 1137e884886fSBarry Smith 1138e884886fSBarry Smith Input Parameters: 1139e884886fSBarry Smith + J - the MatMFFD matrix 1140e884886fSBarry Smith . fun - the function that checks h 1141e884886fSBarry Smith - ctx - any context needed by the function 1142e884886fSBarry Smith 1143e884886fSBarry Smith Options Database Keys: 1144e884886fSBarry Smith . -mat_mffd_check_positivity 1145e884886fSBarry Smith 1146e884886fSBarry Smith Level: advanced 1147e884886fSBarry Smith 114895452b02SPatrick Sanan Notes: 114995452b02SPatrick Sanan For example, MatMFFDCheckPositivity() insures that all entries 1150e884886fSBarry Smith of U + h*a are non-negative 1151e884886fSBarry Smith 1152755b7f64SBarry Smith The function you provide is called after the default h has been computed and allows you to 1153755b7f64SBarry Smith modify it. 1154755b7f64SBarry Smith 1155755b7f64SBarry Smith .seealso: MatMFFDCheckPositivity() 1156e884886fSBarry Smith @*/ 11577087cfbeSBarry Smith PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx) 1158e884886fSBarry Smith { 11594ac538c5SBarry Smith PetscErrorCode ierr; 1160e884886fSBarry Smith 1161e884886fSBarry Smith PetscFunctionBegin; 11620700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 11634ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr); 1164e884886fSBarry Smith PetscFunctionReturn(0); 1165e884886fSBarry Smith } 1166e884886fSBarry Smith 1167e884886fSBarry Smith /*@ 1168e884886fSBarry Smith MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1169e884886fSBarry Smith zero, decreases h until this is satisfied. 1170e884886fSBarry Smith 11713f9fe445SBarry Smith Logically Collective on Vec 1172e884886fSBarry Smith 1173e884886fSBarry Smith Input Parameters: 1174e884886fSBarry Smith + U - base vector that is added to 1175e884886fSBarry Smith . a - vector that is added 1176e884886fSBarry Smith . h - scaling factor on a 1177e884886fSBarry Smith - dummy - context variable (unused) 1178e884886fSBarry Smith 1179e884886fSBarry Smith Options Database Keys: 1180e884886fSBarry Smith . -mat_mffd_check_positivity 1181e884886fSBarry Smith 1182e884886fSBarry Smith Level: advanced 1183e884886fSBarry Smith 118495452b02SPatrick Sanan Notes: 118595452b02SPatrick Sanan This is rarely used directly, rather it is passed as an argument to 1186e884886fSBarry Smith MatMFFDSetCheckh() 1187e884886fSBarry Smith 1188e884886fSBarry Smith .seealso: MatMFFDSetCheckh() 1189e884886fSBarry Smith @*/ 11907087cfbeSBarry Smith PetscErrorCode MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h) 1191e884886fSBarry Smith { 1192e884886fSBarry Smith PetscReal val, minval; 1193e884886fSBarry Smith PetscScalar *u_vec, *a_vec; 1194e884886fSBarry Smith PetscErrorCode ierr; 1195e884886fSBarry Smith PetscInt i,n; 1196e884886fSBarry Smith MPI_Comm comm; 1197e884886fSBarry Smith 1198e884886fSBarry Smith PetscFunctionBegin; 119988b4c220SStefano Zampini PetscValidHeaderSpecific(U,VEC_CLASSID,2); 120088b4c220SStefano Zampini PetscValidHeaderSpecific(a,VEC_CLASSID,3); 120188b4c220SStefano Zampini PetscValidPointer(h,4); 1202e884886fSBarry Smith ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); 1203e884886fSBarry Smith ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); 1204e884886fSBarry Smith ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); 1205e884886fSBarry Smith ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); 1206c068d9bbSLisandro Dalcin minval = PetscAbsScalar(*h)*PetscRealConstant(1.01); 1207e884886fSBarry Smith for (i=0; i<n; i++) { 1208e884886fSBarry Smith if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1209e884886fSBarry Smith val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1210e884886fSBarry Smith if (val < minval) minval = val; 1211e884886fSBarry Smith } 1212e884886fSBarry Smith } 1213e884886fSBarry Smith ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); 1214e884886fSBarry Smith ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); 1215b2566f29SBarry Smith ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr); 1216e884886fSBarry Smith if (val <= PetscAbsScalar(*h)) { 12176712e2f1SBarry Smith ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr); 1218e884886fSBarry Smith if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1219e884886fSBarry Smith else *h = -0.99*val; 1220e884886fSBarry Smith } 1221e884886fSBarry Smith PetscFunctionReturn(0); 1222e884886fSBarry Smith } 1223