xref: /petsc/src/snes/mf/snesmfj.c (revision 3eee16b04de4467482e6ce56cbcc0108d577ea39)
1 
2 #ifndef lint
3 static char vcid[] = "$Id: snesmfj.c,v 1.52 1997/07/02 11:21:45 bsmith Exp bsmith $";
4 #endif
5 
6 #include "src/snes/snesimpl.h"   /*I  "snes.h"   I*/
7 #include "pinclude/pviewer.h"
8 #include <math.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" /* ADIC Ignore */
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" /* ADIC Ignore */
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 extern int VecDot_Seq(Vec,Vec,Scalar *);
59 extern int VecNorm_Seq(Vec,NormType,double *);
60 
61 #undef __FUNC__
62 #define __FUNC__ "SNESMatrixFreeMult_Private"
63 /*
64   SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector
65   product, y = F'(u)*a:
66         y = ( F(u + ha) - F(u) ) /h,
67   where F = nonlinear function, as set by SNESSetFunction()
68         u = current iterate
69         h = difference interval
70 */
71 int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y)
72 {
73   MFCtx_Private *ctx;
74   SNES          snes;
75   double        ovalues[3],values[3],norm, sum, umin;
76   Scalar        h, dot, mone = -1.0;
77   Vec           w,U,F;
78   int           ierr, (*eval_fct)(SNES,Vec,Vec);
79   MPI_Comm      comm;
80 
81   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
82 
83   PetscObjectGetComm((PetscObject)mat,&comm);
84   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
85   snes = ctx->snes;
86   w    = ctx->w;
87   umin = ctx->umin;
88 
89   /* We log matrix-free matrix-vector products separately, so that we can
90      separate the performance monitoring from the cases that use conventional
91      storage.  We may eventually modify event logging to associate events
92      with particular objects, hence alleviating the more general problem. */
93 
94   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
95   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
96     eval_fct = SNESComputeFunction;
97     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
98   }
99   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
100     eval_fct = SNESComputeGradient;
101     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
102   }
103   else SETERRQ(1,0,"Invalid method class");
104 
105   /* Determine a "good" step size, h */
106 
107   /*
108     ierr = VecDot(U,a,&dot); CHKERRQ(ierr);
109     ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr);
110     ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr);
111   */
112 
113   /*
114      Call the Seq Vector routines and then do a single reduction
115      to reduce the number of communications required
116   */
117 
118 #if !defined(PETSC_COMPLEX)
119   PLogEventBegin(VEC_Dot,U,a,0,0);
120   ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr);
121   PLogEventEnd(VEC_Dot,U,a,0,0);
122   PLogEventBegin(VEC_Norm,a,0,0,0);
123   ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
124   ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
125   ovalues[2] = ovalues[2]*ovalues[2];
126   MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );
127   dot = values[0]; sum = values[1]; norm = sqrt(values[2]);
128   PLogEventEnd(VEC_Norm,a,0,0,0);
129 #else
130   {
131     Scalar cvalues[3],covalues[3];
132 
133     PLogEventBegin(VEC_Dot,U,a,0,0);
134     ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr);
135     PLogEventEnd(VEC_Dot,U,a,0,0);
136     PLogEventBegin(VEC_Norm,a,0,0,0);
137     ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
138     ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
139     covalues[1] = ovalues[1];
140     covalues[2] = ovalues[2]*ovalues[2];
141     MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm );
142     dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2]));
143     PLogEventEnd(VEC_Norm,a,0,0,0);
144   }
145 #endif
146 
147 
148   /* Safeguard for step sizes too small */
149   if (sum == 0.0) {dot = 1.0; norm = 1.0;}
150 #if defined(PETSC_COMPLEX)
151   else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum;
152   else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum;
153 #else
154   else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
155   else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
156 #endif
157   h = ctx->error_rel*dot/(norm*norm);
158 
159   /* Evaluate function at F(u + ha) */
160   ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr);
161   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
162   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
163   h = 1.0/h;
164   ierr = VecScale(&h,y); CHKERRQ(ierr);
165   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
166 
167   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
168   return 0;
169 }
170 
171 #undef __FUNC__
172 #define __FUNC__ "SNESDefaultMatrixFreeMatCreate"
173 /*@C
174    SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix
175    context for use with a SNES solver.  This matrix can be used as
176    the Jacobian argument for the routine SNESSetJacobian().
177 
178    Input Parameters:
179 .  snes - the SNES context
180 .  x - vector where SNES solution is to be stored.
181 
182    Output Parameter:
183 .  J - the matrix-free matrix
184 
185    Notes:
186    The matrix-free matrix context merely contains the function pointers
187    and work space for performing finite difference approximations of
188    Jacobian-vector products, J(u)*a, via
189 
190 $       J(u)*a = [J(u+h*a) - J(u)]/h where
191 $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
192 $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
193 $   where
194 $        error_rel = square root of relative error in
195 $                    function evaluation
196 $        umin = minimum iterate parameter
197 
198    The user can set these parameters via SNESSetMatrixFreeParameters().
199    See the nonlinear solvers chapter of the users manual for details.
200 
201    The user should call MatDestroy() when finished with the matrix-free
202    matrix context.
203 
204    Options Database Keys:
205 $  -snes_mf_err <error_rel>
206 $  -snes_mf_unim <umin>
207 
208 .keywords: SNES, default, matrix-free, create, matrix
209 
210 .seealso: MatDestroy(), SNESSetMatrixFreeParameters()
211 @*/
212 int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J)
213 {
214   MPI_Comm      comm;
215   MFCtx_Private *mfctx;
216   int           n, nloc, ierr, flg;
217   char          p[64];
218 
219   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
220   PLogObjectMemory(snes,sizeof(MFCtx_Private));
221   mfctx->sp   = 0;
222   mfctx->snes = snes;
223   mfctx->error_rel = 1.e-8; /* assumes double precision */
224   mfctx->umin      = 1.e-6;
225   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr);
226   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr);
227   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
228   PetscStrcpy(p,"-");
229   if (snes->prefix) PetscStrcat(p,snes->prefix);
230   if (flg) {
231     PetscPrintf(snes->comm,"   %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel);
232     PetscPrintf(snes->comm,"   %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin);
233   }
234   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
235   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
236   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
237   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
238   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr);
239   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr);
240   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr);
241   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr);
242   PLogObjectParent(*J,mfctx->w);
243   PLogObjectParent(snes,*J);
244   return 0;
245 }
246 
247 #undef __FUNC__
248 #define __FUNC__ "SNESSetMatrixFreeParameters"
249 /*@
250    SNESSetMatrixFreeParameters - Sets the parameters for the approximation of
251    matrix-vector products using finite differences.
252 
253 $       J(u)*a = [J(u+h*a) - J(u)]/h where
254 $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
255 $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
256 $
257    Input Parameters:
258 .  snes - the SNES context
259 .  error_rel - relative error (should be set to the square root of
260      the relative error in the function evaluations)
261 .  umin - minimum allowable u-value
262 
263 .keywords: SNES, matrix-free, parameters
264 
265 .seealso: SNESDefaultMatrixFreeMatCreate()
266 @*/
267 int SNESSetMatrixFreeParameters(SNES snes,double error,double umin)
268 {
269   MFCtx_Private *ctx;
270   int           ierr;
271   Mat           mat;
272 
273   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
274   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
275   if (ctx) {
276     if (error != PETSC_DEFAULT) ctx->error_rel = error;
277     if (umin != PETSC_DEFAULT)  ctx->umin = umin;
278   }
279   return 0;
280 }
281 
282 #undef __FUNC__
283 #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace"
284 /*@
285    SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that
286    an operator is supposed to have.  Since roundoff will create a
287    small component in the null space, if you know the null space
288    you may have it automatically removed.
289 
290    Input Parameters:
291 .  J - the matrix-free matrix context
292 .  has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants
293 .  n - number of vectors (excluding constant vector) in null space
294 .  vecs - the vectors that span the null space (excluding the constant vector);
295 .         these vectors must be orthonormal
296 
297 .keywords: SNES, matrix-free, null space
298 @*/
299 int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs)
300 {
301   int           ierr;
302   MFCtx_Private *ctx;
303   MPI_Comm      comm;
304 
305   PetscObjectGetComm((PetscObject)J,&comm);
306 
307   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
308   /* no context indicates that it is not the "matrix free" matrix type */
309   if (!ctx) return 0;
310   ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr);
311   return 0;
312 }
313 
314 
315 
316 
317