xref: /petsc/src/snes/mf/snesmfj.c (revision e1311b9049e89cb3452dcd306fde571f4b440ff2)
1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER
2*e1311b90SBarry Smith static char vcid[] = "$Id: snesmfj.c,v 1.61 1998/03/16 18:45:31 balay 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__
18d4bb536fSBarry Smith #define __FUNC__ "SNESMatrixFreeDestroy_Private"
19*e1311b90SBarry 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);
54eb9086c3SLois Curfman McInnes   }
553a40ed3dSBarry Smith   PetscFunctionReturn(0);
56eb9086c3SLois Curfman McInnes }
57eb9086c3SLois Curfman McInnes 
58e9cf3c21SBarry Smith extern int VecDot_Seq(Vec,Vec,Scalar *);
59c481317fSBarry Smith extern int VecNorm_Seq(Vec,NormType,double *);
60c481317fSBarry Smith 
615615d1e5SSatish Balay #undef __FUNC__
625615d1e5SSatish Balay #define __FUNC__ "SNESMatrixFreeMult_Private"
63eb9086c3SLois Curfman McInnes /*
64eb9086c3SLois Curfman McInnes   SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector
65eb9086c3SLois Curfman McInnes   product, y = F'(u)*a:
66eb9086c3SLois Curfman McInnes         y = ( F(u + ha) - F(u) ) /h,
67eb9086c3SLois Curfman McInnes   where F = nonlinear function, as set by SNESSetFunction()
68eb9086c3SLois Curfman McInnes         u = current iterate
69eb9086c3SLois Curfman McInnes         h = difference interval
70eb9086c3SLois Curfman McInnes */
71eb9086c3SLois Curfman McInnes int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y)
7239e2f89bSBarry Smith {
73fae171e0SBarry Smith   MFCtx_Private *ctx;
74fae171e0SBarry Smith   SNES          snes;
75c31ba22aSSatish Balay   double        ovalues[3],norm, sum, umin;
765334005bSBarry Smith   Scalar        h, dot, mone = -1.0;
77fae171e0SBarry Smith   Vec           w,U,F;
7883e56afdSLois Curfman McInnes   int           ierr, (*eval_fct)(SNES,Vec,Vec);
79c481317fSBarry Smith   MPI_Comm      comm;
80c31ba22aSSatish Balay #if !defined(USE_PETSC_COMPLEX)
81c31ba22aSSatish Balay   double        values[3];
82c31ba22aSSatish Balay #endif
8339e2f89bSBarry Smith 
843a40ed3dSBarry Smith   PetscFunctionBegin;
8556cd22aeSBarry Smith   PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0);
8656cd22aeSBarry Smith 
87c481317fSBarry Smith   PetscObjectGetComm((PetscObject)mat,&comm);
886e38a044SBarry Smith   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
89fae171e0SBarry Smith   snes = ctx->snes;
90fae171e0SBarry Smith   w    = ctx->w;
91eb9086c3SLois Curfman McInnes   umin = ctx->umin;
92fae171e0SBarry Smith 
9357c37382SLois Curfman McInnes   /* We log matrix-free matrix-vector products separately, so that we can
9457c37382SLois Curfman McInnes      separate the performance monitoring from the cases that use conventional
9557c37382SLois Curfman McInnes      storage.  We may eventually modify event logging to associate events
9657c37382SLois Curfman McInnes      with particular objects, hence alleviating the more general problem. */
9757c37382SLois Curfman McInnes 
9878b31e54SBarry Smith   ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr);
9983e56afdSLois Curfman McInnes   if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
10083e56afdSLois Curfman McInnes     eval_fct = SNESComputeFunction;
10178b31e54SBarry Smith     ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr);
10283e56afdSLois Curfman McInnes   }
10383e56afdSLois Curfman McInnes   else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
10483e56afdSLois Curfman McInnes     eval_fct = SNESComputeGradient;
10583e56afdSLois Curfman McInnes     ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr);
10683e56afdSLois Curfman McInnes   }
107a8c6a408SBarry Smith   else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class");
10850361f65SLois Curfman McInnes 
109eb9086c3SLois Curfman McInnes   /* Determine a "good" step size, h */
110c481317fSBarry Smith 
111c481317fSBarry Smith   /*
112eb9086c3SLois Curfman McInnes     ierr = VecDot(U,a,&dot); CHKERRQ(ierr);
113eb9086c3SLois Curfman McInnes     ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr);
114eb9086c3SLois Curfman McInnes     ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr);
115c481317fSBarry Smith   */
116c481317fSBarry Smith 
117c481317fSBarry Smith   /*
11886f5173aSBarry Smith      Call the Seq Vector routines and then do a single reduction
11986f5173aSBarry Smith      to reduce the number of communications required
120c481317fSBarry Smith   */
121c481317fSBarry Smith 
1223a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX)
123c481317fSBarry Smith   PLogEventBegin(VEC_Dot,U,a,0,0);
124c481317fSBarry Smith   ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr);
125c481317fSBarry Smith   PLogEventEnd(VEC_Dot,U,a,0,0);
126c481317fSBarry Smith   PLogEventBegin(VEC_Norm,a,0,0,0);
127c481317fSBarry Smith   ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
128c481317fSBarry Smith   ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
129c481317fSBarry Smith   ovalues[2] = ovalues[2]*ovalues[2];
130005c665bSBarry Smith   PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
131ca161407SBarry Smith   ierr = MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );CHKERRQ(ierr);
132005c665bSBarry Smith   PLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm);
133c481317fSBarry Smith   dot = values[0]; sum = values[1]; norm = sqrt(values[2]);
134c481317fSBarry Smith   PLogEventEnd(VEC_Norm,a,0,0,0);
13586f5173aSBarry Smith #else
13686f5173aSBarry Smith   {
13786f5173aSBarry Smith     Scalar cvalues[3],covalues[3];
13886f5173aSBarry Smith 
13986f5173aSBarry Smith     PLogEventBegin(VEC_Dot,U,a,0,0);
14086f5173aSBarry Smith     ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr);
14186f5173aSBarry Smith     PLogEventEnd(VEC_Dot,U,a,0,0);
14286f5173aSBarry Smith     PLogEventBegin(VEC_Norm,a,0,0,0);
14386f5173aSBarry Smith     ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr);
14486f5173aSBarry Smith     ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr);
14586f5173aSBarry Smith     covalues[1] = ovalues[1];
14686f5173aSBarry Smith     covalues[2] = ovalues[2]*ovalues[2];
147005c665bSBarry Smith     PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
148ca161407SBarry Smith     ierr = MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr);
149005c665bSBarry Smith     PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm);
15086f5173aSBarry Smith     dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2]));
15186f5173aSBarry Smith     PLogEventEnd(VEC_Norm,a,0,0,0);
15286f5173aSBarry Smith   }
15386f5173aSBarry Smith #endif
154c481317fSBarry Smith 
155eb9086c3SLois Curfman McInnes 
156eb9086c3SLois Curfman McInnes   /* Safeguard for step sizes too small */
157edd2f0e1SBarry Smith   if (sum == 0.0) {dot = 1.0; norm = 1.0;}
1583a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX)
159eb9086c3SLois Curfman McInnes   else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum;
160eb9086c3SLois Curfman McInnes   else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum;
16119a167f6SBarry Smith #else
162eb9086c3SLois Curfman McInnes   else if (dot < umin*sum && dot >= 0.0) dot = umin*sum;
163eb9086c3SLois Curfman McInnes   else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum;
16419a167f6SBarry Smith #endif
165eb9086c3SLois Curfman McInnes   h = ctx->error_rel*dot/(norm*norm);
16639e2f89bSBarry Smith 
167eb9086c3SLois Curfman McInnes   /* Evaluate function at F(u + ha) */
168eb9086c3SLois Curfman McInnes   ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr);
16983e56afdSLois Curfman McInnes   ierr = eval_fct(snes,w,y); CHKERRQ(ierr);
170195ee392SLois Curfman McInnes   ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr);
1715334005bSBarry Smith   h = 1.0/h;
172195ee392SLois Curfman McInnes   ierr = VecScale(&h,y); CHKERRQ(ierr);
173b4fd4287SBarry Smith   if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);}
17457c37382SLois Curfman McInnes 
175eb9086c3SLois Curfman McInnes   PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0);
1763a40ed3dSBarry Smith   PetscFunctionReturn(0);
17739e2f89bSBarry Smith }
17883e56afdSLois Curfman McInnes 
1795615d1e5SSatish Balay #undef __FUNC__
1805615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatCreate"
1814b828684SBarry Smith /*@C
1825392566eSBarry Smith    SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix
18350361f65SLois Curfman McInnes    context for use with a SNES solver.  This matrix can be used as
18450361f65SLois Curfman McInnes    the Jacobian argument for the routine SNESSetJacobian().
1855392566eSBarry Smith 
1865392566eSBarry Smith    Input Parameters:
187eb9086c3SLois Curfman McInnes .  snes - the SNES context
1885392566eSBarry Smith .  x - vector where SNES solution is to be stored.
1895392566eSBarry Smith 
190eb9086c3SLois Curfman McInnes    Output Parameter:
1915392566eSBarry Smith .  J - the matrix-free matrix
1925392566eSBarry Smith 
19365afa06eSLois Curfman McInnes    Notes:
19465afa06eSLois Curfman McInnes    The matrix-free matrix context merely contains the function pointers
19565afa06eSLois Curfman McInnes    and work space for performing finite difference approximations of
1963d5db8e1SLois Curfman McInnes    Jacobian-vector products, J(u)*a, via
1973d5db8e1SLois Curfman McInnes 
1983d5db8e1SLois Curfman McInnes $       J(u)*a = [J(u+h*a) - J(u)]/h where
1993d5db8e1SLois Curfman McInnes $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
2003d5db8e1SLois Curfman McInnes $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   otherwise
2013d5db8e1SLois Curfman McInnes $   where
2023d5db8e1SLois Curfman McInnes $        error_rel = square root of relative error in
2033d5db8e1SLois Curfman McInnes $                    function evaluation
2043d5db8e1SLois Curfman McInnes $        umin = minimum iterate parameter
2053d5db8e1SLois Curfman McInnes 
2063d5db8e1SLois Curfman McInnes    The user can set these parameters via SNESSetMatrixFreeParameters().
2073d5db8e1SLois Curfman McInnes    See the nonlinear solvers chapter of the users manual for details.
20865afa06eSLois Curfman McInnes 
20965afa06eSLois Curfman McInnes    The user should call MatDestroy() when finished with the matrix-free
21065afa06eSLois Curfman McInnes    matrix context.
21165afa06eSLois Curfman McInnes 
2123d5db8e1SLois Curfman McInnes    Options Database Keys:
2133d5db8e1SLois Curfman McInnes $  -snes_mf_err <error_rel>
2143d5db8e1SLois Curfman McInnes $  -snes_mf_unim <umin>
2153d5db8e1SLois Curfman McInnes 
21665afa06eSLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix
21765afa06eSLois Curfman McInnes 
2183d5db8e1SLois Curfman McInnes .seealso: MatDestroy(), SNESSetMatrixFreeParameters()
2195392566eSBarry Smith @*/
2205392566eSBarry Smith int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J)
2215392566eSBarry Smith {
2225392566eSBarry Smith   MPI_Comm      comm;
2235392566eSBarry Smith   MFCtx_Private *mfctx;
224eb9086c3SLois Curfman McInnes   int           n, nloc, ierr, flg;
2253e966953SLois Curfman McInnes   char          p[64];
2265392566eSBarry Smith 
2273a40ed3dSBarry Smith   PetscFunctionBegin;
2280452661fSBarry Smith   mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx);
229464493b3SBarry Smith   PLogObjectMemory(snes,sizeof(MFCtx_Private));
230b4fd4287SBarry Smith   mfctx->sp   = 0;
2315392566eSBarry Smith   mfctx->snes = snes;
232eb9086c3SLois Curfman McInnes   mfctx->error_rel = 1.e-8; /* assumes double precision */
233639f9d9dSBarry Smith   mfctx->umin      = 1.e-6;
234eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr);
235eb9086c3SLois Curfman McInnes   ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr);
236639f9d9dSBarry Smith   ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr);
2373e966953SLois Curfman McInnes   PetscStrcpy(p,"-");
2383e966953SLois Curfman McInnes   if (snes->prefix) PetscStrcat(p,snes->prefix);
239639f9d9dSBarry Smith   if (flg) {
24076be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel);
24176be9ce4SBarry Smith     (*PetscHelpPrintf)(snes->comm,"   %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin);
242639f9d9dSBarry Smith   }
2435392566eSBarry Smith   ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr);
244195ee392SLois Curfman McInnes   ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr);
245195ee392SLois Curfman McInnes   ierr = VecGetSize(x,&n); CHKERRQ(ierr);
2467ddc982cSLois Curfman McInnes   ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr);
247029af93fSBarry Smith   ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr);
2481c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr);
2491c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr);
2501c1c02c0SLois Curfman McInnes   ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr);
251b9fa9cd0SBarry Smith   PLogObjectParent(*J,mfctx->w);
252d370d78aSBarry Smith   PLogObjectParent(snes,*J);
2533a40ed3dSBarry Smith   PetscFunctionReturn(0);
25481e6777dSBarry Smith }
25581e6777dSBarry Smith 
2565615d1e5SSatish Balay #undef __FUNC__
2575615d1e5SSatish Balay #define __FUNC__ "SNESSetMatrixFreeParameters"
258b4fd4287SBarry Smith /*@
259eb9086c3SLois Curfman McInnes    SNESSetMatrixFreeParameters - Sets the parameters for the approximation of
260eb9086c3SLois Curfman McInnes    matrix-vector products using finite differences.
261eb9086c3SLois Curfman McInnes 
262639f9d9dSBarry Smith $       J(u)*a = [J(u+h*a) - J(u)]/h where
263639f9d9dSBarry Smith $        h = error_rel*u'a/||a||^2                        if  |u'a| > umin*||a||_{1}
264639f9d9dSBarry Smith $          = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2   else
265639f9d9dSBarry Smith $
266eb9086c3SLois Curfman McInnes    Input Parameters:
267eb9086c3SLois Curfman McInnes .  snes - the SNES context
268ccb6204bSLois Curfman McInnes .  error_rel - relative error (should be set to the square root of
269ccb6204bSLois Curfman McInnes      the relative error in the function evaluations)
270eb9086c3SLois Curfman McInnes .  umin - minimum allowable u-value
271eb9086c3SLois Curfman McInnes 
272eb9086c3SLois Curfman McInnes .keywords: SNES, matrix-free, parameters
2733d5db8e1SLois Curfman McInnes 
2743d5db8e1SLois Curfman McInnes .seealso: SNESDefaultMatrixFreeMatCreate()
275eb9086c3SLois Curfman McInnes @*/
276eb9086c3SLois Curfman McInnes int SNESSetMatrixFreeParameters(SNES snes,double error,double umin)
277eb9086c3SLois Curfman McInnes {
278eb9086c3SLois Curfman McInnes   MFCtx_Private *ctx;
279eb9086c3SLois Curfman McInnes   int           ierr;
280eb9086c3SLois Curfman McInnes   Mat           mat;
281eb9086c3SLois Curfman McInnes 
2823a40ed3dSBarry Smith   PetscFunctionBegin;
283eb9086c3SLois Curfman McInnes   ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr);
284eb9086c3SLois Curfman McInnes   ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr);
285eb9086c3SLois Curfman McInnes   if (ctx) {
286eb9086c3SLois Curfman McInnes     if (error != PETSC_DEFAULT) ctx->error_rel = error;
287eb9086c3SLois Curfman McInnes     if (umin != PETSC_DEFAULT)  ctx->umin = umin;
288eb9086c3SLois Curfman McInnes   }
2893a40ed3dSBarry Smith   PetscFunctionReturn(0);
290eb9086c3SLois Curfman McInnes }
291eb9086c3SLois Curfman McInnes 
2925615d1e5SSatish Balay #undef __FUNC__
2935615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace"
294eb9086c3SLois Curfman McInnes /*@
29521a45821SLois Curfman McInnes    SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that
296f525115eSLois Curfman McInnes    an operator is supposed to have.  Since roundoff will create a
297b4fd4287SBarry Smith    small component in the null space, if you know the null space
298b4fd4287SBarry Smith    you may have it automatically removed.
299b4fd4287SBarry Smith 
300b4fd4287SBarry Smith    Input Parameters:
30121a45821SLois Curfman McInnes .  J - the matrix-free matrix context
30221a45821SLois Curfman McInnes .  has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants
303b4fd4287SBarry Smith .  n - number of vectors (excluding constant vector) in null space
30421a45821SLois Curfman McInnes .  vecs - the vectors that span the null space (excluding the constant vector);
305b4fd4287SBarry Smith .         these vectors must be orthonormal
306b4fd4287SBarry Smith 
307b4fd4287SBarry Smith .keywords: SNES, matrix-free, null space
308b4fd4287SBarry Smith @*/
309b4fd4287SBarry Smith int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs)
310b4fd4287SBarry Smith {
311b4fd4287SBarry Smith   int           ierr;
312b4fd4287SBarry Smith   MFCtx_Private *ctx;
313b4fd4287SBarry Smith   MPI_Comm      comm;
314b4fd4287SBarry Smith 
3153a40ed3dSBarry Smith   PetscFunctionBegin;
316b4fd4287SBarry Smith   PetscObjectGetComm((PetscObject)J,&comm);
317b4fd4287SBarry Smith 
318b4fd4287SBarry Smith   ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr);
319b4fd4287SBarry Smith   /* no context indicates that it is not the "matrix free" matrix type */
3203a40ed3dSBarry Smith   if (!ctx) PetscFunctionReturn(0);
321b4fd4287SBarry Smith   ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr);
3223a40ed3dSBarry Smith   PetscFunctionReturn(0);
323b4fd4287SBarry Smith }
324b4fd4287SBarry Smith 
3255334005bSBarry Smith 
3265334005bSBarry Smith 
3275334005bSBarry Smith 
328