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