xref: /petsc/src/snes/mf/snesmfj.c (revision d4bb536f0e426e9a0292bbfd5743770a9b03f0d5)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*d4bb536fSBarry Smith static char vcid[] = "$Id: snesmfj.c,v 1.55 1997/08/13 22:26:08 bsmith Exp bsmith $";
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__
18*d4bb536fSBarry 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 
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);
29b9fa9cd0SBarry Smith   return 0;
30b9fa9cd0SBarry Smith }
3150361f65SLois Curfman McInnes 
325615d1e5SSatish Balay #undef __FUNC__
33*d4bb536fSBarry 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 
45eb9086c3SLois Curfman McInnes   PetscObjectGetComm((PetscObject)J,&comm);
46eb9086c3SLois Curfman McInnes   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
47eb9086c3SLois Curfman McInnes   ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr);
48eb9086c3SLois Curfman McInnes   ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr);
49eb9086c3SLois Curfman McInnes   if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) {
50eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"  SNES matrix-free approximation:\n");
51eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"    err=%g (relative error in function evaluation)\n",ctx->error_rel);
52eb9086c3SLois Curfman McInnes      PetscFPrintf(comm,fd,"    umin=%g (minimum iterate parameter)\n",ctx->umin);
53eb9086c3SLois Curfman McInnes   }
54eb9086c3SLois Curfman McInnes   return 0;
55eb9086c3SLois Curfman McInnes }
56eb9086c3SLois Curfman McInnes 
57e9cf3c21SBarry Smith extern int VecDot_Seq(Vec,Vec,Scalar *);
58c481317fSBarry Smith extern int VecNorm_Seq(Vec,NormType,double *);
59c481317fSBarry Smith 
605615d1e5SSatish Balay #undef __FUNC__
615615d1e5SSatish Balay #define __FUNC__ "SNESMatrixFreeMult_Private"
62eb9086c3SLois Curfman McInnes /*
63eb9086c3SLois Curfman McInnes   SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector
64eb9086c3SLois Curfman McInnes   product, y = F'(u)*a:
65eb9086c3SLois Curfman McInnes         y = ( F(u + ha) - F(u) ) /h,
66eb9086c3SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
67eb9086c3SLois Curfman McInnes         u = current iterate
68eb9086c3SLois Curfman McInnes         h = difference interval
69eb9086c3SLois Curfman McInnes */
70eb9086c3SLois Curfman McInnes int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y)
7139e2f89bSBarry Smith {
72fae171e0SBarry Smith   MFCtx_Private *ctx;
73fae171e0SBarry Smith   SNES          snes;
74c481317fSBarry Smith   double        ovalues[3],values[3],norm, sum, umin;
755334005bSBarry Smith   Scalar        h, dot, mone = -1.0;
76fae171e0SBarry Smith   Vec           w,U,F;
7783e56afdSLois Curfman McInnes   int           ierr, (*eval_fct)(SNES,Vec,Vec);
78c481317fSBarry Smith   MPI_Comm      comm;
7939e2f89bSBarry Smith 
8056cd22aeSBarry Smith   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
8156cd22aeSBarry Smith 
82c481317fSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
836e38a044SBarry Smith   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
84fae171e0SBarry Smith   snes = ctx->snes;
85fae171e0SBarry Smith   w    = ctx->w;
86eb9086c3SLois Curfman McInnes   umin = ctx->umin;
87fae171e0SBarry Smith 
8857c37382SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
8957c37382SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
9057c37382SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
9157c37382SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
9257c37382SLois Curfman McInnes 
9378b31e54SBarry Smith   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
9483e56afdSLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
9583e56afdSLois Curfman McInnes     eval_fct = SNESComputeFunction;
9678b31e54SBarry Smith     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
9783e56afdSLois Curfman McInnes   }
9883e56afdSLois Curfman McInnes   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
9983e56afdSLois Curfman McInnes     eval_fct = SNESComputeGradient;
10083e56afdSLois Curfman McInnes     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
10183e56afdSLois Curfman McInnes   }
102e3372554SBarry Smith   else SETERRQ(1,0,"Invalid method class");
10350361f65SLois Curfman McInnes 
104eb9086c3SLois Curfman McInnes   /* Determine a "good" step size, h */
105c481317fSBarry Smith 
106c481317fSBarry Smith   /*
107eb9086c3SLois Curfman McInnes     ierr = VecDot(U,a,&dot); CHKERRQ(ierr);
108eb9086c3SLois Curfman McInnes     ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr);
109eb9086c3SLois Curfman McInnes     ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr);
110c481317fSBarry Smith   */
111c481317fSBarry Smith 
112c481317fSBarry Smith   /*
11386f5173aSBarry Smith      Call the Seq Vector routines and then do a single reduction
11486f5173aSBarry Smith      to reduce the number of communications required
115c481317fSBarry Smith   */
116c481317fSBarry Smith 
11786f5173aSBarry Smith #if !defined(PETSC_COMPLEX)
118c481317fSBarry Smith   PLogEventBegin(VEC_Dot,U,a,0,0);
119c481317fSBarry Smith   ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr);
120c481317fSBarry Smith   PLogEventEnd(VEC_Dot,U,a,0,0);
121c481317fSBarry Smith   PLogEventBegin(VEC_Norm,a,0,0,0);
122c481317fSBarry Smith   ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
123c481317fSBarry Smith   ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
124c481317fSBarry Smith   ovalues[2] = ovalues[2]*ovalues[2];
125005c665bSBarry Smith   PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
126c481317fSBarry Smith   MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );
127005c665bSBarry Smith   PLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm);
128c481317fSBarry Smith   dot = values[0]; sum = values[1]; norm = sqrt(values[2]);
129c481317fSBarry Smith   PLogEventEnd(VEC_Norm,a,0,0,0);
13086f5173aSBarry Smith #else
13186f5173aSBarry Smith   {
13286f5173aSBarry Smith     Scalar cvalues[3],covalues[3];
13386f5173aSBarry Smith 
13486f5173aSBarry Smith     PLogEventBegin(VEC_Dot,U,a,0,0);
13586f5173aSBarry Smith     ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr);
13686f5173aSBarry Smith     PLogEventEnd(VEC_Dot,U,a,0,0);
13786f5173aSBarry Smith     PLogEventBegin(VEC_Norm,a,0,0,0);
13886f5173aSBarry Smith     ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
13986f5173aSBarry Smith     ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
14086f5173aSBarry Smith     covalues[1] = ovalues[1];
14186f5173aSBarry Smith     covalues[2] = ovalues[2]*ovalues[2];
142005c665bSBarry Smith     PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
14386f5173aSBarry Smith     MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm );
144005c665bSBarry Smith     PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
14586f5173aSBarry Smith     dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2]));
14686f5173aSBarry Smith     PLogEventEnd(VEC_Norm,a,0,0,0);
14786f5173aSBarry Smith   }
14886f5173aSBarry Smith #endif
149c481317fSBarry Smith 
150eb9086c3SLois Curfman McInnes 
151eb9086c3SLois Curfman McInnes   /* Safeguard for step sizes too small */
152edd2f0e1SBarry Smith   if (sum == 0.0) {dot = 1.0; norm = 1.0;}
15319a167f6SBarry Smith #if defined(PETSC_COMPLEX)
154eb9086c3SLois Curfman McInnes   else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum;
155eb9086c3SLois Curfman McInnes   else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum;
15619a167f6SBarry Smith #else
157eb9086c3SLois Curfman McInnes   else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
158eb9086c3SLois Curfman McInnes   else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
15919a167f6SBarry Smith #endif
160eb9086c3SLois Curfman McInnes   h = ctx->error_rel*dot/(norm*norm);
16139e2f89bSBarry Smith 
162eb9086c3SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
163eb9086c3SLois Curfman McInnes   ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr);
16483e56afdSLois Curfman McInnes   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
165195ee392SLois Curfman McInnes   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
1665334005bSBarry Smith   h = 1.0/h;
167195ee392SLois Curfman McInnes   ierr = VecScale(&h,y); CHKERRQ(ierr);
168b4fd4287SBarry Smith   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
16957c37382SLois Curfman McInnes 
170eb9086c3SLois Curfman McInnes   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
17139e2f89bSBarry Smith   return 0;
17239e2f89bSBarry Smith }
17383e56afdSLois Curfman McInnes 
1745615d1e5SSatish Balay #undef __FUNC__
1755615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatCreate"
1764b828684SBarry Smith /*@C
1775392566eSBarry Smith    SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix
17850361f65SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
17950361f65SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
1805392566eSBarry Smith 
1815392566eSBarry Smith    Input Parameters:
182eb9086c3SLois Curfman McInnes .  snes - the SNES context
1835392566eSBarry Smith .  x - vector where SNES solution is to be stored.
1845392566eSBarry Smith 
185eb9086c3SLois Curfman McInnes    Output Parameter:
1865392566eSBarry Smith .  J - the matrix-free matrix
1875392566eSBarry Smith 
18865afa06eSLois Curfman McInnes    Notes:
18965afa06eSLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
19065afa06eSLois Curfman McInnes    and work space for performing finite difference approximations of
1913d5db8e1SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
1923d5db8e1SLois Curfman McInnes 
1933d5db8e1SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
1943d5db8e1SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
1953d5db8e1SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
1963d5db8e1SLois Curfman McInnes $   where
1973d5db8e1SLois Curfman McInnes $        error_rel = square root of relative error in
1983d5db8e1SLois Curfman McInnes $                    function evaluation
1993d5db8e1SLois Curfman McInnes $        umin = minimum iterate parameter
2003d5db8e1SLois Curfman McInnes 
2013d5db8e1SLois Curfman McInnes    The user can set these parameters via SNESSetMatrixFreeParameters().
2023d5db8e1SLois Curfman McInnes    See the nonlinear solvers chapter of the users manual for details.
20365afa06eSLois Curfman McInnes 
20465afa06eSLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
20565afa06eSLois Curfman McInnes    matrix context.
20665afa06eSLois Curfman McInnes 
2073d5db8e1SLois Curfman McInnes    Options Database Keys:
2083d5db8e1SLois Curfman McInnes $  -snes_mf_err <error_rel>
2093d5db8e1SLois Curfman McInnes $  -snes_mf_unim <umin>
2103d5db8e1SLois Curfman McInnes 
21165afa06eSLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
21265afa06eSLois Curfman McInnes 
2133d5db8e1SLois Curfman McInnes .seealso: MatDestroy(), SNESSetMatrixFreeParameters()
2145392566eSBarry Smith @*/
2155392566eSBarry Smith int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J)
2165392566eSBarry Smith {
2175392566eSBarry Smith   MPI_Comm      comm;
2185392566eSBarry Smith   MFCtx_Private *mfctx;
219eb9086c3SLois Curfman McInnes   int           n, nloc, ierr, flg;
2203e966953SLois Curfman McInnes   char          p[64];
2215392566eSBarry Smith 
2220452661fSBarry Smith   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
223464493b3SBarry Smith   PLogObjectMemory(snes,sizeof(MFCtx_Private));
224b4fd4287SBarry Smith   mfctx->sp   = 0;
2255392566eSBarry Smith   mfctx->snes = snes;
226eb9086c3SLois Curfman McInnes   mfctx->error_rel = 1.e-8; /* assumes double precision */
227639f9d9dSBarry Smith   mfctx->umin      = 1.e-6;
228eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr);
229eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr);
230639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
2313e966953SLois Curfman McInnes   PetscStrcpy(p,"-");
2323e966953SLois Curfman McInnes   if (snes->prefix) PetscStrcat(p,snes->prefix);
233639f9d9dSBarry Smith   if (flg) {
2343e966953SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel);
2353e966953SLois Curfman McInnes     PetscPrintf(snes->comm,"   %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin);
236639f9d9dSBarry Smith   }
2375392566eSBarry Smith   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
238195ee392SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
239195ee392SLois Curfman McInnes   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
2407ddc982cSLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
241029af93fSBarry Smith   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr);
2421c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr);
2431c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr);
2441c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr);
245b9fa9cd0SBarry Smith   PLogObjectParent(*J,mfctx->w);
246d370d78aSBarry Smith   PLogObjectParent(snes,*J);
24781e6777dSBarry Smith   return 0;
24881e6777dSBarry Smith }
24981e6777dSBarry Smith 
2505615d1e5SSatish Balay #undef __FUNC__
2515615d1e5SSatish Balay #define __FUNC__ "SNESSetMatrixFreeParameters"
252b4fd4287SBarry Smith /*@
253eb9086c3SLois Curfman McInnes    SNESSetMatrixFreeParameters - Sets the parameters for the approximation of
254eb9086c3SLois Curfman McInnes    matrix-vector products using finite differences.
255eb9086c3SLois Curfman McInnes 
256639f9d9dSBarry Smith $       J(u)*a = [J(u+h*a) - J(u)]/h where
257639f9d9dSBarry Smith $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
258639f9d9dSBarry Smith $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
259639f9d9dSBarry Smith $
260eb9086c3SLois Curfman McInnes    Input Parameters:
261eb9086c3SLois Curfman McInnes .  snes - the SNES context
262ccb6204bSLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
263ccb6204bSLois Curfman McInnes      the relative error in the function evaluations)
264eb9086c3SLois Curfman McInnes .  umin - minimum allowable u-value
265eb9086c3SLois Curfman McInnes 
266eb9086c3SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
2673d5db8e1SLois Curfman McInnes 
2683d5db8e1SLois Curfman McInnes .seealso: SNESDefaultMatrixFreeMatCreate()
269eb9086c3SLois Curfman McInnes @*/
270eb9086c3SLois Curfman McInnes int SNESSetMatrixFreeParameters(SNES snes,double error,double umin)
271eb9086c3SLois Curfman McInnes {
272eb9086c3SLois Curfman McInnes   MFCtx_Private *ctx;
273eb9086c3SLois Curfman McInnes   int           ierr;
274eb9086c3SLois Curfman McInnes   Mat           mat;
275eb9086c3SLois Curfman McInnes 
276eb9086c3SLois Curfman McInnes   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
277eb9086c3SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
278eb9086c3SLois Curfman McInnes   if (ctx) {
279eb9086c3SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
280eb9086c3SLois Curfman McInnes     if (umin != PETSC_DEFAULT)  ctx->umin = umin;
281eb9086c3SLois Curfman McInnes   }
282eb9086c3SLois Curfman McInnes   return 0;
283eb9086c3SLois Curfman McInnes }
284eb9086c3SLois Curfman McInnes 
2855615d1e5SSatish Balay #undef __FUNC__
2865615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace"
287eb9086c3SLois Curfman McInnes /*@
28821a45821SLois Curfman McInnes    SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that
289f525115eSLois Curfman McInnes    an operator is supposed to have.  Since roundoff will create a
290b4fd4287SBarry Smith    small component in the null space, if you know the null space
291b4fd4287SBarry Smith    you may have it automatically removed.
292b4fd4287SBarry Smith 
293b4fd4287SBarry Smith    Input Parameters:
29421a45821SLois Curfman McInnes .  J - the matrix-free matrix context
29521a45821SLois Curfman McInnes .  has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants
296b4fd4287SBarry Smith .  n - number of vectors (excluding constant vector) in null space
29721a45821SLois Curfman McInnes .  vecs - the vectors that span the null space (excluding the constant vector);
298b4fd4287SBarry Smith .         these vectors must be orthonormal
299b4fd4287SBarry Smith 
300b4fd4287SBarry Smith .keywords: SNES, matrix-free, null space
301b4fd4287SBarry Smith @*/
302b4fd4287SBarry Smith int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs)
303b4fd4287SBarry Smith {
304b4fd4287SBarry Smith   int           ierr;
305b4fd4287SBarry Smith   MFCtx_Private *ctx;
306b4fd4287SBarry Smith   MPI_Comm      comm;
307b4fd4287SBarry Smith 
308b4fd4287SBarry Smith   PetscObjectGetComm((PetscObject)J,&comm);
309b4fd4287SBarry Smith 
310b4fd4287SBarry Smith   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
311b4fd4287SBarry Smith   /* no context indicates that it is not the "matrix free" matrix type */
312b4fd4287SBarry Smith   if (!ctx) return 0;
313b4fd4287SBarry Smith   ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr);
314b4fd4287SBarry Smith   return 0;
315b4fd4287SBarry Smith }
316b4fd4287SBarry Smith 
3175334005bSBarry Smith 
3185334005bSBarry Smith 
3195334005bSBarry Smith 
320