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 18755b7f64SBarry Smith .seealso: PetscFinalize(), MatCreateMFFD(), MatCreateSNESMF() 19b022a5c1SBarry Smith @*/ 207087cfbeSBarry Smith PetscErrorCode MatMFFDFinalizePackage(void) 21b022a5c1SBarry Smith { 22a703d84cSPatrick Lacasse PetscErrorCode ierr; 23a703d84cSPatrick Lacasse 24b022a5c1SBarry Smith PetscFunctionBegin; 2537e93019SBarry Smith ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr); 26b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_FALSE; 27b022a5c1SBarry Smith MatMFFDRegisterAllCalled = PETSC_FALSE; 28b022a5c1SBarry Smith PetscFunctionReturn(0); 29b022a5c1SBarry Smith } 30b022a5c1SBarry Smith 313ec795f1SBarry Smith /*@C 323ec795f1SBarry Smith MatMFFDInitializePackage - This function initializes everything in the MatMFFD package. It is called 338a690491SBarry Smith from MatInitializePackage(). 343ec795f1SBarry Smith 353ec795f1SBarry Smith Level: developer 363ec795f1SBarry Smith 373ec795f1SBarry Smith .seealso: PetscInitialize() 383ec795f1SBarry Smith @*/ 39607a6623SBarry Smith PetscErrorCode MatMFFDInitializePackage(void) 403ec795f1SBarry Smith { 413ec795f1SBarry Smith char logList[256]; 428e81d068SLisandro Dalcin PetscBool opt,pkg; 433ec795f1SBarry Smith PetscErrorCode ierr; 443ec795f1SBarry Smith 453ec795f1SBarry Smith PetscFunctionBegin; 46b022a5c1SBarry Smith if (MatMFFDPackageInitialized) PetscFunctionReturn(0); 47b022a5c1SBarry Smith MatMFFDPackageInitialized = PETSC_TRUE; 483ec795f1SBarry Smith /* Register Classes */ 490700a824SBarry Smith ierr = PetscClassIdRegister("MatMFFD",&MATMFFD_CLASSID);CHKERRQ(ierr); 503ec795f1SBarry Smith /* Register Constructors */ 51607a6623SBarry Smith ierr = MatMFFDRegisterAll();CHKERRQ(ierr); 523ec795f1SBarry Smith /* Register Events */ 530700a824SBarry Smith ierr = PetscLogEventRegister("MatMult MF",MATMFFD_CLASSID,&MATMFFD_Mult);CHKERRQ(ierr); 543ec795f1SBarry Smith /* Process info exclusions */ 558e81d068SLisandro Dalcin ierr = PetscOptionsGetString(NULL,NULL,"-info_exclude",logList,sizeof(logList),&opt);CHKERRQ(ierr); 563ec795f1SBarry Smith if (opt) { 578e81d068SLisandro Dalcin ierr = PetscStrInList("matmffd",logList,',',&pkg);CHKERRQ(ierr); 588e81d068SLisandro Dalcin if (pkg) {ierr = PetscInfoDeactivateClass(MATMFFD_CLASSID);CHKERRQ(ierr);} 593ec795f1SBarry Smith } 603ec795f1SBarry Smith /* Process summary exclusions */ 618e81d068SLisandro Dalcin ierr = PetscOptionsGetString(NULL,NULL,"-log_exclude",logList,sizeof(logList),&opt);CHKERRQ(ierr); 623ec795f1SBarry Smith if (opt) { 638e81d068SLisandro Dalcin ierr = PetscStrInList("matmffd",logList,',',&pkg);CHKERRQ(ierr); 64fa2bb9feSLisandro Dalcin if (pkg) {ierr = PetscLogEventExcludeClass(MATMFFD_CLASSID);CHKERRQ(ierr);} 653ec795f1SBarry Smith } 668e81d068SLisandro Dalcin /* Register package finalizer */ 67b022a5c1SBarry Smith ierr = PetscRegisterFinalize(MatMFFDFinalizePackage);CHKERRQ(ierr); 683ec795f1SBarry Smith PetscFunctionReturn(0); 693ec795f1SBarry Smith } 703ec795f1SBarry Smith 71789d8953SBarry Smith static PetscErrorCode MatMFFDSetType_MFFD(Mat mat,MatMFFDType ftype) 72789d8953SBarry Smith { 73789d8953SBarry Smith PetscErrorCode ierr,(*r)(MatMFFD); 74789d8953SBarry Smith MatMFFD ctx; 75789d8953SBarry Smith PetscBool match; 76789d8953SBarry Smith 77789d8953SBarry Smith PetscFunctionBegin; 78789d8953SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 79789d8953SBarry Smith PetscValidCharPointer(ftype,2); 80789d8953SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 81789d8953SBarry Smith 82789d8953SBarry Smith /* already set, so just return */ 83789d8953SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr); 84789d8953SBarry Smith if (match) PetscFunctionReturn(0); 85789d8953SBarry Smith 86789d8953SBarry Smith /* destroy the old one if it exists */ 87789d8953SBarry Smith if (ctx->ops->destroy) { 88789d8953SBarry Smith ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr); 89789d8953SBarry Smith } 90789d8953SBarry Smith 91789d8953SBarry Smith ierr = PetscFunctionListFind(MatMFFDList,ftype,&r);CHKERRQ(ierr); 92789d8953SBarry Smith if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype); 93789d8953SBarry Smith ierr = (*r)(ctx);CHKERRQ(ierr); 94789d8953SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr); 95789d8953SBarry Smith PetscFunctionReturn(0); 96789d8953SBarry Smith } 97789d8953SBarry Smith 98e884886fSBarry Smith /*@C 99e884886fSBarry Smith MatMFFDSetType - Sets the method that is used to compute the 100e884886fSBarry Smith differencing parameter for finite differene matrix-free formulations. 101e884886fSBarry Smith 102e884886fSBarry Smith Input Parameters: 103e884886fSBarry Smith + mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD() 104e884886fSBarry Smith or MatSetType(mat,MATMFFD); 105e884886fSBarry Smith - ftype - the type requested, either MATMFFD_WP or MATMFFD_DS 106e884886fSBarry Smith 107e884886fSBarry Smith Level: advanced 108e884886fSBarry Smith 109e884886fSBarry Smith Notes: 110e884886fSBarry Smith For example, such routines can compute h for use in 111e884886fSBarry Smith Jacobian-vector products of the form 112e884886fSBarry Smith 113e884886fSBarry Smith F(x+ha) - F(x) 114e884886fSBarry Smith F'(u)a ~= ---------------- 115e884886fSBarry Smith h 116e884886fSBarry Smith 117755b7f64SBarry Smith .seealso: MatCreateSNESMF(), MatMFFDRegister(), MatMFFDSetFunction(), MatCreateMFFD() 118e884886fSBarry Smith @*/ 11919fd82e9SBarry Smith PetscErrorCode MatMFFDSetType(Mat mat,MatMFFDType ftype) 120e884886fSBarry Smith { 121789d8953SBarry Smith PetscErrorCode ierr; 122e884886fSBarry Smith 123e884886fSBarry Smith PetscFunctionBegin; 1240700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 125e884886fSBarry Smith PetscValidCharPointer(ftype,2); 126789d8953SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetType_C",(Mat,MatMFFDType),(mat,ftype));CHKERRQ(ierr); 127e884886fSBarry Smith PetscFunctionReturn(0); 128e884886fSBarry Smith } 129e884886fSBarry Smith 1305d21e1f8SStefano Zampini static PetscErrorCode MatGetDiagonal_MFFD(Mat,Vec); 1315d21e1f8SStefano Zampini 132e884886fSBarry Smith typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/ 13340244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func) 134e884886fSBarry Smith { 135789d8953SBarry Smith MatMFFD ctx; 136789d8953SBarry Smith PetscErrorCode ierr; 137e884886fSBarry Smith 138e884886fSBarry Smith PetscFunctionBegin; 139789d8953SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 140e884886fSBarry Smith ctx->funcisetbase = func; 141e884886fSBarry Smith PetscFunctionReturn(0); 142e884886fSBarry Smith } 143e884886fSBarry Smith 144e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 14540244768SBarry Smith static PetscErrorCode MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci) 146e884886fSBarry Smith { 147789d8953SBarry Smith MatMFFD ctx; 148789d8953SBarry Smith PetscErrorCode ierr; 149e884886fSBarry Smith 150e884886fSBarry Smith PetscFunctionBegin; 151789d8953SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 152e884886fSBarry Smith ctx->funci = funci; 153789d8953SBarry Smith ierr = MatShellSetOperation(mat,MATOP_GET_DIAGONAL,(void (*)(void))MatGetDiagonal_MFFD);CHKERRQ(ierr); 154789d8953SBarry Smith PetscFunctionReturn(0); 1555d21e1f8SStefano Zampini } 156789d8953SBarry Smith 157789d8953SBarry Smith static PetscErrorCode MatMFFDGetH_MFFD(Mat mat,PetscScalar *h) 158789d8953SBarry Smith { 159789d8953SBarry Smith MatMFFD ctx; 160789d8953SBarry Smith PetscErrorCode ierr; 161789d8953SBarry Smith 162789d8953SBarry Smith PetscFunctionBegin; 163789d8953SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 164789d8953SBarry Smith *h = ctx->currenth; 165e884886fSBarry Smith PetscFunctionReturn(0); 166e884886fSBarry Smith } 167e884886fSBarry Smith 16840244768SBarry Smith static PetscErrorCode MatMFFDResetHHistory_MFFD(Mat J) 169bc0f08ceSBarry Smith { 170789d8953SBarry Smith MatMFFD ctx; 171789d8953SBarry Smith PetscErrorCode ierr; 172bc0f08ceSBarry Smith 173bc0f08ceSBarry Smith PetscFunctionBegin; 174789d8953SBarry Smith ierr = MatShellGetContext(J,&ctx);CHKERRQ(ierr); 175bc0f08ceSBarry Smith ctx->ncurrenth = 0; 176bc0f08ceSBarry Smith PetscFunctionReturn(0); 177bc0f08ceSBarry Smith } 178e884886fSBarry Smith 1791c84c290SBarry Smith /*@C 1801c84c290SBarry Smith MatMFFDRegister - Adds a method to the MatMFFD registry. 1811c84c290SBarry Smith 1821c84c290SBarry Smith Not Collective 1831c84c290SBarry Smith 1841c84c290SBarry Smith Input Parameters: 1851c84c290SBarry Smith + name_solver - name of a new user-defined compute-h module 1861c84c290SBarry Smith - routine_create - routine to create method context 1871c84c290SBarry Smith 1881c84c290SBarry Smith Level: developer 1891c84c290SBarry Smith 1901c84c290SBarry Smith Notes: 1911c84c290SBarry Smith MatMFFDRegister() may be called multiple times to add several user-defined solvers. 1921c84c290SBarry Smith 1931c84c290SBarry Smith Sample usage: 1941c84c290SBarry Smith .vb 195bdf89e91SBarry Smith MatMFFDRegister("my_h",MyHCreate); 1961c84c290SBarry Smith .ve 1971c84c290SBarry Smith 1981c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 1991c84c290SBarry Smith $ MatMFFDSetType(mfctx,"my_h") 2001c84c290SBarry Smith or at runtime via the option 201be7a6d03SBarry Smith $ -mat_mffd_type my_h 2021c84c290SBarry Smith 2031c84c290SBarry Smith .seealso: MatMFFDRegisterAll(), MatMFFDRegisterDestroy() 2041c84c290SBarry Smith @*/ 205bdf89e91SBarry Smith PetscErrorCode MatMFFDRegister(const char sname[],PetscErrorCode (*function)(MatMFFD)) 206e884886fSBarry Smith { 207e884886fSBarry Smith PetscErrorCode ierr; 208e884886fSBarry Smith 209e884886fSBarry Smith PetscFunctionBegin; 2109cc31a68SJed Brown ierr = MatInitializePackage();CHKERRQ(ierr); 211a240a19fSJed Brown ierr = PetscFunctionListAdd(&MatMFFDList,sname,function);CHKERRQ(ierr); 212e884886fSBarry Smith PetscFunctionReturn(0); 213e884886fSBarry Smith } 214e884886fSBarry Smith 215e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/ 21640244768SBarry Smith static PetscErrorCode MatDestroy_MFFD(Mat mat) 217e884886fSBarry Smith { 218e884886fSBarry Smith PetscErrorCode ierr; 219789d8953SBarry Smith MatMFFD ctx; 220e884886fSBarry Smith 221e884886fSBarry Smith PetscFunctionBegin; 222789d8953SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 2236bf464f9SBarry Smith ierr = VecDestroy(&ctx->w);CHKERRQ(ierr); 224c51ad4d4SStefano Zampini ierr = VecDestroy(&ctx->current_u);CHKERRQ(ierr); 225cfe22f5eSBarry Smith if (ctx->current_f_allocated) { 2266bf464f9SBarry Smith ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr); 227cfe22f5eSBarry Smith } 228e884886fSBarry Smith if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);} 2296bf464f9SBarry Smith ierr = PetscHeaderDestroy(&ctx);CHKERRQ(ierr); 230e884886fSBarry Smith 231bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C",NULL);CHKERRQ(ierr); 232bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",NULL);CHKERRQ(ierr); 233bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",NULL);CHKERRQ(ierr); 234bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C",NULL);CHKERRQ(ierr); 235bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C",NULL);CHKERRQ(ierr); 236bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C",NULL);CHKERRQ(ierr); 237bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C",NULL);CHKERRQ(ierr); 238bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C",NULL);CHKERRQ(ierr); 239789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetHHistory_C",NULL);CHKERRQ(ierr); 240789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetType_C",NULL);CHKERRQ(ierr); 241789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDGetH_C",NULL);CHKERRQ(ierr); 242e884886fSBarry Smith PetscFunctionReturn(0); 243e884886fSBarry Smith } 244e884886fSBarry Smith 245e884886fSBarry Smith /* 246e884886fSBarry Smith MatMFFDView_MFFD - Views matrix-free parameters. 247e884886fSBarry Smith 248e884886fSBarry Smith */ 24940244768SBarry Smith static PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer) 250e884886fSBarry Smith { 251e884886fSBarry Smith PetscErrorCode ierr; 252789d8953SBarry Smith MatMFFD ctx; 253a283635eSDmitry Karpeev PetscBool iascii, viewbase, viewfunction; 254a283635eSDmitry Karpeev const char *prefix; 255e884886fSBarry Smith 256e884886fSBarry Smith PetscFunctionBegin; 257789d8953SBarry Smith ierr = MatShellGetContext(J,&ctx);CHKERRQ(ierr); 258251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 259e884886fSBarry Smith if (iascii) { 260a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr); 261a283635eSDmitry Karpeev ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 26257622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"err=%g (relative error in function evaluation)\n",(double)ctx->error_rel);CHKERRQ(ierr); 2637adad957SLisandro Dalcin if (!((PetscObject)ctx)->type_name) { 264e884886fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr); 265e884886fSBarry Smith } else { 2667adad957SLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr); 267e884886fSBarry Smith } 268c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 269c92e3469SBarry Smith if (ctx->usecomplex) { 270c92e3469SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Using Lyness complex number trick to compute the matrix-vector product\n");CHKERRQ(ierr); 271c92e3469SBarry Smith } 272c92e3469SBarry Smith #endif 273e884886fSBarry Smith if (ctx->ops->view) { 274e884886fSBarry Smith ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr); 275e884886fSBarry Smith } 276a283635eSDmitry Karpeev ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr); 277a283635eSDmitry Karpeev 278c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr); 279a283635eSDmitry Karpeev if (viewbase) { 280a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr); 281a283635eSDmitry Karpeev ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr); 282a283635eSDmitry Karpeev } 283c5929fdfSBarry Smith ierr = PetscOptionsHasName(((PetscObject)J)->options,prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr); 284a283635eSDmitry Karpeev if (viewfunction) { 285a283635eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr); 286a283635eSDmitry Karpeev ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr); 287a283635eSDmitry Karpeev } 288a283635eSDmitry Karpeev ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 28911aeaf0aSBarry Smith } 290e884886fSBarry Smith PetscFunctionReturn(0); 291e884886fSBarry Smith } 292e884886fSBarry Smith 293e884886fSBarry Smith /* 294e884886fSBarry Smith MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This 295e884886fSBarry Smith allows the user to indicate the beginning of a new linear solve by calling 296e884886fSBarry Smith MatAssemblyXXX() on the matrix free matrix. This then allows the 2971d0fab5eSBarry Smith MatCreateMFFD_WP() to properly compute ||U|| only the first time 298e884886fSBarry Smith in the linear solver rather than every time. 2995a576424SJed Brown 300cc2e6a90SBarry Smith This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library hence 301cc2e6a90SBarry Smith it must be labeled as PETSC_EXTERN 302e884886fSBarry Smith */ 3035a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt) 304e884886fSBarry Smith { 305e884886fSBarry Smith PetscErrorCode ierr; 306789d8953SBarry Smith MatMFFD j; 307e884886fSBarry Smith 308e884886fSBarry Smith PetscFunctionBegin; 309789d8953SBarry Smith ierr = MatShellGetContext(J,&j);CHKERRQ(ierr); 310e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 311e884886fSBarry Smith PetscFunctionReturn(0); 312e884886fSBarry Smith } 313e884886fSBarry Smith 314e884886fSBarry Smith /* 315e884886fSBarry Smith MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a: 316e884886fSBarry Smith 317e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 318e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 319e884886fSBarry Smith u = current iterate 320e884886fSBarry Smith h = difference interval 321e884886fSBarry Smith */ 32240244768SBarry Smith static PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y) 323e884886fSBarry Smith { 324789d8953SBarry Smith MatMFFD ctx; 325e884886fSBarry Smith PetscScalar h; 326e884886fSBarry Smith Vec w,U,F; 327e884886fSBarry Smith PetscErrorCode ierr; 328ace3abfcSBarry Smith PetscBool zeroa; 329e884886fSBarry Smith 330e884886fSBarry Smith PetscFunctionBegin; 331789d8953SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 332ce94432eSBarry 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"); 333e884886fSBarry Smith /* We log matrix-free matrix-vector products separately, so that we can 334e884886fSBarry Smith separate the performance monitoring from the cases that use conventional 335e884886fSBarry Smith storage. We may eventually modify event logging to associate events 336e884886fSBarry Smith with particular objects, hence alleviating the more general problem. */ 337e884886fSBarry Smith ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 338e884886fSBarry Smith 339e884886fSBarry Smith w = ctx->w; 340e884886fSBarry Smith U = ctx->current_u; 3413ec795f1SBarry Smith F = ctx->current_f; 342e884886fSBarry Smith /* 343e884886fSBarry Smith Compute differencing parameter 344e884886fSBarry Smith */ 3454a2f8832SBarry Smith if (!((PetscObject)ctx)->type_name) { 346e884886fSBarry Smith ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr); 3479c6ac3b3SBarry Smith ierr = MatSetFromOptions(mat);CHKERRQ(ierr); 348e884886fSBarry Smith } 349e884886fSBarry Smith ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr); 350e884886fSBarry Smith if (zeroa) { 351e884886fSBarry Smith ierr = VecSet(y,0.0);CHKERRQ(ierr); 352e884886fSBarry Smith PetscFunctionReturn(0); 353e884886fSBarry Smith } 354e884886fSBarry Smith 35584d44b13SHong Zhang if (mat->erroriffailure && PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h"); 356e884886fSBarry Smith if (ctx->checkh) { 357e884886fSBarry Smith ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr); 358e884886fSBarry Smith } 359e884886fSBarry Smith 360e884886fSBarry Smith /* keep a record of the current differencing parameter h */ 361e884886fSBarry Smith ctx->currenth = h; 362e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX) 36357622a8eSBarry Smith ierr = PetscInfo2(mat,"Current differencing parameter: %g + %g i\n",(double)PetscRealPart(h),(double)PetscImaginaryPart(h));CHKERRQ(ierr); 364e884886fSBarry Smith #else 365e884886fSBarry Smith ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr); 366e884886fSBarry Smith #endif 367e884886fSBarry Smith if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 368e884886fSBarry Smith ctx->historyh[ctx->ncurrenth] = h; 369e884886fSBarry Smith } 370e884886fSBarry Smith ctx->ncurrenth++; 371e884886fSBarry Smith 372c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 373c92e3469SBarry Smith if (ctx->usecomplex) h = PETSC_i*h; 374c92e3469SBarry Smith #endif 375c92e3469SBarry Smith 376e884886fSBarry Smith /* w = u + ha */ 377e884886fSBarry Smith ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr); 378e884886fSBarry Smith 3791b797266SDmitry Karpeev /* compute func(U) as base for differencing; only needed first time in and not when provided by user */ 3801b797266SDmitry Karpeev if (ctx->ncurrenth == 1 && ctx->current_f_allocated) { 381e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr); 38252121784SBarry Smith } 383e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr); 384e884886fSBarry Smith 385c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 386c92e3469SBarry Smith if (ctx->usecomplex) { 387c92e3469SBarry Smith ierr = VecImaginaryPart(y);CHKERRQ(ierr); 388c92e3469SBarry Smith h = PetscImaginaryPart(h); 389c92e3469SBarry Smith } else { 390e884886fSBarry Smith ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 391c92e3469SBarry Smith } 392c92e3469SBarry Smith #else 393c92e3469SBarry Smith ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr); 394c92e3469SBarry Smith #endif 395e884886fSBarry Smith ierr = VecScale(y,1.0/h);CHKERRQ(ierr); 39639601f4eSBarry Smith if (mat->nullsp) {ierr = MatNullSpaceRemove(mat->nullsp,y);CHKERRQ(ierr);} 397e884886fSBarry Smith 398e884886fSBarry Smith ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr); 399e884886fSBarry Smith PetscFunctionReturn(0); 400e884886fSBarry Smith } 401e884886fSBarry Smith 402e884886fSBarry Smith /* 403e884886fSBarry Smith MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix 404e884886fSBarry Smith 405e884886fSBarry Smith y ~= (F(u + ha) - F(u))/h, 406e884886fSBarry Smith where F = nonlinear function, as set by SNESSetFunction() 407e884886fSBarry Smith u = current iterate 408e884886fSBarry Smith h = difference interval 409e884886fSBarry Smith */ 4105d21e1f8SStefano Zampini PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a) 411e884886fSBarry Smith { 412789d8953SBarry Smith MatMFFD ctx; 413e884886fSBarry Smith PetscScalar h,*aa,*ww,v; 414e884886fSBarry Smith PetscReal epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON; 415e884886fSBarry Smith Vec w,U; 416e884886fSBarry Smith PetscErrorCode ierr; 417e884886fSBarry Smith PetscInt i,rstart,rend; 418e884886fSBarry Smith 419e884886fSBarry Smith PetscFunctionBegin; 420789d8953SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 421486afcceSStefano Zampini if (!ctx->func) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunction() first"); 422486afcceSStefano Zampini if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Requires calling MatMFFDSetFunctioni() first"); 423e884886fSBarry Smith w = ctx->w; 424e884886fSBarry Smith U = ctx->current_u; 425e884886fSBarry Smith ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr); 426789d8953SBarry Smith if (ctx->funcisetbase) { 427e884886fSBarry Smith ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr); 428789d8953SBarry Smith } 429e884886fSBarry Smith ierr = VecCopy(U,w);CHKERRQ(ierr); 430e884886fSBarry Smith 431e884886fSBarry Smith ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr); 432e884886fSBarry Smith ierr = VecGetArray(a,&aa);CHKERRQ(ierr); 433e884886fSBarry Smith for (i=rstart; i<rend; i++) { 434e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 435e884886fSBarry Smith h = ww[i-rstart]; 436e884886fSBarry Smith if (h == 0.0) h = 1.0; 437e884886fSBarry Smith if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0) h = umin; 438e884886fSBarry Smith else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin; 439e884886fSBarry Smith h *= epsilon; 440e884886fSBarry Smith 441e884886fSBarry Smith ww[i-rstart] += h; 442e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 443e884886fSBarry Smith ierr = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr); 444e884886fSBarry Smith aa[i-rstart] = (v - aa[i-rstart])/h; 445e884886fSBarry Smith 446e884886fSBarry Smith ierr = VecGetArray(w,&ww);CHKERRQ(ierr); 447e884886fSBarry Smith ww[i-rstart] -= h; 448e884886fSBarry Smith ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr); 449e884886fSBarry Smith } 450e884886fSBarry Smith ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr); 451e884886fSBarry Smith PetscFunctionReturn(0); 452e884886fSBarry Smith } 453e884886fSBarry Smith 454d8fa6a54SBarry Smith PETSC_EXTERN PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F) 455e884886fSBarry Smith { 456e884886fSBarry Smith PetscErrorCode ierr; 457789d8953SBarry Smith MatMFFD ctx; 458e884886fSBarry Smith 459e884886fSBarry Smith PetscFunctionBegin; 460789d8953SBarry Smith ierr = MatShellGetContext(J,&ctx);CHKERRQ(ierr); 461e884886fSBarry Smith ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr); 462c51ad4d4SStefano Zampini if (!ctx->current_u) { 463c51ad4d4SStefano Zampini ierr = VecDuplicate(U,&ctx->current_u);CHKERRQ(ierr); 4648860a134SJunchao Zhang ierr = VecLockReadPush(ctx->current_u);CHKERRQ(ierr); 465c51ad4d4SStefano Zampini } 4668860a134SJunchao Zhang ierr = VecLockReadPop(ctx->current_u);CHKERRQ(ierr); 467c51ad4d4SStefano Zampini ierr = VecCopy(U,ctx->current_u);CHKERRQ(ierr); 4688860a134SJunchao Zhang ierr = VecLockReadPush(ctx->current_u);CHKERRQ(ierr); 46952121784SBarry Smith if (F) { 4706bf464f9SBarry Smith if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);} 4713ec795f1SBarry Smith ctx->current_f = F; 472cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_FALSE; 473cfe22f5eSBarry Smith } else if (!ctx->current_f_allocated) { 47406c4b6baSBarry Smith ierr = MatCreateVecs(J,NULL,&ctx->current_f);CHKERRQ(ierr); 4752205254eSKarl Rupp 476cfe22f5eSBarry Smith ctx->current_f_allocated = PETSC_TRUE; 47752121784SBarry Smith } 478e884886fSBarry Smith if (!ctx->w) { 479e884886fSBarry Smith ierr = VecDuplicate(ctx->current_u,&ctx->w);CHKERRQ(ierr); 480e884886fSBarry Smith } 481e884886fSBarry Smith J->assembled = PETSC_TRUE; 482e884886fSBarry Smith PetscFunctionReturn(0); 483e884886fSBarry Smith } 4845a576424SJed Brown 485e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/ 486b2573a8aSBarry Smith 48740244768SBarry Smith static PetscErrorCode MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx) 488e884886fSBarry Smith { 489789d8953SBarry Smith MatMFFD ctx; 490789d8953SBarry Smith PetscErrorCode ierr; 491e884886fSBarry Smith 492e884886fSBarry Smith PetscFunctionBegin; 493789d8953SBarry Smith ierr = MatShellGetContext(J,&ctx);CHKERRQ(ierr); 494e884886fSBarry Smith ctx->checkh = fun; 495e884886fSBarry Smith ctx->checkhctx = ectx; 496e884886fSBarry Smith PetscFunctionReturn(0); 497e884886fSBarry Smith } 498e884886fSBarry Smith 4996aa9148fSLisandro Dalcin /*@C 5006aa9148fSLisandro Dalcin MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all 5016aa9148fSLisandro Dalcin MatMFFD options in the database. 5026aa9148fSLisandro Dalcin 5036aa9148fSLisandro Dalcin Collective on Mat 5046aa9148fSLisandro Dalcin 5056aa9148fSLisandro Dalcin Input Parameter: 5066aa9148fSLisandro Dalcin + A - the Mat context 5076aa9148fSLisandro Dalcin - prefix - the prefix to prepend to all option names 5086aa9148fSLisandro Dalcin 5096aa9148fSLisandro Dalcin Notes: 5106aa9148fSLisandro Dalcin A hyphen (-) must NOT be given at the beginning of the prefix name. 5116aa9148fSLisandro Dalcin The first character of all runtime options is AUTOMATICALLY the hyphen. 5126aa9148fSLisandro Dalcin 5136aa9148fSLisandro Dalcin Level: advanced 5146aa9148fSLisandro Dalcin 515755b7f64SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF(), MatCreateMFFD() 5166aa9148fSLisandro Dalcin @*/ 5177087cfbeSBarry Smith PetscErrorCode MatMFFDSetOptionsPrefix(Mat mat,const char prefix[]) 5186aa9148fSLisandro Dalcin 5196aa9148fSLisandro Dalcin { 520789d8953SBarry Smith MatMFFD mfctx; 5216aa9148fSLisandro Dalcin PetscErrorCode ierr; 5225fd66863SKarl Rupp 5236aa9148fSLisandro Dalcin PetscFunctionBegin; 5240700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 525789d8953SBarry Smith ierr = MatShellGetContext(mat,&mfctx);CHKERRQ(ierr); 5260700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 5276aa9148fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr); 5286aa9148fSLisandro Dalcin PetscFunctionReturn(0); 5296aa9148fSLisandro Dalcin } 5306aa9148fSLisandro Dalcin 53140244768SBarry Smith static PetscErrorCode MatSetFromOptions_MFFD(PetscOptionItems *PetscOptionsObject,Mat mat) 532e884886fSBarry Smith { 533789d8953SBarry Smith MatMFFD mfctx; 534e884886fSBarry Smith PetscErrorCode ierr; 535ace3abfcSBarry Smith PetscBool flg; 536e884886fSBarry Smith char ftype[256]; 537e884886fSBarry Smith 538e884886fSBarry Smith PetscFunctionBegin; 5390700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 540789d8953SBarry Smith ierr = MatShellGetContext(mat,&mfctx);CHKERRQ(ierr); 5410700a824SBarry Smith PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1); 5423194b578SJed Brown ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr); 543a264d7a6SBarry Smith ierr = PetscOptionsFList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr); 544e884886fSBarry Smith if (flg) { 545e884886fSBarry Smith ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr); 546e884886fSBarry Smith } 547e884886fSBarry Smith 548e884886fSBarry Smith ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr); 549e884886fSBarry Smith ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr); 550e884886fSBarry Smith 55190d69ab7SBarry Smith flg = PETSC_FALSE; 5520298fd71SBarry Smith ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr); 553e884886fSBarry Smith if (flg) { 554e884886fSBarry Smith ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr); 555e884886fSBarry Smith } 556c92e3469SBarry Smith #if defined(PETSC_USE_COMPLEX) 557c92e3469SBarry 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); 558c92e3469SBarry Smith #endif 559e884886fSBarry Smith if (mfctx->ops->setfromoptions) { 560e55864a3SBarry Smith ierr = (*mfctx->ops->setfromoptions)(PetscOptionsObject,mfctx);CHKERRQ(ierr); 561e884886fSBarry Smith } 562e884886fSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 563e884886fSBarry Smith PetscFunctionReturn(0); 564e884886fSBarry Smith } 565e884886fSBarry Smith 56640244768SBarry Smith static PetscErrorCode MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period) 567bc0f08ceSBarry Smith { 568789d8953SBarry Smith MatMFFD ctx; 569789d8953SBarry Smith PetscErrorCode ierr; 570bc0f08ceSBarry Smith 571bc0f08ceSBarry Smith PetscFunctionBegin; 572789d8953SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 573bc0f08ceSBarry Smith ctx->recomputeperiod = period; 574bc0f08ceSBarry Smith PetscFunctionReturn(0); 575bc0f08ceSBarry Smith } 576bc0f08ceSBarry Smith 57740244768SBarry Smith static PetscErrorCode MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 578bc0f08ceSBarry Smith { 579789d8953SBarry Smith MatMFFD ctx; 580789d8953SBarry Smith PetscErrorCode ierr; 581bc0f08ceSBarry Smith 582bc0f08ceSBarry Smith PetscFunctionBegin; 583789d8953SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 584bc0f08ceSBarry Smith ctx->func = func; 585bc0f08ceSBarry Smith ctx->funcctx = funcctx; 586bc0f08ceSBarry Smith PetscFunctionReturn(0); 587bc0f08ceSBarry Smith } 588bc0f08ceSBarry Smith 58940244768SBarry Smith static PetscErrorCode MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error) 590bc0f08ceSBarry Smith { 591789d8953SBarry Smith MatMFFD ctx; 592789d8953SBarry Smith PetscErrorCode ierr; 593bc0f08ceSBarry Smith 594bc0f08ceSBarry Smith PetscFunctionBegin; 595789d8953SBarry Smith ierr = MatShellGetContext(mat,&ctx);CHKERRQ(ierr); 596bc0f08ceSBarry Smith if (error != PETSC_DEFAULT) ctx->error_rel = error; 597bc0f08ceSBarry Smith PetscFunctionReturn(0); 598bc0f08ceSBarry Smith } 599bc0f08ceSBarry Smith 600789d8953SBarry Smith PetscErrorCode MatMFFDSetHHistory_MFFD(Mat J,PetscScalar history[],PetscInt nhistory) 6013b49f96aSBarry Smith { 602789d8953SBarry Smith MatMFFD ctx; 603789d8953SBarry Smith PetscErrorCode ierr; 604789d8953SBarry Smith 6053b49f96aSBarry Smith PetscFunctionBegin; 606789d8953SBarry Smith ierr = MatShellGetContext(J,&ctx);CHKERRQ(ierr); 607789d8953SBarry Smith ctx->historyh = history; 608789d8953SBarry Smith ctx->maxcurrenth = nhistory; 609789d8953SBarry Smith ctx->currenth = 0.; 6103b49f96aSBarry Smith PetscFunctionReturn(0); 6113b49f96aSBarry Smith } 6123b49f96aSBarry Smith 613e884886fSBarry Smith /*MC 614e884886fSBarry Smith MATMFFD - MATMFFD = "mffd" - A matrix free matrix type. 615e884886fSBarry Smith 616e884886fSBarry Smith Level: advanced 617e884886fSBarry Smith 618789d8953SBarry Smith Developers Note: This is implemented on top of MATSHELL to get support for scaling and shifting without requiring duplicate code 619789d8953SBarry Smith 620755b7f64SBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction(), MatMFFDSetType(), 621755b7f64SBarry Smith MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 622755b7f64SBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 623755b7f64SBarry Smith MatMFFDGetH(), 624e884886fSBarry Smith M*/ 6258cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_MFFD(Mat A) 626e884886fSBarry Smith { 627e884886fSBarry Smith MatMFFD mfctx; 628e884886fSBarry Smith PetscErrorCode ierr; 629e884886fSBarry Smith 630e884886fSBarry Smith PetscFunctionBegin; 631607a6623SBarry Smith ierr = MatMFFDInitializePackage();CHKERRQ(ierr); 632e884886fSBarry Smith 633789d8953SBarry Smith ierr = PetscHeaderCreate(mfctx,MATMFFD_CLASSID,"MatMFFD","Matrix-free Finite Differencing","Mat",PetscObjectComm((PetscObject)A),NULL,NULL);CHKERRQ(ierr); 6342205254eSKarl Rupp 635e884886fSBarry Smith mfctx->error_rel = PETSC_SQRT_MACHINE_EPSILON; 636e884886fSBarry Smith mfctx->recomputeperiod = 1; 637e884886fSBarry Smith mfctx->count = 0; 638e884886fSBarry Smith mfctx->currenth = 0.0; 6390298fd71SBarry Smith mfctx->historyh = NULL; 640e884886fSBarry Smith mfctx->ncurrenth = 0; 641e884886fSBarry Smith mfctx->maxcurrenth = 0; 6427adad957SLisandro Dalcin ((PetscObject)mfctx)->type_name = 0; 643e884886fSBarry Smith 644e884886fSBarry Smith /* 645e884886fSBarry Smith Create the empty data structure to contain compute-h routines. 646e884886fSBarry Smith These will be filled in below from the command line options or 647e884886fSBarry Smith a later call with MatMFFDSetType() or if that is not called 648e884886fSBarry Smith then it will default in the first use of MatMult_MFFD() 649e884886fSBarry Smith */ 650e884886fSBarry Smith mfctx->ops->compute = 0; 651e884886fSBarry Smith mfctx->ops->destroy = 0; 652e884886fSBarry Smith mfctx->ops->view = 0; 653e884886fSBarry Smith mfctx->ops->setfromoptions = 0; 654e884886fSBarry Smith mfctx->hctx = 0; 655e884886fSBarry Smith 656e884886fSBarry Smith mfctx->func = 0; 657e884886fSBarry Smith mfctx->funcctx = 0; 6580298fd71SBarry Smith mfctx->w = NULL; 659789d8953SBarry Smith mfctx->mat = A; 660e884886fSBarry Smith 661789d8953SBarry Smith ierr = MatSetType(A,MATSHELL);CHKERRQ(ierr); 662789d8953SBarry Smith ierr = MatShellSetContext(A,mfctx);CHKERRQ(ierr); 663789d8953SBarry Smith ierr = MatShellSetOperation(A,MATOP_MULT,(void (*)(void))MatMult_MFFD);CHKERRQ(ierr); 664789d8953SBarry Smith ierr = MatShellSetOperation(A,MATOP_DESTROY,(void (*)(void))MatDestroy_MFFD);CHKERRQ(ierr); 665789d8953SBarry Smith ierr = MatShellSetOperation(A,MATOP_VIEW,(void (*)(void))MatView_MFFD);CHKERRQ(ierr); 666789d8953SBarry Smith ierr = MatShellSetOperation(A,MATOP_ASSEMBLY_END,(void (*)(void))MatAssemblyEnd_MFFD);CHKERRQ(ierr); 667789d8953SBarry Smith ierr = MatShellSetOperation(A,MATOP_SET_FROM_OPTIONS,(void (*)(void))MatSetFromOptions_MFFD);CHKERRQ(ierr); 668e884886fSBarry Smith 669bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C",MatMFFDSetBase_MFFD);CHKERRQ(ierr); 670bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr); 671bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr); 672bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C",MatMFFDSetFunction_MFFD);CHKERRQ(ierr); 673bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr); 674bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr); 675bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr); 676bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr); 677789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetHHistory_C",MatMFFDSetHHistory_MFFD);CHKERRQ(ierr); 678789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetType_C",MatMFFDSetType_MFFD);CHKERRQ(ierr); 679789d8953SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDGetH_C",MatMFFDGetH_MFFD);CHKERRQ(ierr); 680e884886fSBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr); 681e884886fSBarry Smith PetscFunctionReturn(0); 682e884886fSBarry Smith } 683e884886fSBarry Smith 684e884886fSBarry Smith /*@ 685e884886fSBarry Smith MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF() 686e884886fSBarry Smith 687e884886fSBarry Smith Collective on Vec 688e884886fSBarry Smith 689e884886fSBarry Smith Input Parameters: 690fef1beadSBarry Smith + comm - MPI communicator 691fef1beadSBarry Smith . m - number of local rows (or PETSC_DECIDE to have calculated if M is given) 692fef1beadSBarry Smith This value should be the same as the local size used in creating the 693fef1beadSBarry Smith y vector for the matrix-vector product y = Ax. 694fef1beadSBarry Smith . n - This value should be the same as the local size used in creating the 695fef1beadSBarry Smith x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have 696fef1beadSBarry Smith calculated if N is given) For square matrices n is almost always m. 697fef1beadSBarry Smith . M - number of global rows (or PETSC_DETERMINE to have calculated if m is given) 698fef1beadSBarry Smith - N - number of global columns (or PETSC_DETERMINE to have calculated if n is given) 699fef1beadSBarry Smith 700e884886fSBarry Smith 701e884886fSBarry Smith Output Parameter: 702e884886fSBarry Smith . J - the matrix-free matrix 703e884886fSBarry Smith 7049c6ac3b3SBarry Smith Options Database Keys: call MatSetFromOptions() to trigger these 7059c6ac3b3SBarry Smith + -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS) 706c92e3469SBarry Smith . -mat_mffd_err - square root of estimated relative error in function evaluation 707c92e3469SBarry Smith . -mat_mffd_period - how often h is recomputed, defaults to 1, everytime 708c92e3469SBarry Smith . -mat_mffd_check_positivity - possibly decrease h until U + h*a has only positive values 709c92e3469SBarry Smith - -mat_mffd_complex - use the Lyness trick with complex numbers to compute the matrix-vector product instead of differencing 710c92e3469SBarry Smith (requires real valued functions but that PETSc be configured for complex numbers) 7119c6ac3b3SBarry Smith 7129c6ac3b3SBarry Smith 713e884886fSBarry Smith Level: advanced 714e884886fSBarry Smith 715e884886fSBarry Smith Notes: 716e884886fSBarry Smith The matrix-free matrix context merely contains the function pointers 717e884886fSBarry Smith and work space for performing finite difference approximations of 718e884886fSBarry Smith Jacobian-vector products, F'(u)*a, 719e884886fSBarry Smith 720e884886fSBarry Smith The default code uses the following approach to compute h 721e884886fSBarry Smith 722e884886fSBarry Smith .vb 723e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 724e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 725e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 726e884886fSBarry Smith where 727e884886fSBarry Smith error_rel = square root of relative error in function evaluation 728e884886fSBarry Smith umin = minimum iterate parameter 729e884886fSBarry Smith .ve 730e884886fSBarry Smith 731ca93e954SBarry Smith You can call SNESSetJacobian() with MatMFFDComputeJacobian() if you are using matrix and not a different 732ca93e954SBarry Smith preconditioner matrix 733ca93e954SBarry Smith 734e884886fSBarry Smith The user can set the error_rel via MatMFFDSetFunctionError() and 735a7f22e61SSatish Balay umin via MatMFFDDSSetUmin(); see Users-Manual: ch_snes for details. 736e884886fSBarry Smith 737e884886fSBarry Smith The user should call MatDestroy() when finished with the matrix-free 738e884886fSBarry Smith matrix context. 739e884886fSBarry Smith 740e884886fSBarry Smith Options Database Keys: 741e884886fSBarry Smith + -mat_mffd_err <error_rel> - Sets error_rel 742e884886fSBarry Smith . -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only) 743e884886fSBarry Smith - -mat_mffd_check_positivity 744e884886fSBarry Smith 7451d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction() 746e884886fSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(), 74781242352SJed Brown MatMFFDGetH(), MatMFFDRegister(), MatMFFDComputeJacobian() 748e884886fSBarry Smith 749e884886fSBarry Smith @*/ 7507087cfbeSBarry Smith PetscErrorCode MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J) 751e884886fSBarry Smith { 752e884886fSBarry Smith PetscErrorCode ierr; 753e884886fSBarry Smith 754e884886fSBarry Smith PetscFunctionBegin; 755e884886fSBarry Smith ierr = MatCreate(comm,J);CHKERRQ(ierr); 756fef1beadSBarry Smith ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr); 757e884886fSBarry Smith ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr); 758be7c243fSJed Brown ierr = MatSetUp(*J);CHKERRQ(ierr); 759e884886fSBarry Smith PetscFunctionReturn(0); 760e884886fSBarry Smith } 761e884886fSBarry Smith 762e884886fSBarry Smith /*@ 763e884886fSBarry Smith MatMFFDGetH - Gets the last value that was used as the differencing 764e884886fSBarry Smith parameter. 765e884886fSBarry Smith 766e884886fSBarry Smith Not Collective 767e884886fSBarry Smith 768e884886fSBarry Smith Input Parameters: 769e884886fSBarry Smith . mat - the matrix obtained with MatCreateSNESMF() 770e884886fSBarry Smith 771*fd292e60Sprj- Output Parameter: 772e884886fSBarry Smith . h - the differencing step size 773e884886fSBarry Smith 774e884886fSBarry Smith Level: advanced 775e884886fSBarry Smith 7761d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory() 777e884886fSBarry Smith @*/ 7787087cfbeSBarry Smith PetscErrorCode MatMFFDGetH(Mat mat,PetscScalar *h) 779e884886fSBarry Smith { 780bc0f08ceSBarry Smith PetscErrorCode ierr; 781e884886fSBarry Smith 782e884886fSBarry Smith PetscFunctionBegin; 78388b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 78488b4c220SStefano Zampini PetscValidPointer(h,2); 785789d8953SBarry Smith ierr = PetscUseMethod(mat,"MatMFFDGetH_C",(Mat,PetscScalar*),(mat,h));CHKERRQ(ierr); 786e884886fSBarry Smith PetscFunctionReturn(0); 787e884886fSBarry Smith } 788e884886fSBarry Smith 789e884886fSBarry Smith /*@C 790e884886fSBarry Smith MatMFFDSetFunction - Sets the function used in applying the matrix free. 791e884886fSBarry Smith 7923f9fe445SBarry Smith Logically Collective on Mat 793e884886fSBarry Smith 794e884886fSBarry Smith Input Parameters: 79514f633e2SBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() or MatCreateMFFD() 796e884886fSBarry Smith . func - the function to use 797e884886fSBarry Smith - funcctx - optional function context passed to function 798e884886fSBarry Smith 79914f633e2SBarry Smith Calling Sequence of func: 80014f633e2SBarry Smith $ func (void *funcctx, Vec x, Vec f) 80114f633e2SBarry Smith 80214f633e2SBarry Smith + funcctx - user provided context 80314f633e2SBarry Smith . x - input vector 80414f633e2SBarry Smith - f - computed output function 80514f633e2SBarry Smith 806e884886fSBarry Smith Level: advanced 807e884886fSBarry Smith 808e884886fSBarry Smith Notes: 809e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 810e884886fSBarry Smith matrix inside your compute Jacobian routine 811e884886fSBarry Smith 8123ec795f1SBarry Smith If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used. 813e884886fSBarry Smith 8141d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD, 8151d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction() 816e884886fSBarry Smith @*/ 8177087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx) 818e884886fSBarry Smith { 819bc0f08ceSBarry Smith PetscErrorCode ierr; 820e884886fSBarry Smith 821e884886fSBarry Smith PetscFunctionBegin; 82288b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 823bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr); 824e884886fSBarry Smith PetscFunctionReturn(0); 825e884886fSBarry Smith } 826e884886fSBarry Smith 827e884886fSBarry Smith /*@C 828e884886fSBarry Smith MatMFFDSetFunctioni - Sets the function for a single component 829e884886fSBarry Smith 8303f9fe445SBarry Smith Logically Collective on Mat 831e884886fSBarry Smith 832e884886fSBarry Smith Input Parameters: 833e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 834e884886fSBarry Smith - funci - the function to use 835e884886fSBarry Smith 836e884886fSBarry Smith Level: advanced 837e884886fSBarry Smith 838e884886fSBarry Smith Notes: 839e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 84094eb4328SStefano Zampini matrix inside your compute Jacobian routine. 84194eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 842486afcceSStefano Zampini funci must not contain any MPI call as it is called inside a loop on the local portion of the vector. 843e884886fSBarry Smith 84494eb4328SStefano Zampini .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 8451d0fab5eSBarry Smith 846e884886fSBarry Smith @*/ 8477087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*)) 848e884886fSBarry Smith { 8494ac538c5SBarry Smith PetscErrorCode ierr; 850e884886fSBarry Smith 851e884886fSBarry Smith PetscFunctionBegin; 8520700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 8534ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr); 854e884886fSBarry Smith PetscFunctionReturn(0); 855e884886fSBarry Smith } 856e884886fSBarry Smith 857e884886fSBarry Smith /*@C 858e884886fSBarry Smith MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation 859e884886fSBarry Smith 8603f9fe445SBarry Smith Logically Collective on Mat 861e884886fSBarry Smith 862e884886fSBarry Smith Input Parameters: 863e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 864e884886fSBarry Smith - func - the function to use 865e884886fSBarry Smith 866e884886fSBarry Smith Level: advanced 867e884886fSBarry Smith 868e884886fSBarry Smith Notes: 869e884886fSBarry Smith If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free 87094eb4328SStefano Zampini matrix inside your compute Jacobian routine. 87194eb4328SStefano Zampini This function is necessary to compute the diagonal of the matrix. 872e884886fSBarry Smith 873e884886fSBarry Smith 874e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 87594eb4328SStefano Zampini MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction(), MatGetDiagonal() 876e884886fSBarry Smith @*/ 8777087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec)) 878e884886fSBarry Smith { 8794ac538c5SBarry Smith PetscErrorCode ierr; 880e884886fSBarry Smith 881e884886fSBarry Smith PetscFunctionBegin; 8820700a824SBarry Smith PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 8834ac538c5SBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr); 884e884886fSBarry Smith PetscFunctionReturn(0); 885e884886fSBarry Smith } 886e884886fSBarry Smith 887e884886fSBarry Smith /*@ 888e884886fSBarry Smith MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime 889e884886fSBarry Smith 8903f9fe445SBarry Smith Logically Collective on Mat 891e884886fSBarry Smith 892e884886fSBarry Smith Input Parameters: 893e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateSNESMF() 894e884886fSBarry Smith - period - 1 for everytime, 2 for every second etc 895e884886fSBarry Smith 896e884886fSBarry Smith Options Database Keys: 897a2b725a8SWilliam Gropp . -mat_mffd_period <period> 898e884886fSBarry Smith 899e884886fSBarry Smith Level: advanced 900e884886fSBarry Smith 901e884886fSBarry Smith 902e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), 9031d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 904e884886fSBarry Smith @*/ 9057087cfbeSBarry Smith PetscErrorCode MatMFFDSetPeriod(Mat mat,PetscInt period) 906e884886fSBarry Smith { 907bc0f08ceSBarry Smith PetscErrorCode ierr; 908e884886fSBarry Smith 909e884886fSBarry Smith PetscFunctionBegin; 91088b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 91188b4c220SStefano Zampini PetscValidLogicalCollectiveInt(mat,period,2); 912bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr); 913e884886fSBarry Smith PetscFunctionReturn(0); 914e884886fSBarry Smith } 915e884886fSBarry Smith 916e884886fSBarry Smith /*@ 917e884886fSBarry Smith MatMFFDSetFunctionError - Sets the error_rel for the approximation of 918e884886fSBarry Smith matrix-vector products using finite differences. 919e884886fSBarry Smith 9203f9fe445SBarry Smith Logically Collective on Mat 921e884886fSBarry Smith 922e884886fSBarry Smith Input Parameters: 923e884886fSBarry Smith + mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF() 924e884886fSBarry Smith - error_rel - relative error (should be set to the square root of 925e884886fSBarry Smith the relative error in the function evaluations) 926e884886fSBarry Smith 927e884886fSBarry Smith Options Database Keys: 928a2b725a8SWilliam Gropp . -mat_mffd_err <error_rel> - Sets error_rel 929e884886fSBarry Smith 930e884886fSBarry Smith Level: advanced 931e884886fSBarry Smith 932e884886fSBarry Smith Notes: 933e884886fSBarry Smith The default matrix-free matrix-vector product routine computes 934e884886fSBarry Smith .vb 935e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 936e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 937e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 938e884886fSBarry Smith .ve 939e884886fSBarry Smith 940e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD 9411d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDResetHHistory() 942e884886fSBarry Smith @*/ 9437087cfbeSBarry Smith PetscErrorCode MatMFFDSetFunctionError(Mat mat,PetscReal error) 944e884886fSBarry Smith { 945bc0f08ceSBarry Smith PetscErrorCode ierr; 946e884886fSBarry Smith 947e884886fSBarry Smith PetscFunctionBegin; 94888b4c220SStefano Zampini PetscValidHeaderSpecific(mat,MAT_CLASSID,1); 94988b4c220SStefano Zampini PetscValidLogicalCollectiveReal(mat,error,2); 950bc0f08ceSBarry Smith ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr); 951e884886fSBarry Smith PetscFunctionReturn(0); 952e884886fSBarry Smith } 953e884886fSBarry Smith 954e884886fSBarry Smith /*@ 955e884886fSBarry Smith MatMFFDSetHHistory - Sets an array to collect a history of the 956e884886fSBarry Smith differencing values (h) computed for the matrix-free product. 957e884886fSBarry Smith 9583f9fe445SBarry Smith Logically Collective on Mat 959e884886fSBarry Smith 960e884886fSBarry Smith Input Parameters: 961e884886fSBarry Smith + J - the matrix-free matrix context 962e884886fSBarry Smith . histroy - space to hold the history 963e884886fSBarry Smith - nhistory - number of entries in history, if more entries are generated than 964e884886fSBarry Smith nhistory, then the later ones are discarded 965e884886fSBarry Smith 966e884886fSBarry Smith Level: advanced 967e884886fSBarry Smith 968e884886fSBarry Smith Notes: 969e884886fSBarry Smith Use MatMFFDResetHHistory() to reset the history counter and collect 970e884886fSBarry Smith a new batch of differencing parameters, h. 971e884886fSBarry Smith 972e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 9731d0fab5eSBarry Smith MatMFFDResetHHistory(), MatMFFDSetFunctionError() 974e884886fSBarry Smith 975e884886fSBarry Smith @*/ 9767087cfbeSBarry Smith PetscErrorCode MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory) 977e884886fSBarry Smith { 978bc0f08ceSBarry Smith PetscErrorCode ierr; 979e884886fSBarry Smith 980e884886fSBarry Smith PetscFunctionBegin; 98188b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 98288b4c220SStefano Zampini if (history) PetscValidPointer(history,2); 98388b4c220SStefano Zampini PetscValidLogicalCollectiveInt(J,nhistory,3); 984789d8953SBarry Smith ierr = PetscUseMethod(J,"MatMFFDSetHHistory_C",(Mat,PetscScalar[],PetscInt),(J,history,nhistory));CHKERRQ(ierr); 985e884886fSBarry Smith PetscFunctionReturn(0); 986e884886fSBarry Smith } 987e884886fSBarry Smith 988e884886fSBarry Smith /*@ 989e884886fSBarry Smith MatMFFDResetHHistory - Resets the counter to zero to begin 990e884886fSBarry Smith collecting a new set of differencing histories. 991e884886fSBarry Smith 9923f9fe445SBarry Smith Logically Collective on Mat 993e884886fSBarry Smith 994e884886fSBarry Smith Input Parameters: 995e884886fSBarry Smith . J - the matrix-free matrix context 996e884886fSBarry Smith 997e884886fSBarry Smith Level: advanced 998e884886fSBarry Smith 999e884886fSBarry Smith Notes: 1000e884886fSBarry Smith Use MatMFFDSetHHistory() to create the original history counter. 1001e884886fSBarry Smith 1002e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(), 10031d0fab5eSBarry Smith MatMFFDSetHHistory(), MatMFFDSetFunctionError() 1004e884886fSBarry Smith 1005e884886fSBarry Smith @*/ 10067087cfbeSBarry Smith PetscErrorCode MatMFFDResetHHistory(Mat J) 1007e884886fSBarry Smith { 1008bc0f08ceSBarry Smith PetscErrorCode ierr; 1009e884886fSBarry Smith 1010e884886fSBarry Smith PetscFunctionBegin; 101188b4c220SStefano Zampini PetscValidHeaderSpecific(J,MAT_CLASSID,1); 1012bc0f08ceSBarry Smith ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr); 1013e884886fSBarry Smith PetscFunctionReturn(0); 1014e884886fSBarry Smith } 1015e884886fSBarry Smith 1016487a658cSBarry Smith /*@ 1017e884886fSBarry Smith MatMFFDSetBase - Sets the vector U at which matrix vector products of the 1018e884886fSBarry Smith Jacobian are computed 1019e884886fSBarry Smith 10203f9fe445SBarry Smith Logically Collective on Mat 1021e884886fSBarry Smith 1022e884886fSBarry Smith Input Parameters: 1023e884886fSBarry Smith + J - the MatMFFD matrix 10243ec795f1SBarry Smith . U - the vector 1025bcddec3dSBarry Smith - F - (optional) vector that contains F(u) if it has been already computed 1026e884886fSBarry Smith 102795452b02SPatrick Sanan Notes: 102895452b02SPatrick Sanan This is rarely used directly 1029e884886fSBarry Smith 10308af5ae88SBarry Smith If F is provided then it is not recomputed. Otherwise the function is evaluated at the base 10318af5ae88SBarry Smith point during the first MatMult() after each call to MatMFFDSetBase(). 1032dff2f722SBarry Smith 1033e884886fSBarry Smith Level: advanced 1034e884886fSBarry Smith 1035e884886fSBarry Smith @*/ 10367087cfbeSBarry Smith PetscErrorCode MatMFFDSetBase(Mat J,Vec U,Vec F) 1037e884886fSBarry Smith { 10384ac538c5SBarry Smith PetscErrorCode ierr; 1039e884886fSBarry Smith 1040e884886fSBarry Smith PetscFunctionBegin; 10410700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 10420700a824SBarry Smith PetscValidHeaderSpecific(U,VEC_CLASSID,2); 10430700a824SBarry Smith if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3); 10444ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr); 1045e884886fSBarry Smith PetscFunctionReturn(0); 1046e884886fSBarry Smith } 1047e884886fSBarry Smith 1048e884886fSBarry Smith /*@C 1049e884886fSBarry Smith MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts 1050e884886fSBarry Smith it to satisfy some criteria 1051e884886fSBarry Smith 10523f9fe445SBarry Smith Logically Collective on Mat 1053e884886fSBarry Smith 1054e884886fSBarry Smith Input Parameters: 1055e884886fSBarry Smith + J - the MatMFFD matrix 1056e884886fSBarry Smith . fun - the function that checks h 1057e884886fSBarry Smith - ctx - any context needed by the function 1058e884886fSBarry Smith 1059e884886fSBarry Smith Options Database Keys: 1060e884886fSBarry Smith . -mat_mffd_check_positivity 1061e884886fSBarry Smith 1062e884886fSBarry Smith Level: advanced 1063e884886fSBarry Smith 106495452b02SPatrick Sanan Notes: 106595452b02SPatrick Sanan For example, MatMFFDCheckPositivity() insures that all entries 1066e884886fSBarry Smith of U + h*a are non-negative 1067e884886fSBarry Smith 1068755b7f64SBarry Smith The function you provide is called after the default h has been computed and allows you to 1069755b7f64SBarry Smith modify it. 1070755b7f64SBarry Smith 1071755b7f64SBarry Smith .seealso: MatMFFDCheckPositivity() 1072e884886fSBarry Smith @*/ 10737087cfbeSBarry Smith PetscErrorCode MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx) 1074e884886fSBarry Smith { 10754ac538c5SBarry Smith PetscErrorCode ierr; 1076e884886fSBarry Smith 1077e884886fSBarry Smith PetscFunctionBegin; 10780700a824SBarry Smith PetscValidHeaderSpecific(J,MAT_CLASSID,1); 10794ac538c5SBarry Smith ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr); 1080e884886fSBarry Smith PetscFunctionReturn(0); 1081e884886fSBarry Smith } 1082e884886fSBarry Smith 1083e884886fSBarry Smith /*@ 1084e884886fSBarry Smith MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or 1085e884886fSBarry Smith zero, decreases h until this is satisfied. 1086e884886fSBarry Smith 10873f9fe445SBarry Smith Logically Collective on Vec 1088e884886fSBarry Smith 1089e884886fSBarry Smith Input Parameters: 1090e884886fSBarry Smith + U - base vector that is added to 1091e884886fSBarry Smith . a - vector that is added 1092e884886fSBarry Smith . h - scaling factor on a 1093e884886fSBarry Smith - dummy - context variable (unused) 1094e884886fSBarry Smith 1095e884886fSBarry Smith Options Database Keys: 1096e884886fSBarry Smith . -mat_mffd_check_positivity 1097e884886fSBarry Smith 1098e884886fSBarry Smith Level: advanced 1099e884886fSBarry Smith 110095452b02SPatrick Sanan Notes: 110195452b02SPatrick Sanan This is rarely used directly, rather it is passed as an argument to 1102e884886fSBarry Smith MatMFFDSetCheckh() 1103e884886fSBarry Smith 1104e884886fSBarry Smith .seealso: MatMFFDSetCheckh() 1105e884886fSBarry Smith @*/ 11067087cfbeSBarry Smith PetscErrorCode MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h) 1107e884886fSBarry Smith { 1108e884886fSBarry Smith PetscReal val, minval; 1109e884886fSBarry Smith PetscScalar *u_vec, *a_vec; 1110e884886fSBarry Smith PetscErrorCode ierr; 1111e884886fSBarry Smith PetscInt i,n; 1112e884886fSBarry Smith MPI_Comm comm; 1113e884886fSBarry Smith 1114e884886fSBarry Smith PetscFunctionBegin; 111588b4c220SStefano Zampini PetscValidHeaderSpecific(U,VEC_CLASSID,2); 111688b4c220SStefano Zampini PetscValidHeaderSpecific(a,VEC_CLASSID,3); 111788b4c220SStefano Zampini PetscValidPointer(h,4); 1118e884886fSBarry Smith ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr); 1119e884886fSBarry Smith ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr); 1120e884886fSBarry Smith ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr); 1121e884886fSBarry Smith ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr); 1122c068d9bbSLisandro Dalcin minval = PetscAbsScalar(*h)*PetscRealConstant(1.01); 1123e884886fSBarry Smith for (i=0; i<n; i++) { 1124e884886fSBarry Smith if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) { 1125e884886fSBarry Smith val = PetscAbsScalar(u_vec[i]/a_vec[i]); 1126e884886fSBarry Smith if (val < minval) minval = val; 1127e884886fSBarry Smith } 1128e884886fSBarry Smith } 1129e884886fSBarry Smith ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr); 1130e884886fSBarry Smith ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr); 1131b2566f29SBarry Smith ierr = MPIU_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr); 1132e884886fSBarry Smith if (val <= PetscAbsScalar(*h)) { 11336712e2f1SBarry Smith ierr = PetscInfo2(U,"Scaling back h from %g to %g\n",(double)PetscRealPart(*h),(double)(.99*val));CHKERRQ(ierr); 1134e884886fSBarry Smith if (PetscRealPart(*h) > 0.0) *h = 0.99*val; 1135e884886fSBarry Smith else *h = -0.99*val; 1136e884886fSBarry Smith } 1137e884886fSBarry Smith PetscFunctionReturn(0); 1138e884886fSBarry Smith } 1139