xref: /petsc/src/snes/mf/snesmfj.c (revision 4fb1253bef8d5f66fbaff399b342ce3a5bb8e686)
1 #ifdef PETSC_RCS_HEADER
2 static char vcid[] = "$Id: snesmfj.c,v 1.70 1998/10/31 23:40:15 bsmith Exp bsmith $";
3 #endif
4 
5 #include "src/snes/snesimpl.h"   /*I  "snes.h"   I*/
6 #include "pinclude/pviewer.h"
7 #include <math.h>
8 
9 typedef struct {  /* default context for matrix-free SNES */
10   SNES        snes;      /* SNES context */
11   Vec         w;         /* work vector */
12   PCNullSpace sp;        /* null space context */
13   double      error_rel; /* square root of relative error in computing function */
14   double      umin;      /* minimum allowable u'a value relative to |u|_1 */
15   Scalar      currenth;  /* last differencing parameter used */
16   Scalar      *historyh; /* history of h */
17   int         ncurrenth,maxcurrenth;
18 } MFCtx_Private;
19 
20 #undef __FUNC__
21 #define __FUNC__ "SNESMatrixFreeDestroy_Private"
22 int SNESMatrixFreeDestroy_Private(Mat mat)
23 {
24   int           ierr;
25   MFCtx_Private *ctx;
26 
27   PetscFunctionBegin;
28   ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr);
29   ierr = VecDestroy(ctx->w); CHKERRQ(ierr);
30   if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);}
31   PetscFree(ctx);
32   PetscFunctionReturn(0);
33 }
34 
35 #undef __FUNC__
36 #define __FUNC__ "SNESMatrixFreeView_Private"
37 /*
38    SNESMatrixFreeView_Private - Views matrix-free parameters.
39 
40 */
41 int SNESMatrixFreeView_Private(Mat J,Viewer viewer)
42 {
43   int           ierr;
44   MFCtx_Private *ctx;
45   MPI_Comm      comm;
46   FILE          *fd;
47   ViewerType    vtype;
48 
49   PetscFunctionBegin;
50   PetscObjectGetComm((PetscObject)J,&comm);
51   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
52   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
53   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
54   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
55      PetscFPrintf(comm,fd,"  SNES matrix-free approximation:\n");
56      PetscFPrintf(comm,fd,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);
57      PetscFPrintf(comm,fd,"    umin=%g (minimum iterate parameter)\n",ctx->umin);
58   } else {
59     SETERRQ(1,1,"Viewer type not supported for this object");
60   }
61   PetscFunctionReturn(0);
62 }
63 
64 extern int VecDot_Seq(Vec,Vec,Scalar *);
65 extern int VecNorm_Seq(Vec,NormType,double *);
66 
67 #undef __FUNC__
68 #define __FUNC__ "SNESMatrixFreeMult_Private"
69 /*
70   SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector
71   product, y = F'(u)*a:
72         y = ( F(u + ha) - F(u) ) /h,
73   where F = nonlinear function, as set by SNESSetFunction()
74         u = current iterate
75         h = difference interval
76 */
77 int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y)
78 {
79   MFCtx_Private *ctx;
80   SNES          snes;
81   double        ovalues[3],norm, sum, umin;
82   Scalar        h, dot, mone = -1.0;
83   Vec           w,U,F;
84   int           ierr, (*eval_fct)(SNES,Vec,Vec);
85   MPI_Comm      comm;
86 #if !defined(USE_PETSC_COMPLEX)
87   double        values[3];
88 #endif
89 
90   PetscFunctionBegin;
91   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
92 
93   PetscObjectGetComm((PetscObject)mat,&comm);
94   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
95   snes = ctx->snes;
96   w    = ctx->w;
97   umin = ctx->umin;
98 
99   /* We log matrix-free matrix-vector products separately, so that we can
100      separate the performance monitoring from the cases that use conventional
101      storage.  We may eventually modify event logging to associate events
102      with particular objects, hence alleviating the more general problem. */
103 
104   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
105   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
106     eval_fct = SNESComputeFunction;
107     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
108   }
109   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
110     eval_fct = SNESComputeGradient;
111     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
112   }
113   else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class");
114 
115   /* Determine a "good" step size, h */
116 
117   /*
118     ierr = VecDot(U,a,&dot); CHKERRQ(ierr);
119     ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr);
120     ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr);
121   */
122 
123   /*
124      Call the Seq Vector routines and then do a single reduction
125      to reduce the number of communications required
126   */
127 
128 #if !defined(USE_PETSC_COMPLEX)
129   PLogEventBegin(VEC_Dot,U,a,0,0);
130   ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr);
131   PLogEventEnd(VEC_Dot,U,a,0,0);
132   PLogEventBegin(VEC_Norm,a,0,0,0);
133   ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
134   ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
135   ovalues[2] = ovalues[2]*ovalues[2];
136   PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
137   ierr = MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );CHKERRQ(ierr);
138   PLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm);
139   dot = values[0]; sum = values[1]; norm = sqrt(values[2]);
140   PLogEventEnd(VEC_Norm,a,0,0,0);
141 #else
142   {
143     Scalar cvalues[3],covalues[3];
144 
145     PLogEventBegin(VEC_Dot,U,a,0,0);
146     ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr);
147     PLogEventEnd(VEC_Dot,U,a,0,0);
148     PLogEventBegin(VEC_Norm,a,0,0,0);
149     ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
150     ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
151     covalues[1] = ovalues[1];
152     covalues[2] = ovalues[2]*ovalues[2];
153     PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
154     ierr = MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
155     PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
156     dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2]));
157     PLogEventEnd(VEC_Norm,a,0,0,0);
158   }
159 #endif
160 
161   /* Safeguard for step sizes too small */
162   if (sum == 0.0) {dot = 1.0; norm = 1.0;}
163 #if defined(USE_PETSC_COMPLEX)
164   else if (PetscAbsScalar(dot) < umin*sum && PetscReal(dot) >= 0.0) dot = umin*sum;
165   else if (PetscAbsScalar(dot) < 0.0 && PetscReal(dot) > -umin*sum) dot = -umin*sum;
166 #else
167   else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
168   else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
169 #endif
170   h = ctx->error_rel*dot/(norm*norm);
171 
172   /* keep a record of the current differencing parameter h */
173   ctx->currenth = h;
174 #if defined(USE_PETSC_COMPLEX)
175   PLogInfo(mat,"Current differencing parameter: %g + %g i\n",PetscReal(h),PetscImaginary(h));
176 #else
177   PLogInfo(mat,"Current differencing parameter: %g\n",h);
178 #endif
179   if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) {
180     ctx->historyh[ctx->ncurrenth++] = h;
181   }
182 
183   /* Evaluate function at F(u + ha) */
184   ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr);
185   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
186 
187   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
188   h = 1.0/h;
189   ierr = VecScale(&h,y); CHKERRQ(ierr);
190   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
191 
192   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
193   PetscFunctionReturn(0);
194 }
195 
196 #undef __FUNC__
197 #define __FUNC__ "SNESDefaultMatrixFreeCreate"
198 /*@C
199    SNESDefaultMatrixFreeCreate - Creates a matrix-free matrix
200    context for use with a SNES solver.  This matrix can be used as
201    the Jacobian argument for the routine SNESSetJacobian().
202 
203    Collective on SNES and Vec
204 
205    Input Parameters:
206 +  snes - the SNES context
207 -  x - vector where SNES solution is to be stored.
208 
209    Output Parameter:
210 .  J - the matrix-free matrix
211 
212    Notes:
213    The matrix-free matrix context merely contains the function pointers
214    and work space for performing finite difference approximations of
215    Jacobian-vector products, J(u)*a, via
216 
217 .vb
218      J(u)*a = [J(u+h*a) - J(u)]/h where
219      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
220        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
221  where
222      error_rel = square root of relative error in function evaluation
223      umin = minimum iterate parameter
224 .ve
225 
226    The user can set these parameters via SNESDefaultMatrixFreeSetParameters().
227    See the nonlinear solvers chapter of the users manual for details.
228 
229    The user should call MatDestroy() when finished with the matrix-free
230    matrix context.
231 
232    Options Database Keys:
233 +  -snes_mf_err <error_rel> - Sets error_rel
234 .  -snes_mf_unim <umin> - Sets umin
235 -  -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h
236 
237 .keywords: SNES, default, matrix-free, create, matrix
238 
239 .seealso: MatDestroy(), SNESDefaultMatrixFreeSetParameters(),
240           SNESDefaultMatrixFreeSetHHistory(), SNESDefaultMatrixFreeResetHHistory(),
241           SNESDefaultMatrixFreeGetH(),SNESDefaultMatrixFreeKSPMonitor()
242 
243 @*/
244 int SNESDefaultMatrixFreeCreate(SNES snes,Vec x, Mat *J)
245 {
246   MPI_Comm      comm;
247   MFCtx_Private *mfctx;
248   int           n, nloc, ierr, flg;
249   char          p[64];
250 
251   PetscFunctionBegin;
252   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
253   PLogObjectMemory(snes,sizeof(MFCtx_Private));
254   mfctx->sp          = 0;
255   mfctx->snes        = snes;
256   mfctx->error_rel   = 1.e-8; /* assumes double precision */
257   mfctx->umin        = 1.e-6;
258   mfctx->currenth    = 0.0;
259   mfctx->historyh    = PETSC_NULL;
260   mfctx->ncurrenth   = 0;
261   mfctx->maxcurrenth = 0;
262 
263   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr);
264   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr);
265   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
266   PetscStrcpy(p,"-");
267   if (snes->prefix) PetscStrcat(p,snes->prefix);
268   if (flg) {
269     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel);
270     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin);
271   }
272   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
273   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
274   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
275   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
276   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr);
277   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr);
278   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr);
279   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr);
280   PLogObjectParent(*J,mfctx->w);
281   PLogObjectParent(snes,*J);
282   PetscFunctionReturn(0);
283 }
284 
285 #undef __FUNC__
286 #define __FUNC__ "SNESDefaultMatrixFreeGetH"
287 /*@
288    SNESDefaultMatrixFreeGetH - Gets the last h that was used as the differencing
289      parameter.
290 
291    Not Collective
292 
293    Input Parameters:
294 .   mat - the matrix obtained with SNESDefaultMatrixFreeCreate()
295 
296    Output Paramter:
297 .  h - the differencing step size
298 
299 .keywords: SNES, matrix-free, parameters
300 
301 .seealso: SNESDefaultMatrixFreeCreate(),SNESDefaultMatrixFreeSetHHistory(),
302           SNESDefaultMatrixFreeResetHHistory(),SNESDefaultMatrixFreeKSPMonitor()
303 @*/
304 int SNESDefaultMatrixFreeGetH(Mat mat,Scalar *h)
305 {
306   MFCtx_Private *ctx;
307   int           ierr;
308 
309   PetscFunctionBegin;
310   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
311   if (ctx) {
312     *h = ctx->currenth;
313   }
314   PetscFunctionReturn(0);
315 }
316 
317 #undef __FUNC__
318 #define __FUNC__ "SNESDefaultMatrixFreeKSPMonitor"
319 /*
320    SNESDefaultMatrixFreeKSPMonitor - A KSP monitor for use with the default PETSc
321       SNES matrix free routines. Prints the h differencing parameter used at each
322       timestep.
323 
324 */
325 int SNESDefaultMatrixFreeKSPMonitor(KSP ksp,int n,double rnorm,void *dummy)
326 {
327   PC            pc;
328   MFCtx_Private *ctx;
329   int           ierr;
330   Mat           mat;
331   MPI_Comm      comm;
332   PetscTruth    nonzeroinitialguess;
333 
334   PetscFunctionBegin;
335   ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr);
336   ierr = KSPGetPC(ksp,&pc); CHKERRQ(ierr);
337   ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr);
338   ierr = PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
339   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
340   if (n > 0 || nonzeroinitialguess) {
341 #if defined(USE_PETSC_COMPLEX)
342     PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm,
343                 PetscReal(ctx->currenth),PetscImaginary(ctx->currenth));
344 #else
345     PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth);
346 #endif
347   } else {
348     PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm);
349   }
350   PetscFunctionReturn(0);
351 }
352 
353 #undef __FUNC__
354 #define __FUNC__ "SNESDefaultMatrixFreeSetParameters"
355 /*@
356    SNESDefaultMatrixFreeSetParameters - Sets the parameters for the approximation of
357    matrix-vector products using finite differences.
358 
359    Collective on Mat
360 
361    Input Parameters:
362 +  mat - the matrix free matrix created via SNESDefaultMatrixFreeCreate()
363 .  error_rel - relative error (should be set to the square root of
364                the relative error in the function evaluations)
365 -  umin - minimum allowable u-value
366 
367    Notes:
368    The default matrix-free matrix-vector product routine computes
369 .vb
370      J(u)*a = [J(u+h*a) - J(u)]/h where
371      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
372        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
373 .ve
374 
375    Options Database Keys:
376 +  -snes_mf_err <error_rel> - Sets error_rel
377 -  -snes_mf_unim <umin> - Sets umin
378 
379 .keywords: SNES, matrix-free, parameters
380 
381 .seealso: SNESDefaultMatrixFreeCreate(),SNESDefaultMatrixFreeGetH(),
382           SNESDefaultMatrixFreeSetHHistory(), SNESDefaultMatrixFreeResetHHistory(),
383           SNESDefaultMatrixFreeKSPMonitor()
384 @*/
385 int SNESDefaultMatrixFreeSetParameters(Mat mat,double error,double umin)
386 {
387   MFCtx_Private *ctx;
388   int           ierr;
389 
390   PetscFunctionBegin;
391   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
392   if (ctx) {
393     if (error != PETSC_DEFAULT) ctx->error_rel = error;
394     if (umin != PETSC_DEFAULT)  ctx->umin = umin;
395   }
396   PetscFunctionReturn(0);
397 }
398 
399 #undef __FUNC__
400 #define __FUNC__ "SNESDefaultMatrixFreeAddNullSpace"
401 /*@
402    SNESDefaultMatrixFreeAddNullSpace - Provides a null space that
403    an operator is supposed to have.  Since roundoff will create a
404    small component in the null space, if you know the null space
405    you may have it automatically removed.
406 
407    Collective on Mat
408 
409    Input Parameters:
410 +  J - the matrix-free matrix context
411 .  has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants
412 .  n - number of vectors (excluding constant vector) in null space
413 -  vecs - the vectors that span the null space (excluding the constant vector);
414           these vectors must be orthonormal
415 
416 .keywords: SNES, matrix-free, null space
417 
418 .seealso: SNESDefaultMatrixFreeGetH(), SNESDefaultMatrixFreeCreate(),
419           SNESDefaultMatrixFreeSetHHistory(), SNESDefaultMatrixFreeResetHHistory(),
420           SNESDefaultMatrixFreeKSPMonitor(), SNESDefaultMatrixFreeSetParameters()
421 
422 @*/
423 int SNESDefaultMatrixFreeAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs)
424 {
425   int           ierr;
426   MFCtx_Private *ctx;
427   MPI_Comm      comm;
428 
429   PetscFunctionBegin;
430   PetscObjectGetComm((PetscObject)J,&comm);
431 
432   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
433   /* no context indicates that it is not the "matrix free" matrix type */
434   if (!ctx) PetscFunctionReturn(0);
435   ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr);
436   PetscFunctionReturn(0);
437 }
438 
439 #undef __FUNC__
440 #define __FUNC__ "SNESDefaultMatrixFreeSetHHistory"
441 /*@
442    SNESDefaultMatrixFreeSetHHistory - Sets an array to collect a history
443       of the differencing values h computed for the matrix free product
444 
445    Collective on Mat
446 
447    Input Parameters:
448 +  J - the matrix-free matrix context
449 .  histroy - space to hold the h history
450 -  nhistory - number of entries in history, if more h are generated than
451               nhistory the later ones are discarded
452 
453    Notes:
454     Use SNESDefaultMatrixFreeResetHHistory() to reset the history counter
455     and collect a new batch of h.
456 
457 .keywords: SNES, matrix-free, h history, differencing history
458 
459 .seealso: SNESDefaultMatrixFreeGetH(), SNESDefaultMatrixFreeCreate(),
460           SNESDefaultMatrixFreeResetHHistory(),
461           SNESDefaultMatrixFreeKSPMonitor(), SNESDefaultMatrixFreeSetParameters()
462 
463 @*/
464 int SNESDefaultMatrixFreeSetHHistory(Mat J,Scalar *history,int nhistory)
465 {
466   int           ierr;
467   MFCtx_Private *ctx;
468 
469   PetscFunctionBegin;
470 
471   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
472   /* no context indicates that it is not the "matrix free" matrix type */
473   if (!ctx) PetscFunctionReturn(0);
474   ctx->historyh    = history;
475   ctx->maxcurrenth = nhistory;
476   ctx->currenth    = 0;
477 
478   PetscFunctionReturn(0);
479 }
480 
481 /*@
482    SNESDefaultMatrixFreeResetHHistory - Resets the counter to zero to begin
483       collecting a new set of differencing histories.
484 
485    Collective on Mat
486 
487    Input Parameters:
488 .  J - the matrix-free matrix context
489 
490    Notes:
491     Use SNESDefaultMatrixFreeSetHHistory() to create the original history counter
492 
493 .keywords: SNES, matrix-free, h history, differencing history
494 
495 .seealso: SNESDefaultMatrixFreeGetH(), SNESDefaultMatrixFreeCreate(),
496           SNESDefaultMatrixFreeSetHHistory(),
497           SNESDefaultMatrixFreeKSPMonitor(), SNESDefaultMatrixFreeSetParameters()
498 
499 @*/
500 int SNESDefaultMatrixFreeResetHHistory(Mat J)
501 {
502   int           ierr;
503   MFCtx_Private *ctx;
504 
505   PetscFunctionBegin;
506 
507   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
508   /* no context indicates that it is not the "matrix free" matrix type */
509   if (!ctx) PetscFunctionReturn(0);
510   ctx->currenth    = 0;
511 
512   PetscFunctionReturn(0);
513 }
514 
515 
516 
517 
518 
519 
520