xref: /petsc/src/mat/impls/mffd/mffd.c (revision e884886f040fd140b63a22c479f657cf835a1983)
1 #define PETSCMAT_DLL
2 
3 #include "include/private/matimpl.h"
4 #include "src/mat/impls/mffd/mffdimpl.h"   /*I  "petscmat.h"   I*/
5 
6 PetscFList MatMFFDPetscFList         = 0;
7 PetscTruth MatMFFDRegisterAllCalled = PETSC_FALSE;
8 
9 PetscCookie PETSCMAT_DLLEXPORT MATMFFD_COOKIE = 0;
10 PetscEvent  MATMFFD_Mult = 0;
11 
12 #undef __FUNCT__
13 #define __FUNCT__ "MatMFFDSetType"
14 /*@C
15     MatMFFDSetType - Sets the method that is used to compute the
16     differencing parameter for finite differene matrix-free formulations.
17 
18     Input Parameters:
19 +   mat - the "matrix-free" matrix created via MatCreateSNESMF(), or MatCreateMFFD()
20           or MatSetType(mat,MATMFFD);
21 -   ftype - the type requested, either MATMFFD_WP or MATMFFD_DS
22 
23     Level: advanced
24 
25     Notes:
26     For example, such routines can compute h for use in
27     Jacobian-vector products of the form
28 
29                         F(x+ha) - F(x)
30           F'(u)a  ~=  ----------------
31                               h
32 
33 .seealso: MatCreateSNESMF(), MatMFFDRegisterDynamic)
34 @*/
35 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetType(Mat mat,MatMFFDType ftype)
36 {
37   PetscErrorCode ierr,(*r)(MatMFFD);
38   MatMFFD        ctx = (MatMFFD)mat->data;
39   PetscTruth     match;
40 
41   PetscFunctionBegin;
42   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
43   PetscValidCharPointer(ftype,2);
44 
45   /* already set, so just return */
46   ierr = PetscTypeCompare((PetscObject)ctx,ftype,&match);CHKERRQ(ierr);
47   if (match) PetscFunctionReturn(0);
48 
49   /* destroy the old one if it exists */
50   if (ctx->ops->destroy) {
51     ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);
52   }
53 
54   ierr =  PetscFListFind(ctx->comm,MatMFFDPetscFList,ftype,(void (**)(void)) &r);CHKERRQ(ierr);
55   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown MatMFFD type %s given",ftype);
56   ierr = (*r)(ctx);CHKERRQ(ierr);
57   ierr = PetscObjectChangeTypeName((PetscObject)ctx,ftype);CHKERRQ(ierr);
58   PetscFunctionReturn(0);
59 }
60 
61 typedef PetscErrorCode (*FCN1)(void*,Vec); /* force argument to next function to not be extern C*/
62 EXTERN_C_BEGIN
63 #undef __FUNCT__
64 #define __FUNCT__ "MatMFFDSetFunctioniBase_FD"
65 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase_FD(Mat mat,FCN1 func)
66 {
67   MatMFFD ctx = (MatMFFD)mat->data;
68 
69   PetscFunctionBegin;
70   ctx->funcisetbase = func;
71   PetscFunctionReturn(0);
72 }
73 EXTERN_C_END
74 
75 typedef PetscErrorCode (*FCN2)(void*,PetscInt,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
76 EXTERN_C_BEGIN
77 #undef __FUNCT__
78 #define __FUNCT__ "MatMFFDSetFunctioni_FD"
79 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni_FD(Mat mat,FCN2 funci)
80 {
81   MatMFFD ctx = (MatMFFD)mat->data;
82 
83   PetscFunctionBegin;
84   ctx->funci = funci;
85   PetscFunctionReturn(0);
86 }
87 EXTERN_C_END
88 
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "MatMFFDRegister"
92 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(MatMFFD))
93 {
94   PetscErrorCode ierr;
95   char           fullname[PETSC_MAX_PATH_LEN];
96 
97   PetscFunctionBegin;
98   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
99   ierr = PetscFListAdd(&MatMFFDPetscFList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
100   PetscFunctionReturn(0);
101 }
102 
103 
104 #undef __FUNCT__
105 #define __FUNCT__ "MatMFFDRegisterDestroy"
106 /*@C
107    MatMFFDRegisterDestroy - Frees the list of MatMFFD methods that were
108    registered by MatMFFDRegisterDynamic).
109 
110    Not Collective
111 
112    Level: developer
113 
114 .keywords: MatMFFD, register, destroy
115 
116 .seealso: MatMFFDRegisterDynamic), MatMFFDRegisterAll()
117 @*/
118 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDRegisterDestroy(void)
119 {
120   PetscErrorCode ierr;
121 
122   PetscFunctionBegin;
123   ierr = PetscFListDestroy(&MatMFFDPetscFList);CHKERRQ(ierr);
124   MatMFFDRegisterAllCalled = PETSC_FALSE;
125   PetscFunctionReturn(0);
126 }
127 
128 /* ----------------------------------------------------------------------------------------*/
129 #undef __FUNCT__
130 #define __FUNCT__ "MatDestroy_MFFD"
131 PetscErrorCode MatDestroy_MFFD(Mat mat)
132 {
133   PetscErrorCode ierr;
134   MatMFFD        ctx = (MatMFFD)mat->data;
135 
136   PetscFunctionBegin;
137   if (ctx->w) {
138     ierr = VecDestroy(ctx->w);CHKERRQ(ierr);
139   }
140   if (ctx->ops->destroy) {ierr = (*ctx->ops->destroy)(ctx);CHKERRQ(ierr);}
141   if (ctx->sp) {ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);}
142   ierr = PetscHeaderDestroy(ctx);CHKERRQ(ierr);
143 
144   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetBase_C","",PETSC_NULL);CHKERRQ(ierr);
145   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C","",PETSC_NULL);CHKERRQ(ierr);
146   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetFunctioni_C","",PETSC_NULL);CHKERRQ(ierr);
147   ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMFFDSetCheckh_C","",PETSC_NULL);CHKERRQ(ierr);
148 
149   PetscFunctionReturn(0);
150 }
151 
152 #undef __FUNCT__
153 #define __FUNCT__ "MatView_MFFD"
154 /*
155    MatMFFDView_MFFD - Views matrix-free parameters.
156 
157 */
158 PetscErrorCode MatView_MFFD(Mat J,PetscViewer viewer)
159 {
160   PetscErrorCode ierr;
161   MatMFFD        ctx = (MatMFFD)J->data;
162   PetscTruth     iascii;
163 
164   PetscFunctionBegin;
165   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
166   if (iascii) {
167      ierr = PetscViewerASCIIPrintf(viewer,"  matrix-free approximation:\n");CHKERRQ(ierr);
168      ierr = PetscViewerASCIIPrintf(viewer,"    err=%G (relative error in function evaluation)\n",ctx->error_rel);CHKERRQ(ierr);
169      if (!ctx->type_name) {
170        ierr = PetscViewerASCIIPrintf(viewer,"    The compute h routine has not yet been set\n");CHKERRQ(ierr);
171      } else {
172        ierr = PetscViewerASCIIPrintf(viewer,"    Using %s compute h routine\n",ctx->type_name);CHKERRQ(ierr);
173      }
174      if (ctx->ops->view) {
175        ierr = (*ctx->ops->view)(ctx,viewer);CHKERRQ(ierr);
176      }
177   } else {
178     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for matrix-free matrix",((PetscObject)viewer)->type_name);
179   }
180   PetscFunctionReturn(0);
181 }
182 
183 #undef __FUNCT__
184 #define __FUNCT__ "MatAssemblyEnd_MFFD"
185 /*
186    MatAssemblyEnd_MFFD - Resets the ctx->ncurrenth to zero. This
187    allows the user to indicate the beginning of a new linear solve by calling
188    MatAssemblyXXX() on the matrix free matrix. This then allows the
189    MatMFFDCreate_WP() to properly compute ||U|| only the first time
190    in the linear solver rather than every time.
191 */
192 PetscErrorCode MatAssemblyEnd_MFFD(Mat J,MatAssemblyType mt)
193 {
194   PetscErrorCode ierr;
195   MatMFFD        j = (MatMFFD)J->data;
196 
197   PetscFunctionBegin;
198   ierr      = MatMFFDResetHHistory(J);CHKERRQ(ierr);
199   j->vshift = 0.0;
200   j->vscale = 1.0;
201   PetscFunctionReturn(0);
202 }
203 
204 #undef __FUNCT__
205 #define __FUNCT__ "MatMult_MFFD"
206 /*
207   MatMult_MFFD - Default matrix-free form for Jacobian-vector product, y = F'(u)*a:
208 
209         y ~= (F(u + ha) - F(u))/h,
210   where F = nonlinear function, as set by SNESSetFunction()
211         u = current iterate
212         h = difference interval
213 */
214 PetscErrorCode MatMult_MFFD(Mat mat,Vec a,Vec y)
215 {
216   MatMFFD        ctx = (MatMFFD)mat->data;
217   PetscScalar    h;
218   Vec            w,U,F;
219   PetscErrorCode ierr;
220   PetscTruth     zeroa;
221 
222   PetscFunctionBegin;
223   /* We log matrix-free matrix-vector products separately, so that we can
224      separate the performance monitoring from the cases that use conventional
225      storage.  We may eventually modify event logging to associate events
226      with particular objects, hence alleviating the more general problem. */
227   ierr = PetscLogEventBegin(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
228 
229   w    = ctx->w;
230   U    = ctx->current_u;
231 
232   /*
233       Compute differencing parameter
234   */
235   if (!ctx->ops->compute) {
236     ierr = MatMFFDSetType(mat,MATMFFD_WP);CHKERRQ(ierr);
237     ierr = MatMFFDSetFromOptions(mat);CHKERRQ(ierr);
238   }
239   ierr = (*ctx->ops->compute)(ctx,U,a,&h,&zeroa);CHKERRQ(ierr);
240   if (zeroa) {
241     ierr = VecSet(y,0.0);CHKERRQ(ierr);
242     PetscFunctionReturn(0);
243   }
244 
245   if (ctx->checkh) {
246     ierr = (*ctx->checkh)(ctx->checkhctx,U,a,&h);CHKERRQ(ierr);
247   }
248 
249   /* keep a record of the current differencing parameter h */
250   ctx->currenth = h;
251 #if defined(PETSC_USE_COMPLEX)
252   ierr = PetscInfo2(mat,"Current differencing parameter: %G + %G i\n",PetscRealPart(h),PetscImaginaryPart(h));CHKERRQ(ierr);
253 #else
254   ierr = PetscInfo1(mat,"Current differencing parameter: %15.12e\n",h);CHKERRQ(ierr);
255 #endif
256   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
257     ctx->historyh[ctx->ncurrenth] = h;
258   }
259   ctx->ncurrenth++;
260 
261   /* w = u + ha */
262   ierr = VecWAXPY(w,h,a,U);CHKERRQ(ierr);
263 
264   F = ctx->funcvec;
265   /* compute func(U) as base for differencing */
266   if (ctx->ncurrenth == 1) {
267     ierr = (*ctx->func)(ctx->funcctx,U,F);CHKERRQ(ierr);
268   }
269   ierr = (*ctx->func)(ctx->funcctx,w,y);CHKERRQ(ierr);
270 
271   ierr = VecAXPY(y,-1.0,F);CHKERRQ(ierr);
272   ierr = VecScale(y,1.0/h);CHKERRQ(ierr);
273 
274   ierr = VecAXPBY(y,ctx->vshift,ctx->vscale,a);CHKERRQ(ierr);
275 
276   if (ctx->sp) {ierr = MatNullSpaceRemove(ctx->sp,y,PETSC_NULL);CHKERRQ(ierr);}
277 
278   ierr = PetscLogEventEnd(MATMFFD_Mult,a,y,0,0);CHKERRQ(ierr);
279   PetscFunctionReturn(0);
280 }
281 
282 #undef __FUNCT__
283 #define __FUNCT__ "MatGetDiagonal_MFFD"
284 /*
285   MatGetDiagonal_MFFD - Gets the diagonal for a matrix free matrix
286 
287         y ~= (F(u + ha) - F(u))/h,
288   where F = nonlinear function, as set by SNESSetFunction()
289         u = current iterate
290         h = difference interval
291 */
292 PetscErrorCode MatGetDiagonal_MFFD(Mat mat,Vec a)
293 {
294   MatMFFD        ctx = (MatMFFD)mat->data;
295   PetscScalar    h,*aa,*ww,v;
296   PetscReal      epsilon = PETSC_SQRT_MACHINE_EPSILON,umin = 100.0*PETSC_SQRT_MACHINE_EPSILON;
297   Vec            w,U;
298   PetscErrorCode ierr;
299   PetscInt       i,rstart,rend;
300 
301   PetscFunctionBegin;
302   if (!ctx->funci) {
303     SETERRQ(PETSC_ERR_ORDER,"Requires calling MatMFFDSetFunctioni() first");
304   }
305 
306   w    = ctx->w;
307   U    = ctx->current_u;
308   ierr = (*ctx->func)(ctx->funcctx,U,a);CHKERRQ(ierr);
309   ierr = (*ctx->funcisetbase)(ctx->funcctx,U);CHKERRQ(ierr);
310   ierr = VecCopy(U,w);CHKERRQ(ierr);
311 
312   ierr = VecGetOwnershipRange(a,&rstart,&rend);CHKERRQ(ierr);
313   ierr = VecGetArray(a,&aa);CHKERRQ(ierr);
314   for (i=rstart; i<rend; i++) {
315     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
316     h  = ww[i-rstart];
317     if (h == 0.0) h = 1.0;
318 #if !defined(PETSC_USE_COMPLEX)
319     if (h < umin && h >= 0.0)      h = umin;
320     else if (h < 0.0 && h > -umin) h = -umin;
321 #else
322     if (PetscAbsScalar(h) < umin && PetscRealPart(h) >= 0.0)     h = umin;
323     else if (PetscRealPart(h) < 0.0 && PetscAbsScalar(h) < umin) h = -umin;
324 #endif
325     h     *= epsilon;
326 
327     ww[i-rstart] += h;
328     ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr);
329     ierr          = (*ctx->funci)(ctx->funcctx,i,w,&v);CHKERRQ(ierr);
330     aa[i-rstart]  = (v - aa[i-rstart])/h;
331 
332     /* possibly shift and scale result */
333     aa[i - rstart] = ctx->vshift + ctx->vscale*aa[i-rstart];
334 
335     ierr = VecGetArray(w,&ww);CHKERRQ(ierr);
336     ww[i-rstart] -= h;
337     ierr = VecRestoreArray(w,&ww);CHKERRQ(ierr);
338   }
339   ierr = VecRestoreArray(a,&aa);CHKERRQ(ierr);
340   PetscFunctionReturn(0);
341 }
342 
343 #undef __FUNCT__
344 #define __FUNCT__ "MatShift_MFFD"
345 PetscErrorCode MatShift_MFFD(Mat Y,PetscScalar a)
346 {
347   MatMFFD shell = (MatMFFD)Y->data;
348   PetscFunctionBegin;
349   shell->vshift += a;
350   PetscFunctionReturn(0);
351 }
352 
353 #undef __FUNCT__
354 #define __FUNCT__ "MatScale_MFFD"
355 PetscErrorCode MatScale_MFFD(Mat Y,PetscScalar a)
356 {
357   MatMFFD shell = (MatMFFD)Y->data;
358   PetscFunctionBegin;
359   shell->vscale *= a;
360   PetscFunctionReturn(0);
361 }
362 
363 EXTERN_C_BEGIN
364 #undef __FUNCT__
365 #define __FUNCT__ "MatMFFDSetBase_FD"
366 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase_FD(Mat J,Vec U)
367 {
368   PetscErrorCode ierr;
369   MatMFFD        ctx = (MatMFFD)J->data;
370 
371   PetscFunctionBegin;
372   ierr = MatMFFDResetHHistory(J);CHKERRQ(ierr);
373   ctx->current_u = U;
374   if (!ctx->w) {
375     ierr = VecDuplicate(ctx->current_u, &ctx->w);CHKERRQ(ierr);
376   }
377   J->assembled = PETSC_TRUE;
378   PetscFunctionReturn(0);
379 }
380 EXTERN_C_END
381 typedef PetscErrorCode (*FCN3)(void*,Vec,Vec,PetscScalar*); /* force argument to next function to not be extern C*/
382 EXTERN_C_BEGIN
383 #undef __FUNCT__
384 #define __FUNCT__ "MatMFFDSetCheckh_FD"
385 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh_FD(Mat J,FCN3 fun,void*ectx)
386 {
387   MatMFFD ctx = (MatMFFD)J->data;
388 
389   PetscFunctionBegin;
390   ctx->checkh    = fun;
391   ctx->checkhctx = ectx;
392   PetscFunctionReturn(0);
393 }
394 EXTERN_C_END
395 
396 #undef __FUNCT__
397 #define __FUNCT__ "MatMFFDSetFromOptions"
398 /*@
399    MatMFFDSetFromOptions - Sets the MatMFFD options from the command line
400    parameter.
401 
402    Collective on Mat
403 
404    Input Parameters:
405 .  mat - the matrix obtained with MatCreateMFFD() or MatCreateSNESMF()
406 
407    Options Database Keys:
408 +  -mat_mffd_type - wp or ds (see MATMFFD_WP or MATMFFD_DS)
409 -  -mat_mffd_err - square root of estimated relative error in function evaluation
410 -  -mat_mffd_period - how often h is recomputed, defaults to 1, everytime
411 
412    Level: advanced
413 
414 .keywords: SNES, matrix-free, parameters
415 
416 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(),
417           MatMFFDResetHHistory(), MatMFFDKSPMonitor()
418 @*/
419 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFromOptions(Mat mat)
420 {
421   MatMFFD        mfctx = (MatMFFD)mat->data;
422   PetscErrorCode ierr;
423   PetscTruth     flg;
424   char           ftype[256];
425 
426   PetscFunctionBegin;
427   if (!MatMFFDRegisterAllCalled) {ierr = MatMFFDRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
428 
429   ierr = PetscOptionsBegin(mfctx->comm,mfctx->prefix,"Set matrix free computation parameters","MatMFFD");CHKERRQ(ierr);
430   ierr = PetscOptionsList("-mat_mffd_type","Matrix free type","MatMFFDSetType",MatSNESMPetscFList,mfctx->type_name,ftype,256,&flg);CHKERRQ(ierr);
431   if (flg) {
432     ierr = MatMFFDSetType(mat,ftype);CHKERRQ(ierr);
433   }
434 
435   ierr = PetscOptionsReal("-mat_mffd_err","set sqrt relative error in function","MatMFFDSetFunctionError",mfctx->error_rel,&mfctx->error_rel,0);CHKERRQ(ierr);
436   ierr = PetscOptionsInt("-mat_mffd_period","how often h is recomputed","MatMFFDSetPeriod",mfctx->recomputeperiod,&mfctx->recomputeperiod,0);CHKERRQ(ierr);
437 
438   ierr = PetscOptionsName("-mat_mffd_check_positivity","Insure that U + h*a is nonnegative","MatMFFDSetCheckh",&flg);CHKERRQ(ierr);
439   if (flg) {
440     ierr = MatMFFDSetCheckh(mat,MatMFFDCheckPositivity,0);CHKERRQ(ierr);
441   }
442   if (mfctx->ops->setfromoptions) {
443     ierr = (*mfctx->ops->setfromoptions)(mfctx);CHKERRQ(ierr);
444   }
445   ierr = PetscOptionsEnd();CHKERRQ(ierr);
446   PetscFunctionReturn(0);
447 }
448 
449 /*MC
450   MATMFFD - MATMFFD = "mffd" - A matrix free matrix type.
451 
452   Level: advanced
453 
454 .seealso: MatCreateMFFD(), MatCreateSNESMF()
455 M*/
456 EXTERN_C_BEGIN
457 #undef __FUNCT__
458 #define __FUNCT__ "MatCreate_MFFD"
459 PetscErrorCode PETSCMAT_DLLEXPORT MatCreate_MFFD(Mat A)
460 {
461   MatMFFD         mfctx;
462   PetscErrorCode  ierr;
463 
464   PetscFunctionBegin;
465 
466   ierr = PetscHeaderCreate(mfctx,_p_MatMFFD,struct _MFOps,MATMFFD_COOKIE,0,"MatMFFD",A->comm,MatDestroy_MFFD,MatView_MFFD);CHKERRQ(ierr);
467   mfctx->sp              = 0;
468   mfctx->error_rel       = PETSC_SQRT_MACHINE_EPSILON;
469   mfctx->recomputeperiod = 1;
470   mfctx->count           = 0;
471   mfctx->currenth        = 0.0;
472   mfctx->historyh        = PETSC_NULL;
473   mfctx->ncurrenth       = 0;
474   mfctx->maxcurrenth     = 0;
475   mfctx->type_name       = 0;
476 
477   mfctx->vshift          = 0.0;
478   mfctx->vscale          = 1.0;
479 
480   /*
481      Create the empty data structure to contain compute-h routines.
482      These will be filled in below from the command line options or
483      a later call with MatMFFDSetType() or if that is not called
484      then it will default in the first use of MatMult_MFFD()
485   */
486   mfctx->ops->compute        = 0;
487   mfctx->ops->destroy        = 0;
488   mfctx->ops->view           = 0;
489   mfctx->ops->setfromoptions = 0;
490   mfctx->hctx                = 0;
491 
492   mfctx->func                = 0;
493   mfctx->funcctx             = 0;
494   mfctx->funcvec             = 0;
495   mfctx->w                   = PETSC_NULL;
496 
497   A->data                = mfctx;
498 
499   A->ops->mult           = MatMult_MFFD;
500   A->ops->destroy        = MatDestroy_MFFD;
501   A->ops->view           = MatView_MFFD;
502   A->ops->assemblyend    = MatAssemblyEnd_MFFD;
503   A->ops->getdiagonal    = MatGetDiagonal_MFFD;
504   A->ops->scale          = MatScale_MFFD;
505   A->ops->shift          = MatShift_MFFD;
506   A->ops->setfromoptions = MatMFFDSetFromOptions;
507   A->assembled = PETSC_TRUE;
508 
509   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetBase_C","MatMFFDSetBase_FD",MatMFFDSetBase_FD);CHKERRQ(ierr);
510   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioniBase_C","MatMFFDSetFunctioniBase_FD",MatMFFDSetFunctioniBase_FD);CHKERRQ(ierr);
511   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetFunctioni_C","MatMFFDSetFunctioni_FD",MatMFFDSetFunctioni_FD);CHKERRQ(ierr);
512   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMFFDSetCheckh_C","MatMFFDSetCheckh_FD",MatMFFDSetCheckh_FD);CHKERRQ(ierr);
513   mfctx->mat = A;
514   ierr = PetscObjectChangeTypeName((PetscObject)A,MATMFFD);CHKERRQ(ierr);
515   PetscFunctionReturn(0);
516 }
517 EXTERN_C_END
518 
519 #undef __FUNCT__
520 #define __FUNCT__ "MatCreateMFFD"
521 /*@
522    MatCreateMFFD - Creates a matrix-free matrix. See also MatCreateSNESMF()
523 
524    Collective on Vec
525 
526    Input Parameters:
527 .  x - vector that defines layout of the vectors and matrices
528 
529    Output Parameter:
530 .  J - the matrix-free matrix
531 
532    Level: advanced
533 
534    Notes:
535    The matrix-free matrix context merely contains the function pointers
536    and work space for performing finite difference approximations of
537    Jacobian-vector products, F'(u)*a,
538 
539    The default code uses the following approach to compute h
540 
541 .vb
542      F'(u)*a = [F(u+h*a) - F(u)]/h where
543      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
544        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
545  where
546      error_rel = square root of relative error in function evaluation
547      umin = minimum iterate parameter
548 .ve
549 
550    The user can set the error_rel via MatMFFDSetFunctionError() and
551    umin via MatMFFDDefaultSetUmin(); see the nonlinear solvers chapter
552    of the users manual for details.
553 
554    The user should call MatDestroy() when finished with the matrix-free
555    matrix context.
556 
557    Options Database Keys:
558 +  -mat_mffd_err <error_rel> - Sets error_rel
559 .  -mat_mffd_unim <umin> - Sets umin (for default PETSc routine that computes h only)
560 .  -mat_mffd_ksp_monitor - KSP monitor routine that prints differencing h
561 -  -mat_mffd_check_positivity
562 
563 .keywords: default, matrix-free, create, matrix
564 
565 .seealso: MatDestroy(), MatMFFDSetFunctionError(), MatMFFDDefaultSetUmin(), MatMFFDSetFunction()
566           MatMFFDSetHHistory(), MatMFFDResetHHistory(), MatCreateSNESMF(),
567           MatMFFDGetH(),MatMFFDKSPMonitor(), MatMFFDRegisterDynamic),, MatMFFDComputeJacobian()
568 
569 @*/
570 PetscErrorCode PETSCMAT_DLLEXPORT MatCreateMFFD(Vec x,Mat *J)
571 {
572   MPI_Comm       comm;
573   PetscErrorCode ierr;
574   PetscInt       n,nloc;
575 
576   PetscFunctionBegin;
577   ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr);
578   ierr = VecGetSize(x,&n);CHKERRQ(ierr);
579   ierr = VecGetLocalSize(x,&nloc);CHKERRQ(ierr);
580   ierr = MatCreate(comm,J);CHKERRQ(ierr);
581   ierr = MatSetSizes(*J,nloc,nloc,n,n);CHKERRQ(ierr);
582   ierr = MatSetType(*J,MATMFFD);CHKERRQ(ierr);
583   PetscFunctionReturn(0);
584 }
585 
586 
587 #undef __FUNCT__
588 #define __FUNCT__ "MatMFFDGetH"
589 /*@
590    MatMFFDGetH - Gets the last value that was used as the differencing
591    parameter.
592 
593    Not Collective
594 
595    Input Parameters:
596 .  mat - the matrix obtained with MatCreateSNESMF()
597 
598    Output Paramter:
599 .  h - the differencing step size
600 
601    Level: advanced
602 
603 .keywords: SNES, matrix-free, parameters
604 
605 .seealso: MatCreateSNESMF(),MatMFFDSetHHistory(), MatCreateMFFD(), MATMFFD
606           MatMFFDResetHHistory(),MatMFFDKSPMonitor()
607 @*/
608 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDGetH(Mat mat,PetscScalar *h)
609 {
610   MatMFFD ctx = (MatMFFD)mat->data;
611 
612   PetscFunctionBegin;
613   *h = ctx->currenth;
614   PetscFunctionReturn(0);
615 }
616 
617 
618 #undef __FUNCT__
619 #define __FUNCT__ "MatMFFDSetFunction"
620 /*@C
621    MatMFFDSetFunction - Sets the function used in applying the matrix free.
622 
623    Collective on Mat
624 
625    Input Parameters:
626 +  mat - the matrix free matrix created via MatCreateSNESMF()
627 .  v   - workspace vector
628 .  func - the function to use
629 -  funcctx - optional function context passed to function
630 
631    Level: advanced
632 
633    Notes:
634     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
635     matrix inside your compute Jacobian routine
636 
637     If this is not set then it will use the function set with SNESSetFunction()
638 
639 .keywords: SNES, matrix-free, function
640 
641 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
642           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
643           MatMFFDKSPMonitor(), SNESetFunction()
644 @*/
645 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunction(Mat mat,Vec v,PetscErrorCode (*func)(void*,Vec,Vec),void *funcctx)
646 {
647   MatMFFD ctx = (MatMFFD)mat->data;
648 
649   PetscFunctionBegin;
650   ctx->func    = func;
651   ctx->funcctx = funcctx;
652   ctx->funcvec = v;
653   PetscFunctionReturn(0);
654 }
655 
656 #undef __FUNCT__
657 #define __FUNCT__ "MatMFFDSetFunctioni"
658 /*@C
659    MatMFFDSetFunctioni - Sets the function for a single component
660 
661    Collective on Mat
662 
663    Input Parameters:
664 +  mat - the matrix free matrix created via MatCreateSNESMF()
665 -  funci - the function to use
666 
667    Level: advanced
668 
669    Notes:
670     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
671     matrix inside your compute Jacobian routine
672 
673 
674 .keywords: SNES, matrix-free, function
675 
676 .seealso: MatCreateSNESMF(),MatMFFDGetH(),
677           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
678           MatMFFDKSPMonitor(), SNESetFunction()
679 @*/
680 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioni(Mat mat,PetscErrorCode (*funci)(void*,PetscInt,Vec,PetscScalar*))
681 {
682   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(void*,PetscInt,Vec,PetscScalar*));
683 
684   PetscFunctionBegin;
685   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
686   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMFFDSetFunctioni_C",(void (**)(void))&f);CHKERRQ(ierr);
687   if (f) {
688     ierr = (*f)(mat,funci);CHKERRQ(ierr);
689   }
690   PetscFunctionReturn(0);
691 }
692 
693 
694 #undef __FUNCT__
695 #define __FUNCT__ "MatMFFDSetFunctioniBase"
696 /*@C
697    MatMFFDSetFunctioniBase - Sets the base vector for a single component function evaluation
698 
699    Collective on Mat
700 
701    Input Parameters:
702 +  mat - the matrix free matrix created via MatCreateSNESMF()
703 -  func - the function to use
704 
705    Level: advanced
706 
707    Notes:
708     If you use this you MUST call MatAssemblyBegin()/MatAssemblyEnd() on the matrix free
709     matrix inside your compute Jacobian routine
710 
711 
712 .keywords: SNES, matrix-free, function
713 
714 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
715           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
716           MatMFFDKSPMonitor(), SNESetFunction()
717 @*/
718 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctioniBase(Mat mat,PetscErrorCode (*func)(void*,Vec))
719 {
720   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(void*,Vec));
721 
722   PetscFunctionBegin;
723   PetscValidHeaderSpecific(mat,MAT_COOKIE,1);
724   ierr = PetscObjectQueryFunction((PetscObject)mat,"MatMFFDSetFunctioniBase_C",(void (**)(void))&f);CHKERRQ(ierr);
725   if (f) {
726     ierr = (*f)(mat,func);CHKERRQ(ierr);
727   }
728   PetscFunctionReturn(0);
729 }
730 
731 
732 #undef __FUNCT__
733 #define __FUNCT__ "MatMFFDSetPeriod"
734 /*@
735    MatMFFDSetPeriod - Sets how often h is recomputed, by default it is everytime
736 
737    Collective on Mat
738 
739    Input Parameters:
740 +  mat - the matrix free matrix created via MatCreateSNESMF()
741 -  period - 1 for everytime, 2 for every second etc
742 
743    Options Database Keys:
744 +  -mat_mffd_period <period>
745 
746    Level: advanced
747 
748 
749 .keywords: SNES, matrix-free, parameters
750 
751 .seealso: MatCreateSNESMF(),MatMFFDGetH(),
752           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
753           MatMFFDKSPMonitor()
754 @*/
755 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetPeriod(Mat mat,PetscInt period)
756 {
757   MatMFFD ctx = (MatMFFD)mat->data;
758 
759   PetscFunctionBegin;
760   ctx->recomputeperiod = period;
761   PetscFunctionReturn(0);
762 }
763 
764 #undef __FUNCT__
765 #define __FUNCT__ "MatMFFDSetFunctionError"
766 /*@
767    MatMFFDSetFunctionError - Sets the error_rel for the approximation of
768    matrix-vector products using finite differences.
769 
770    Collective on Mat
771 
772    Input Parameters:
773 +  mat - the matrix free matrix created via MatCreateMFFD() or MatCreateSNESMF()
774 -  error_rel - relative error (should be set to the square root of
775                the relative error in the function evaluations)
776 
777    Options Database Keys:
778 +  -mat_mffd_err <error_rel> - Sets error_rel
779 
780    Level: advanced
781 
782    Notes:
783    The default matrix-free matrix-vector product routine computes
784 .vb
785      F'(u)*a = [F(u+h*a) - F(u)]/h where
786      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
787        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
788 .ve
789 
790 .keywords: SNES, matrix-free, parameters
791 
792 .seealso: MatCreateSNESMF(),MatMFFDGetH(), MatCreateMFFD(), MATMFFD
793           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
794           MatMFFDKSPMonitor()
795 @*/
796 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetFunctionError(Mat mat,PetscReal error)
797 {
798   MatMFFD ctx = (MatMFFD)mat->data;
799 
800   PetscFunctionBegin;
801   if (error != PETSC_DEFAULT) ctx->error_rel = error;
802   PetscFunctionReturn(0);
803 }
804 
805 #undef __FUNCT__
806 #define __FUNCT__ "MatMFFDAddNullSpace"
807 /*@
808    MatMFFDAddNullSpace - Provides a null space that an operator is
809    supposed to have.  Since roundoff will create a small component in
810    the null space, if you know the null space you may have it
811    automatically removed.
812 
813    Collective on Mat
814 
815    Input Parameters:
816 +  J - the matrix-free matrix context
817 -  nullsp - object created with MatNullSpaceCreate()
818 
819    Level: advanced
820 
821 .keywords: SNES, matrix-free, null space
822 
823 .seealso: MatNullSpaceCreate(), MatMFFDGetH(), MatCreateSNESMF(), MatCreateMFFD(), MATMFFD
824           MatMFFDSetHHistory(), MatMFFDResetHHistory(),
825           MatMFFDKSPMonitor(), MatMFFDErrorRel()
826 @*/
827 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDAddNullSpace(Mat J,MatNullSpace nullsp)
828 {
829   PetscErrorCode ierr;
830   MatMFFD      ctx = (MatMFFD)J->data;
831 
832   PetscFunctionBegin;
833   ierr = PetscObjectReference((PetscObject)nullsp);CHKERRQ(ierr);
834   if (ctx->sp) { ierr = MatNullSpaceDestroy(ctx->sp);CHKERRQ(ierr); }
835   ctx->sp = nullsp;
836   PetscFunctionReturn(0);
837 }
838 
839 #undef __FUNCT__
840 #define __FUNCT__ "MatMFFDSetHHistory"
841 /*@
842    MatMFFDSetHHistory - Sets an array to collect a history of the
843    differencing values (h) computed for the matrix-free product.
844 
845    Collective on Mat
846 
847    Input Parameters:
848 +  J - the matrix-free matrix context
849 .  histroy - space to hold the history
850 -  nhistory - number of entries in history, if more entries are generated than
851               nhistory, then the later ones are discarded
852 
853    Level: advanced
854 
855    Notes:
856    Use MatMFFDResetHHistory() to reset the history counter and collect
857    a new batch of differencing parameters, h.
858 
859 .keywords: SNES, matrix-free, h history, differencing history
860 
861 .seealso: MatMFFDGetH(), MatCreateSNESMF(),
862           MatMFFDResetHHistory(),
863           MatMFFDKSPMonitor(), MatMFFDSetFunctionError()
864 
865 @*/
866 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetHHistory(Mat J,PetscScalar history[],PetscInt nhistory)
867 {
868   MatMFFD ctx = (MatMFFD)J->data;
869 
870   PetscFunctionBegin;
871   ctx->historyh    = history;
872   ctx->maxcurrenth = nhistory;
873   ctx->currenth    = 0;
874   PetscFunctionReturn(0);
875 }
876 
877 #undef __FUNCT__
878 #define __FUNCT__ "MatMFFDResetHHistory"
879 /*@
880    MatMFFDResetHHistory - Resets the counter to zero to begin
881    collecting a new set of differencing histories.
882 
883    Collective on Mat
884 
885    Input Parameters:
886 .  J - the matrix-free matrix context
887 
888    Level: advanced
889 
890    Notes:
891    Use MatMFFDSetHHistory() to create the original history counter.
892 
893 .keywords: SNES, matrix-free, h history, differencing history
894 
895 .seealso: MatMFFDGetH(), MatCreateSNESMF(),
896           MatMFFDSetHHistory(),
897           MatMFFDKSPMonitor(), MatMFFDSetFunctionError()
898 
899 @*/
900 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDResetHHistory(Mat J)
901 {
902   MatMFFD ctx = (MatMFFD)J->data;
903 
904   PetscFunctionBegin;
905   ctx->ncurrenth    = 0;
906   PetscFunctionReturn(0);
907 }
908 
909 
910 #undef __FUNCT__
911 #define __FUNCT__ "MatMFFDSetBase"
912 /*@
913     MatMFFDSetBase - Sets the vector U at which matrix vector products of the
914         Jacobian are computed
915 
916     Collective on Mat
917 
918     Input Parameters:
919 +   J - the MatMFFD matrix
920 -   U - the vector
921 
922     Notes: This is rarely used directly
923 
924     Level: advanced
925 
926 @*/
927 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetBase(Mat J,Vec U)
928 {
929   PetscErrorCode ierr,(*f)(Mat,Vec);
930 
931   PetscFunctionBegin;
932   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
933   PetscValidHeaderSpecific(U,VEC_COOKIE,2);
934   ierr = PetscObjectQueryFunction((PetscObject)J,"MatMFFDSetBase_C",(void (**)(void))&f);CHKERRQ(ierr);
935   if (f) {
936     ierr = (*f)(J,U);CHKERRQ(ierr);
937   }
938   PetscFunctionReturn(0);
939 }
940 
941 #undef __FUNCT__
942 #define __FUNCT__ "MatMFFDSetCheckh"
943 /*@C
944     MatMFFDSetCheckh - Sets a function that checks the computed h and adjusts
945         it to satisfy some criteria
946 
947     Collective on Mat
948 
949     Input Parameters:
950 +   J - the MatMFFD matrix
951 .   fun - the function that checks h
952 -   ctx - any context needed by the function
953 
954     Options Database Keys:
955 .   -mat_mffd_check_positivity
956 
957     Level: advanced
958 
959     Notes: For example, MatMFFDSetCheckPositivity() insures that all entries
960        of U + h*a are non-negative
961 
962 .seealso:  MatMFFDSetCheckPositivity()
963 @*/
964 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDSetCheckh(Mat J,PetscErrorCode (*fun)(void*,Vec,Vec,PetscScalar*),void* ctx)
965 {
966   PetscErrorCode ierr,(*f)(Mat,PetscErrorCode (*)(void*,Vec,Vec,PetscScalar*),void*);
967 
968   PetscFunctionBegin;
969   PetscValidHeaderSpecific(J,MAT_COOKIE,1);
970   ierr = PetscObjectQueryFunction((PetscObject)J,"MatMFFDSetCheckh_C",(void (**)(void))&f);CHKERRQ(ierr);
971   if (f) {
972     ierr = (*f)(J,fun,ctx);CHKERRQ(ierr);
973   }
974   PetscFunctionReturn(0);
975 }
976 
977 #undef __FUNCT__
978 #define __FUNCT__ "MatMFFDSetCheckPositivity"
979 /*@
980     MatMFFDCheckPositivity - Checks that all entries in U + h*a are positive or
981         zero, decreases h until this is satisfied.
982 
983     Collective on Vec
984 
985     Input Parameters:
986 +   U - base vector that is added to
987 .   a - vector that is added
988 .   h - scaling factor on a
989 -   dummy - context variable (unused)
990 
991     Options Database Keys:
992 .   -mat_mffd_check_positivity
993 
994     Level: advanced
995 
996     Notes: This is rarely used directly, rather it is passed as an argument to
997            MatMFFDSetCheckh()
998 
999 .seealso:  MatMFFDSetCheckh()
1000 @*/
1001 PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCheckPositivity(void* dummy,Vec U,Vec a,PetscScalar *h)
1002 {
1003   PetscReal      val, minval;
1004   PetscScalar    *u_vec, *a_vec;
1005   PetscErrorCode ierr;
1006   PetscInt       i,n;
1007   MPI_Comm       comm;
1008 
1009   PetscFunctionBegin;
1010   ierr = PetscObjectGetComm((PetscObject)U,&comm);CHKERRQ(ierr);
1011   ierr = VecGetArray(U,&u_vec);CHKERRQ(ierr);
1012   ierr = VecGetArray(a,&a_vec);CHKERRQ(ierr);
1013   ierr = VecGetLocalSize(U,&n);CHKERRQ(ierr);
1014   minval = PetscAbsScalar(*h*1.01);
1015   for(i=0;i<n;i++) {
1016     if (PetscRealPart(u_vec[i] + *h*a_vec[i]) <= 0.0) {
1017       val = PetscAbsScalar(u_vec[i]/a_vec[i]);
1018       if (val < minval) minval = val;
1019     }
1020   }
1021   ierr = VecRestoreArray(U,&u_vec);CHKERRQ(ierr);
1022   ierr = VecRestoreArray(a,&a_vec);CHKERRQ(ierr);
1023   ierr = PetscGlobalMin(&minval,&val,comm);CHKERRQ(ierr);
1024   if (val <= PetscAbsScalar(*h)) {
1025     ierr = PetscInfo2(U,"Scaling back h from %G to %G\n",PetscRealPart(*h),.99*val);CHKERRQ(ierr);
1026     if (PetscRealPart(*h) > 0.0) *h =  0.99*val;
1027     else                         *h = -0.99*val;
1028   }
1029   PetscFunctionReturn(0);
1030 }
1031 
1032 
1033 
1034 
1035 
1036 
1037 
1038