xref: /petsc/src/snes/mf/snesmfj.c (revision c31ba22a230fe2c23596c66917081dfec3ccf152)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*c31ba22aSSatish Balay static char vcid[] = "$Id: snesmfj.c,v 1.60 1998/01/14 02:44:45 bsmith 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"
19fae171e0SBarry Smith int SNESMatrixFreeDestroy_Private(PetscObject obj)
20b9fa9cd0SBarry Smith {
21b9fa9cd0SBarry Smith   int           ierr;
22fae171e0SBarry Smith   Mat           mat = (Mat) obj;
23fae171e0SBarry Smith   MFCtx_Private *ctx;
24fae171e0SBarry Smith 
253a40ed3dSBarry Smith   PetscFunctionBegin;
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);
303a40ed3dSBarry Smith   PetscFunctionReturn(0);
31b9fa9cd0SBarry Smith }
3250361f65SLois Curfman McInnes 
335615d1e5SSatish Balay #undef __FUNC__
34d4bb536fSBarry Smith #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 
463a40ed3dSBarry Smith   PetscFunctionBegin;
47eb9086c3SLois Curfman McInnes   PetscObjectGetComm((PetscObject)J,&comm);
48eb9086c3SLois Curfman McInnes   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
49eb9086c3SLois Curfman McInnes   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
50eb9086c3SLois Curfman McInnes   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
51eb9086c3SLois Curfman McInnes   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
52eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"  SNES matrix-free approximation:\n");
53eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);
54eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"    umin=%g (minimum iterate parameter)\n",ctx->umin);
55eb9086c3SLois Curfman McInnes   }
563a40ed3dSBarry Smith   PetscFunctionReturn(0);
57eb9086c3SLois Curfman McInnes }
58eb9086c3SLois Curfman McInnes 
59e9cf3c21SBarry Smith extern int VecDot_Seq(Vec,Vec,Scalar *);
60c481317fSBarry Smith extern int VecNorm_Seq(Vec,NormType,double *);
61c481317fSBarry Smith 
625615d1e5SSatish Balay #undef __FUNC__
635615d1e5SSatish Balay #define __FUNC__ "SNESMatrixFreeMult_Private"
64eb9086c3SLois Curfman McInnes /*
65eb9086c3SLois Curfman McInnes   SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector
66eb9086c3SLois Curfman McInnes   product, y = F'(u)*a:
67eb9086c3SLois Curfman McInnes         y = ( F(u + ha) - F(u) ) /h,
68eb9086c3SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
69eb9086c3SLois Curfman McInnes         u = current iterate
70eb9086c3SLois Curfman McInnes         h = difference interval
71eb9086c3SLois Curfman McInnes */
72eb9086c3SLois Curfman McInnes int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y)
7339e2f89bSBarry Smith {
74fae171e0SBarry Smith   MFCtx_Private *ctx;
75fae171e0SBarry Smith   SNES          snes;
76*c31ba22aSSatish Balay   double        ovalues[3],norm, sum, umin;
775334005bSBarry Smith   Scalar        h, dot, mone = -1.0;
78fae171e0SBarry Smith   Vec           w,U,F;
7983e56afdSLois Curfman McInnes   int           ierr, (*eval_fct)(SNES,Vec,Vec);
80c481317fSBarry Smith   MPI_Comm      comm;
81*c31ba22aSSatish Balay #if !defined(USE_PETSC_COMPLEX)
82*c31ba22aSSatish Balay   double        values[3];
83*c31ba22aSSatish Balay #endif
8439e2f89bSBarry Smith 
853a40ed3dSBarry Smith   PetscFunctionBegin;
8656cd22aeSBarry Smith   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
8756cd22aeSBarry Smith 
88c481317fSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
896e38a044SBarry Smith   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
90fae171e0SBarry Smith   snes = ctx->snes;
91fae171e0SBarry Smith   w    = ctx->w;
92eb9086c3SLois Curfman McInnes   umin = ctx->umin;
93fae171e0SBarry Smith 
9457c37382SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
9557c37382SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
9657c37382SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
9757c37382SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
9857c37382SLois Curfman McInnes 
9978b31e54SBarry Smith   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
10083e56afdSLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
10183e56afdSLois Curfman McInnes     eval_fct = SNESComputeFunction;
10278b31e54SBarry Smith     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
10383e56afdSLois Curfman McInnes   }
10483e56afdSLois Curfman McInnes   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
10583e56afdSLois Curfman McInnes     eval_fct = SNESComputeGradient;
10683e56afdSLois Curfman McInnes     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
10783e56afdSLois Curfman McInnes   }
108a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class");
10950361f65SLois Curfman McInnes 
110eb9086c3SLois Curfman McInnes   /* Determine a "good" step size, h */
111c481317fSBarry Smith 
112c481317fSBarry Smith   /*
113eb9086c3SLois Curfman McInnes     ierr = VecDot(U,a,&dot); CHKERRQ(ierr);
114eb9086c3SLois Curfman McInnes     ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr);
115eb9086c3SLois Curfman McInnes     ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr);
116c481317fSBarry Smith   */
117c481317fSBarry Smith 
118c481317fSBarry Smith   /*
11986f5173aSBarry Smith      Call the Seq Vector routines and then do a single reduction
12086f5173aSBarry Smith      to reduce the number of communications required
121c481317fSBarry Smith   */
122c481317fSBarry Smith 
1233a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
124c481317fSBarry Smith   PLogEventBegin(VEC_Dot,U,a,0,0);
125c481317fSBarry Smith   ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr);
126c481317fSBarry Smith   PLogEventEnd(VEC_Dot,U,a,0,0);
127c481317fSBarry Smith   PLogEventBegin(VEC_Norm,a,0,0,0);
128c481317fSBarry Smith   ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
129c481317fSBarry Smith   ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
130c481317fSBarry Smith   ovalues[2] = ovalues[2]*ovalues[2];
131005c665bSBarry Smith   PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
132ca161407SBarry Smith   ierr = MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );CHKERRQ(ierr);
133005c665bSBarry Smith   PLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm);
134c481317fSBarry Smith   dot = values[0]; sum = values[1]; norm = sqrt(values[2]);
135c481317fSBarry Smith   PLogEventEnd(VEC_Norm,a,0,0,0);
13686f5173aSBarry Smith #else
13786f5173aSBarry Smith   {
13886f5173aSBarry Smith     Scalar cvalues[3],covalues[3];
13986f5173aSBarry Smith 
14086f5173aSBarry Smith     PLogEventBegin(VEC_Dot,U,a,0,0);
14186f5173aSBarry Smith     ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr);
14286f5173aSBarry Smith     PLogEventEnd(VEC_Dot,U,a,0,0);
14386f5173aSBarry Smith     PLogEventBegin(VEC_Norm,a,0,0,0);
14486f5173aSBarry Smith     ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
14586f5173aSBarry Smith     ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
14686f5173aSBarry Smith     covalues[1] = ovalues[1];
14786f5173aSBarry Smith     covalues[2] = ovalues[2]*ovalues[2];
148005c665bSBarry Smith     PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
149ca161407SBarry Smith     ierr = MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
150005c665bSBarry Smith     PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
15186f5173aSBarry Smith     dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2]));
15286f5173aSBarry Smith     PLogEventEnd(VEC_Norm,a,0,0,0);
15386f5173aSBarry Smith   }
15486f5173aSBarry Smith #endif
155c481317fSBarry Smith 
156eb9086c3SLois Curfman McInnes 
157eb9086c3SLois Curfman McInnes   /* Safeguard for step sizes too small */
158edd2f0e1SBarry Smith   if (sum == 0.0) {dot = 1.0; norm = 1.0;}
1593a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
160eb9086c3SLois Curfman McInnes   else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum;
161eb9086c3SLois Curfman McInnes   else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum;
16219a167f6SBarry Smith #else
163eb9086c3SLois Curfman McInnes   else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
164eb9086c3SLois Curfman McInnes   else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
16519a167f6SBarry Smith #endif
166eb9086c3SLois Curfman McInnes   h = ctx->error_rel*dot/(norm*norm);
16739e2f89bSBarry Smith 
168eb9086c3SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
169eb9086c3SLois Curfman McInnes   ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr);
17083e56afdSLois Curfman McInnes   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
171195ee392SLois Curfman McInnes   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
1725334005bSBarry Smith   h = 1.0/h;
173195ee392SLois Curfman McInnes   ierr = VecScale(&h,y); CHKERRQ(ierr);
174b4fd4287SBarry Smith   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
17557c37382SLois Curfman McInnes 
176eb9086c3SLois Curfman McInnes   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
1773a40ed3dSBarry Smith   PetscFunctionReturn(0);
17839e2f89bSBarry Smith }
17983e56afdSLois Curfman McInnes 
1805615d1e5SSatish Balay #undef __FUNC__
1815615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatCreate"
1824b828684SBarry Smith /*@C
1835392566eSBarry Smith    SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix
18450361f65SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
18550361f65SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
1865392566eSBarry Smith 
1875392566eSBarry Smith    Input Parameters:
188eb9086c3SLois Curfman McInnes .  snes - the SNES context
1895392566eSBarry Smith .  x - vector where SNES solution is to be stored.
1905392566eSBarry Smith 
191eb9086c3SLois Curfman McInnes    Output Parameter:
1925392566eSBarry Smith .  J - the matrix-free matrix
1935392566eSBarry Smith 
19465afa06eSLois Curfman McInnes    Notes:
19565afa06eSLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
19665afa06eSLois Curfman McInnes    and work space for performing finite difference approximations of
1973d5db8e1SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
1983d5db8e1SLois Curfman McInnes 
1993d5db8e1SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
2003d5db8e1SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
2013d5db8e1SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
2023d5db8e1SLois Curfman McInnes $   where
2033d5db8e1SLois Curfman McInnes $        error_rel = square root of relative error in
2043d5db8e1SLois Curfman McInnes $                    function evaluation
2053d5db8e1SLois Curfman McInnes $        umin = minimum iterate parameter
2063d5db8e1SLois Curfman McInnes 
2073d5db8e1SLois Curfman McInnes    The user can set these parameters via SNESSetMatrixFreeParameters().
2083d5db8e1SLois Curfman McInnes    See the nonlinear solvers chapter of the users manual for details.
20965afa06eSLois Curfman McInnes 
21065afa06eSLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
21165afa06eSLois Curfman McInnes    matrix context.
21265afa06eSLois Curfman McInnes 
2133d5db8e1SLois Curfman McInnes    Options Database Keys:
2143d5db8e1SLois Curfman McInnes $  -snes_mf_err <error_rel>
2153d5db8e1SLois Curfman McInnes $  -snes_mf_unim <umin>
2163d5db8e1SLois Curfman McInnes 
21765afa06eSLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
21865afa06eSLois Curfman McInnes 
2193d5db8e1SLois Curfman McInnes .seealso: MatDestroy(), SNESSetMatrixFreeParameters()
2205392566eSBarry Smith @*/
2215392566eSBarry Smith int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J)
2225392566eSBarry Smith {
2235392566eSBarry Smith   MPI_Comm      comm;
2245392566eSBarry Smith   MFCtx_Private *mfctx;
225eb9086c3SLois Curfman McInnes   int           n, nloc, ierr, flg;
2263e966953SLois Curfman McInnes   char          p[64];
2275392566eSBarry Smith 
2283a40ed3dSBarry Smith   PetscFunctionBegin;
2290452661fSBarry Smith   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
230464493b3SBarry Smith   PLogObjectMemory(snes,sizeof(MFCtx_Private));
231b4fd4287SBarry Smith   mfctx->sp   = 0;
2325392566eSBarry Smith   mfctx->snes = snes;
233eb9086c3SLois Curfman McInnes   mfctx->error_rel = 1.e-8; /* assumes double precision */
234639f9d9dSBarry Smith   mfctx->umin      = 1.e-6;
235eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr);
236eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr);
237639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
2383e966953SLois Curfman McInnes   PetscStrcpy(p,"-");
2393e966953SLois Curfman McInnes   if (snes->prefix) PetscStrcat(p,snes->prefix);
240639f9d9dSBarry Smith   if (flg) {
24176be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel);
24276be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin);
243639f9d9dSBarry Smith   }
2445392566eSBarry Smith   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
245195ee392SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
246195ee392SLois Curfman McInnes   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
2477ddc982cSLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
248029af93fSBarry Smith   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr);
2491c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr);
2501c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr);
2511c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr);
252b9fa9cd0SBarry Smith   PLogObjectParent(*J,mfctx->w);
253d370d78aSBarry Smith   PLogObjectParent(snes,*J);
2543a40ed3dSBarry Smith   PetscFunctionReturn(0);
25581e6777dSBarry Smith }
25681e6777dSBarry Smith 
2575615d1e5SSatish Balay #undef __FUNC__
2585615d1e5SSatish Balay #define __FUNC__ "SNESSetMatrixFreeParameters"
259b4fd4287SBarry Smith /*@
260eb9086c3SLois Curfman McInnes    SNESSetMatrixFreeParameters - Sets the parameters for the approximation of
261eb9086c3SLois Curfman McInnes    matrix-vector products using finite differences.
262eb9086c3SLois Curfman McInnes 
263639f9d9dSBarry Smith $       J(u)*a = [J(u+h*a) - J(u)]/h where
264639f9d9dSBarry Smith $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
265639f9d9dSBarry Smith $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
266639f9d9dSBarry Smith $
267eb9086c3SLois Curfman McInnes    Input Parameters:
268eb9086c3SLois Curfman McInnes .  snes - the SNES context
269ccb6204bSLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
270ccb6204bSLois Curfman McInnes      the relative error in the function evaluations)
271eb9086c3SLois Curfman McInnes .  umin - minimum allowable u-value
272eb9086c3SLois Curfman McInnes 
273eb9086c3SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
2743d5db8e1SLois Curfman McInnes 
2753d5db8e1SLois Curfman McInnes .seealso: SNESDefaultMatrixFreeMatCreate()
276eb9086c3SLois Curfman McInnes @*/
277eb9086c3SLois Curfman McInnes int SNESSetMatrixFreeParameters(SNES snes,double error,double umin)
278eb9086c3SLois Curfman McInnes {
279eb9086c3SLois Curfman McInnes   MFCtx_Private *ctx;
280eb9086c3SLois Curfman McInnes   int           ierr;
281eb9086c3SLois Curfman McInnes   Mat           mat;
282eb9086c3SLois Curfman McInnes 
2833a40ed3dSBarry Smith   PetscFunctionBegin;
284eb9086c3SLois Curfman McInnes   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
285eb9086c3SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
286eb9086c3SLois Curfman McInnes   if (ctx) {
287eb9086c3SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
288eb9086c3SLois Curfman McInnes     if (umin != PETSC_DEFAULT)  ctx->umin = umin;
289eb9086c3SLois Curfman McInnes   }
2903a40ed3dSBarry Smith   PetscFunctionReturn(0);
291eb9086c3SLois Curfman McInnes }
292eb9086c3SLois Curfman McInnes 
2935615d1e5SSatish Balay #undef __FUNC__
2945615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace"
295eb9086c3SLois Curfman McInnes /*@
29621a45821SLois Curfman McInnes    SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that
297f525115eSLois Curfman McInnes    an operator is supposed to have.  Since roundoff will create a
298b4fd4287SBarry Smith    small component in the null space, if you know the null space
299b4fd4287SBarry Smith    you may have it automatically removed.
300b4fd4287SBarry Smith 
301b4fd4287SBarry Smith    Input Parameters:
30221a45821SLois Curfman McInnes .  J - the matrix-free matrix context
30321a45821SLois Curfman McInnes .  has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants
304b4fd4287SBarry Smith .  n - number of vectors (excluding constant vector) in null space
30521a45821SLois Curfman McInnes .  vecs - the vectors that span the null space (excluding the constant vector);
306b4fd4287SBarry Smith .         these vectors must be orthonormal
307b4fd4287SBarry Smith 
308b4fd4287SBarry Smith .keywords: SNES, matrix-free, null space
309b4fd4287SBarry Smith @*/
310b4fd4287SBarry Smith int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs)
311b4fd4287SBarry Smith {
312b4fd4287SBarry Smith   int           ierr;
313b4fd4287SBarry Smith   MFCtx_Private *ctx;
314b4fd4287SBarry Smith   MPI_Comm      comm;
315b4fd4287SBarry Smith 
3163a40ed3dSBarry Smith   PetscFunctionBegin;
317b4fd4287SBarry Smith   PetscObjectGetComm((PetscObject)J,&comm);
318b4fd4287SBarry Smith 
319b4fd4287SBarry Smith   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
320b4fd4287SBarry Smith   /* no context indicates that it is not the "matrix free" matrix type */
3213a40ed3dSBarry Smith   if (!ctx) PetscFunctionReturn(0);
322b4fd4287SBarry Smith   ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr);
3233a40ed3dSBarry Smith   PetscFunctionReturn(0);
324b4fd4287SBarry Smith }
325b4fd4287SBarry Smith 
3265334005bSBarry Smith 
3275334005bSBarry Smith 
3285334005bSBarry Smith 
329