1e884886fSBarry Smith 2b45d2f2cSJed Brown #include <petsc-private/matimpl.h> 3c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h> /*I "petscmat.h" I*/ 4e884886fSBarry Smith 5140e18c1SBarry Smith PetscFunctionList MatMFFDList = 0; 6ace3abfcSBarry Smith PetscBool MatMFFDRegisterAllCalled = PETSC_FALSE; 7e884886fSBarry Smith 87087cfbeSBarry Smith PetscClassId MATMFFD_CLASSID; 9166c7f25SBarry Smith PetscLogEvent MATMFFD_Mult; 10e884886fSBarry Smith 11ace3abfcSBarry Smith static PetscBool MatMFFDPackageInitialized = PETSC_FALSE; 12b022a5c1SBarry Smith #undef __FUNCT__ 13b022a5c1SBarry Smith #define __FUNCT__ "MatMFFDFinalizePackage" 14b022a5c1SBarry Smith /*@C 152390153bSJed Brown MatMFFDFinalizePackage - This function destroys everything in the MatMFFD package. It is 16b022a5c1SBarry Smith called from PetscFinalize(). 17b022a5c1SBarry Smith 18b022a5c1SBarry Smith Level: developer 19b022a5c1SBarry Smith 202390153bSJed Brown .keywords: Petsc, destroy, package 21b022a5c1SBarry Smith .seealso: PetscFinalize() 22b022a5c1SBarry Smith @*/ 237087cfbeSBarry Smith PetscErrorCode MatMFFDFinalizePackage(void) 24b022a5c1SBarry Smith { 25a703d84cSPatrick Lacasse PetscErrorCode ierr; 26a703d84cSPatrick Lacasse 27b022a5c1SBarry Smith PetscFunctionBegin; 28b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_FALSE; 29b022a5c1SBarry Smith MatMFFDRegisterAllCalled = PETSC_FALSE; 30a703d84cSPatrick Lacasse ierr = MatMFFDRegisterDestroy();CHKERRQ(ierr); 310298fd71SBarry Smith MatMFFDList = NULL; 32b022a5c1SBarry Smith PetscFunctionReturn(0); 33b022a5c1SBarry Smith } 34b022a5c1SBarry Smith 35e884886fSBarry Smith #undef __FUNCT__ 363ec795f1SBarry Smith #define __FUNCT__ "MatMFFDInitializePackage" 373ec795f1SBarry Smith /*@C 383ec795f1SBarry Smith MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called 393ec795f1SBarry Smith from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to MatCreate_MFFD() 403ec795f1SBarry Smith when using static libraries. 413ec795f1SBarry Smith 423ec795f1SBarry Smith Input Parameter: 430298fd71SBarry Smith . path - The dynamic library path, or NULL 443ec795f1SBarry Smith 453ec795f1SBarry Smith Level: developer 463ec795f1SBarry Smith 473ec795f1SBarry Smith .keywords: Vec, initialize, package 483ec795f1SBarry Smith .seealso: PetscInitialize() 493ec795f1SBarry Smith @*/ 507087cfbeSBarry Smith PetscErrorCode MatMFFDInitializePackage(const char path[]) 513ec795f1SBarry Smith { 523ec795f1SBarry Smith char logList[256]; 533ec795f1SBarry Smith char *className; 54ace3abfcSBarry Smith PetscBool opt; 553ec795f1SBarry Smith PetscErrorCode ierr; 563ec795f1SBarry Smith 573ec795f1SBarry Smith PetscFunctionBegin; 58b022a5c1SBarry Smith if (MatMFFDPackageInitialized) PetscFunctionReturn(0); 59b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_TRUE; 603ec795f1SBarry Smith /* Register Classes */ 610700a824SBarry Smith ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr); 623ec795f1SBarry Smith /* Register Constructors */ 633ec795f1SBarry Smith ierr = MatMFFDRegisterAll(path);CHKERRQ(ierr); 643ec795f1SBarry Smith /* Register Events */ 650700a824SBarry Smith ierr = PetscLogEventRegister("MatMult MF", MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr); 663ec795f1SBarry Smith 673ec795f1SBarry Smith /* Process info exclusions */ 680298fd71SBarry Smith ierr = PetscOptionsGetString(NULL, "-info_exclude", logList, 256, &opt);CHKERRQ(ierr); 693ec795f1SBarry Smith if (opt) { 703ec795f1SBarry Smith ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr); 713ec795f1SBarry Smith if (className) { 720700a824SBarry Smith ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr); 733ec795f1SBarry Smith } 743ec795f1SBarry Smith } 753ec795f1SBarry Smith /* Process summary exclusions */ 760298fd71SBarry Smith ierr = PetscOptionsGetString(NULL, "-log_summary_exclude", logList, 256, &opt);CHKERRQ(ierr); 773ec795f1SBarry Smith if (opt) { 783ec795f1SBarry Smith ierr = PetscStrstr(logList, "matmffd", &className);CHKERRQ(ierr); 793ec795f1SBarry Smith if (className) { 800700a824SBarry Smith ierr = PetscLogEventDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr); 813ec795f1SBarry Smith } 823ec795f1SBarry Smith } 83b022a5c1SBarry Smith ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr); 843ec795f1SBarry Smith PetscFunctionReturn(0); 853ec795f1SBarry Smith } 863ec795f1SBarry Smith 873ec795f1SBarry Smith #undef __FUNCT__ 88e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetType" 89e884886fSBarry Smith /*@C 90e884886fSBarry Smith MatMFFDSetType - Sets the method that is used to compute the 91e884886fSBarry Smith differencing parameter for finite differene matrix-free formulations. 92e884886fSBarry Smith 93e884886fSBarry Smith Input Parameters: 94e884886fSBarry Smith + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD() 95e884886fSBarry Smith or MatSetType(mat,MATMFFD); 96e884886fSBarry Smith - ftype - the type requested, either MATMFFD_WP or MATMFFD_DS 97e884886fSBarry Smith 98e884886fSBarry Smith Level: advanced 99e884886fSBarry Smith 100e884886fSBarry Smith Notes: 101e884886fSBarry Smith For example, such routines can compute h for use in 102e884886fSBarry Smith Jacobian-vector products of the form 103e884886fSBarry Smith 104e884886fSBarry Smith F(x+ha) - F(x) 105e884886fSBarry Smith F'(u)a ~= ---------------- 106e884886fSBarry Smith h 107e884886fSBarry Smith 108d2d6cebeSBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegisterDynamic(), MatMFFDSetFunction() 109e884886fSBarry Smith @*/ 11019fd82e9SBarry Smith PetscErrorCode MatMFFDSetType(Mat mat,MatMFFDType ftype) 111e884886fSBarry Smith { 112e884886fSBarry Smith PetscErrorCode ierr,(*r)(MatMFFD); 113e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 114ace3abfcSBarry Smith PetscBool match; 115e884886fSBarry Smith 116e884886fSBarry Smith PetscFunctionBegin; 1170700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 118e884886fSBarry Smith PetscValidCharPointer(ftype,2); 119e884886fSBarry Smith 120251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 1219a13eae7SJed Brown if (!match) PetscFunctionReturn(0); 1229a13eae7SJed Brown 123e884886fSBarry Smith /* already set, so just return */ 124251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr); 125e884886fSBarry Smith if (match) PetscFunctionReturn(0); 126e884886fSBarry Smith 127e884886fSBarry Smith /* destroy the old one if it exists */ 128e884886fSBarry Smith if (ctx->ops->destroy) { 129e884886fSBarry Smith ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr); 130e884886fSBarry Smith } 131e884886fSBarry Smith 132ce94432eSBarry Smith ierr = PetscFunctionListFind(PetscObjectComm((PetscObject)ctx),MatMFFDList,ftype,PETSC_TRUE,(void (**)(void)) &r);CHKERRQ(ierr); 133e32f2f54SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype); 134e884886fSBarry Smith ierr = (*r)(ctx);CHKERRQ(ierr); 135e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr); 136e884886fSBarry Smith PetscFunctionReturn(0); 137e884886fSBarry Smith } 138e884886fSBarry Smith 139e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/ 140e884886fSBarry Smith #undef __FUNCT__ 141c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase_MFFD" 1427087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func) 143e884886fSBarry Smith { 144e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 145e884886fSBarry Smith 146e884886fSBarry Smith PetscFunctionBegin; 147e884886fSBarry Smith ctx->funcisetbase = func; 148e884886fSBarry Smith PetscFunctionReturn(0); 149e884886fSBarry Smith } 150e884886fSBarry Smith 151e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 152e884886fSBarry Smith #undef __FUNCT__ 153c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni_MFFD" 1547087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci) 155e884886fSBarry Smith { 156e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 157e884886fSBarry Smith 158e884886fSBarry Smith PetscFunctionBegin; 159e884886fSBarry Smith ctx->funci = funci; 160e884886fSBarry Smith PetscFunctionReturn(0); 161e884886fSBarry Smith } 162e884886fSBarry Smith 163bc0f08ceSBarry Smith #undef __FUNCT__ 164bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory_MFFD" 165bc0f08ceSBarry Smith PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) 166bc0f08ceSBarry Smith { 167bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 168bc0f08ceSBarry Smith 169bc0f08ceSBarry Smith PetscFunctionBegin; 170bc0f08ceSBarry Smith ctx->ncurrenth = 0; 171bc0f08ceSBarry Smith PetscFunctionReturn(0); 172bc0f08ceSBarry Smith } 173e884886fSBarry Smith 174e884886fSBarry Smith #undef __FUNCT__ 175e884886fSBarry Smith #define __FUNCT__ "MatMFFDRegister" 1767087cfbeSBarry Smith PetscErrorCode MatMFFDRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(MatMFFD)) 177e884886fSBarry Smith { 178e884886fSBarry Smith PetscErrorCode ierr; 179e884886fSBarry Smith char fullname[PETSC_MAX_PATH_LEN]; 180e884886fSBarry Smith 181e884886fSBarry Smith PetscFunctionBegin; 182140e18c1SBarry Smith ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr); 183140e18c1SBarry Smith ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&MatMFFDList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 184e884886fSBarry Smith PetscFunctionReturn(0); 185e884886fSBarry Smith } 186e884886fSBarry Smith 187e884886fSBarry Smith 188e884886fSBarry Smith #undef __FUNCT__ 189e884886fSBarry Smith #define __FUNCT__ "MatMFFDRegisterDestroy" 190e884886fSBarry Smith /*@C 191e884886fSBarry Smith MatMFFDRegisterDestroy - Frees the list of MatMFFD methods that were 192e884886fSBarry Smith registered by MatMFFDRegisterDynamic). 193e884886fSBarry Smith 194e884886fSBarry Smith Not Collective 195e884886fSBarry Smith 196e884886fSBarry Smith Level: developer 197e884886fSBarry Smith 198e884886fSBarry Smith .keywords: MatMFFD, register, destroy 199e884886fSBarry Smith 200e884886fSBarry Smith .seealso: MatMFFDRegisterDynamic), MatMFFDRegisterAll() 201e884886fSBarry Smith @*/ 2027087cfbeSBarry Smith PetscErrorCode MatMFFDRegisterDestroy(void) 203e884886fSBarry Smith { 204e884886fSBarry Smith PetscErrorCode ierr; 205e884886fSBarry Smith 206e884886fSBarry Smith PetscFunctionBegin; 207140e18c1SBarry Smith ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr); 2082205254eSKarl Rupp 209e884886fSBarry Smith MatMFFDRegisterAllCalled = PETSC_FALSE; 210e884886fSBarry Smith PetscFunctionReturn(0); 211e884886fSBarry Smith } 212e884886fSBarry Smith 213bc0f08ceSBarry Smith #undef __FUNCT__ 214bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDAddNullSpace_MFFD" 215bc0f08ceSBarry Smith PetscErrorCode MatMFFDAddNullSpace_MFFD(Mat J,MatNullSpace nullsp) 216bc0f08ceSBarry Smith { 217bc0f08ceSBarry Smith PetscErrorCode ierr; 218bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 219bc0f08ceSBarry Smith 220bc0f08ceSBarry Smith PetscFunctionBegin; 221bc0f08ceSBarry Smith ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr); 222bc0f08ceSBarry Smith if (ctx->sp) { ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr); } 223bc0f08ceSBarry Smith ctx->sp = nullsp; 224bc0f08ceSBarry Smith PetscFunctionReturn(0); 225bc0f08ceSBarry Smith } 226bc0f08ceSBarry Smith 227e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/ 228e884886fSBarry Smith #undef __FUNCT__ 229e884886fSBarry Smith #define __FUNCT__ "MatDestroy_MFFD" 230e884886fSBarry Smith PetscErrorCode MatDestroy_MFFD(Mat mat) 231e884886fSBarry Smith { 232e884886fSBarry Smith PetscErrorCode ierr; 233e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 234e884886fSBarry Smith 235e884886fSBarry Smith PetscFunctionBegin; 2366bf464f9SBarry Smith ierr = VecDestroy(&ctx->w);CHKERRQ(ierr); 2376bf464f9SBarry Smith ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr); 2386bf464f9SBarry Smith ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr); 2396bf464f9SBarry Smith ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr); 240cfe22f5eSBarry Smith if (ctx->current_f_allocated) { 2416bf464f9SBarry Smith ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr); 242cfe22f5eSBarry Smith } 243e884886fSBarry Smith if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);} 2446bf464f9SBarry Smith ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr); 2456bf464f9SBarry Smith ierr = PetscHeaderDestroy(&ctx);CHKERRQ(ierr); 246bf0cc555SLisandro Dalcin mat->data = 0; 247e884886fSBarry Smith 2480298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C","",NULL);CHKERRQ(ierr); 2490298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C","",NULL);CHKERRQ(ierr); 2500298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C","",NULL);CHKERRQ(ierr); 2510298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C","",NULL);CHKERRQ(ierr); 2520298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C","",NULL);CHKERRQ(ierr); 2530298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C","",NULL);CHKERRQ(ierr); 2540298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C","",NULL);CHKERRQ(ierr); 2550298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C","",NULL);CHKERRQ(ierr); 2560298fd71SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDAddNullSpace_C","",NULL);CHKERRQ(ierr); 257e884886fSBarry Smith PetscFunctionReturn(0); 258e884886fSBarry Smith } 259e884886fSBarry Smith 260e884886fSBarry Smith #undef __FUNCT__ 261e884886fSBarry Smith #define __FUNCT__ "MatView_MFFD" 262e884886fSBarry Smith /* 263e884886fSBarry Smith MatMFFDView_MFFD - Views matrix-free parameters. 264e884886fSBarry Smith 265e884886fSBarry Smith */ 266e884886fSBarry Smith PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer) 267e884886fSBarry Smith { 268e884886fSBarry Smith PetscErrorCode ierr; 269e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 270a283635eSDmitry Karpeev PetscBool iascii, viewbase, viewfunction; 271a283635eSDmitry Karpeev const char *prefix; 272e884886fSBarry Smith 273e884886fSBarry Smith PetscFunctionBegin; 274251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 275e884886fSBarry Smith if (iascii) { 276a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr); 277a283635eSDmitry Karpeev ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 278e884886fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr); 2797adad957SLisandro Dalcin if (!((PetscObject)ctx)->type_name) { 280e884886fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr); 281e884886fSBarry Smith } else { 2827adad957SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr); 283e884886fSBarry Smith } 284e884886fSBarry Smith if (ctx->ops->view) { 285e884886fSBarry Smith ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr); 286e884886fSBarry Smith } 287a283635eSDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr); 288a283635eSDmitry Karpeev 289a283635eSDmitry Karpeev ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr); 290a283635eSDmitry Karpeev if (viewbase) { 291a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr); 292a283635eSDmitry Karpeev ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr); 293a283635eSDmitry Karpeev } 294a283635eSDmitry Karpeev ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr); 295a283635eSDmitry Karpeev if (viewfunction) { 296a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr); 297a283635eSDmitry Karpeev ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr); 298a283635eSDmitry Karpeev } 299a283635eSDmitry Karpeev ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 30011aeaf0aSBarry Smith } 301e884886fSBarry Smith PetscFunctionReturn(0); 302e884886fSBarry Smith } 303e884886fSBarry Smith 304e884886fSBarry Smith #undef __FUNCT__ 305e884886fSBarry Smith #define __FUNCT__ "MatAssemblyEnd_MFFD" 306e884886fSBarry Smith /* 307e884886fSBarry Smith MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 308e884886fSBarry Smith allows the user to indicate the beginning of a new linear solve by calling 309e884886fSBarry Smith MatAssemblyXXX() on the matrix free matrix. This then allows the 3101d0fab5eSBarry Smith MatCreateMFFD_WP() to properly compute ||U|| only the first time 311e884886fSBarry Smith in the linear solver rather than every time. 3125a576424SJed Brown 3135a576424SJed Brown This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library. 314e884886fSBarry Smith */ 3155a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt) 316e884886fSBarry Smith { 317e884886fSBarry Smith PetscErrorCode ierr; 318e884886fSBarry Smith MatMFFD j = (MatMFFD)J->data; 319e884886fSBarry Smith 320e884886fSBarry Smith PetscFunctionBegin; 321e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 322e884886fSBarry Smith j->vshift = 0.0; 323e884886fSBarry Smith j->vscale = 1.0; 324e884886fSBarry Smith PetscFunctionReturn(0); 325e884886fSBarry Smith } 326e884886fSBarry Smith 327e884886fSBarry Smith #undef __FUNCT__ 328e884886fSBarry Smith #define __FUNCT__ "MatMult_MFFD" 329e884886fSBarry Smith /* 330e884886fSBarry Smith MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 331e884886fSBarry Smith 332e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 333e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 334e884886fSBarry Smith u = current iterate 335e884886fSBarry Smith h = difference interval 336e884886fSBarry Smith */ 337e884886fSBarry Smith PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y) 338e884886fSBarry Smith { 339e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 340e884886fSBarry Smith PetscScalar h; 341e884886fSBarry Smith Vec w,U,F; 342e884886fSBarry Smith PetscErrorCode ierr; 343ace3abfcSBarry Smith PetscBool zeroa; 344e884886fSBarry Smith 345e884886fSBarry Smith PetscFunctionBegin; 346ce94432eSBarry 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"); 347e884886fSBarry Smith /* We log matrix-free matrix-vector products separately, so that we can 348e884886fSBarry Smith separate the performance monitoring from the cases that use conventional 349e884886fSBarry Smith storage. We may eventually modify event logging to associate events 350e884886fSBarry Smith with particular objects, hence alleviating the more general problem. */ 351e884886fSBarry Smith ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 352e884886fSBarry Smith 353e884886fSBarry Smith w = ctx->w; 354e884886fSBarry Smith U = ctx->current_u; 3553ec795f1SBarry Smith F = ctx->current_f; 356e884886fSBarry Smith /* 357e884886fSBarry Smith Compute differencing parameter 358e884886fSBarry Smith */ 359e884886fSBarry Smith if (!ctx->ops->compute) { 360e884886fSBarry Smith ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr); 3619c6ac3b3SBarry Smith ierr = MatSetFromOptions(mat);CHKERRQ(ierr); 362e884886fSBarry Smith } 363e884886fSBarry Smith ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr); 364e884886fSBarry Smith if (zeroa) { 365e884886fSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 366e884886fSBarry Smith PetscFunctionReturn(0); 367e884886fSBarry Smith } 368e884886fSBarry Smith 369e32f2f54SBarry Smith if (PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h"); 370e884886fSBarry Smith if (ctx->checkh) { 371e884886fSBarry Smith ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr); 372e884886fSBarry Smith } 373e884886fSBarry Smith 374e884886fSBarry Smith /* keep a record of the current differencing parameter h */ 375e884886fSBarry Smith ctx->currenth = h; 376e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX) 377e884886fSBarry Smith ierr = PetscInfo2(mat,"Current differencing parameter: %G + %G i\n",PetscRealPart(h),PetscImaginaryPart(h));CHKERRQ(ierr); 378e884886fSBarry Smith #else 379e884886fSBarry Smith ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr); 380e884886fSBarry Smith #endif 381e884886fSBarry Smith if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 382e884886fSBarry Smith ctx->historyh[ctx->ncurrenth] = h; 383e884886fSBarry Smith } 384e884886fSBarry Smith ctx->ncurrenth++; 385e884886fSBarry Smith 386e884886fSBarry Smith /* w = u + ha */ 387c3bb7e23SBarry Smith if (ctx->drscale) { 388c3bb7e23SBarry Smith ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr); 389c3bb7e23SBarry Smith ierr = VecAYPX(U,h,w);CHKERRQ(ierr); 390c3bb7e23SBarry Smith } else { 391e884886fSBarry Smith ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr); 392c3bb7e23SBarry Smith } 393e884886fSBarry Smith 3941b797266SDmitry Karpeev /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 3951b797266SDmitry Karpeev if (ctx->ncurrenth == 1 && ctx->current_f_allocated) { 396e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr); 39752121784SBarry Smith } 398e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr); 399e884886fSBarry Smith 400e884886fSBarry Smith ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 401e884886fSBarry Smith ierr = VecScale(y,1.0/h);CHKERRQ(ierr); 402e884886fSBarry Smith 403c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 404e884886fSBarry Smith ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr); 405c3bb7e23SBarry Smith } 406c3bb7e23SBarry Smith if (ctx->dlscale) { 407c3bb7e23SBarry Smith ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr); 408c3bb7e23SBarry Smith } 409c3bb7e23SBarry Smith if (ctx->dshift) { 410c3bb7e23SBarry Smith ierr = VecPointwiseMult(ctx->dshift,a,U);CHKERRQ(ierr); 411c3bb7e23SBarry Smith ierr = VecAXPY(y,1.0,U);CHKERRQ(ierr); 412c3bb7e23SBarry Smith } 413e884886fSBarry Smith 4140298fd71SBarry Smith if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,NULL);CHKERRQ(ierr);} 415e884886fSBarry Smith 416e884886fSBarry Smith ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 417e884886fSBarry Smith PetscFunctionReturn(0); 418e884886fSBarry Smith } 419e884886fSBarry Smith 420e884886fSBarry Smith #undef __FUNCT__ 421e884886fSBarry Smith #define __FUNCT__ "MatGetDiagonal_MFFD" 422e884886fSBarry Smith /* 423e884886fSBarry Smith MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 424e884886fSBarry Smith 425e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 426e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 427e884886fSBarry Smith u = current iterate 428e884886fSBarry Smith h = difference interval 429e884886fSBarry Smith */ 430e884886fSBarry Smith PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 431e884886fSBarry Smith { 432e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 433e884886fSBarry Smith PetscScalar h,*aa,*ww,v; 434e884886fSBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 435e884886fSBarry Smith Vec w,U; 436e884886fSBarry Smith PetscErrorCode ierr; 437e884886fSBarry Smith PetscInt i,rstart,rend; 438e884886fSBarry Smith 439e884886fSBarry Smith PetscFunctionBegin; 440e7e72b3dSBarry Smith if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first"); 441e884886fSBarry Smith 442e884886fSBarry Smith w = ctx->w; 443e884886fSBarry Smith U = ctx->current_u; 444e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr); 445e884886fSBarry Smith ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr); 446e884886fSBarry Smith ierr = VecCopy(U,w);CHKERRQ(ierr); 447e884886fSBarry Smith 448e884886fSBarry Smith ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr); 449e884886fSBarry Smith ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 450e884886fSBarry Smith for (i=rstart; i<rend; i++) { 451e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 452e884886fSBarry Smith h = ww[i-rstart]; 453e884886fSBarry Smith if (h == 0.0) h = 1.0; 454e884886fSBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 455e884886fSBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 456e884886fSBarry Smith h *= epsilon; 457e884886fSBarry Smith 458e884886fSBarry Smith ww[i-rstart] += h; 459e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 460e884886fSBarry Smith ierr = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr); 461e884886fSBarry Smith aa[i-rstart] = (v - aa[i-rstart])/h; 462e884886fSBarry Smith 463e884886fSBarry Smith /* possibly shift and scale result */ 464c35ec66cSMatthew G Knepley if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) { 465e884886fSBarry Smith aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart]; 466c3bb7e23SBarry Smith } 467e884886fSBarry Smith 468e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 469e884886fSBarry Smith ww[i-rstart] -= h; 470e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 471e884886fSBarry Smith } 472e884886fSBarry Smith ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 473e884886fSBarry Smith PetscFunctionReturn(0); 474e884886fSBarry Smith } 475e884886fSBarry Smith 476e884886fSBarry Smith #undef __FUNCT__ 477c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalScale_MFFD" 478c3bb7e23SBarry Smith PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr) 479c3bb7e23SBarry Smith { 480c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 481c3bb7e23SBarry Smith PetscErrorCode ierr; 482c3bb7e23SBarry Smith 483c3bb7e23SBarry Smith PetscFunctionBegin; 484c3bb7e23SBarry Smith if (ll && !aij->dlscale) { 485c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr); 486c3bb7e23SBarry Smith } 487c3bb7e23SBarry Smith if (rr && !aij->drscale) { 488c3bb7e23SBarry Smith ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr); 489c3bb7e23SBarry Smith } 490c3bb7e23SBarry Smith if (ll) { 491c3bb7e23SBarry Smith ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr); 492c3bb7e23SBarry Smith } 493c3bb7e23SBarry Smith if (rr) { 494c3bb7e23SBarry Smith ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr); 495c3bb7e23SBarry Smith } 496c3bb7e23SBarry Smith PetscFunctionReturn(0); 497c3bb7e23SBarry Smith } 498c3bb7e23SBarry Smith 499c3bb7e23SBarry Smith #undef __FUNCT__ 500c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalSet_MFFD" 501c3bb7e23SBarry Smith PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode) 502c3bb7e23SBarry Smith { 503c3bb7e23SBarry Smith MatMFFD aij = (MatMFFD)mat->data; 504c3bb7e23SBarry Smith PetscErrorCode ierr; 505c3bb7e23SBarry Smith 506c3bb7e23SBarry Smith PetscFunctionBegin; 507ce94432eSBarry Smith if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES"); 508c3bb7e23SBarry Smith if (!aij->dshift) { 509c3bb7e23SBarry Smith ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr); 510c3bb7e23SBarry Smith } 511c3bb7e23SBarry Smith ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr); 512c3bb7e23SBarry Smith PetscFunctionReturn(0); 513c3bb7e23SBarry Smith } 514c3bb7e23SBarry Smith 515c3bb7e23SBarry Smith #undef __FUNCT__ 516e884886fSBarry Smith #define __FUNCT__ "MatShift_MFFD" 517e884886fSBarry Smith PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a) 518e884886fSBarry Smith { 519e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 5205fd66863SKarl Rupp 521e884886fSBarry Smith PetscFunctionBegin; 522e884886fSBarry Smith shell->vshift += a; 523e884886fSBarry Smith PetscFunctionReturn(0); 524e884886fSBarry Smith } 525e884886fSBarry Smith 526e884886fSBarry Smith #undef __FUNCT__ 527e884886fSBarry Smith #define __FUNCT__ "MatScale_MFFD" 528e884886fSBarry Smith PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a) 529e884886fSBarry Smith { 530e884886fSBarry Smith MatMFFD shell = (MatMFFD)Y->data; 5315fd66863SKarl Rupp 532e884886fSBarry Smith PetscFunctionBegin; 533e884886fSBarry Smith shell->vscale *= a; 534e884886fSBarry Smith PetscFunctionReturn(0); 535e884886fSBarry Smith } 536e884886fSBarry Smith 537e884886fSBarry Smith #undef __FUNCT__ 538c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetBase_MFFD" 539d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F) 540e884886fSBarry Smith { 541e884886fSBarry Smith PetscErrorCode ierr; 542e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 543e884886fSBarry Smith 544e884886fSBarry Smith PetscFunctionBegin; 545e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 5462205254eSKarl Rupp 547e884886fSBarry Smith ctx->current_u = U; 54852121784SBarry Smith if (F) { 5496bf464f9SBarry Smith if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);} 5503ec795f1SBarry Smith ctx->current_f = F; 551cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_FALSE; 552cfe22f5eSBarry Smith } else if (!ctx->current_f_allocated) { 55352121784SBarry Smith ierr = VecDuplicate(ctx->current_u, &ctx->current_f);CHKERRQ(ierr); 5542205254eSKarl Rupp 555cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_TRUE; 55652121784SBarry Smith } 557e884886fSBarry Smith if (!ctx->w) { 558e884886fSBarry Smith ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr); 559e884886fSBarry Smith } 560e884886fSBarry Smith J->assembled = PETSC_TRUE; 561e884886fSBarry Smith PetscFunctionReturn(0); 562e884886fSBarry Smith } 5635a576424SJed Brown 564e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 565b2573a8aSBarry Smith 566e884886fSBarry Smith #undef __FUNCT__ 567c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh_MFFD" 5687087cfbeSBarry Smith PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx) 569e884886fSBarry Smith { 570e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 571e884886fSBarry Smith 572e884886fSBarry Smith PetscFunctionBegin; 573e884886fSBarry Smith ctx->checkh = fun; 574e884886fSBarry Smith ctx->checkhctx = ectx; 575e884886fSBarry Smith PetscFunctionReturn(0); 576e884886fSBarry Smith } 577e884886fSBarry Smith 578e884886fSBarry Smith #undef __FUNCT__ 5796aa9148fSLisandro Dalcin #define __FUNCT__ "MatMFFDSetOptionsPrefix" 5806aa9148fSLisandro Dalcin /*@C 5816aa9148fSLisandro Dalcin MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 5826aa9148fSLisandro Dalcin MatMFFD options in the database. 5836aa9148fSLisandro Dalcin 5846aa9148fSLisandro Dalcin Collective on Mat 5856aa9148fSLisandro Dalcin 5866aa9148fSLisandro Dalcin Input Parameter: 5876aa9148fSLisandro Dalcin + A - the Mat context 5886aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names 5896aa9148fSLisandro Dalcin 5906aa9148fSLisandro Dalcin Notes: 5916aa9148fSLisandro Dalcin A hyphen (-) must NOT be given at the beginning of the prefix name. 5926aa9148fSLisandro Dalcin The first character of all runtime options is AUTOMATICALLY the hyphen. 5936aa9148fSLisandro Dalcin 5946aa9148fSLisandro Dalcin Level: advanced 5956aa9148fSLisandro Dalcin 5966aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters 5976aa9148fSLisandro Dalcin 5989c6ac3b3SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF() 5996aa9148fSLisandro Dalcin @*/ 6007087cfbeSBarry Smith PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[]) 6016aa9148fSLisandro Dalcin 6026aa9148fSLisandro Dalcin { 6030298fd71SBarry Smith MatMFFD mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL; 6046aa9148fSLisandro Dalcin PetscErrorCode ierr; 6055fd66863SKarl Rupp 6066aa9148fSLisandro Dalcin PetscFunctionBegin; 6070700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 6080700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 6096aa9148fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr); 6106aa9148fSLisandro Dalcin PetscFunctionReturn(0); 6116aa9148fSLisandro Dalcin } 6126aa9148fSLisandro Dalcin 6136aa9148fSLisandro Dalcin #undef __FUNCT__ 6149c6ac3b3SBarry Smith #define __FUNCT__ "MatSetFromOptions_MFFD" 6159c6ac3b3SBarry Smith PetscErrorCode MatSetFromOptions_MFFD(Mat mat) 616e884886fSBarry Smith { 61771f08665SBarry Smith MatMFFD mfctx = (MatMFFD)mat->data; 618e884886fSBarry Smith PetscErrorCode ierr; 619ace3abfcSBarry Smith PetscBool flg; 620e884886fSBarry Smith char ftype[256]; 621e884886fSBarry Smith 622e884886fSBarry Smith PetscFunctionBegin; 6230700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 6240700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 6253194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr); 626b022a5c1SBarry Smith ierr = PetscOptionsList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr); 627e884886fSBarry Smith if (flg) { 628e884886fSBarry Smith ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr); 629e884886fSBarry Smith } 630e884886fSBarry Smith 631e884886fSBarry Smith ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 632e884886fSBarry Smith ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 633e884886fSBarry Smith 63490d69ab7SBarry Smith flg = PETSC_FALSE; 6350298fd71SBarry Smith ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr); 636e884886fSBarry Smith if (flg) { 637e884886fSBarry Smith ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr); 638e884886fSBarry Smith } 639e884886fSBarry Smith if (mfctx->ops->setfromoptions) { 640e884886fSBarry Smith ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr); 641e884886fSBarry Smith } 642e884886fSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 643e884886fSBarry Smith PetscFunctionReturn(0); 644e884886fSBarry Smith } 645e884886fSBarry Smith 646bc0f08ceSBarry Smith #undef __FUNCT__ 647bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod_MFFD" 648bc0f08ceSBarry Smith PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period) 649bc0f08ceSBarry Smith { 650bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 651bc0f08ceSBarry Smith 652bc0f08ceSBarry Smith PetscFunctionBegin; 653bc0f08ceSBarry Smith PetscValidLogicalCollectiveInt(mat,period,2); 654bc0f08ceSBarry Smith ctx->recomputeperiod = period; 655bc0f08ceSBarry Smith PetscFunctionReturn(0); 656bc0f08ceSBarry Smith } 657bc0f08ceSBarry Smith 658bc0f08ceSBarry Smith #undef __FUNCT__ 659bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunction_MFFD" 660bc0f08ceSBarry Smith PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 661bc0f08ceSBarry Smith { 662bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 663bc0f08ceSBarry Smith 664bc0f08ceSBarry Smith PetscFunctionBegin; 665bc0f08ceSBarry Smith ctx->func = func; 666bc0f08ceSBarry Smith ctx->funcctx = funcctx; 667bc0f08ceSBarry Smith PetscFunctionReturn(0); 668bc0f08ceSBarry Smith } 669bc0f08ceSBarry Smith 670bc0f08ceSBarry Smith #undef __FUNCT__ 671bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError_MFFD" 672bc0f08ceSBarry Smith PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error) 673bc0f08ceSBarry Smith { 674bc0f08ceSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 675bc0f08ceSBarry Smith 676bc0f08ceSBarry Smith PetscFunctionBegin; 677bc0f08ceSBarry Smith PetscValidLogicalCollectiveReal(mat,error,2); 678bc0f08ceSBarry Smith if (error != PETSC_DEFAULT) ctx->error_rel = error; 679bc0f08ceSBarry Smith PetscFunctionReturn(0); 680bc0f08ceSBarry Smith } 681bc0f08ceSBarry Smith 682e884886fSBarry Smith /*MC 683e884886fSBarry Smith MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 684e884886fSBarry Smith 685e884886fSBarry Smith Level: advanced 686e884886fSBarry Smith 687d2d6cebeSBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction() 688e884886fSBarry Smith M*/ 689e884886fSBarry Smith #undef __FUNCT__ 690e884886fSBarry Smith #define __FUNCT__ "MatCreate_MFFD" 691*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 692e884886fSBarry Smith { 693e884886fSBarry Smith MatMFFD mfctx; 694e884886fSBarry Smith PetscErrorCode ierr; 695e884886fSBarry Smith 696e884886fSBarry Smith PetscFunctionBegin; 697519f805aSKarl Rupp #if !defined(PETSC_USE_DYNAMIC_LIBRARIES) 6980298fd71SBarry Smith ierr = MatMFFDInitializePackage(NULL);CHKERRQ(ierr); 6993ec795f1SBarry Smith #endif 700e884886fSBarry Smith 701ce94432eSBarry Smith ierr = PetscHeaderCreate(mfctx,_p_MatMFFD,struct _MFOps,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr); 7022205254eSKarl Rupp 703e884886fSBarry Smith mfctx->sp = 0; 704e884886fSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 705e884886fSBarry Smith mfctx->recomputeperiod = 1; 706e884886fSBarry Smith mfctx->count = 0; 707e884886fSBarry Smith mfctx->currenth = 0.0; 7080298fd71SBarry Smith mfctx->historyh = NULL; 709e884886fSBarry Smith mfctx->ncurrenth = 0; 710e884886fSBarry Smith mfctx->maxcurrenth = 0; 7117adad957SLisandro Dalcin ((PetscObject)mfctx)->type_name = 0; 712e884886fSBarry Smith 713e884886fSBarry Smith mfctx->vshift = 0.0; 714e884886fSBarry Smith mfctx->vscale = 1.0; 715e884886fSBarry Smith 716e884886fSBarry Smith /* 717e884886fSBarry Smith Create the empty data structure to contain compute-h routines. 718e884886fSBarry Smith These will be filled in below from the command line options or 719e884886fSBarry Smith a later call with MatMFFDSetType() or if that is not called 720e884886fSBarry Smith then it will default in the first use of MatMult_MFFD() 721e884886fSBarry Smith */ 722e884886fSBarry Smith mfctx->ops->compute = 0; 723e884886fSBarry Smith mfctx->ops->destroy = 0; 724e884886fSBarry Smith mfctx->ops->view = 0; 725e884886fSBarry Smith mfctx->ops->setfromoptions = 0; 726e884886fSBarry Smith mfctx->hctx = 0; 727e884886fSBarry Smith 728e884886fSBarry Smith mfctx->func = 0; 729e884886fSBarry Smith mfctx->funcctx = 0; 7300298fd71SBarry Smith mfctx->w = NULL; 731e884886fSBarry Smith 732e884886fSBarry Smith A->data = mfctx; 733e884886fSBarry Smith 734e884886fSBarry Smith A->ops->mult = MatMult_MFFD; 735e884886fSBarry Smith A->ops->destroy = MatDestroy_MFFD; 736e884886fSBarry Smith A->ops->view = MatView_MFFD; 737e884886fSBarry Smith A->ops->assemblyend = MatAssemblyEnd_MFFD; 738e884886fSBarry Smith A->ops->getdiagonal = MatGetDiagonal_MFFD; 739e884886fSBarry Smith A->ops->scale = MatScale_MFFD; 740e884886fSBarry Smith A->ops->shift = MatShift_MFFD; 741c3bb7e23SBarry Smith A->ops->diagonalscale = MatDiagonalScale_MFFD; 742c3bb7e23SBarry Smith A->ops->diagonalset = MatDiagonalSet_MFFD; 7439c6ac3b3SBarry Smith A->ops->setfromoptions = MatSetFromOptions_MFFD; 744e884886fSBarry Smith A->assembled = PETSC_TRUE; 745e884886fSBarry Smith 74626283091SBarry Smith ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 74726283091SBarry Smith ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 748ee1cef2cSJed Brown 74900de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C","MatMFFDSetBase_MFFD",MatMFFDSetBase_MFFD);CHKERRQ(ierr); 75000de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C","MatMFFDSetFunctioniBase_MFFD",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr); 75100de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C","MatMFFDSetFunctioni_MFFD",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr); 75200de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C","MatMFFDSetFunction_MFFD",MatMFFDSetFunction_MFFD);CHKERRQ(ierr); 75300de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C","MatMFFDSetCheckh_MFFD",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr); 75400de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C","MatMFFDSetPeriod_MFFD",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr); 75500de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C","MatMFFDSetFunctionError_MFFD",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr); 75600de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C","MatMFFDResetHHistory_MFFD",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr); 75700de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDAddNullSpace_C","MatMFFDAddNullSpace_MFFD",MatMFFDAddNullSpace_MFFD);CHKERRQ(ierr); 7582205254eSKarl Rupp 759e884886fSBarry Smith mfctx->mat = A; 7602205254eSKarl Rupp 761e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr); 762e884886fSBarry Smith PetscFunctionReturn(0); 763e884886fSBarry Smith } 764e884886fSBarry Smith 765e884886fSBarry Smith #undef __FUNCT__ 766e884886fSBarry Smith #define __FUNCT__ "MatCreateMFFD" 767e884886fSBarry Smith /*@ 768e884886fSBarry Smith MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 769e884886fSBarry Smith 770e884886fSBarry Smith Collective on Vec 771e884886fSBarry Smith 772e884886fSBarry Smith Input Parameters: 773fef1beadSBarry Smith + comm - MPI communicator 774fef1beadSBarry Smith . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 775fef1beadSBarry Smith This value should be the same as the local size used in creating the 776fef1beadSBarry Smith y vector for the matrix-vector product y = Ax. 777fef1beadSBarry Smith . n - This value should be the same as the local size used in creating the 778fef1beadSBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 779fef1beadSBarry Smith calculated if N is given) For square matrices n is almost always m. 780fef1beadSBarry Smith . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 781fef1beadSBarry Smith - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 782fef1beadSBarry Smith 783e884886fSBarry Smith 784e884886fSBarry Smith Output Parameter: 785e884886fSBarry Smith . J - the matrix-free matrix 786e884886fSBarry Smith 7879c6ac3b3SBarry Smith Options Database Keys: call MatSetFromOptions() to trigger these 7889c6ac3b3SBarry Smith + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 7899c6ac3b3SBarry Smith - -mat_mffd_err - square root of estimated relative error in function evaluation 7909c6ac3b3SBarry Smith - -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 7919c6ac3b3SBarry Smith 7929c6ac3b3SBarry Smith 793e884886fSBarry Smith Level: advanced 794e884886fSBarry Smith 795e884886fSBarry Smith Notes: 796e884886fSBarry Smith The matrix-free matrix context merely contains the function pointers 797e884886fSBarry Smith and work space for performing finite difference approximations of 798e884886fSBarry Smith Jacobian-vector products, F'(u)*a, 799e884886fSBarry Smith 800e884886fSBarry Smith The default code uses the following approach to compute h 801e884886fSBarry Smith 802e884886fSBarry Smith .vb 803e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 804e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 805e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 806e884886fSBarry Smith where 807e884886fSBarry Smith error_rel = square root of relative error in function evaluation 808e884886fSBarry Smith umin = minimum iterate parameter 809e884886fSBarry Smith .ve 810e884886fSBarry Smith 811e884886fSBarry Smith The user can set the error_rel via MatMFFDSetFunctionError() and 8120598bfebSBarry Smith umin via MatMFFDDSSetUmin(); see the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A> for details. 813e884886fSBarry Smith 814e884886fSBarry Smith The user should call MatDestroy() when finished with the matrix-free 815e884886fSBarry Smith matrix context. 816e884886fSBarry Smith 817e884886fSBarry Smith Options Database Keys: 818e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 819e884886fSBarry Smith . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only) 820e884886fSBarry Smith - -mat_mffd_check_positivity 821e884886fSBarry Smith 822e884886fSBarry Smith .keywords: default, matrix-free, create, matrix 823e884886fSBarry Smith 8241d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 825e884886fSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 8261d0fab5eSBarry Smith MatMFFDGetH(), MatMFFDRegisterDynamic), MatMFFDComputeJacobian() 827e884886fSBarry Smith 828e884886fSBarry Smith @*/ 8297087cfbeSBarry Smith PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 830e884886fSBarry Smith { 831e884886fSBarry Smith PetscErrorCode ierr; 832e884886fSBarry Smith 833e884886fSBarry Smith PetscFunctionBegin; 834e884886fSBarry Smith ierr = MatCreate(comm,J);CHKERRQ(ierr); 835fef1beadSBarry Smith ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr); 836e884886fSBarry Smith ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 837be7c243fSJed Brown ierr = MatSetUp(*J);CHKERRQ(ierr); 838e884886fSBarry Smith PetscFunctionReturn(0); 839e884886fSBarry Smith } 840e884886fSBarry Smith 841e884886fSBarry Smith 842e884886fSBarry Smith #undef __FUNCT__ 843e884886fSBarry Smith #define __FUNCT__ "MatMFFDGetH" 844e884886fSBarry Smith /*@ 845e884886fSBarry Smith MatMFFDGetH - Gets the last value that was used as the differencing 846e884886fSBarry Smith parameter. 847e884886fSBarry Smith 848e884886fSBarry Smith Not Collective 849e884886fSBarry Smith 850e884886fSBarry Smith Input Parameters: 851e884886fSBarry Smith . mat - the matrix obtained with MatCreateSNESMF() 852e884886fSBarry Smith 853e884886fSBarry Smith Output Paramter: 854e884886fSBarry Smith . h - the differencing step size 855e884886fSBarry Smith 856e884886fSBarry Smith Level: advanced 857e884886fSBarry Smith 858e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 859e884886fSBarry Smith 8601d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory() 861e884886fSBarry Smith @*/ 8627087cfbeSBarry Smith PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 863e884886fSBarry Smith { 864e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 865bc0f08ceSBarry Smith PetscErrorCode ierr; 866bc0f08ceSBarry Smith PetscBool match; 867e884886fSBarry Smith 868e884886fSBarry Smith PetscFunctionBegin; 869251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr); 870ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 871bc0f08ceSBarry Smith 872e884886fSBarry Smith *h = ctx->currenth; 873e884886fSBarry Smith PetscFunctionReturn(0); 874e884886fSBarry Smith } 875e884886fSBarry Smith 876e884886fSBarry Smith #undef __FUNCT__ 877e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunction" 878e884886fSBarry Smith /*@C 879e884886fSBarry Smith MatMFFDSetFunction - Sets the function used in applying the matrix free. 880e884886fSBarry Smith 8813f9fe445SBarry Smith Logically Collective on Mat 882e884886fSBarry Smith 883e884886fSBarry Smith Input Parameters: 884e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 885e884886fSBarry Smith . func - the function to use 886e884886fSBarry Smith - funcctx - optional function context passed to function 887e884886fSBarry Smith 888e884886fSBarry Smith Level: advanced 889e884886fSBarry Smith 890e884886fSBarry Smith Notes: 891e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 892e884886fSBarry Smith matrix inside your compute Jacobian routine 893e884886fSBarry Smith 8943ec795f1SBarry Smith If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 895e884886fSBarry Smith 896e884886fSBarry Smith .keywords: SNES, matrix-free, function 897e884886fSBarry Smith 8981d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD, 8991d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 900e884886fSBarry Smith @*/ 9017087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 902e884886fSBarry Smith { 903bc0f08ceSBarry Smith PetscErrorCode ierr; 904e884886fSBarry Smith 905e884886fSBarry Smith PetscFunctionBegin; 906bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr); 907e884886fSBarry Smith PetscFunctionReturn(0); 908e884886fSBarry Smith } 909e884886fSBarry Smith 910e884886fSBarry Smith #undef __FUNCT__ 911e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni" 912e884886fSBarry Smith /*@C 913e884886fSBarry Smith MatMFFDSetFunctioni - Sets the function for a single component 914e884886fSBarry Smith 9153f9fe445SBarry Smith Logically Collective on Mat 916e884886fSBarry Smith 917e884886fSBarry Smith Input Parameters: 918e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 919e884886fSBarry Smith - funci - the function to use 920e884886fSBarry Smith 921e884886fSBarry Smith Level: advanced 922e884886fSBarry Smith 923e884886fSBarry Smith Notes: 924e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 925e884886fSBarry Smith matrix inside your compute Jacobian routine 926e884886fSBarry Smith 927e884886fSBarry Smith 928e884886fSBarry Smith .keywords: SNES, matrix-free, function 929e884886fSBarry Smith 9301d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 9311d0fab5eSBarry Smith 932e884886fSBarry Smith @*/ 9337087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 934e884886fSBarry Smith { 9354ac538c5SBarry Smith PetscErrorCode ierr; 936e884886fSBarry Smith 937e884886fSBarry Smith PetscFunctionBegin; 9380700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 9394ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr); 940e884886fSBarry Smith PetscFunctionReturn(0); 941e884886fSBarry Smith } 942e884886fSBarry Smith 943e884886fSBarry Smith 944e884886fSBarry Smith #undef __FUNCT__ 945e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase" 946e884886fSBarry Smith /*@C 947e884886fSBarry Smith MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 948e884886fSBarry Smith 9493f9fe445SBarry Smith Logically Collective on Mat 950e884886fSBarry Smith 951e884886fSBarry Smith Input Parameters: 952e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 953e884886fSBarry Smith - func - the function to use 954e884886fSBarry Smith 955e884886fSBarry Smith Level: advanced 956e884886fSBarry Smith 957e884886fSBarry Smith Notes: 958e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 959e884886fSBarry Smith matrix inside your compute Jacobian routine 960e884886fSBarry Smith 961e884886fSBarry Smith 962e884886fSBarry Smith .keywords: SNES, matrix-free, function 963e884886fSBarry Smith 964e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 9651d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 966e884886fSBarry Smith @*/ 9677087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 968e884886fSBarry Smith { 9694ac538c5SBarry Smith PetscErrorCode ierr; 970e884886fSBarry Smith 971e884886fSBarry Smith PetscFunctionBegin; 9720700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 9734ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr); 974e884886fSBarry Smith PetscFunctionReturn(0); 975e884886fSBarry Smith } 976e884886fSBarry Smith 977e884886fSBarry Smith #undef __FUNCT__ 978e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod" 979e884886fSBarry Smith /*@ 980e884886fSBarry Smith MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 981e884886fSBarry Smith 9823f9fe445SBarry Smith Logically Collective on Mat 983e884886fSBarry Smith 984e884886fSBarry Smith Input Parameters: 985e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 986e884886fSBarry Smith - period - 1 for everytime, 2 for every second etc 987e884886fSBarry Smith 988e884886fSBarry Smith Options Database Keys: 989e884886fSBarry Smith + -mat_mffd_period <period> 990e884886fSBarry Smith 991e884886fSBarry Smith Level: advanced 992e884886fSBarry Smith 993e884886fSBarry Smith 994e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 995e884886fSBarry Smith 996e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), 9971d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 998e884886fSBarry Smith @*/ 9997087cfbeSBarry Smith PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 1000e884886fSBarry Smith { 1001bc0f08ceSBarry Smith PetscErrorCode ierr; 1002e884886fSBarry Smith 1003e884886fSBarry Smith PetscFunctionBegin; 1004bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr); 1005e884886fSBarry Smith PetscFunctionReturn(0); 1006e884886fSBarry Smith } 1007e884886fSBarry Smith 1008e884886fSBarry Smith #undef __FUNCT__ 1009e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError" 1010e884886fSBarry Smith /*@ 1011e884886fSBarry Smith MatMFFDSetFunctionError - Sets the error_rel for the approximation of 1012e884886fSBarry Smith matrix-vector products using finite differences. 1013e884886fSBarry Smith 10143f9fe445SBarry Smith Logically Collective on Mat 1015e884886fSBarry Smith 1016e884886fSBarry Smith Input Parameters: 1017e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 1018e884886fSBarry Smith - error_rel - relative error (should be set to the square root of 1019e884886fSBarry Smith the relative error in the function evaluations) 1020e884886fSBarry Smith 1021e884886fSBarry Smith Options Database Keys: 1022e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 1023e884886fSBarry Smith 1024e884886fSBarry Smith Level: advanced 1025e884886fSBarry Smith 1026e884886fSBarry Smith Notes: 1027e884886fSBarry Smith The default matrix-free matrix-vector product routine computes 1028e884886fSBarry Smith .vb 1029e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 1030e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 1031e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 1032e884886fSBarry Smith .ve 1033e884886fSBarry Smith 1034e884886fSBarry Smith .keywords: SNES, matrix-free, parameters 1035e884886fSBarry Smith 1036e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 10371d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 1038e884886fSBarry Smith @*/ 10397087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 1040e884886fSBarry Smith { 1041bc0f08ceSBarry Smith PetscErrorCode ierr; 1042e884886fSBarry Smith 1043e884886fSBarry Smith PetscFunctionBegin; 1044bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr); 1045e884886fSBarry Smith PetscFunctionReturn(0); 1046e884886fSBarry Smith } 1047e884886fSBarry Smith 1048e884886fSBarry Smith #undef __FUNCT__ 1049e884886fSBarry Smith #define __FUNCT__ "MatMFFDAddNullSpace" 1050e884886fSBarry Smith /*@ 1051e884886fSBarry Smith MatMFFDAddNullSpace - Provides a null space that an operator is 1052e884886fSBarry Smith supposed to have. Since roundoff will create a small component in 1053e884886fSBarry Smith the null space, if you know the null space you may have it 1054e884886fSBarry Smith automatically removed. 1055e884886fSBarry Smith 10563f9fe445SBarry Smith Logically Collective on Mat 1057e884886fSBarry Smith 1058e884886fSBarry Smith Input Parameters: 1059e884886fSBarry Smith + J - the matrix-free matrix context 1060e884886fSBarry Smith - nullsp - object created with MatNullSpaceCreate() 1061e884886fSBarry Smith 1062e884886fSBarry Smith Level: advanced 1063e884886fSBarry Smith 1064e884886fSBarry Smith .keywords: SNES, matrix-free, null space 1065e884886fSBarry Smith 1066e884886fSBarry Smith .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD 10671d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 1068e884886fSBarry Smith @*/ 10697087cfbeSBarry Smith PetscErrorCode MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp) 1070e884886fSBarry Smith { 1071e884886fSBarry Smith PetscErrorCode ierr; 1072e884886fSBarry Smith 1073e884886fSBarry Smith PetscFunctionBegin; 1074bc0f08ceSBarry Smith ierr = PetscTryMethod(J,"MatMFFDAddNullSpace_C",(Mat,MatNullSpace),(J,nullsp));CHKERRQ(ierr); 1075e884886fSBarry Smith PetscFunctionReturn(0); 1076e884886fSBarry Smith } 1077e884886fSBarry Smith 1078e884886fSBarry Smith #undef __FUNCT__ 1079e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetHHistory" 1080e884886fSBarry Smith /*@ 1081e884886fSBarry Smith MatMFFDSetHHistory - Sets an array to collect a history of the 1082e884886fSBarry Smith differencing values (h) computed for the matrix-free product. 1083e884886fSBarry Smith 10843f9fe445SBarry Smith Logically Collective on Mat 1085e884886fSBarry Smith 1086e884886fSBarry Smith Input Parameters: 1087e884886fSBarry Smith + J - the matrix-free matrix context 1088e884886fSBarry Smith . histroy - space to hold the history 1089e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than 1090e884886fSBarry Smith nhistory, then the later ones are discarded 1091e884886fSBarry Smith 1092e884886fSBarry Smith Level: advanced 1093e884886fSBarry Smith 1094e884886fSBarry Smith Notes: 1095e884886fSBarry Smith Use MatMFFDResetHHistory() to reset the history counter and collect 1096e884886fSBarry Smith a new batch of differencing parameters, h. 1097e884886fSBarry Smith 1098e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1099e884886fSBarry Smith 1100e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 11011d0fab5eSBarry Smith MatMFFDResetHHistory(), MatMFFDSetFunctionError() 1102e884886fSBarry Smith 1103e884886fSBarry Smith @*/ 11047087cfbeSBarry Smith PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 1105e884886fSBarry Smith { 1106e884886fSBarry Smith MatMFFD ctx = (MatMFFD)J->data; 1107bc0f08ceSBarry Smith PetscErrorCode ierr; 1108bc0f08ceSBarry Smith PetscBool match; 1109e884886fSBarry Smith 1110e884886fSBarry Smith PetscFunctionBegin; 1111251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr); 1112ce94432eSBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix"); 1113e884886fSBarry Smith ctx->historyh = history; 1114e884886fSBarry Smith ctx->maxcurrenth = nhistory; 111575567043SBarry Smith ctx->currenth = 0.; 1116e884886fSBarry Smith PetscFunctionReturn(0); 1117e884886fSBarry Smith } 1118e884886fSBarry Smith 1119bc0f08ceSBarry Smith 1120e884886fSBarry Smith #undef __FUNCT__ 1121e884886fSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory" 1122e884886fSBarry Smith /*@ 1123e884886fSBarry Smith MatMFFDResetHHistory - Resets the counter to zero to begin 1124e884886fSBarry Smith collecting a new set of differencing histories. 1125e884886fSBarry Smith 11263f9fe445SBarry Smith Logically Collective on Mat 1127e884886fSBarry Smith 1128e884886fSBarry Smith Input Parameters: 1129e884886fSBarry Smith . J - the matrix-free matrix context 1130e884886fSBarry Smith 1131e884886fSBarry Smith Level: advanced 1132e884886fSBarry Smith 1133e884886fSBarry Smith Notes: 1134e884886fSBarry Smith Use MatMFFDSetHHistory() to create the original history counter. 1135e884886fSBarry Smith 1136e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history 1137e884886fSBarry Smith 1138e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 11391d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDSetFunctionError() 1140e884886fSBarry Smith 1141e884886fSBarry Smith @*/ 11427087cfbeSBarry Smith PetscErrorCode MatMFFDResetHHistory(Mat J) 1143e884886fSBarry Smith { 1144bc0f08ceSBarry Smith PetscErrorCode ierr; 1145e884886fSBarry Smith 1146e884886fSBarry Smith PetscFunctionBegin; 1147bc0f08ceSBarry Smith ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr); 1148e884886fSBarry Smith PetscFunctionReturn(0); 1149e884886fSBarry Smith } 1150e884886fSBarry Smith 1151e884886fSBarry Smith 1152e884886fSBarry Smith #undef __FUNCT__ 1153e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetBase" 1154e884886fSBarry Smith /*@ 1155e884886fSBarry Smith MatMFFDSetBase - Sets the vector U at which matrix vector products of the 1156e884886fSBarry Smith Jacobian are computed 1157e884886fSBarry Smith 11583f9fe445SBarry Smith Logically Collective on Mat 1159e884886fSBarry Smith 1160e884886fSBarry Smith Input Parameters: 1161e884886fSBarry Smith + J - the MatMFFD matrix 11623ec795f1SBarry Smith . U - the vector 1163bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed 1164e884886fSBarry Smith 1165e884886fSBarry Smith Notes: This is rarely used directly 1166e884886fSBarry Smith 11678af5ae88SBarry Smith If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 11688af5ae88SBarry Smith point during the first MatMult() after each call to MatMFFDSetBase(). 1169dff2f722SBarry Smith 1170e884886fSBarry Smith Level: advanced 1171e884886fSBarry Smith 1172e884886fSBarry Smith @*/ 11737087cfbeSBarry Smith PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 1174e884886fSBarry Smith { 11754ac538c5SBarry Smith PetscErrorCode ierr; 1176e884886fSBarry Smith 1177e884886fSBarry Smith PetscFunctionBegin; 11780700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 11790700a824SBarry Smith PetscValidHeaderSpecific(U,VEC_CLASSID,2); 11800700a824SBarry Smith if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 11814ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr); 1182e884886fSBarry Smith PetscFunctionReturn(0); 1183e884886fSBarry Smith } 1184e884886fSBarry Smith 1185e884886fSBarry Smith #undef __FUNCT__ 1186e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh" 1187e884886fSBarry Smith /*@C 1188e884886fSBarry Smith MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 1189e884886fSBarry Smith it to satisfy some criteria 1190e884886fSBarry Smith 11913f9fe445SBarry Smith Logically Collective on Mat 1192e884886fSBarry Smith 1193e884886fSBarry Smith Input Parameters: 1194e884886fSBarry Smith + J - the MatMFFD matrix 1195e884886fSBarry Smith . fun - the function that checks h 1196e884886fSBarry Smith - ctx - any context needed by the function 1197e884886fSBarry Smith 1198e884886fSBarry Smith Options Database Keys: 1199e884886fSBarry Smith . -mat_mffd_check_positivity 1200e884886fSBarry Smith 1201e884886fSBarry Smith Level: advanced 1202e884886fSBarry Smith 1203e884886fSBarry Smith Notes: For example, MatMFFDSetCheckPositivity() insures that all entries 1204e884886fSBarry Smith of U + h*a are non-negative 1205e884886fSBarry Smith 1206e884886fSBarry Smith .seealso: MatMFFDSetCheckPositivity() 1207e884886fSBarry Smith @*/ 12087087cfbeSBarry Smith PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx) 1209e884886fSBarry Smith { 12104ac538c5SBarry Smith PetscErrorCode ierr; 1211e884886fSBarry Smith 1212e884886fSBarry Smith PetscFunctionBegin; 12130700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 12144ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr); 1215e884886fSBarry Smith PetscFunctionReturn(0); 1216e884886fSBarry Smith } 1217e884886fSBarry Smith 1218e884886fSBarry Smith #undef __FUNCT__ 1219e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckPositivity" 1220e884886fSBarry Smith /*@ 1221e884886fSBarry Smith MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1222e884886fSBarry Smith zero, decreases h until this is satisfied. 1223e884886fSBarry Smith 12243f9fe445SBarry Smith Logically Collective on Vec 1225e884886fSBarry Smith 1226e884886fSBarry Smith Input Parameters: 1227e884886fSBarry Smith + U - base vector that is added to 1228e884886fSBarry Smith . a - vector that is added 1229e884886fSBarry Smith . h - scaling factor on a 1230e884886fSBarry Smith - dummy - context variable (unused) 1231e884886fSBarry Smith 1232e884886fSBarry Smith Options Database Keys: 1233e884886fSBarry Smith . -mat_mffd_check_positivity 1234e884886fSBarry Smith 1235e884886fSBarry Smith Level: advanced 1236e884886fSBarry Smith 1237e884886fSBarry Smith Notes: This is rarely used directly, rather it is passed as an argument to 1238e884886fSBarry Smith MatMFFDSetCheckh() 1239e884886fSBarry Smith 1240e884886fSBarry Smith .seealso: MatMFFDSetCheckh() 1241e884886fSBarry Smith @*/ 12427087cfbeSBarry Smith PetscErrorCode MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h) 1243e884886fSBarry Smith { 1244e884886fSBarry Smith PetscReal val, minval; 1245e884886fSBarry Smith PetscScalar *u_vec, *a_vec; 1246e884886fSBarry Smith PetscErrorCode ierr; 1247e884886fSBarry Smith PetscInt i,n; 1248e884886fSBarry Smith MPI_Comm comm; 1249e884886fSBarry Smith 1250e884886fSBarry Smith PetscFunctionBegin; 1251e884886fSBarry Smith ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); 1252e884886fSBarry Smith ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); 1253e884886fSBarry Smith ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); 1254e884886fSBarry Smith ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); 1255e884886fSBarry Smith minval = PetscAbsScalar(*h*1.01); 1256e884886fSBarry Smith for (i=0; i<n; i++) { 1257e884886fSBarry Smith if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1258e884886fSBarry Smith val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1259e884886fSBarry Smith if (val < minval) minval = val; 1260e884886fSBarry Smith } 1261e884886fSBarry Smith } 1262e884886fSBarry Smith ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); 1263e884886fSBarry Smith ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); 1264d9822059SBarry Smith ierr = MPI_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr); 1265e884886fSBarry Smith if (val <= PetscAbsScalar(*h)) { 1266e884886fSBarry Smith ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr); 1267e884886fSBarry Smith if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1268e884886fSBarry Smith else *h = -0.99*val; 1269e884886fSBarry Smith } 1270e884886fSBarry Smith PetscFunctionReturn(0); 1271e884886fSBarry Smith } 1272e884886fSBarry Smith 1273e884886fSBarry Smith 1274e884886fSBarry Smith 1275e884886fSBarry Smith 1276e884886fSBarry Smith 1277e884886fSBarry Smith 1278e884886fSBarry Smith 1279