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