xref: /petsc/src/mat/impls/mffd/mffd.c (revision 5a5764247e8e42b2972d2371d5c619b20bcb2a28)
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 EXTERN_C_BEGIN
141e884886fSBarry Smith #undef __FUNCT__
142c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase_MFFD"
1437087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase_MFFD(Mat mat,FCN1 func)
144e884886fSBarry Smith {
145e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
146e884886fSBarry Smith 
147e884886fSBarry Smith   PetscFunctionBegin;
148e884886fSBarry Smith   ctx->funcisetbase = func;
149e884886fSBarry Smith   PetscFunctionReturn(0);
150e884886fSBarry Smith }
151e884886fSBarry Smith EXTERN_C_END
152e884886fSBarry Smith 
153e884886fSBarry Smith typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
154e884886fSBarry Smith EXTERN_C_BEGIN
155e884886fSBarry Smith #undef __FUNCT__
156c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni_MFFD"
1577087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni_MFFD(Mat mat,FCN2 funci)
158e884886fSBarry Smith {
159e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
160e884886fSBarry Smith 
161e884886fSBarry Smith   PetscFunctionBegin;
162e884886fSBarry Smith   ctx->funci = funci;
163e884886fSBarry Smith   PetscFunctionReturn(0);
164e884886fSBarry Smith }
165e884886fSBarry Smith EXTERN_C_END
166e884886fSBarry Smith 
167bc0f08ceSBarry Smith EXTERN_C_BEGIN
168bc0f08ceSBarry Smith #undef __FUNCT__
169bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory_MFFD"
170bc0f08ceSBarry Smith PetscErrorCode  MatMFFDResetHHistory_MFFD(Mat J)
171bc0f08ceSBarry Smith {
172bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
173bc0f08ceSBarry Smith 
174bc0f08ceSBarry Smith   PetscFunctionBegin;
175bc0f08ceSBarry Smith   ctx->ncurrenth = 0;
176bc0f08ceSBarry Smith   PetscFunctionReturn(0);
177bc0f08ceSBarry Smith }
178bc0f08ceSBarry Smith EXTERN_C_END
179e884886fSBarry Smith 
180e884886fSBarry Smith #undef __FUNCT__
181e884886fSBarry Smith #define __FUNCT__ "MatMFFDRegister"
1827087cfbeSBarry Smith PetscErrorCode  MatMFFDRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(MatMFFD))
183e884886fSBarry Smith {
184e884886fSBarry Smith   PetscErrorCode ierr;
185e884886fSBarry Smith   char           fullname[PETSC_MAX_PATH_LEN];
186e884886fSBarry Smith 
187e884886fSBarry Smith   PetscFunctionBegin;
188140e18c1SBarry Smith   ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr);
189140e18c1SBarry Smith   ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&MatMFFDList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
190e884886fSBarry Smith   PetscFunctionReturn(0);
191e884886fSBarry Smith }
192e884886fSBarry Smith 
193e884886fSBarry Smith 
194e884886fSBarry Smith #undef __FUNCT__
195e884886fSBarry Smith #define __FUNCT__ "MatMFFDRegisterDestroy"
196e884886fSBarry Smith /*@C
197e884886fSBarry Smith    MatMFFDRegisterDestroy - Frees the list of MatMFFD methods that were
198e884886fSBarry Smith    registered by MatMFFDRegisterDynamic).
199e884886fSBarry Smith 
200e884886fSBarry Smith    Not Collective
201e884886fSBarry Smith 
202e884886fSBarry Smith    Level: developer
203e884886fSBarry Smith 
204e884886fSBarry Smith .keywords: MatMFFD, register, destroy
205e884886fSBarry Smith 
206e884886fSBarry Smith .seealso: MatMFFDRegisterDynamic), MatMFFDRegisterAll()
207e884886fSBarry Smith @*/
2087087cfbeSBarry Smith PetscErrorCode  MatMFFDRegisterDestroy(void)
209e884886fSBarry Smith {
210e884886fSBarry Smith   PetscErrorCode ierr;
211e884886fSBarry Smith 
212e884886fSBarry Smith   PetscFunctionBegin;
213140e18c1SBarry Smith   ierr = PetscFunctionListDestroy(&MatMFFDList);CHKERRQ(ierr);
2142205254eSKarl Rupp 
215e884886fSBarry Smith   MatMFFDRegisterAllCalled = PETSC_FALSE;
216e884886fSBarry Smith   PetscFunctionReturn(0);
217e884886fSBarry Smith }
218e884886fSBarry Smith 
219bc0f08ceSBarry Smith EXTERN_C_BEGIN
220bc0f08ceSBarry Smith #undef __FUNCT__
221bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDAddNullSpace_MFFD"
222bc0f08ceSBarry Smith PetscErrorCode  MatMFFDAddNullSpace_MFFD(Mat J,MatNullSpace nullsp)
223bc0f08ceSBarry Smith {
224bc0f08ceSBarry Smith   PetscErrorCode ierr;
225bc0f08ceSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
226bc0f08ceSBarry Smith 
227bc0f08ceSBarry Smith   PetscFunctionBegin;
228bc0f08ceSBarry Smith   ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr);
229bc0f08ceSBarry Smith   if (ctx->sp) { ierr = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr); }
230bc0f08ceSBarry Smith   ctx->sp = nullsp;
231bc0f08ceSBarry Smith   PetscFunctionReturn(0);
232bc0f08ceSBarry Smith }
233bc0f08ceSBarry Smith EXTERN_C_END
234bc0f08ceSBarry Smith 
235e884886fSBarry Smith /* ----------------------------------------------------------------------------------------*/
236e884886fSBarry Smith #undef __FUNCT__
237e884886fSBarry Smith #define __FUNCT__ "MatDestroy_MFFD"
238e884886fSBarry Smith PetscErrorCode MatDestroy_MFFD(Mat mat)
239e884886fSBarry Smith {
240e884886fSBarry Smith   PetscErrorCode ierr;
241e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
242e884886fSBarry Smith 
243e884886fSBarry Smith   PetscFunctionBegin;
2446bf464f9SBarry Smith   ierr = VecDestroy(&ctx->w);CHKERRQ(ierr);
2456bf464f9SBarry Smith   ierr = VecDestroy(&ctx->drscale);CHKERRQ(ierr);
2466bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dlscale);CHKERRQ(ierr);
2476bf464f9SBarry Smith   ierr = VecDestroy(&ctx->dshift);CHKERRQ(ierr);
248cfe22f5eSBarry Smith   if (ctx->current_f_allocated) {
2496bf464f9SBarry Smith     ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);
250cfe22f5eSBarry Smith   }
251e884886fSBarry Smith   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
2526bf464f9SBarry Smith   ierr      = MatNullSpaceDestroy(&ctx->sp);CHKERRQ(ierr);
2536bf464f9SBarry Smith   ierr      = PetscHeaderDestroy(&ctx);CHKERRQ(ierr);
254bf0cc555SLisandro Dalcin   mat->data = 0;
255e884886fSBarry Smith 
2560298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C","",NULL);CHKERRQ(ierr);
2570298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C","",NULL);CHKERRQ(ierr);
2580298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C","",NULL);CHKERRQ(ierr);
2590298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunction_C","",NULL);CHKERRQ(ierr);
2600298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctionError_C","",NULL);CHKERRQ(ierr);
2610298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C","",NULL);CHKERRQ(ierr);
2620298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetPeriod_C","",NULL);CHKERRQ(ierr);
2630298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDResetHHistory_C","",NULL);CHKERRQ(ierr);
2640298fd71SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDAddNullSpace_C","",NULL);CHKERRQ(ierr);
265e884886fSBarry Smith   PetscFunctionReturn(0);
266e884886fSBarry Smith }
267e884886fSBarry Smith 
268e884886fSBarry Smith #undef __FUNCT__
269e884886fSBarry Smith #define __FUNCT__ "MatView_MFFD"
270e884886fSBarry Smith /*
271e884886fSBarry Smith    MatMFFDView_MFFD - Views matrix-free parameters.
272e884886fSBarry Smith 
273e884886fSBarry Smith */
274e884886fSBarry Smith PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
275e884886fSBarry Smith {
276e884886fSBarry Smith   PetscErrorCode ierr;
277e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
278a283635eSDmitry Karpeev   PetscBool      iascii, viewbase, viewfunction;
279a283635eSDmitry Karpeev   const char     *prefix;
280e884886fSBarry Smith 
281e884886fSBarry Smith   PetscFunctionBegin;
282251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
283e884886fSBarry Smith   if (iascii) {
284a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPrintf(viewer,"Matrix-free approximation:\n");CHKERRQ(ierr);
285a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
286e884886fSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
2877adad957SLisandro Dalcin     if (!((PetscObject)ctx)->type_name) {
288e884886fSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"The compute h routine has not yet been set\n");CHKERRQ(ierr);
289e884886fSBarry Smith     } else {
2907adad957SLisandro Dalcin       ierr = PetscViewerASCIIPrintf(viewer,"Using %s compute h routine\n",((PetscObject)ctx)->type_name);CHKERRQ(ierr);
291e884886fSBarry Smith     }
292e884886fSBarry Smith     if (ctx->ops->view) {
293e884886fSBarry Smith       ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
294e884886fSBarry Smith     }
295a283635eSDmitry Karpeev     ierr = PetscObjectGetOptionsPrefix((PetscObject)J, &prefix);CHKERRQ(ierr);
296a283635eSDmitry Karpeev 
297a283635eSDmitry Karpeev     ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_base", &viewbase);CHKERRQ(ierr);
298a283635eSDmitry Karpeev     if (viewbase) {
299a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Base:\n");CHKERRQ(ierr);
300a283635eSDmitry Karpeev       ierr = VecView(ctx->current_u, viewer);CHKERRQ(ierr);
301a283635eSDmitry Karpeev     }
302a283635eSDmitry Karpeev     ierr = PetscOptionsHasName(prefix, "-mat_mffd_view_function", &viewfunction);CHKERRQ(ierr);
303a283635eSDmitry Karpeev     if (viewfunction) {
304a283635eSDmitry Karpeev       ierr = PetscViewerASCIIPrintf(viewer, "Function:\n");CHKERRQ(ierr);
305a283635eSDmitry Karpeev       ierr = VecView(ctx->current_f, viewer);CHKERRQ(ierr);
306a283635eSDmitry Karpeev     }
307a283635eSDmitry Karpeev     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
30811aeaf0aSBarry Smith   }
309e884886fSBarry Smith   PetscFunctionReturn(0);
310e884886fSBarry Smith }
311e884886fSBarry Smith 
312e884886fSBarry Smith #undef __FUNCT__
313e884886fSBarry Smith #define __FUNCT__ "MatAssemblyEnd_MFFD"
314e884886fSBarry Smith /*
315e884886fSBarry Smith    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
316e884886fSBarry Smith    allows the user to indicate the beginning of a new linear solve by calling
317e884886fSBarry Smith    MatAssemblyXXX() on the matrix free matrix. This then allows the
3181d0fab5eSBarry Smith    MatCreateMFFD_WP() to properly compute ||U|| only the first time
319e884886fSBarry Smith    in the linear solver rather than every time.
320*5a576424SJed Brown 
321*5a576424SJed Brown    This function is referenced directly from MatAssemblyEnd_SNESMF(), which may be in a different shared library.
322e884886fSBarry Smith */
323*5a576424SJed Brown PETSC_EXTERN PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
324e884886fSBarry Smith {
325e884886fSBarry Smith   PetscErrorCode ierr;
326e884886fSBarry Smith   MatMFFD        j = (MatMFFD)J->data;
327e884886fSBarry Smith 
328e884886fSBarry Smith   PetscFunctionBegin;
329e884886fSBarry Smith   ierr      = MatMFFDResetHHistory(J);CHKERRQ(ierr);
330e884886fSBarry Smith   j->vshift = 0.0;
331e884886fSBarry Smith   j->vscale = 1.0;
332e884886fSBarry Smith   PetscFunctionReturn(0);
333e884886fSBarry Smith }
334e884886fSBarry Smith 
335e884886fSBarry Smith #undef __FUNCT__
336e884886fSBarry Smith #define __FUNCT__ "MatMult_MFFD"
337e884886fSBarry Smith /*
338e884886fSBarry Smith   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
339e884886fSBarry Smith 
340e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
341e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
342e884886fSBarry Smith         u = current iterate
343e884886fSBarry Smith         h = difference interval
344e884886fSBarry Smith */
345e884886fSBarry Smith PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
346e884886fSBarry Smith {
347e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
348e884886fSBarry Smith   PetscScalar    h;
349e884886fSBarry Smith   Vec            w,U,F;
350e884886fSBarry Smith   PetscErrorCode ierr;
351ace3abfcSBarry Smith   PetscBool      zeroa;
352e884886fSBarry Smith 
353e884886fSBarry Smith   PetscFunctionBegin;
354ce94432eSBarry 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");
355e884886fSBarry Smith   /* We log matrix-free matrix-vector products separately, so that we can
356e884886fSBarry Smith      separate the performance monitoring from the cases that use conventional
357e884886fSBarry Smith      storage.  We may eventually modify event logging to associate events
358e884886fSBarry Smith      with particular objects, hence alleviating the more general problem. */
359e884886fSBarry Smith   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
360e884886fSBarry Smith 
361e884886fSBarry Smith   w = ctx->w;
362e884886fSBarry Smith   U = ctx->current_u;
3633ec795f1SBarry Smith   F = ctx->current_f;
364e884886fSBarry Smith   /*
365e884886fSBarry Smith       Compute differencing parameter
366e884886fSBarry Smith   */
367e884886fSBarry Smith   if (!ctx->ops->compute) {
368e884886fSBarry Smith     ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr);
3699c6ac3b3SBarry Smith     ierr = MatSetFromOptions(mat);CHKERRQ(ierr);
370e884886fSBarry Smith   }
371e884886fSBarry Smith   ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr);
372e884886fSBarry Smith   if (zeroa) {
373e884886fSBarry Smith     ierr = VecSet(y,0.0);CHKERRQ(ierr);
374e884886fSBarry Smith     PetscFunctionReturn(0);
375e884886fSBarry Smith   }
376e884886fSBarry Smith 
377e32f2f54SBarry Smith   if (PetscIsInfOrNanScalar(h)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Computed Nan differencing parameter h");
378e884886fSBarry Smith   if (ctx->checkh) {
379e884886fSBarry Smith     ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr);
380e884886fSBarry Smith   }
381e884886fSBarry Smith 
382e884886fSBarry Smith   /* keep a record of the current differencing parameter h */
383e884886fSBarry Smith   ctx->currenth = h;
384e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX)
385e884886fSBarry Smith   ierr = PetscInfo2(mat,"Current differencing parameter: %G + %G i\n",PetscRealPart(h),PetscImaginaryPart(h));CHKERRQ(ierr);
386e884886fSBarry Smith #else
387e884886fSBarry Smith   ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr);
388e884886fSBarry Smith #endif
389e884886fSBarry Smith   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
390e884886fSBarry Smith     ctx->historyh[ctx->ncurrenth] = h;
391e884886fSBarry Smith   }
392e884886fSBarry Smith   ctx->ncurrenth++;
393e884886fSBarry Smith 
394e884886fSBarry Smith   /* w = u + ha */
395c3bb7e23SBarry Smith   if (ctx->drscale) {
396c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->drscale,a,U);CHKERRQ(ierr);
397c3bb7e23SBarry Smith     ierr = VecAYPX(U,h,w);CHKERRQ(ierr);
398c3bb7e23SBarry Smith   } else {
399e884886fSBarry Smith     ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr);
400c3bb7e23SBarry Smith   }
401e884886fSBarry Smith 
4021b797266SDmitry Karpeev   /* compute func(U) as base for differencing; only needed first time in and not when provided by user */
4031b797266SDmitry Karpeev   if (ctx->ncurrenth == 1 && ctx->current_f_allocated) {
404e884886fSBarry Smith     ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr);
40552121784SBarry Smith   }
406e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr);
407e884886fSBarry Smith 
408e884886fSBarry Smith   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
409e884886fSBarry Smith   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
410e884886fSBarry Smith 
411c35ec66cSMatthew G Knepley   if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
412e884886fSBarry Smith     ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
413c3bb7e23SBarry Smith   }
414c3bb7e23SBarry Smith   if (ctx->dlscale) {
415c3bb7e23SBarry Smith     ierr = VecPointwiseMult(y,ctx->dlscale,y);CHKERRQ(ierr);
416c3bb7e23SBarry Smith   }
417c3bb7e23SBarry Smith   if (ctx->dshift) {
418c3bb7e23SBarry Smith     ierr = VecPointwiseMult(ctx->dshift,a,U);CHKERRQ(ierr);
419c3bb7e23SBarry Smith     ierr = VecAXPY(y,1.0,U);CHKERRQ(ierr);
420c3bb7e23SBarry Smith   }
421e884886fSBarry Smith 
4220298fd71SBarry Smith   if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,NULL);CHKERRQ(ierr);}
423e884886fSBarry Smith 
424e884886fSBarry Smith   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
425e884886fSBarry Smith   PetscFunctionReturn(0);
426e884886fSBarry Smith }
427e884886fSBarry Smith 
428e884886fSBarry Smith #undef __FUNCT__
429e884886fSBarry Smith #define __FUNCT__ "MatGetDiagonal_MFFD"
430e884886fSBarry Smith /*
431e884886fSBarry Smith   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
432e884886fSBarry Smith 
433e884886fSBarry Smith         y ~= (F(u + ha) - F(u))/h,
434e884886fSBarry Smith   where F = nonlinear function, as set by SNESSetFunction()
435e884886fSBarry Smith         u = current iterate
436e884886fSBarry Smith         h = difference interval
437e884886fSBarry Smith */
438e884886fSBarry Smith PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
439e884886fSBarry Smith {
440e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
441e884886fSBarry Smith   PetscScalar    h,*aa,*ww,v;
442e884886fSBarry Smith   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
443e884886fSBarry Smith   Vec            w,U;
444e884886fSBarry Smith   PetscErrorCode ierr;
445e884886fSBarry Smith   PetscInt       i,rstart,rend;
446e884886fSBarry Smith 
447e884886fSBarry Smith   PetscFunctionBegin;
448e7e72b3dSBarry Smith   if (!ctx->funci) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first");
449e884886fSBarry Smith 
450e884886fSBarry Smith   w    = ctx->w;
451e884886fSBarry Smith   U    = ctx->current_u;
452e884886fSBarry Smith   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
453e884886fSBarry Smith   ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
454e884886fSBarry Smith   ierr = VecCopy(U,w);CHKERRQ(ierr);
455e884886fSBarry Smith 
456e884886fSBarry Smith   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
457e884886fSBarry Smith   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
458e884886fSBarry Smith   for (i=rstart; i<rend; i++) {
459e884886fSBarry Smith     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
460e884886fSBarry Smith     h    = ww[i-rstart];
461e884886fSBarry Smith     if (h == 0.0) h = 1.0;
462e884886fSBarry Smith     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
463e884886fSBarry Smith     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
464e884886fSBarry Smith     h *= epsilon;
465e884886fSBarry Smith 
466e884886fSBarry Smith     ww[i-rstart] += h;
467e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
468e884886fSBarry Smith     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
469e884886fSBarry Smith     aa[i-rstart]  = (v - aa[i-rstart])/h;
470e884886fSBarry Smith 
471e884886fSBarry Smith     /* possibly shift and scale result */
472c35ec66cSMatthew G Knepley     if ((ctx->vshift != 0.0) || (ctx->vscale != 1.0)) {
473e884886fSBarry Smith       aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
474c3bb7e23SBarry Smith     }
475e884886fSBarry Smith 
476e884886fSBarry Smith     ierr          = VecGetArray(w,&ww);CHKERRQ(ierr);
477e884886fSBarry Smith     ww[i-rstart] -= h;
478e884886fSBarry Smith     ierr          = VecRestoreArray(w,&ww);CHKERRQ(ierr);
479e884886fSBarry Smith   }
480e884886fSBarry Smith   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
481e884886fSBarry Smith   PetscFunctionReturn(0);
482e884886fSBarry Smith }
483e884886fSBarry Smith 
484e884886fSBarry Smith #undef __FUNCT__
485c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalScale_MFFD"
486c3bb7e23SBarry Smith PetscErrorCode MatDiagonalScale_MFFD(Mat mat,Vec ll,Vec rr)
487c3bb7e23SBarry Smith {
488c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
489c3bb7e23SBarry Smith   PetscErrorCode ierr;
490c3bb7e23SBarry Smith 
491c3bb7e23SBarry Smith   PetscFunctionBegin;
492c3bb7e23SBarry Smith   if (ll && !aij->dlscale) {
493c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dlscale);CHKERRQ(ierr);
494c3bb7e23SBarry Smith   }
495c3bb7e23SBarry Smith   if (rr && !aij->drscale) {
496c3bb7e23SBarry Smith     ierr = VecDuplicate(rr,&aij->drscale);CHKERRQ(ierr);
497c3bb7e23SBarry Smith   }
498c3bb7e23SBarry Smith   if (ll) {
499c3bb7e23SBarry Smith     ierr = VecCopy(ll,aij->dlscale);CHKERRQ(ierr);
500c3bb7e23SBarry Smith   }
501c3bb7e23SBarry Smith   if (rr) {
502c3bb7e23SBarry Smith     ierr = VecCopy(rr,aij->drscale);CHKERRQ(ierr);
503c3bb7e23SBarry Smith   }
504c3bb7e23SBarry Smith   PetscFunctionReturn(0);
505c3bb7e23SBarry Smith }
506c3bb7e23SBarry Smith 
507c3bb7e23SBarry Smith #undef __FUNCT__
508c3bb7e23SBarry Smith #define __FUNCT__ "MatDiagonalSet_MFFD"
509c3bb7e23SBarry Smith PetscErrorCode MatDiagonalSet_MFFD(Mat mat,Vec ll,InsertMode mode)
510c3bb7e23SBarry Smith {
511c3bb7e23SBarry Smith   MatMFFD        aij = (MatMFFD)mat->data;
512c3bb7e23SBarry Smith   PetscErrorCode ierr;
513c3bb7e23SBarry Smith 
514c3bb7e23SBarry Smith   PetscFunctionBegin;
515ce94432eSBarry Smith   if (mode == INSERT_VALUES) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No diagonal set with INSERT_VALUES");
516c3bb7e23SBarry Smith   if (!aij->dshift) {
517c3bb7e23SBarry Smith     ierr = VecDuplicate(ll,&aij->dshift);CHKERRQ(ierr);
518c3bb7e23SBarry Smith   }
519c3bb7e23SBarry Smith   ierr = VecAXPY(aij->dshift,1.0,ll);CHKERRQ(ierr);
520c3bb7e23SBarry Smith   PetscFunctionReturn(0);
521c3bb7e23SBarry Smith }
522c3bb7e23SBarry Smith 
523c3bb7e23SBarry Smith #undef __FUNCT__
524e884886fSBarry Smith #define __FUNCT__ "MatShift_MFFD"
525e884886fSBarry Smith PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
526e884886fSBarry Smith {
527e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
5285fd66863SKarl Rupp 
529e884886fSBarry Smith   PetscFunctionBegin;
530e884886fSBarry Smith   shell->vshift += a;
531e884886fSBarry Smith   PetscFunctionReturn(0);
532e884886fSBarry Smith }
533e884886fSBarry Smith 
534e884886fSBarry Smith #undef __FUNCT__
535e884886fSBarry Smith #define __FUNCT__ "MatScale_MFFD"
536e884886fSBarry Smith PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
537e884886fSBarry Smith {
538e884886fSBarry Smith   MatMFFD shell = (MatMFFD)Y->data;
5395fd66863SKarl Rupp 
540e884886fSBarry Smith   PetscFunctionBegin;
541e884886fSBarry Smith   shell->vscale *= a;
542e884886fSBarry Smith   PetscFunctionReturn(0);
543e884886fSBarry Smith }
544e884886fSBarry Smith 
545e884886fSBarry Smith #undef __FUNCT__
546c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetBase_MFFD"
547*5a576424SJed Brown /* PETSC_EXTERN_C because this function is referenced directly from MatMFFDSetBase_SNESMF(). */
548*5a576424SJed Brown PETSC_EXTERN_C PetscErrorCode MatMFFDSetBase_MFFD(Mat J,Vec U,Vec F)
549e884886fSBarry Smith {
550e884886fSBarry Smith   PetscErrorCode ierr;
551e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
552e884886fSBarry Smith 
553e884886fSBarry Smith   PetscFunctionBegin;
554e884886fSBarry Smith   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
5552205254eSKarl Rupp 
556e884886fSBarry Smith   ctx->current_u = U;
55752121784SBarry Smith   if (F) {
5586bf464f9SBarry Smith     if (ctx->current_f_allocated) {ierr = VecDestroy(&ctx->current_f);CHKERRQ(ierr);}
5593ec795f1SBarry Smith     ctx->current_f           = F;
560cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_FALSE;
561cfe22f5eSBarry Smith   } else if (!ctx->current_f_allocated) {
56252121784SBarry Smith     ierr = VecDuplicate(ctx->current_u, &ctx->current_f);CHKERRQ(ierr);
5632205254eSKarl Rupp 
564cfe22f5eSBarry Smith     ctx->current_f_allocated = PETSC_TRUE;
56552121784SBarry Smith   }
566e884886fSBarry Smith   if (!ctx->w) {
567e884886fSBarry Smith     ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr);
568e884886fSBarry Smith   }
569e884886fSBarry Smith   J->assembled = PETSC_TRUE;
570e884886fSBarry Smith   PetscFunctionReturn(0);
571e884886fSBarry Smith }
572*5a576424SJed Brown 
573e884886fSBarry Smith typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
574e884886fSBarry Smith EXTERN_C_BEGIN
575e884886fSBarry Smith #undef __FUNCT__
576c879e56bSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh_MFFD"
5777087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh_MFFD(Mat J,FCN3 fun,void *ectx)
578e884886fSBarry Smith {
579e884886fSBarry Smith   MatMFFD ctx = (MatMFFD)J->data;
580e884886fSBarry Smith 
581e884886fSBarry Smith   PetscFunctionBegin;
582e884886fSBarry Smith   ctx->checkh    = fun;
583e884886fSBarry Smith   ctx->checkhctx = ectx;
584e884886fSBarry Smith   PetscFunctionReturn(0);
585e884886fSBarry Smith }
586e884886fSBarry Smith EXTERN_C_END
587e884886fSBarry Smith 
588e884886fSBarry Smith #undef __FUNCT__
5896aa9148fSLisandro Dalcin #define __FUNCT__ "MatMFFDSetOptionsPrefix"
5906aa9148fSLisandro Dalcin /*@C
5916aa9148fSLisandro Dalcin    MatMFFDSetOptionsPrefix - Sets the prefix used for searching for all
5926aa9148fSLisandro Dalcin    MatMFFD options in the database.
5936aa9148fSLisandro Dalcin 
5946aa9148fSLisandro Dalcin    Collective on Mat
5956aa9148fSLisandro Dalcin 
5966aa9148fSLisandro Dalcin    Input Parameter:
5976aa9148fSLisandro Dalcin +  A - the Mat context
5986aa9148fSLisandro Dalcin -  prefix - the prefix to prepend to all option names
5996aa9148fSLisandro Dalcin 
6006aa9148fSLisandro Dalcin    Notes:
6016aa9148fSLisandro Dalcin    A hyphen (-) must NOT be given at the beginning of the prefix name.
6026aa9148fSLisandro Dalcin    The first character of all runtime options is AUTOMATICALLY the hyphen.
6036aa9148fSLisandro Dalcin 
6046aa9148fSLisandro Dalcin    Level: advanced
6056aa9148fSLisandro Dalcin 
6066aa9148fSLisandro Dalcin .keywords: SNES, matrix-free, parameters
6076aa9148fSLisandro Dalcin 
6089c6ac3b3SBarry Smith .seealso: MatSetFromOptions(), MatCreateSNESMF()
6096aa9148fSLisandro Dalcin @*/
6107087cfbeSBarry Smith PetscErrorCode  MatMFFDSetOptionsPrefix(Mat mat,const char prefix[])
6116aa9148fSLisandro Dalcin 
6126aa9148fSLisandro Dalcin {
6130298fd71SBarry Smith   MatMFFD        mfctx = mat ? (MatMFFD)mat->data : (MatMFFD)NULL;
6146aa9148fSLisandro Dalcin   PetscErrorCode ierr;
6155fd66863SKarl Rupp 
6166aa9148fSLisandro Dalcin   PetscFunctionBegin;
6170700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6180700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
6196aa9148fSLisandro Dalcin   ierr = PetscObjectSetOptionsPrefix((PetscObject)mfctx,prefix);CHKERRQ(ierr);
6206aa9148fSLisandro Dalcin   PetscFunctionReturn(0);
6216aa9148fSLisandro Dalcin }
6226aa9148fSLisandro Dalcin 
6236aa9148fSLisandro Dalcin #undef __FUNCT__
6249c6ac3b3SBarry Smith #define __FUNCT__ "MatSetFromOptions_MFFD"
6259c6ac3b3SBarry Smith PetscErrorCode  MatSetFromOptions_MFFD(Mat mat)
626e884886fSBarry Smith {
62771f08665SBarry Smith   MatMFFD        mfctx = (MatMFFD)mat->data;
628e884886fSBarry Smith   PetscErrorCode ierr;
629ace3abfcSBarry Smith   PetscBool      flg;
630e884886fSBarry Smith   char           ftype[256];
631e884886fSBarry Smith 
632e884886fSBarry Smith   PetscFunctionBegin;
6330700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
6340700a824SBarry Smith   PetscValidHeaderSpecific(mfctx,MATMFFD_CLASSID,1);
6353194b578SJed Brown   ierr = PetscObjectOptionsBegin((PetscObject)mfctx);CHKERRQ(ierr);
636b022a5c1SBarry Smith   ierr = PetscOptionsList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatMFFDList,((PetscObject)mfctx)->type_name,ftype,256,&flg);CHKERRQ(ierr);
637e884886fSBarry Smith   if (flg) {
638e884886fSBarry Smith     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
639e884886fSBarry Smith   }
640e884886fSBarry Smith 
641e884886fSBarry Smith   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
642e884886fSBarry Smith   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
643e884886fSBarry Smith 
64490d69ab7SBarry Smith   flg  = PETSC_FALSE;
6450298fd71SBarry Smith   ierr = PetscOptionsBool("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",flg,&flg,NULL);CHKERRQ(ierr);
646e884886fSBarry Smith   if (flg) {
647e884886fSBarry Smith     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
648e884886fSBarry Smith   }
649e884886fSBarry Smith   if (mfctx->ops->setfromoptions) {
650e884886fSBarry Smith     ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr);
651e884886fSBarry Smith   }
652e884886fSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
653e884886fSBarry Smith   PetscFunctionReturn(0);
654e884886fSBarry Smith }
655e884886fSBarry Smith 
656bc0f08ceSBarry Smith EXTERN_C_BEGIN
657bc0f08ceSBarry Smith #undef __FUNCT__
658bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod_MFFD"
659bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetPeriod_MFFD(Mat mat,PetscInt period)
660bc0f08ceSBarry Smith {
661bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
662bc0f08ceSBarry Smith 
663bc0f08ceSBarry Smith   PetscFunctionBegin;
664bc0f08ceSBarry Smith   PetscValidLogicalCollectiveInt(mat,period,2);
665bc0f08ceSBarry Smith   ctx->recomputeperiod = period;
666bc0f08ceSBarry Smith   PetscFunctionReturn(0);
667bc0f08ceSBarry Smith }
668bc0f08ceSBarry Smith EXTERN_C_END
669bc0f08ceSBarry Smith 
670bc0f08ceSBarry Smith EXTERN_C_BEGIN
671bc0f08ceSBarry Smith #undef __FUNCT__
672bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunction_MFFD"
673bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetFunction_MFFD(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
674bc0f08ceSBarry Smith {
675bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
676bc0f08ceSBarry Smith 
677bc0f08ceSBarry Smith   PetscFunctionBegin;
678bc0f08ceSBarry Smith   ctx->func    = func;
679bc0f08ceSBarry Smith   ctx->funcctx = funcctx;
680bc0f08ceSBarry Smith   PetscFunctionReturn(0);
681bc0f08ceSBarry Smith }
682bc0f08ceSBarry Smith EXTERN_C_END
683bc0f08ceSBarry Smith 
684bc0f08ceSBarry Smith EXTERN_C_BEGIN
685bc0f08ceSBarry Smith #undef __FUNCT__
686bc0f08ceSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError_MFFD"
687bc0f08ceSBarry Smith PetscErrorCode  MatMFFDSetFunctionError_MFFD(Mat mat,PetscReal error)
688bc0f08ceSBarry Smith {
689bc0f08ceSBarry Smith   MatMFFD ctx = (MatMFFD)mat->data;
690bc0f08ceSBarry Smith 
691bc0f08ceSBarry Smith   PetscFunctionBegin;
692bc0f08ceSBarry Smith   PetscValidLogicalCollectiveReal(mat,error,2);
693bc0f08ceSBarry Smith   if (error != PETSC_DEFAULT) ctx->error_rel = error;
694bc0f08ceSBarry Smith   PetscFunctionReturn(0);
695bc0f08ceSBarry Smith }
696bc0f08ceSBarry Smith EXTERN_C_END
697bc0f08ceSBarry Smith 
698e884886fSBarry Smith /*MC
699e884886fSBarry Smith   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
700e884886fSBarry Smith 
701e884886fSBarry Smith   Level: advanced
702e884886fSBarry Smith 
703d2d6cebeSBarry Smith .seealso: MatCreateMFFD(), MatCreateSNESMF(), MatMFFDSetFunction()
704e884886fSBarry Smith M*/
705e884886fSBarry Smith EXTERN_C_BEGIN
706e884886fSBarry Smith #undef __FUNCT__
707e884886fSBarry Smith #define __FUNCT__ "MatCreate_MFFD"
7087087cfbeSBarry Smith PetscErrorCode  MatCreate_MFFD(Mat A)
709e884886fSBarry Smith {
710e884886fSBarry Smith   MatMFFD        mfctx;
711e884886fSBarry Smith   PetscErrorCode ierr;
712e884886fSBarry Smith 
713e884886fSBarry Smith   PetscFunctionBegin;
714519f805aSKarl Rupp #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
7150298fd71SBarry Smith   ierr = MatMFFDInitializePackage(NULL);CHKERRQ(ierr);
7163ec795f1SBarry Smith #endif
717e884886fSBarry Smith 
718ce94432eSBarry 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);
7192205254eSKarl Rupp 
720e884886fSBarry Smith   mfctx->sp                       = 0;
721e884886fSBarry Smith   mfctx->error_rel                = PETSC_SQRT_MACHINE_EPSILON;
722e884886fSBarry Smith   mfctx->recomputeperiod          = 1;
723e884886fSBarry Smith   mfctx->count                    = 0;
724e884886fSBarry Smith   mfctx->currenth                 = 0.0;
7250298fd71SBarry Smith   mfctx->historyh                 = NULL;
726e884886fSBarry Smith   mfctx->ncurrenth                = 0;
727e884886fSBarry Smith   mfctx->maxcurrenth              = 0;
7287adad957SLisandro Dalcin   ((PetscObject)mfctx)->type_name = 0;
729e884886fSBarry Smith 
730e884886fSBarry Smith   mfctx->vshift = 0.0;
731e884886fSBarry Smith   mfctx->vscale = 1.0;
732e884886fSBarry Smith 
733e884886fSBarry Smith   /*
734e884886fSBarry Smith      Create the empty data structure to contain compute-h routines.
735e884886fSBarry Smith      These will be filled in below from the command line options or
736e884886fSBarry Smith      a later call with MatMFFDSetType() or if that is not called
737e884886fSBarry Smith      then it will default in the first use of MatMult_MFFD()
738e884886fSBarry Smith   */
739e884886fSBarry Smith   mfctx->ops->compute        = 0;
740e884886fSBarry Smith   mfctx->ops->destroy        = 0;
741e884886fSBarry Smith   mfctx->ops->view           = 0;
742e884886fSBarry Smith   mfctx->ops->setfromoptions = 0;
743e884886fSBarry Smith   mfctx->hctx                = 0;
744e884886fSBarry Smith 
745e884886fSBarry Smith   mfctx->func    = 0;
746e884886fSBarry Smith   mfctx->funcctx = 0;
7470298fd71SBarry Smith   mfctx->w       = NULL;
748e884886fSBarry Smith 
749e884886fSBarry Smith   A->data = mfctx;
750e884886fSBarry Smith 
751e884886fSBarry Smith   A->ops->mult           = MatMult_MFFD;
752e884886fSBarry Smith   A->ops->destroy        = MatDestroy_MFFD;
753e884886fSBarry Smith   A->ops->view           = MatView_MFFD;
754e884886fSBarry Smith   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
755e884886fSBarry Smith   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
756e884886fSBarry Smith   A->ops->scale          = MatScale_MFFD;
757e884886fSBarry Smith   A->ops->shift          = MatShift_MFFD;
758c3bb7e23SBarry Smith   A->ops->diagonalscale  = MatDiagonalScale_MFFD;
759c3bb7e23SBarry Smith   A->ops->diagonalset    = MatDiagonalSet_MFFD;
7609c6ac3b3SBarry Smith   A->ops->setfromoptions = MatSetFromOptions_MFFD;
761e884886fSBarry Smith   A->assembled           = PETSC_TRUE;
762e884886fSBarry Smith 
76326283091SBarry Smith   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
76426283091SBarry Smith   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
765ee1cef2cSJed Brown 
76600de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetBase_C","MatMFFDSetBase_MFFD",MatMFFDSetBase_MFFD);CHKERRQ(ierr);
76700de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioniBase_C","MatMFFDSetFunctioniBase_MFFD",MatMFFDSetFunctioniBase_MFFD);CHKERRQ(ierr);
76800de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctioni_C","MatMFFDSetFunctioni_MFFD",MatMFFDSetFunctioni_MFFD);CHKERRQ(ierr);
76900de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunction_C","MatMFFDSetFunction_MFFD",MatMFFDSetFunction_MFFD);CHKERRQ(ierr);
77000de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetCheckh_C","MatMFFDSetCheckh_MFFD",MatMFFDSetCheckh_MFFD);CHKERRQ(ierr);
77100de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetPeriod_C","MatMFFDSetPeriod_MFFD",MatMFFDSetPeriod_MFFD);CHKERRQ(ierr);
77200de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDSetFunctionError_C","MatMFFDSetFunctionError_MFFD",MatMFFDSetFunctionError_MFFD);CHKERRQ(ierr);
77300de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDResetHHistory_C","MatMFFDResetHHistory_MFFD",MatMFFDResetHHistory_MFFD);CHKERRQ(ierr);
77400de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatMFFDAddNullSpace_C","MatMFFDAddNullSpace_MFFD",MatMFFDAddNullSpace_MFFD);CHKERRQ(ierr);
7752205254eSKarl Rupp 
776e884886fSBarry Smith   mfctx->mat = A;
7772205254eSKarl Rupp 
778e884886fSBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
779e884886fSBarry Smith   PetscFunctionReturn(0);
780e884886fSBarry Smith }
781e884886fSBarry Smith EXTERN_C_END
782e884886fSBarry Smith 
783e884886fSBarry Smith #undef __FUNCT__
784e884886fSBarry Smith #define __FUNCT__ "MatCreateMFFD"
785e884886fSBarry Smith /*@
786e884886fSBarry Smith    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
787e884886fSBarry Smith 
788e884886fSBarry Smith    Collective on Vec
789e884886fSBarry Smith 
790e884886fSBarry Smith    Input Parameters:
791fef1beadSBarry Smith +  comm - MPI communicator
792fef1beadSBarry Smith .  m - number of local rows (or PETSC_DECIDE to have calculated if M is given)
793fef1beadSBarry Smith            This value should be the same as the local size used in creating the
794fef1beadSBarry Smith            y vector for the matrix-vector product y = Ax.
795fef1beadSBarry Smith .  n - This value should be the same as the local size used in creating the
796fef1beadSBarry Smith        x vector for the matrix-vector product y = Ax. (or PETSC_DECIDE to have
797fef1beadSBarry Smith        calculated if N is given) For square matrices n is almost always m.
798fef1beadSBarry Smith .  M - number of global rows (or PETSC_DETERMINE to have calculated if m is given)
799fef1beadSBarry Smith -  N - number of global columns (or PETSC_DETERMINE to have calculated if n is given)
800fef1beadSBarry Smith 
801e884886fSBarry Smith 
802e884886fSBarry Smith    Output Parameter:
803e884886fSBarry Smith .  J - the matrix-free matrix
804e884886fSBarry Smith 
8059c6ac3b3SBarry Smith    Options Database Keys: call MatSetFromOptions() to trigger these
8069c6ac3b3SBarry Smith +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
8079c6ac3b3SBarry Smith -  -mat_mffd_err - square root of estimated relative error in function evaluation
8089c6ac3b3SBarry Smith -  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
8099c6ac3b3SBarry Smith 
8109c6ac3b3SBarry Smith 
811e884886fSBarry Smith    Level: advanced
812e884886fSBarry Smith 
813e884886fSBarry Smith    Notes:
814e884886fSBarry Smith    The matrix-free matrix context merely contains the function pointers
815e884886fSBarry Smith    and work space for performing finite difference approximations of
816e884886fSBarry Smith    Jacobian-vector products, F'(u)*a,
817e884886fSBarry Smith 
818e884886fSBarry Smith    The default code uses the following approach to compute h
819e884886fSBarry Smith 
820e884886fSBarry Smith .vb
821e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
822e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
823e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
824e884886fSBarry Smith  where
825e884886fSBarry Smith      error_rel = square root of relative error in function evaluation
826e884886fSBarry Smith      umin = minimum iterate parameter
827e884886fSBarry Smith .ve
828e884886fSBarry Smith 
829e884886fSBarry Smith    The user can set the error_rel via MatMFFDSetFunctionError() and
8300598bfebSBarry Smith    umin via MatMFFDDSSetUmin(); see the <A href="../../docs/manual.pdf#nameddest=ch_snes">SNES chapter of the users manual</A> for details.
831e884886fSBarry Smith 
832e884886fSBarry Smith    The user should call MatDestroy() when finished with the matrix-free
833e884886fSBarry Smith    matrix context.
834e884886fSBarry Smith 
835e884886fSBarry Smith    Options Database Keys:
836e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
837e884886fSBarry Smith .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
838e884886fSBarry Smith -  -mat_mffd_check_positivity
839e884886fSBarry Smith 
840e884886fSBarry Smith .keywords: default, matrix-free, create, matrix
841e884886fSBarry Smith 
8421d0fab5eSBarry Smith .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDSSetUmin(), MatMFFDSetFunction()
843e884886fSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
8441d0fab5eSBarry Smith           MatMFFDGetH(), MatMFFDRegisterDynamic), MatMFFDComputeJacobian()
845e884886fSBarry Smith 
846e884886fSBarry Smith @*/
8477087cfbeSBarry Smith PetscErrorCode  MatCreateMFFD(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,Mat *J)
848e884886fSBarry Smith {
849e884886fSBarry Smith   PetscErrorCode ierr;
850e884886fSBarry Smith 
851e884886fSBarry Smith   PetscFunctionBegin;
852e884886fSBarry Smith   ierr = MatCreate(comm,J);CHKERRQ(ierr);
853fef1beadSBarry Smith   ierr = MatSetSizes(*J,m,n,M,N);CHKERRQ(ierr);
854e884886fSBarry Smith   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
855be7c243fSJed Brown   ierr = MatSetUp(*J);CHKERRQ(ierr);
856e884886fSBarry Smith   PetscFunctionReturn(0);
857e884886fSBarry Smith }
858e884886fSBarry Smith 
859e884886fSBarry Smith 
860e884886fSBarry Smith #undef __FUNCT__
861e884886fSBarry Smith #define __FUNCT__ "MatMFFDGetH"
862e884886fSBarry Smith /*@
863e884886fSBarry Smith    MatMFFDGetH - Gets the last value that was used as the differencing
864e884886fSBarry Smith    parameter.
865e884886fSBarry Smith 
866e884886fSBarry Smith    Not Collective
867e884886fSBarry Smith 
868e884886fSBarry Smith    Input Parameters:
869e884886fSBarry Smith .  mat - the matrix obtained with MatCreateSNESMF()
870e884886fSBarry Smith 
871e884886fSBarry Smith    Output Paramter:
872e884886fSBarry Smith .  h - the differencing step size
873e884886fSBarry Smith 
874e884886fSBarry Smith    Level: advanced
875e884886fSBarry Smith 
876e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
877e884886fSBarry Smith 
8781d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD, MatMFFDResetHHistory()
879e884886fSBarry Smith @*/
8807087cfbeSBarry Smith PetscErrorCode  MatMFFDGetH(Mat mat,PetscScalar *h)
881e884886fSBarry Smith {
882e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)mat->data;
883bc0f08ceSBarry Smith   PetscErrorCode ierr;
884bc0f08ceSBarry Smith   PetscBool      match;
885e884886fSBarry Smith 
886e884886fSBarry Smith   PetscFunctionBegin;
887251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mat,MATMFFD,&match);CHKERRQ(ierr);
888ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
889bc0f08ceSBarry Smith 
890e884886fSBarry Smith   *h = ctx->currenth;
891e884886fSBarry Smith   PetscFunctionReturn(0);
892e884886fSBarry Smith }
893e884886fSBarry Smith 
894e884886fSBarry Smith #undef __FUNCT__
895e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunction"
896e884886fSBarry Smith /*@C
897e884886fSBarry Smith    MatMFFDSetFunction - Sets the function used in applying the matrix free.
898e884886fSBarry Smith 
8993f9fe445SBarry Smith    Logically Collective on Mat
900e884886fSBarry Smith 
901e884886fSBarry Smith    Input Parameters:
902e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
903e884886fSBarry Smith .  func - the function to use
904e884886fSBarry Smith -  funcctx - optional function context passed to function
905e884886fSBarry Smith 
906e884886fSBarry Smith    Level: advanced
907e884886fSBarry Smith 
908e884886fSBarry Smith    Notes:
909e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
910e884886fSBarry Smith     matrix inside your compute Jacobian routine
911e884886fSBarry Smith 
9123ec795f1SBarry Smith     If this is not set then it will use the function set with SNESSetFunction() if MatCreateSNESMF() was used.
913e884886fSBarry Smith 
914e884886fSBarry Smith .keywords: SNES, matrix-free, function
915e884886fSBarry Smith 
9161d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD,
9171d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
918e884886fSBarry Smith @*/
9197087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunction(Mat mat,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
920e884886fSBarry Smith {
921bc0f08ceSBarry Smith   PetscErrorCode ierr;
922e884886fSBarry Smith 
923e884886fSBarry Smith   PetscFunctionBegin;
924bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunction_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec),void*),(mat,func,funcctx));CHKERRQ(ierr);
925e884886fSBarry Smith   PetscFunctionReturn(0);
926e884886fSBarry Smith }
927e884886fSBarry Smith 
928e884886fSBarry Smith #undef __FUNCT__
929e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioni"
930e884886fSBarry Smith /*@C
931e884886fSBarry Smith    MatMFFDSetFunctioni - Sets the function for a single component
932e884886fSBarry Smith 
9333f9fe445SBarry Smith    Logically Collective on Mat
934e884886fSBarry Smith 
935e884886fSBarry Smith    Input Parameters:
936e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
937e884886fSBarry Smith -  funci - the function to use
938e884886fSBarry Smith 
939e884886fSBarry Smith    Level: advanced
940e884886fSBarry Smith 
941e884886fSBarry Smith    Notes:
942e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
943e884886fSBarry Smith     matrix inside your compute Jacobian routine
944e884886fSBarry Smith 
945e884886fSBarry Smith 
946e884886fSBarry Smith .keywords: SNES, matrix-free, function
947e884886fSBarry Smith 
9481d0fab5eSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
9491d0fab5eSBarry Smith 
950e884886fSBarry Smith @*/
9517087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
952e884886fSBarry Smith {
9534ac538c5SBarry Smith   PetscErrorCode ierr;
954e884886fSBarry Smith 
955e884886fSBarry Smith   PetscFunctionBegin;
9560700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9574ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioni_C",(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*)),(mat,funci));CHKERRQ(ierr);
958e884886fSBarry Smith   PetscFunctionReturn(0);
959e884886fSBarry Smith }
960e884886fSBarry Smith 
961e884886fSBarry Smith 
962e884886fSBarry Smith #undef __FUNCT__
963e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctioniBase"
964e884886fSBarry Smith /*@C
965e884886fSBarry Smith    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
966e884886fSBarry Smith 
9673f9fe445SBarry Smith    Logically Collective on Mat
968e884886fSBarry Smith 
969e884886fSBarry Smith    Input Parameters:
970e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
971e884886fSBarry Smith -  func - the function to use
972e884886fSBarry Smith 
973e884886fSBarry Smith    Level: advanced
974e884886fSBarry Smith 
975e884886fSBarry Smith    Notes:
976e884886fSBarry Smith     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
977e884886fSBarry Smith     matrix inside your compute Jacobian routine
978e884886fSBarry Smith 
979e884886fSBarry Smith 
980e884886fSBarry Smith .keywords: SNES, matrix-free, function
981e884886fSBarry Smith 
982e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
9831d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory(), SNESetFunction()
984e884886fSBarry Smith @*/
9857087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
986e884886fSBarry Smith {
9874ac538c5SBarry Smith   PetscErrorCode ierr;
988e884886fSBarry Smith 
989e884886fSBarry Smith   PetscFunctionBegin;
9900700a824SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
9914ac538c5SBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctioniBase_C",(Mat,PetscErrorCode (*)(void*,Vec)),(mat,func));CHKERRQ(ierr);
992e884886fSBarry Smith   PetscFunctionReturn(0);
993e884886fSBarry Smith }
994e884886fSBarry Smith 
995e884886fSBarry Smith #undef __FUNCT__
996e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetPeriod"
997e884886fSBarry Smith /*@
998e884886fSBarry Smith    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
999e884886fSBarry Smith 
10003f9fe445SBarry Smith    Logically Collective on Mat
1001e884886fSBarry Smith 
1002e884886fSBarry Smith    Input Parameters:
1003e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateSNESMF()
1004e884886fSBarry Smith -  period - 1 for everytime, 2 for every second etc
1005e884886fSBarry Smith 
1006e884886fSBarry Smith    Options Database Keys:
1007e884886fSBarry Smith +  -mat_mffd_period <period>
1008e884886fSBarry Smith 
1009e884886fSBarry Smith    Level: advanced
1010e884886fSBarry Smith 
1011e884886fSBarry Smith 
1012e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
1013e884886fSBarry Smith 
1014e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(),
10151d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1016e884886fSBarry Smith @*/
10177087cfbeSBarry Smith PetscErrorCode  MatMFFDSetPeriod(Mat mat,PetscInt period)
1018e884886fSBarry Smith {
1019bc0f08ceSBarry Smith   PetscErrorCode ierr;
1020e884886fSBarry Smith 
1021e884886fSBarry Smith   PetscFunctionBegin;
1022bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetPeriod_C",(Mat,PetscInt),(mat,period));CHKERRQ(ierr);
1023e884886fSBarry Smith   PetscFunctionReturn(0);
1024e884886fSBarry Smith }
1025e884886fSBarry Smith 
1026e884886fSBarry Smith #undef __FUNCT__
1027e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFunctionError"
1028e884886fSBarry Smith /*@
1029e884886fSBarry Smith    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
1030e884886fSBarry Smith    matrix-vector products using finite differences.
1031e884886fSBarry Smith 
10323f9fe445SBarry Smith    Logically Collective on Mat
1033e884886fSBarry Smith 
1034e884886fSBarry Smith    Input Parameters:
1035e884886fSBarry Smith +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
1036e884886fSBarry Smith -  error_rel - relative error (should be set to the square root of
1037e884886fSBarry Smith                the relative error in the function evaluations)
1038e884886fSBarry Smith 
1039e884886fSBarry Smith    Options Database Keys:
1040e884886fSBarry Smith +  -mat_mffd_err <error_rel> - Sets error_rel
1041e884886fSBarry Smith 
1042e884886fSBarry Smith    Level: advanced
1043e884886fSBarry Smith 
1044e884886fSBarry Smith    Notes:
1045e884886fSBarry Smith    The default matrix-free matrix-vector product routine computes
1046e884886fSBarry Smith .vb
1047e884886fSBarry Smith      F'(u)*a = [F(u+h*a) - F(u)]/h where
1048e884886fSBarry Smith      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
1049e884886fSBarry Smith        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
1050e884886fSBarry Smith .ve
1051e884886fSBarry Smith 
1052e884886fSBarry Smith .keywords: SNES, matrix-free, parameters
1053e884886fSBarry Smith 
1054e884886fSBarry Smith .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
10551d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1056e884886fSBarry Smith @*/
10577087cfbeSBarry Smith PetscErrorCode  MatMFFDSetFunctionError(Mat mat,PetscReal error)
1058e884886fSBarry Smith {
1059bc0f08ceSBarry Smith   PetscErrorCode ierr;
1060e884886fSBarry Smith 
1061e884886fSBarry Smith   PetscFunctionBegin;
1062bc0f08ceSBarry Smith   ierr = PetscTryMethod(mat,"MatMFFDSetFunctionError_C",(Mat,PetscReal),(mat,error));CHKERRQ(ierr);
1063e884886fSBarry Smith   PetscFunctionReturn(0);
1064e884886fSBarry Smith }
1065e884886fSBarry Smith 
1066e884886fSBarry Smith #undef __FUNCT__
1067e884886fSBarry Smith #define __FUNCT__ "MatMFFDAddNullSpace"
1068e884886fSBarry Smith /*@
1069e884886fSBarry Smith    MatMFFDAddNullSpace - Provides a null space that an operator is
1070e884886fSBarry Smith    supposed to have.  Since roundoff will create a small component in
1071e884886fSBarry Smith    the null space, if you know the null space you may have it
1072e884886fSBarry Smith    automatically removed.
1073e884886fSBarry Smith 
10743f9fe445SBarry Smith    Logically Collective on Mat
1075e884886fSBarry Smith 
1076e884886fSBarry Smith    Input Parameters:
1077e884886fSBarry Smith +  J - the matrix-free matrix context
1078e884886fSBarry Smith -  nullsp - object created with MatNullSpaceCreate()
1079e884886fSBarry Smith 
1080e884886fSBarry Smith    Level: advanced
1081e884886fSBarry Smith 
1082e884886fSBarry Smith .keywords: SNES, matrix-free, null space
1083e884886fSBarry Smith 
1084e884886fSBarry Smith .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD
10851d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDResetHHistory()
1086e884886fSBarry Smith @*/
10877087cfbeSBarry Smith PetscErrorCode  MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp)
1088e884886fSBarry Smith {
1089e884886fSBarry Smith   PetscErrorCode ierr;
1090e884886fSBarry Smith 
1091e884886fSBarry Smith   PetscFunctionBegin;
1092bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDAddNullSpace_C",(Mat,MatNullSpace),(J,nullsp));CHKERRQ(ierr);
1093e884886fSBarry Smith   PetscFunctionReturn(0);
1094e884886fSBarry Smith }
1095e884886fSBarry Smith 
1096e884886fSBarry Smith #undef __FUNCT__
1097e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetHHistory"
1098e884886fSBarry Smith /*@
1099e884886fSBarry Smith    MatMFFDSetHHistory - Sets an array to collect a history of the
1100e884886fSBarry Smith    differencing values (h) computed for the matrix-free product.
1101e884886fSBarry Smith 
11023f9fe445SBarry Smith    Logically Collective on Mat
1103e884886fSBarry Smith 
1104e884886fSBarry Smith    Input Parameters:
1105e884886fSBarry Smith +  J - the matrix-free matrix context
1106e884886fSBarry Smith .  histroy - space to hold the history
1107e884886fSBarry Smith -  nhistory - number of entries in history, if more entries are generated than
1108e884886fSBarry Smith               nhistory, then the later ones are discarded
1109e884886fSBarry Smith 
1110e884886fSBarry Smith    Level: advanced
1111e884886fSBarry Smith 
1112e884886fSBarry Smith    Notes:
1113e884886fSBarry Smith    Use MatMFFDResetHHistory() to reset the history counter and collect
1114e884886fSBarry Smith    a new batch of differencing parameters, h.
1115e884886fSBarry Smith 
1116e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1117e884886fSBarry Smith 
1118e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
11191d0fab5eSBarry Smith           MatMFFDResetHHistory(), MatMFFDSetFunctionError()
1120e884886fSBarry Smith 
1121e884886fSBarry Smith @*/
11227087cfbeSBarry Smith PetscErrorCode  MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
1123e884886fSBarry Smith {
1124e884886fSBarry Smith   MatMFFD        ctx = (MatMFFD)J->data;
1125bc0f08ceSBarry Smith   PetscErrorCode ierr;
1126bc0f08ceSBarry Smith   PetscBool      match;
1127e884886fSBarry Smith 
1128e884886fSBarry Smith   PetscFunctionBegin;
1129251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)J,MATMFFD,&match);CHKERRQ(ierr);
1130ce94432eSBarry Smith   if (!match) SETERRQ(PetscObjectComm((PetscObject)J),PETSC_ERR_ARG_WRONG,"Not a MFFD matrix");
1131e884886fSBarry Smith   ctx->historyh    = history;
1132e884886fSBarry Smith   ctx->maxcurrenth = nhistory;
113375567043SBarry Smith   ctx->currenth    = 0.;
1134e884886fSBarry Smith   PetscFunctionReturn(0);
1135e884886fSBarry Smith }
1136e884886fSBarry Smith 
1137bc0f08ceSBarry Smith 
1138e884886fSBarry Smith #undef __FUNCT__
1139e884886fSBarry Smith #define __FUNCT__ "MatMFFDResetHHistory"
1140e884886fSBarry Smith /*@
1141e884886fSBarry Smith    MatMFFDResetHHistory - Resets the counter to zero to begin
1142e884886fSBarry Smith    collecting a new set of differencing histories.
1143e884886fSBarry Smith 
11443f9fe445SBarry Smith    Logically Collective on Mat
1145e884886fSBarry Smith 
1146e884886fSBarry Smith    Input Parameters:
1147e884886fSBarry Smith .  J - the matrix-free matrix context
1148e884886fSBarry Smith 
1149e884886fSBarry Smith    Level: advanced
1150e884886fSBarry Smith 
1151e884886fSBarry Smith    Notes:
1152e884886fSBarry Smith    Use MatMFFDSetHHistory() to create the original history counter.
1153e884886fSBarry Smith 
1154e884886fSBarry Smith .keywords: SNES, matrix-free, h history, differencing history
1155e884886fSBarry Smith 
1156e884886fSBarry Smith .seealso: MatMFFDGetH(), MatCreateSNESMF(),
11571d0fab5eSBarry Smith           MatMFFDSetHHistory(), MatMFFDSetFunctionError()
1158e884886fSBarry Smith 
1159e884886fSBarry Smith @*/
11607087cfbeSBarry Smith PetscErrorCode  MatMFFDResetHHistory(Mat J)
1161e884886fSBarry Smith {
1162bc0f08ceSBarry Smith   PetscErrorCode ierr;
1163e884886fSBarry Smith 
1164e884886fSBarry Smith   PetscFunctionBegin;
1165bc0f08ceSBarry Smith   ierr = PetscTryMethod(J,"MatMFFDResetHHistory_C",(Mat),(J));CHKERRQ(ierr);
1166e884886fSBarry Smith   PetscFunctionReturn(0);
1167e884886fSBarry Smith }
1168e884886fSBarry Smith 
1169e884886fSBarry Smith 
1170e884886fSBarry Smith #undef __FUNCT__
1171e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetBase"
1172e884886fSBarry Smith /*@
1173e884886fSBarry Smith     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
1174e884886fSBarry Smith         Jacobian are computed
1175e884886fSBarry Smith 
11763f9fe445SBarry Smith     Logically Collective on Mat
1177e884886fSBarry Smith 
1178e884886fSBarry Smith     Input Parameters:
1179e884886fSBarry Smith +   J - the MatMFFD matrix
11803ec795f1SBarry Smith .   U - the vector
1181bcddec3dSBarry Smith -   F - (optional) vector that contains F(u) if it has been already computed
1182e884886fSBarry Smith 
1183e884886fSBarry Smith     Notes: This is rarely used directly
1184e884886fSBarry Smith 
11858af5ae88SBarry Smith     If F is provided then it is not recomputed. Otherwise the function is evaluated at the base
11868af5ae88SBarry Smith     point during the first MatMult() after each call to MatMFFDSetBase().
1187dff2f722SBarry Smith 
1188e884886fSBarry Smith     Level: advanced
1189e884886fSBarry Smith 
1190e884886fSBarry Smith @*/
11917087cfbeSBarry Smith PetscErrorCode  MatMFFDSetBase(Mat J,Vec U,Vec F)
1192e884886fSBarry Smith {
11934ac538c5SBarry Smith   PetscErrorCode ierr;
1194e884886fSBarry Smith 
1195e884886fSBarry Smith   PetscFunctionBegin;
11960700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
11970700a824SBarry Smith   PetscValidHeaderSpecific(U,VEC_CLASSID,2);
11980700a824SBarry Smith   if (F) PetscValidHeaderSpecific(F,VEC_CLASSID,3);
11994ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetBase_C",(Mat,Vec,Vec),(J,U,F));CHKERRQ(ierr);
1200e884886fSBarry Smith   PetscFunctionReturn(0);
1201e884886fSBarry Smith }
1202e884886fSBarry Smith 
1203e884886fSBarry Smith #undef __FUNCT__
1204e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckh"
1205e884886fSBarry Smith /*@C
1206e884886fSBarry Smith     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
1207e884886fSBarry Smith         it to satisfy some criteria
1208e884886fSBarry Smith 
12093f9fe445SBarry Smith     Logically Collective on Mat
1210e884886fSBarry Smith 
1211e884886fSBarry Smith     Input Parameters:
1212e884886fSBarry Smith +   J - the MatMFFD matrix
1213e884886fSBarry Smith .   fun - the function that checks h
1214e884886fSBarry Smith -   ctx - any context needed by the function
1215e884886fSBarry Smith 
1216e884886fSBarry Smith     Options Database Keys:
1217e884886fSBarry Smith .   -mat_mffd_check_positivity
1218e884886fSBarry Smith 
1219e884886fSBarry Smith     Level: advanced
1220e884886fSBarry Smith 
1221e884886fSBarry Smith     Notes: For example, MatMFFDSetCheckPositivity() insures that all entries
1222e884886fSBarry Smith        of U + h*a are non-negative
1223e884886fSBarry Smith 
1224e884886fSBarry Smith .seealso:  MatMFFDSetCheckPositivity()
1225e884886fSBarry Smith @*/
12267087cfbeSBarry Smith PetscErrorCode  MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void *ctx)
1227e884886fSBarry Smith {
12284ac538c5SBarry Smith   PetscErrorCode ierr;
1229e884886fSBarry Smith 
1230e884886fSBarry Smith   PetscFunctionBegin;
12310700a824SBarry Smith   PetscValidHeaderSpecific(J,MAT_CLASSID,1);
12324ac538c5SBarry Smith   ierr = PetscTryMethod(J,"MatMFFDSetCheckh_C",(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*),(J,fun,ctx));CHKERRQ(ierr);
1233e884886fSBarry Smith   PetscFunctionReturn(0);
1234e884886fSBarry Smith }
1235e884886fSBarry Smith 
1236e884886fSBarry Smith #undef __FUNCT__
1237e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetCheckPositivity"
1238e884886fSBarry Smith /*@
1239e884886fSBarry Smith     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
1240e884886fSBarry Smith         zero, decreases h until this is satisfied.
1241e884886fSBarry Smith 
12423f9fe445SBarry Smith     Logically Collective on Vec
1243e884886fSBarry Smith 
1244e884886fSBarry Smith     Input Parameters:
1245e884886fSBarry Smith +   U - base vector that is added to
1246e884886fSBarry Smith .   a - vector that is added
1247e884886fSBarry Smith .   h - scaling factor on a
1248e884886fSBarry Smith -   dummy - context variable (unused)
1249e884886fSBarry Smith 
1250e884886fSBarry Smith     Options Database Keys:
1251e884886fSBarry Smith .   -mat_mffd_check_positivity
1252e884886fSBarry Smith 
1253e884886fSBarry Smith     Level: advanced
1254e884886fSBarry Smith 
1255e884886fSBarry Smith     Notes: This is rarely used directly, rather it is passed as an argument to
1256e884886fSBarry Smith            MatMFFDSetCheckh()
1257e884886fSBarry Smith 
1258e884886fSBarry Smith .seealso:  MatMFFDSetCheckh()
1259e884886fSBarry Smith @*/
12607087cfbeSBarry Smith PetscErrorCode  MatMFFDCheckPositivity(void *dummy,Vec U,Vec a,PetscScalar *h)
1261e884886fSBarry Smith {
1262e884886fSBarry Smith   PetscReal      val, minval;
1263e884886fSBarry Smith   PetscScalar    *u_vec, *a_vec;
1264e884886fSBarry Smith   PetscErrorCode ierr;
1265e884886fSBarry Smith   PetscInt       i,n;
1266e884886fSBarry Smith   MPI_Comm       comm;
1267e884886fSBarry Smith 
1268e884886fSBarry Smith   PetscFunctionBegin;
1269e884886fSBarry Smith   ierr   = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1270e884886fSBarry Smith   ierr   = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1271e884886fSBarry Smith   ierr   = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1272e884886fSBarry Smith   ierr   = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1273e884886fSBarry Smith   minval = PetscAbsScalar(*h*1.01);
1274e884886fSBarry Smith   for (i=0; i<n; i++) {
1275e884886fSBarry Smith     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1276e884886fSBarry Smith       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1277e884886fSBarry Smith       if (val < minval) minval = val;
1278e884886fSBarry Smith     }
1279e884886fSBarry Smith   }
1280e884886fSBarry Smith   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1281e884886fSBarry Smith   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1282d9822059SBarry Smith   ierr = MPI_Allreduce(&minval,&val,1,MPIU_REAL,MPIU_MIN,comm);CHKERRQ(ierr);
1283e884886fSBarry Smith   if (val <= PetscAbsScalar(*h)) {
1284e884886fSBarry Smith     ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr);
1285e884886fSBarry Smith     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1286e884886fSBarry Smith     else                         *h = -0.99*val;
1287e884886fSBarry Smith   }
1288e884886fSBarry Smith   PetscFunctionReturn(0);
1289e884886fSBarry Smith }
1290e884886fSBarry Smith 
1291e884886fSBarry Smith 
1292e884886fSBarry Smith 
1293e884886fSBarry Smith 
1294e884886fSBarry Smith 
1295e884886fSBarry Smith 
1296e884886fSBarry Smith 
1297