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