xref: /petsc/src/snes/mf/snesmfj.c (revision 8ba3a721c9105e9240dbd90dceddbf8833cabcd1)
1 
2 #ifndef lint
3 static char vcid[] = "$Id: snesmfj.c,v 1.44 1997/01/22 01:42:26 curfman Exp curfman $";
4 #endif
5 
6 #include "draw.h"       /*I  "draw.h"   I*/
7 #include "src/snes/snesimpl.h"   /*I  "snes.h"   I*/
8 #include "pinclude/pviewer.h"
9 
10 typedef struct {  /* default context for matrix-free SNES */
11   SNES        snes;      /* SNES context */
12   Vec         w;         /* work vector */
13   PCNullSpace sp;        /* null space context */
14   double      error_rel; /* square root of relative error in computing function */
15   double      umin;      /* minimum allowable u'a value relative to |u|_1 */
16 } MFCtx_Private;
17 
18 #undef __FUNC__
19 #define __FUNC__ "SNESMatrixFreeDestroy_Private"
20 int SNESMatrixFreeDestroy_Private(PetscObject obj)
21 {
22   int           ierr;
23   Mat           mat = (Mat) obj;
24   MFCtx_Private *ctx;
25 
26   ierr = MatShellGetContext(mat,(void **)&ctx);
27   ierr = VecDestroy(ctx->w); CHKERRQ(ierr);
28   if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp); CHKERRQ(ierr);}
29   PetscFree(ctx);
30   return 0;
31 }
32 
33 #undef __FUNC__
34 #define __FUNC__ "SNESMatrixFreeView_Private"
35 /*
36    SNESMatrixFreeView_Private - Views matrix-free parameters.
37  */
38 int SNESMatrixFreeView_Private(Mat J,Viewer viewer)
39 {
40   int           ierr;
41   MFCtx_Private *ctx;
42   MPI_Comm      comm;
43   FILE          *fd;
44   ViewerType    vtype;
45 
46   PetscObjectGetComm((PetscObject)J,&comm);
47   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
48   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
49   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
50   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
51      PetscFPrintf(comm,fd,"  SNES matrix-free approximation:\n");
52      PetscFPrintf(comm,fd,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);
53      PetscFPrintf(comm,fd,"    umin=%g (minimum iterate parameter)\n",ctx->umin);
54   }
55   return 0;
56 }
57 
58 #undef __FUNC__
59 #define __FUNC__ "SNESMatrixFreeMult_Private"
60 /*
61   SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector
62   product, y = F'(u)*a:
63         y = ( F(u + ha) - F(u) ) /h,
64   where F = nonlinear function, as set by SNESSetFunction()
65         u = current iterate
66         h = difference interval
67 */
68 int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y)
69 {
70   MFCtx_Private *ctx;
71   SNES          snes;
72   double        norm, sum, umin;
73   Scalar        h, dot, mone = -1.0;
74   Vec           w,U,F;
75   int           ierr, (*eval_fct)(SNES,Vec,Vec);
76 
77   MatShellGetContext(mat,(void **)&ctx);
78   snes = ctx->snes;
79   w    = ctx->w;
80   umin = ctx->umin;
81 
82   /* We log matrix-free matrix-vector products separately, so that we can
83      separate the performance monitoring from the cases that use conventional
84      storage.  We may eventually modify event logging to associate events
85      with particular objects, hence alleviating the more general problem. */
86   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
87 
88   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
89   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
90     eval_fct = SNESComputeFunction;
91     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
92   }
93   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
94     eval_fct = SNESComputeGradient;
95     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
96   }
97   else SETERRQ(1,0,"Invalid method class");
98 
99   /* Determine a "good" step size, h */
100   ierr = VecDot(U,a,&dot); CHKERRQ(ierr);
101   ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr);
102   ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr);
103 
104   /* Safeguard for step sizes too small */
105   if (sum == 0.0) {dot = 1.0; norm = 1.0;}
106 #if defined(PETSC_COMPLEX)
107   else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum;
108   else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum;
109 #else
110   else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
111   else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
112 #endif
113   h = ctx->error_rel*dot/(norm*norm);
114 
115   /* Evaluate function at F(u + ha) */
116   ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr);
117   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
118   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
119   h = 1.0/h;
120   ierr = VecScale(&h,y); CHKERRQ(ierr);
121   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
122 
123   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
124   return 0;
125 }
126 
127 #undef __FUNC__
128 #define __FUNC__ "SNESDefaultMatrixFreeMatCreate"
129 /*@C
130    SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix
131    context for use with a SNES solver.  This matrix can be used as
132    the Jacobian argument for the routine SNESSetJacobian().
133 
134    Input Parameters:
135 .  snes - the SNES context
136 .  x - vector where SNES solution is to be stored.
137 
138    Output Parameter:
139 .  J - the matrix-free matrix
140 
141    Notes:
142    The matrix-free matrix context merely contains the function pointers
143    and work space for performing finite difference approximations of
144    Jacobian-vector products, J(u)*a, via
145 
146 $       J(u)*a = [J(u+h*a) - J(u)]/h where
147 $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
148 $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
149 $   where
150 $        error_rel = square root of relative error in
151 $                    function evaluation
152 $        umin = minimum iterate parameter
153 
154    The user can set these parameters via SNESSetMatrixFreeParameters().
155    See the nonlinear solvers chapter of the users manual for details.
156 
157    The user should call MatDestroy() when finished with the matrix-free
158    matrix context.
159 
160    Options Database Keys:
161 $  -snes_mf_err <error_rel>
162 $  -snes_mf_unim <umin>
163 
164 .keywords: SNES, default, matrix-free, create, matrix
165 
166 .seealso: MatDestroy(), SNESSetMatrixFreeParameters()
167 @*/
168 int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J)
169 {
170   MPI_Comm      comm;
171   MFCtx_Private *mfctx;
172   int           n, nloc, ierr, flg;
173   char          p[64];
174 
175   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
176   PLogObjectMemory(snes,sizeof(MFCtx_Private));
177   mfctx->sp   = 0;
178   mfctx->snes = snes;
179   mfctx->error_rel = 1.e-8; /* assumes double precision */
180   mfctx->umin      = 1.e-6;
181   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr);
182   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr);
183   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
184   PetscStrcpy(p,"-");
185   if (snes->prefix) PetscStrcat(p,snes->prefix);
186   if (flg) {
187     PetscPrintf(snes->comm,"   %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel);
188     PetscPrintf(snes->comm,"   %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin);
189   }
190   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
191   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
192   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
193   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
194   ierr = MatCreateShell(comm,nloc,n,n,n,(void*)mfctx,J); CHKERRQ(ierr);
195   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private); CHKERRQ(ierr);
196   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private); CHKERRQ(ierr);
197   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr);
198   PLogObjectParent(*J,mfctx->w);
199   PLogObjectParent(snes,*J);
200   return 0;
201 }
202 
203 #undef __FUNC__
204 #define __FUNC__ "SNESSetMatrixFreeParameters"
205 /*@
206    SNESSetMatrixFreeParameters - Sets the parameters for the approximation of
207    matrix-vector products using finite differences.
208 
209 $       J(u)*a = [J(u+h*a) - J(u)]/h where
210 $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
211 $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
212 $
213    Input Parameters:
214 .  snes - the SNES context
215 .  error_rel - relative error (should be set to the square root of
216      the relative error in the function evaluations)
217 .  umin - minimum allowable u-value
218 
219 .keywords: SNES, matrix-free, parameters
220 
221 .seealso: SNESDefaultMatrixFreeMatCreate()
222 @*/
223 int SNESSetMatrixFreeParameters(SNES snes,double error,double umin)
224 {
225   MFCtx_Private *ctx;
226   int           ierr;
227   Mat           mat;
228 
229   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
230   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
231   if (ctx) {
232     if (error != PETSC_DEFAULT) ctx->error_rel = error;
233     if (umin != PETSC_DEFAULT)  ctx->umin = umin;
234   }
235   return 0;
236 }
237 
238 #undef __FUNC__
239 #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace"
240 /*@
241    SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that
242    an operator is supposed to have.  Since roundoff will create a
243    small component in the null space, if you know the null space
244    you may have it automatically removed.
245 
246    Input Parameters:
247 .  J - the matrix-free matrix context
248 .  has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants
249 .  n - number of vectors (excluding constant vector) in null space
250 .  vecs - the vectors that span the null space (excluding the constant vector);
251 .         these vectors must be orthonormal
252 
253 .keywords: SNES, matrix-free, null space
254 @*/
255 int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs)
256 {
257   int           ierr;
258   MFCtx_Private *ctx;
259   MPI_Comm      comm;
260 
261   PetscObjectGetComm((PetscObject)J,&comm);
262 
263   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
264   /* no context indicates that it is not the "matrix free" matrix type */
265   if (!ctx) return 0;
266   ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr);
267   return 0;
268 }
269 
270 
271 
272 
273