xref: /petsc/src/snes/mf/snesmfj.c (revision ccb6204bdbefb9b674fc9e66648c3cafaef5560d)
181e6777dSBarry Smith 
281e6777dSBarry Smith #ifndef lint
3*ccb6204bSLois Curfman McInnes static char vcid[] = "$Id: snesmfj.c,v 1.42 1997/01/06 20:29:45 balay Exp curfman $";
481e6777dSBarry Smith #endif
581e6777dSBarry Smith 
6b1f0a012SBarry Smith #include "draw.h"       /*I  "draw.h"   I*/
770f55243SBarry Smith #include "src/snes/snesimpl.h"   /*I  "snes.h"   I*/
8eb9086c3SLois Curfman McInnes #include "pinclude/pviewer.h"
981e6777dSBarry Smith 
1050361f65SLois Curfman McInnes typedef struct {  /* default context for matrix-free SNES */
11eb9086c3SLois Curfman McInnes   SNES        snes;      /* SNES context */
12eb9086c3SLois Curfman McInnes   Vec         w;         /* work vector */
13eb9086c3SLois Curfman McInnes   PCNullSpace sp;        /* null space context */
14eb9086c3SLois Curfman McInnes   double      error_rel; /* square root of relative error in computing function */
15639f9d9dSBarry Smith   double      umin;      /* minimum allowable u'a value relative to |u|_1 */
1639e2f89bSBarry Smith } MFCtx_Private;
1739e2f89bSBarry Smith 
185615d1e5SSatish Balay #undef __FUNC__
195615d1e5SSatish Balay #define __FUNC__ "SNESMatrixFreeDestroy_Private"
20fae171e0SBarry Smith int SNESMatrixFreeDestroy_Private(PetscObject obj)
21b9fa9cd0SBarry Smith {
22b9fa9cd0SBarry Smith   int           ierr;
23fae171e0SBarry Smith   Mat           mat = (Mat) obj;
24fae171e0SBarry Smith   MFCtx_Private *ctx;
25fae171e0SBarry Smith 
26fae171e0SBarry Smith   ierr = MatShellGetContext(mat,(void **)&ctx);
27b9fa9cd0SBarry Smith   ierr = VecDestroy(ctx->w); CHKERRQ(ierr);
28b4fd4287SBarry Smith   if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp); CHKERRQ(ierr);}
29fae171e0SBarry Smith   PetscFree(ctx);
30b9fa9cd0SBarry Smith   return 0;
31b9fa9cd0SBarry Smith }
3250361f65SLois Curfman McInnes 
335615d1e5SSatish Balay #undef __FUNC__
345615d1e5SSatish Balay #define __FUNC__ "SNESMatrixFreeView_Private"
3539e2f89bSBarry Smith /*
361d1bbde2SLois Curfman McInnes    SNESMatrixFreeView_Private - Views matrix-free parameters.
3739e2f89bSBarry Smith  */
38eb9086c3SLois Curfman McInnes int SNESMatrixFreeView_Private(Mat J,Viewer viewer)
39eb9086c3SLois Curfman McInnes {
40eb9086c3SLois Curfman McInnes   int           ierr;
41eb9086c3SLois Curfman McInnes   MFCtx_Private *ctx;
42eb9086c3SLois Curfman McInnes   MPI_Comm      comm;
43eb9086c3SLois Curfman McInnes   FILE          *fd;
44eb9086c3SLois Curfman McInnes   ViewerType    vtype;
45eb9086c3SLois Curfman McInnes 
46eb9086c3SLois Curfman McInnes   PetscObjectGetComm((PetscObject)J,&comm);
47eb9086c3SLois Curfman McInnes   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
48eb9086c3SLois Curfman McInnes   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
49eb9086c3SLois Curfman McInnes   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
50eb9086c3SLois Curfman McInnes   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
51eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"  SNES matrix-free approximation:\n");
52eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);
53eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"    umin=%g (minimum iterate parameter)\n",ctx->umin);
54eb9086c3SLois Curfman McInnes   }
55eb9086c3SLois Curfman McInnes   return 0;
56eb9086c3SLois Curfman McInnes }
57eb9086c3SLois Curfman McInnes 
585615d1e5SSatish Balay #undef __FUNC__
595615d1e5SSatish Balay #define __FUNC__ "SNESMatrixFreeMult_Private"
60eb9086c3SLois Curfman McInnes /*
61eb9086c3SLois Curfman McInnes   SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector
62eb9086c3SLois Curfman McInnes   product, y = F'(u)*a:
63eb9086c3SLois Curfman McInnes         y = ( F(u + ha) - F(u) ) /h,
64eb9086c3SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
65eb9086c3SLois Curfman McInnes         u = current iterate
66eb9086c3SLois Curfman McInnes         h = difference interval
67eb9086c3SLois Curfman McInnes */
68eb9086c3SLois Curfman McInnes int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y)
6939e2f89bSBarry Smith {
70fae171e0SBarry Smith   MFCtx_Private *ctx;
71fae171e0SBarry Smith   SNES          snes;
72eb9086c3SLois Curfman McInnes   double        norm, sum, umin;
735334005bSBarry Smith   Scalar        h, dot, mone = -1.0;
74fae171e0SBarry Smith   Vec           w,U,F;
7583e56afdSLois Curfman McInnes   int           ierr, (*eval_fct)(SNES,Vec,Vec);
7639e2f89bSBarry Smith 
77fae171e0SBarry Smith   MatShellGetContext(mat,(void **)&ctx);
78fae171e0SBarry Smith   snes = ctx->snes;
79fae171e0SBarry Smith   w    = ctx->w;
80eb9086c3SLois Curfman McInnes   umin = ctx->umin;
81fae171e0SBarry Smith 
8257c37382SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
8357c37382SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
8457c37382SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
8557c37382SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
86eb9086c3SLois Curfman McInnes   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
8757c37382SLois Curfman McInnes 
8878b31e54SBarry Smith   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
8983e56afdSLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
9083e56afdSLois Curfman McInnes     eval_fct = SNESComputeFunction;
9178b31e54SBarry Smith     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
9283e56afdSLois Curfman McInnes   }
9383e56afdSLois Curfman McInnes   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
9483e56afdSLois Curfman McInnes     eval_fct = SNESComputeGradient;
9583e56afdSLois Curfman McInnes     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
9683e56afdSLois Curfman McInnes   }
97e3372554SBarry Smith   else SETERRQ(1,0,"Invalid method class");
9850361f65SLois Curfman McInnes 
99eb9086c3SLois Curfman McInnes   /* Determine a "good" step size, h */
100eb9086c3SLois Curfman McInnes   ierr = VecDot(U,a,&dot); CHKERRQ(ierr);
101eb9086c3SLois Curfman McInnes   ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr);
102eb9086c3SLois Curfman McInnes   ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr);
103eb9086c3SLois Curfman McInnes 
104eb9086c3SLois Curfman McInnes   /* Safeguard for step sizes too small */
105edd2f0e1SBarry Smith   if (sum == 0.0) {dot = 1.0; norm = 1.0;}
10619a167f6SBarry Smith #if defined(PETSC_COMPLEX)
107eb9086c3SLois Curfman McInnes   else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum;
108eb9086c3SLois Curfman McInnes   else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum;
10919a167f6SBarry Smith #else
110eb9086c3SLois Curfman McInnes   else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
111eb9086c3SLois Curfman McInnes   else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
11219a167f6SBarry Smith #endif
113eb9086c3SLois Curfman McInnes   h = ctx->error_rel*dot/(norm*norm);
11439e2f89bSBarry Smith 
115eb9086c3SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
116eb9086c3SLois Curfman McInnes   ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr);
11783e56afdSLois Curfman McInnes   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
118195ee392SLois Curfman McInnes   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
1195334005bSBarry Smith   h = 1.0/h;
120195ee392SLois Curfman McInnes   ierr = VecScale(&h,y); CHKERRQ(ierr);
121b4fd4287SBarry Smith   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
12257c37382SLois Curfman McInnes 
123eb9086c3SLois Curfman McInnes   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
12439e2f89bSBarry Smith   return 0;
12539e2f89bSBarry Smith }
12683e56afdSLois Curfman McInnes 
1275615d1e5SSatish Balay #undef __FUNC__
1285615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatCreate"
1294b828684SBarry Smith /*@C
1305392566eSBarry Smith    SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix
13150361f65SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
13250361f65SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
1335392566eSBarry Smith 
1345392566eSBarry Smith    Input Parameters:
135eb9086c3SLois Curfman McInnes .  snes - the SNES context
1365392566eSBarry Smith .  x - vector where SNES solution is to be stored.
1375392566eSBarry Smith 
138eb9086c3SLois Curfman McInnes    Output Parameter:
1395392566eSBarry Smith .  J - the matrix-free matrix
1405392566eSBarry Smith 
14165afa06eSLois Curfman McInnes    Notes:
14265afa06eSLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
14365afa06eSLois Curfman McInnes    and work space for performing finite difference approximations of
14465afa06eSLois Curfman McInnes    matrix operations such as matrix-vector products.
14565afa06eSLois Curfman McInnes 
14665afa06eSLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
14765afa06eSLois Curfman McInnes    matrix context.
14865afa06eSLois Curfman McInnes 
14965afa06eSLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
15065afa06eSLois Curfman McInnes 
15165afa06eSLois Curfman McInnes .seealso: MatDestroy()
1525392566eSBarry Smith @*/
1535392566eSBarry Smith int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J)
1545392566eSBarry Smith {
1555392566eSBarry Smith   MPI_Comm      comm;
1565392566eSBarry Smith   MFCtx_Private *mfctx;
157eb9086c3SLois Curfman McInnes   int           n, nloc, ierr, flg;
1585392566eSBarry Smith 
1590452661fSBarry Smith   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
160464493b3SBarry Smith   PLogObjectMemory(snes,sizeof(MFCtx_Private));
161b4fd4287SBarry Smith   mfctx->sp   = 0;
1625392566eSBarry Smith   mfctx->snes = snes;
163eb9086c3SLois Curfman McInnes   mfctx->error_rel = 1.e-8; /* assumes double precision */
164639f9d9dSBarry Smith   mfctx->umin      = 1.e-6;
165eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr);
166eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr);
167639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
168639f9d9dSBarry Smith   if (flg) {
169639f9d9dSBarry Smith     PetscPrintf(snes->comm,"-snes_mf_err <err> set sqrt rel tol in function. Default 1.e-8\n");
170639f9d9dSBarry Smith     PetscPrintf(snes->comm,"-snes_mf_umin <umin> see users manual. Default 1.e-8\n");
171639f9d9dSBarry Smith   }
1725392566eSBarry Smith   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
173195ee392SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
174195ee392SLois Curfman McInnes   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
1757ddc982cSLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
1767ddc982cSLois Curfman McInnes   ierr = MatCreateShell(comm,nloc,n,n,n,(void*)mfctx,J); CHKERRQ(ierr);
1771c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private); CHKERRQ(ierr);
1781c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private); CHKERRQ(ierr);
1791c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr);
180b9fa9cd0SBarry Smith   PLogObjectParent(*J,mfctx->w);
181d370d78aSBarry Smith   PLogObjectParent(snes,*J);
18281e6777dSBarry Smith   return 0;
18381e6777dSBarry Smith }
18481e6777dSBarry Smith 
1855615d1e5SSatish Balay #undef __FUNC__
1865615d1e5SSatish Balay #define __FUNC__ "SNESSetMatrixFreeParameters"
187b4fd4287SBarry Smith /*@
188eb9086c3SLois Curfman McInnes    SNESSetMatrixFreeParameters - Sets the parameters for the approximation of
189eb9086c3SLois Curfman McInnes    matrix-vector products using finite differences.
190eb9086c3SLois Curfman McInnes 
191639f9d9dSBarry Smith $       J(u)*a = [J(u+h*a) - J(u)]/h where
192639f9d9dSBarry Smith $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
193639f9d9dSBarry Smith $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
194639f9d9dSBarry Smith $
195eb9086c3SLois Curfman McInnes    Input Parameters:
196eb9086c3SLois Curfman McInnes .  snes - the SNES context
197*ccb6204bSLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
198*ccb6204bSLois Curfman McInnes      the relative error in the function evaluations)
199eb9086c3SLois Curfman McInnes .  umin - minimum allowable u-value
200eb9086c3SLois Curfman McInnes 
201eb9086c3SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
202eb9086c3SLois Curfman McInnes @*/
203eb9086c3SLois Curfman McInnes int SNESSetMatrixFreeParameters(SNES snes,double error,double umin)
204eb9086c3SLois Curfman McInnes {
205eb9086c3SLois Curfman McInnes   MFCtx_Private *ctx;
206eb9086c3SLois Curfman McInnes   int           ierr;
207eb9086c3SLois Curfman McInnes   Mat           mat;
208eb9086c3SLois Curfman McInnes 
209eb9086c3SLois Curfman McInnes   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
210eb9086c3SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
211eb9086c3SLois Curfman McInnes   if (ctx) {
212eb9086c3SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
213eb9086c3SLois Curfman McInnes     if (umin != PETSC_DEFAULT)  ctx->umin = umin;
214eb9086c3SLois Curfman McInnes   }
215eb9086c3SLois Curfman McInnes   return 0;
216eb9086c3SLois Curfman McInnes }
217eb9086c3SLois Curfman McInnes 
2185615d1e5SSatish Balay #undef __FUNC__
2195615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace"
220eb9086c3SLois Curfman McInnes /*@
22121a45821SLois Curfman McInnes    SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that
222f525115eSLois Curfman McInnes    an operator is supposed to have.  Since roundoff will create a
223b4fd4287SBarry Smith    small component in the null space, if you know the null space
224b4fd4287SBarry Smith    you may have it automatically removed.
225b4fd4287SBarry Smith 
226b4fd4287SBarry Smith    Input Parameters:
22721a45821SLois Curfman McInnes .  J - the matrix-free matrix context
22821a45821SLois Curfman McInnes .  has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants
229b4fd4287SBarry Smith .  n - number of vectors (excluding constant vector) in null space
23021a45821SLois Curfman McInnes .  vecs - the vectors that span the null space (excluding the constant vector);
231b4fd4287SBarry Smith .         these vectors must be orthonormal
232b4fd4287SBarry Smith 
233b4fd4287SBarry Smith .keywords: SNES, matrix-free, null space
234b4fd4287SBarry Smith @*/
235b4fd4287SBarry Smith int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs)
236b4fd4287SBarry Smith {
237b4fd4287SBarry Smith   int           ierr;
238b4fd4287SBarry Smith   MFCtx_Private *ctx;
239b4fd4287SBarry Smith   MPI_Comm      comm;
240b4fd4287SBarry Smith 
241b4fd4287SBarry Smith   PetscObjectGetComm((PetscObject)J,&comm);
242b4fd4287SBarry Smith 
243b4fd4287SBarry Smith   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
244b4fd4287SBarry Smith   /* no context indicates that it is not the "matrix free" matrix type */
245b4fd4287SBarry Smith   if (!ctx) return 0;
246b4fd4287SBarry Smith   ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr);
247b4fd4287SBarry Smith   return 0;
248b4fd4287SBarry Smith }
249b4fd4287SBarry Smith 
2505334005bSBarry Smith 
2515334005bSBarry Smith 
2525334005bSBarry Smith 
253