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