xref: /petsc/src/snes/mf/snesmfj.c (revision e20fef11c1a1457fc77d0de77d14a2dd6b3370e4)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*e20fef11SSatish Balay static char vcid[] = "$Id: snesmfj.c,v 1.65 1998/04/24 04:52:27 curfman Exp balay $";
381e6777dSBarry Smith #endif
481e6777dSBarry Smith 
570f55243SBarry Smith #include "src/snes/snesimpl.h"   /*I  "snes.h"   I*/
6eb9086c3SLois Curfman McInnes #include "pinclude/pviewer.h"
7c481317fSBarry Smith #include <math.h>
881e6777dSBarry Smith 
950361f65SLois Curfman McInnes typedef struct {  /* default context for matrix-free SNES */
10eb9086c3SLois Curfman McInnes   SNES        snes;      /* SNES context */
11eb9086c3SLois Curfman McInnes   Vec         w;         /* work vector */
12eb9086c3SLois Curfman McInnes   PCNullSpace sp;        /* null space context */
13eb9086c3SLois Curfman McInnes   double      error_rel; /* square root of relative error in computing function */
14639f9d9dSBarry Smith   double      umin;      /* minimum allowable u'a value relative to |u|_1 */
1539e2f89bSBarry Smith } MFCtx_Private;
1639e2f89bSBarry Smith 
175615d1e5SSatish Balay #undef __FUNC__
18d4bb536fSBarry Smith #define __FUNC__ "SNESMatrixFreeDestroy_Private"
19e1311b90SBarry Smith int SNESMatrixFreeDestroy_Private(Mat mat)
20b9fa9cd0SBarry Smith {
21b9fa9cd0SBarry Smith   int           ierr;
22fae171e0SBarry Smith   MFCtx_Private *ctx;
23fae171e0SBarry Smith 
243a40ed3dSBarry Smith   PetscFunctionBegin;
25fae171e0SBarry Smith   ierr = MatShellGetContext(mat,(void **)&ctx);
26b9fa9cd0SBarry Smith   ierr = VecDestroy(ctx->w); CHKERRQ(ierr);
27b4fd4287SBarry Smith   if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp); CHKERRQ(ierr);}
28fae171e0SBarry Smith   PetscFree(ctx);
293a40ed3dSBarry Smith   PetscFunctionReturn(0);
30b9fa9cd0SBarry Smith }
3150361f65SLois Curfman McInnes 
325615d1e5SSatish Balay #undef __FUNC__
33d4bb536fSBarry Smith #define __FUNC__ "SNESMatrixFreeView_Private"
3439e2f89bSBarry Smith /*
351d1bbde2SLois Curfman McInnes    SNESMatrixFreeView_Private - Views matrix-free parameters.
3639e2f89bSBarry Smith  */
37eb9086c3SLois Curfman McInnes int SNESMatrixFreeView_Private(Mat J,Viewer viewer)
38eb9086c3SLois Curfman McInnes {
39eb9086c3SLois Curfman McInnes   int           ierr;
40eb9086c3SLois Curfman McInnes   MFCtx_Private *ctx;
41eb9086c3SLois Curfman McInnes   MPI_Comm      comm;
42eb9086c3SLois Curfman McInnes   FILE          *fd;
43eb9086c3SLois Curfman McInnes   ViewerType    vtype;
44eb9086c3SLois Curfman McInnes 
453a40ed3dSBarry Smith   PetscFunctionBegin;
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);
545cd90555SBarry Smith   } else {
555cd90555SBarry Smith     SETERRQ(1,1,"Viewer type not supported for this object");
56eb9086c3SLois Curfman McInnes   }
573a40ed3dSBarry Smith   PetscFunctionReturn(0);
58eb9086c3SLois Curfman McInnes }
59eb9086c3SLois Curfman McInnes 
60e9cf3c21SBarry Smith extern int VecDot_Seq(Vec,Vec,Scalar *);
61c481317fSBarry Smith extern int VecNorm_Seq(Vec,NormType,double *);
62c481317fSBarry Smith 
635615d1e5SSatish Balay #undef __FUNC__
645615d1e5SSatish Balay #define __FUNC__ "SNESMatrixFreeMult_Private"
65eb9086c3SLois Curfman McInnes /*
66eb9086c3SLois Curfman McInnes   SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector
67eb9086c3SLois Curfman McInnes   product, y = F'(u)*a:
68eb9086c3SLois Curfman McInnes         y = ( F(u + ha) - F(u) ) /h,
69eb9086c3SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
70eb9086c3SLois Curfman McInnes         u = current iterate
71eb9086c3SLois Curfman McInnes         h = difference interval
72eb9086c3SLois Curfman McInnes */
73eb9086c3SLois Curfman McInnes int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y)
7439e2f89bSBarry Smith {
75fae171e0SBarry Smith   MFCtx_Private *ctx;
76fae171e0SBarry Smith   SNES          snes;
77c31ba22aSSatish Balay   double        ovalues[3],norm, sum, umin;
785334005bSBarry Smith   Scalar        h, dot, mone = -1.0;
79fae171e0SBarry Smith   Vec           w,U,F;
8083e56afdSLois Curfman McInnes   int           ierr, (*eval_fct)(SNES,Vec,Vec);
81c481317fSBarry Smith   MPI_Comm      comm;
82c31ba22aSSatish Balay #if !defined(USE_PETSC_COMPLEX)
83c31ba22aSSatish Balay   double        values[3];
84c31ba22aSSatish Balay #endif
8539e2f89bSBarry Smith 
863a40ed3dSBarry Smith   PetscFunctionBegin;
8756cd22aeSBarry Smith   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
8856cd22aeSBarry Smith 
89c481317fSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
906e38a044SBarry Smith   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
91fae171e0SBarry Smith   snes = ctx->snes;
92fae171e0SBarry Smith   w    = ctx->w;
93eb9086c3SLois Curfman McInnes   umin = ctx->umin;
94fae171e0SBarry Smith 
9557c37382SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
9657c37382SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
9757c37382SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
9857c37382SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
9957c37382SLois Curfman McInnes 
10078b31e54SBarry Smith   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
10183e56afdSLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
10283e56afdSLois Curfman McInnes     eval_fct = SNESComputeFunction;
10378b31e54SBarry Smith     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
10483e56afdSLois Curfman McInnes   }
10583e56afdSLois Curfman McInnes   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
10683e56afdSLois Curfman McInnes     eval_fct = SNESComputeGradient;
10783e56afdSLois Curfman McInnes     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
10883e56afdSLois Curfman McInnes   }
109a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class");
11050361f65SLois Curfman McInnes 
111eb9086c3SLois Curfman McInnes   /* Determine a "good" step size, h */
112c481317fSBarry Smith 
113c481317fSBarry Smith   /*
114eb9086c3SLois Curfman McInnes     ierr = VecDot(U,a,&dot); CHKERRQ(ierr);
115eb9086c3SLois Curfman McInnes     ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr);
116eb9086c3SLois Curfman McInnes     ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr);
117c481317fSBarry Smith   */
118c481317fSBarry Smith 
119c481317fSBarry Smith   /*
12086f5173aSBarry Smith      Call the Seq Vector routines and then do a single reduction
12186f5173aSBarry Smith      to reduce the number of communications required
122c481317fSBarry Smith   */
123c481317fSBarry Smith 
1243a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
125c481317fSBarry Smith   PLogEventBegin(VEC_Dot,U,a,0,0);
126c481317fSBarry Smith   ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr);
127c481317fSBarry Smith   PLogEventEnd(VEC_Dot,U,a,0,0);
128c481317fSBarry Smith   PLogEventBegin(VEC_Norm,a,0,0,0);
129c481317fSBarry Smith   ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
130c481317fSBarry Smith   ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
131c481317fSBarry Smith   ovalues[2] = ovalues[2]*ovalues[2];
132005c665bSBarry Smith   PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
133ca161407SBarry Smith   ierr = MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );CHKERRQ(ierr);
134005c665bSBarry Smith   PLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm);
135c481317fSBarry Smith   dot = values[0]; sum = values[1]; norm = sqrt(values[2]);
136c481317fSBarry Smith   PLogEventEnd(VEC_Norm,a,0,0,0);
13786f5173aSBarry Smith #else
13886f5173aSBarry Smith   {
13986f5173aSBarry Smith     Scalar cvalues[3],covalues[3];
14086f5173aSBarry Smith 
14186f5173aSBarry Smith     PLogEventBegin(VEC_Dot,U,a,0,0);
14286f5173aSBarry Smith     ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr);
14386f5173aSBarry Smith     PLogEventEnd(VEC_Dot,U,a,0,0);
14486f5173aSBarry Smith     PLogEventBegin(VEC_Norm,a,0,0,0);
14586f5173aSBarry Smith     ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
14686f5173aSBarry Smith     ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
14786f5173aSBarry Smith     covalues[1] = ovalues[1];
14886f5173aSBarry Smith     covalues[2] = ovalues[2]*ovalues[2];
149005c665bSBarry Smith     PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
150ca161407SBarry Smith     ierr = MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
151005c665bSBarry Smith     PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
15286f5173aSBarry Smith     dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2]));
15386f5173aSBarry Smith     PLogEventEnd(VEC_Norm,a,0,0,0);
15486f5173aSBarry Smith   }
15586f5173aSBarry Smith #endif
156c481317fSBarry Smith 
157eb9086c3SLois Curfman McInnes 
158eb9086c3SLois Curfman McInnes   /* Safeguard for step sizes too small */
159edd2f0e1SBarry Smith   if (sum == 0.0) {dot = 1.0; norm = 1.0;}
1603a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
161*e20fef11SSatish Balay   else if (PetscAbsScalar(dot) < umin*sum && PetscReal(dot) >= 0.0) dot = umin*sum;
162*e20fef11SSatish Balay   else if (PetscAbsScalar(dot) < 0.0 && PetscReal(dot) > -umin*sum) dot = -umin*sum;
16319a167f6SBarry Smith #else
164eb9086c3SLois Curfman McInnes   else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
165eb9086c3SLois Curfman McInnes   else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
16619a167f6SBarry Smith #endif
167eb9086c3SLois Curfman McInnes   h = ctx->error_rel*dot/(norm*norm);
16839e2f89bSBarry Smith 
169eb9086c3SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
170eb9086c3SLois Curfman McInnes   ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr);
17183e56afdSLois Curfman McInnes   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
172195ee392SLois Curfman McInnes   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
1735334005bSBarry Smith   h = 1.0/h;
174195ee392SLois Curfman McInnes   ierr = VecScale(&h,y); CHKERRQ(ierr);
175b4fd4287SBarry Smith   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
17657c37382SLois Curfman McInnes 
177eb9086c3SLois Curfman McInnes   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
1783a40ed3dSBarry Smith   PetscFunctionReturn(0);
17939e2f89bSBarry Smith }
18083e56afdSLois Curfman McInnes 
1815615d1e5SSatish Balay #undef __FUNC__
1825615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatCreate"
1834b828684SBarry Smith /*@C
1845392566eSBarry Smith    SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix
18550361f65SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
18650361f65SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
1875392566eSBarry Smith 
188c7afd0dbSLois Curfman McInnes    Collective on SNES and Vec
189c7afd0dbSLois Curfman McInnes 
1905392566eSBarry Smith    Input Parameters:
191c7afd0dbSLois Curfman McInnes +  snes - the SNES context
192c7afd0dbSLois Curfman McInnes -  x - vector where SNES solution is to be stored.
1935392566eSBarry Smith 
194eb9086c3SLois Curfman McInnes    Output Parameter:
1955392566eSBarry Smith .  J - the matrix-free matrix
1965392566eSBarry Smith 
19765afa06eSLois Curfman McInnes    Notes:
19865afa06eSLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
19965afa06eSLois Curfman McInnes    and work space for performing finite difference approximations of
2003d5db8e1SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
2013d5db8e1SLois Curfman McInnes 
202c7afd0dbSLois Curfman McInnes .vb
203c7afd0dbSLois Curfman McInnes      J(u)*a = [J(u+h*a) - J(u)]/h where
204c7afd0dbSLois Curfman McInnes      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
205c7afd0dbSLois Curfman McInnes        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
206c7afd0dbSLois Curfman McInnes  where
207c7afd0dbSLois Curfman McInnes      error_rel = square root of relative error in function evaluation
208c7afd0dbSLois Curfman McInnes      umin = minimum iterate parameter
209c7afd0dbSLois Curfman McInnes .ve
2103d5db8e1SLois Curfman McInnes 
2113d5db8e1SLois Curfman McInnes    The user can set these parameters via SNESSetMatrixFreeParameters().
2123d5db8e1SLois Curfman McInnes    See the nonlinear solvers chapter of the users manual for details.
21365afa06eSLois Curfman McInnes 
21465afa06eSLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
21565afa06eSLois Curfman McInnes    matrix context.
21665afa06eSLois Curfman McInnes 
2173d5db8e1SLois Curfman McInnes    Options Database Keys:
218c7afd0dbSLois Curfman McInnes +  -snes_mf_err <error_rel> - Sets error_rel
219c7afd0dbSLois Curfman McInnes -  -snes_mf_unim <umin> - Sets umin
2203d5db8e1SLois Curfman McInnes 
22165afa06eSLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
22265afa06eSLois Curfman McInnes 
2233d5db8e1SLois Curfman McInnes .seealso: MatDestroy(), SNESSetMatrixFreeParameters()
2245392566eSBarry Smith @*/
2255392566eSBarry Smith int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J)
2265392566eSBarry Smith {
2275392566eSBarry Smith   MPI_Comm      comm;
2285392566eSBarry Smith   MFCtx_Private *mfctx;
229eb9086c3SLois Curfman McInnes   int           n, nloc, ierr, flg;
2303e966953SLois Curfman McInnes   char          p[64];
2315392566eSBarry Smith 
2323a40ed3dSBarry Smith   PetscFunctionBegin;
2330452661fSBarry Smith   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
234464493b3SBarry Smith   PLogObjectMemory(snes,sizeof(MFCtx_Private));
235b4fd4287SBarry Smith   mfctx->sp   = 0;
2365392566eSBarry Smith   mfctx->snes = snes;
237eb9086c3SLois Curfman McInnes   mfctx->error_rel = 1.e-8; /* assumes double precision */
238639f9d9dSBarry Smith   mfctx->umin      = 1.e-6;
239eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr);
240eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr);
241639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
2423e966953SLois Curfman McInnes   PetscStrcpy(p,"-");
2433e966953SLois Curfman McInnes   if (snes->prefix) PetscStrcat(p,snes->prefix);
244639f9d9dSBarry Smith   if (flg) {
24576be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel);
24676be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin);
247639f9d9dSBarry Smith   }
2485392566eSBarry Smith   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
249195ee392SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
250195ee392SLois Curfman McInnes   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
2517ddc982cSLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
252029af93fSBarry Smith   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr);
2531c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr);
2541c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr);
2551c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr);
256b9fa9cd0SBarry Smith   PLogObjectParent(*J,mfctx->w);
257d370d78aSBarry Smith   PLogObjectParent(snes,*J);
2583a40ed3dSBarry Smith   PetscFunctionReturn(0);
25981e6777dSBarry Smith }
26081e6777dSBarry Smith 
2615615d1e5SSatish Balay #undef __FUNC__
2625615d1e5SSatish Balay #define __FUNC__ "SNESSetMatrixFreeParameters"
263b4fd4287SBarry Smith /*@
264eb9086c3SLois Curfman McInnes    SNESSetMatrixFreeParameters - Sets the parameters for the approximation of
265eb9086c3SLois Curfman McInnes    matrix-vector products using finite differences.
266eb9086c3SLois Curfman McInnes 
267c7afd0dbSLois Curfman McInnes    Collective on SNES
268c7afd0dbSLois Curfman McInnes 
269eb9086c3SLois Curfman McInnes    Input Parameters:
270c7afd0dbSLois Curfman McInnes +  snes - the SNES context
271ccb6204bSLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
272ccb6204bSLois Curfman McInnes                the relative error in the function evaluations)
273c7afd0dbSLois Curfman McInnes -  umin - minimum allowable u-value
274eb9086c3SLois Curfman McInnes 
275c7afd0dbSLois Curfman McInnes    Notes:
276c7afd0dbSLois Curfman McInnes    The default matrix-free matrix-vector product routine computes
277c7afd0dbSLois Curfman McInnes .vb
278c7afd0dbSLois Curfman McInnes      J(u)*a = [J(u+h*a) - J(u)]/h where
279c7afd0dbSLois Curfman McInnes      h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
280c7afd0dbSLois Curfman McInnes        = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
281c7afd0dbSLois Curfman McInnes .ve
282c7afd0dbSLois Curfman McInnes 
283c7afd0dbSLois Curfman McInnes    Options Database Keys:
284c7afd0dbSLois Curfman McInnes +  -snes_mf_err <error_rel> - Sets error_rel
285c7afd0dbSLois Curfman McInnes -  -snes_mf_unim <umin> - Sets umin
286fee21e36SBarry Smith 
287eb9086c3SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
2883d5db8e1SLois Curfman McInnes 
2893d5db8e1SLois Curfman McInnes .seealso: SNESDefaultMatrixFreeMatCreate()
290eb9086c3SLois Curfman McInnes @*/
291eb9086c3SLois Curfman McInnes int SNESSetMatrixFreeParameters(SNES snes,double error,double umin)
292eb9086c3SLois Curfman McInnes {
293eb9086c3SLois Curfman McInnes   MFCtx_Private *ctx;
294eb9086c3SLois Curfman McInnes   int           ierr;
295eb9086c3SLois Curfman McInnes   Mat           mat;
296eb9086c3SLois Curfman McInnes 
2973a40ed3dSBarry Smith   PetscFunctionBegin;
298eb9086c3SLois Curfman McInnes   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
299eb9086c3SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
300eb9086c3SLois Curfman McInnes   if (ctx) {
301eb9086c3SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
302eb9086c3SLois Curfman McInnes     if (umin != PETSC_DEFAULT)  ctx->umin = umin;
303eb9086c3SLois Curfman McInnes   }
3043a40ed3dSBarry Smith   PetscFunctionReturn(0);
305eb9086c3SLois Curfman McInnes }
306eb9086c3SLois Curfman McInnes 
3075615d1e5SSatish Balay #undef __FUNC__
3085615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace"
309eb9086c3SLois Curfman McInnes /*@
31021a45821SLois Curfman McInnes    SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that
311f525115eSLois Curfman McInnes    an operator is supposed to have.  Since roundoff will create a
312b4fd4287SBarry Smith    small component in the null space, if you know the null space
313b4fd4287SBarry Smith    you may have it automatically removed.
314b4fd4287SBarry Smith 
315c7afd0dbSLois Curfman McInnes    Collective on Mat
316c7afd0dbSLois Curfman McInnes 
317b4fd4287SBarry Smith    Input Parameters:
318c7afd0dbSLois Curfman McInnes +  J - the matrix-free matrix context
31921a45821SLois Curfman McInnes .  has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants
320b4fd4287SBarry Smith .  n - number of vectors (excluding constant vector) in null space
321c7afd0dbSLois Curfman McInnes -  vecs - the vectors that span the null space (excluding the constant vector);
322c7afd0dbSLois Curfman McInnes           these vectors must be orthonormal
323fee21e36SBarry Smith 
324b4fd4287SBarry Smith .keywords: SNES, matrix-free, null space
325b4fd4287SBarry Smith @*/
326b4fd4287SBarry Smith int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs)
327b4fd4287SBarry Smith {
328b4fd4287SBarry Smith   int           ierr;
329b4fd4287SBarry Smith   MFCtx_Private *ctx;
330b4fd4287SBarry Smith   MPI_Comm      comm;
331b4fd4287SBarry Smith 
3323a40ed3dSBarry Smith   PetscFunctionBegin;
333b4fd4287SBarry Smith   PetscObjectGetComm((PetscObject)J,&comm);
334b4fd4287SBarry Smith 
335b4fd4287SBarry Smith   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
336b4fd4287SBarry Smith   /* no context indicates that it is not the "matrix free" matrix type */
3373a40ed3dSBarry Smith   if (!ctx) PetscFunctionReturn(0);
338b4fd4287SBarry Smith   ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr);
3393a40ed3dSBarry Smith   PetscFunctionReturn(0);
340b4fd4287SBarry Smith }
341b4fd4287SBarry Smith 
3425334005bSBarry Smith 
3435334005bSBarry Smith 
3445334005bSBarry Smith 
345