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