181e6777dSBarry Smith 281e6777dSBarry Smith #ifndef lint 3*56cd22aeSBarry Smith static char vcid[] = "$Id: snesmfj.c,v 1.52 1997/07/02 11:21:45 bsmith Exp bsmith $"; 481e6777dSBarry Smith #endif 581e6777dSBarry Smith 670f55243SBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 7eb9086c3SLois Curfman McInnes #include "pinclude/pviewer.h" 8c481317fSBarry Smith #include <math.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__ 195eea60f9SBarry Smith #define __FUNC__ "SNESMatrixFreeDestroy_Private" /* ADIC Ignore */ 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__ 345eea60f9SBarry Smith #define __FUNC__ "SNESMatrixFreeView_Private" /* ADIC Ignore */ 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 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; 75c481317fSBarry Smith double ovalues[3],values[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; 8039e2f89bSBarry Smith 81*56cd22aeSBarry Smith PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0); 82*56cd22aeSBarry Smith 83c481317fSBarry Smith PetscObjectGetComm((PetscObject)mat,&comm); 846e38a044SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 85fae171e0SBarry Smith snes = ctx->snes; 86fae171e0SBarry Smith w = ctx->w; 87eb9086c3SLois Curfman McInnes umin = ctx->umin; 88fae171e0SBarry Smith 8957c37382SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 9057c37382SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 9157c37382SLois Curfman McInnes storage. We may eventually modify event logging to associate events 9257c37382SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 9357c37382SLois Curfman McInnes 9478b31e54SBarry Smith ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 9583e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 9683e56afdSLois Curfman McInnes eval_fct = SNESComputeFunction; 9778b31e54SBarry Smith ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 9883e56afdSLois Curfman McInnes } 9983e56afdSLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 10083e56afdSLois Curfman McInnes eval_fct = SNESComputeGradient; 10183e56afdSLois Curfman McInnes ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr); 10283e56afdSLois Curfman McInnes } 103e3372554SBarry Smith else SETERRQ(1,0,"Invalid method class"); 10450361f65SLois Curfman McInnes 105eb9086c3SLois Curfman McInnes /* Determine a "good" step size, h */ 106c481317fSBarry Smith 107c481317fSBarry Smith /* 108eb9086c3SLois Curfman McInnes ierr = VecDot(U,a,&dot); CHKERRQ(ierr); 109eb9086c3SLois Curfman McInnes ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr); 110eb9086c3SLois Curfman McInnes ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr); 111c481317fSBarry Smith */ 112c481317fSBarry Smith 113c481317fSBarry Smith /* 11486f5173aSBarry Smith Call the Seq Vector routines and then do a single reduction 11586f5173aSBarry Smith to reduce the number of communications required 116c481317fSBarry Smith */ 117c481317fSBarry Smith 11886f5173aSBarry Smith #if !defined(PETSC_COMPLEX) 119c481317fSBarry Smith PLogEventBegin(VEC_Dot,U,a,0,0); 120c481317fSBarry Smith ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr); 121c481317fSBarry Smith PLogEventEnd(VEC_Dot,U,a,0,0); 122c481317fSBarry Smith PLogEventBegin(VEC_Norm,a,0,0,0); 123c481317fSBarry Smith ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 124c481317fSBarry Smith ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 125c481317fSBarry Smith ovalues[2] = ovalues[2]*ovalues[2]; 126c481317fSBarry Smith MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm ); 127c481317fSBarry Smith dot = values[0]; sum = values[1]; norm = sqrt(values[2]); 128c481317fSBarry Smith PLogEventEnd(VEC_Norm,a,0,0,0); 12986f5173aSBarry Smith #else 13086f5173aSBarry Smith { 13186f5173aSBarry Smith Scalar cvalues[3],covalues[3]; 13286f5173aSBarry Smith 13386f5173aSBarry Smith PLogEventBegin(VEC_Dot,U,a,0,0); 13486f5173aSBarry Smith ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr); 13586f5173aSBarry Smith PLogEventEnd(VEC_Dot,U,a,0,0); 13686f5173aSBarry Smith PLogEventBegin(VEC_Norm,a,0,0,0); 13786f5173aSBarry Smith ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 13886f5173aSBarry Smith ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 13986f5173aSBarry Smith covalues[1] = ovalues[1]; 14086f5173aSBarry Smith covalues[2] = ovalues[2]*ovalues[2]; 14186f5173aSBarry Smith MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm ); 14286f5173aSBarry Smith dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2])); 14386f5173aSBarry Smith PLogEventEnd(VEC_Norm,a,0,0,0); 14486f5173aSBarry Smith } 14586f5173aSBarry Smith #endif 146c481317fSBarry Smith 147eb9086c3SLois Curfman McInnes 148eb9086c3SLois Curfman McInnes /* Safeguard for step sizes too small */ 149edd2f0e1SBarry Smith if (sum == 0.0) {dot = 1.0; norm = 1.0;} 15019a167f6SBarry Smith #if defined(PETSC_COMPLEX) 151eb9086c3SLois Curfman McInnes else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum; 152eb9086c3SLois Curfman McInnes else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum; 15319a167f6SBarry Smith #else 154eb9086c3SLois Curfman McInnes else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 155eb9086c3SLois Curfman McInnes else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 15619a167f6SBarry Smith #endif 157eb9086c3SLois Curfman McInnes h = ctx->error_rel*dot/(norm*norm); 15839e2f89bSBarry Smith 159eb9086c3SLois Curfman McInnes /* Evaluate function at F(u + ha) */ 160eb9086c3SLois Curfman McInnes ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr); 16183e56afdSLois Curfman McInnes ierr = eval_fct(snes,w,y); CHKERRQ(ierr); 162195ee392SLois Curfman McInnes ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 1635334005bSBarry Smith h = 1.0/h; 164195ee392SLois Curfman McInnes ierr = VecScale(&h,y); CHKERRQ(ierr); 165b4fd4287SBarry Smith if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);} 16657c37382SLois Curfman McInnes 167eb9086c3SLois Curfman McInnes PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0); 16839e2f89bSBarry Smith return 0; 16939e2f89bSBarry Smith } 17083e56afdSLois Curfman McInnes 1715615d1e5SSatish Balay #undef __FUNC__ 1725615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatCreate" 1734b828684SBarry Smith /*@C 1745392566eSBarry Smith SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix 17550361f65SLois Curfman McInnes context for use with a SNES solver. This matrix can be used as 17650361f65SLois Curfman McInnes the Jacobian argument for the routine SNESSetJacobian(). 1775392566eSBarry Smith 1785392566eSBarry Smith Input Parameters: 179eb9086c3SLois Curfman McInnes . snes - the SNES context 1805392566eSBarry Smith . x - vector where SNES solution is to be stored. 1815392566eSBarry Smith 182eb9086c3SLois Curfman McInnes Output Parameter: 1835392566eSBarry Smith . J - the matrix-free matrix 1845392566eSBarry Smith 18565afa06eSLois Curfman McInnes Notes: 18665afa06eSLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 18765afa06eSLois Curfman McInnes and work space for performing finite difference approximations of 1883d5db8e1SLois Curfman McInnes Jacobian-vector products, J(u)*a, via 1893d5db8e1SLois Curfman McInnes 1903d5db8e1SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h where 1913d5db8e1SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 1923d5db8e1SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 1933d5db8e1SLois Curfman McInnes $ where 1943d5db8e1SLois Curfman McInnes $ error_rel = square root of relative error in 1953d5db8e1SLois Curfman McInnes $ function evaluation 1963d5db8e1SLois Curfman McInnes $ umin = minimum iterate parameter 1973d5db8e1SLois Curfman McInnes 1983d5db8e1SLois Curfman McInnes The user can set these parameters via SNESSetMatrixFreeParameters(). 1993d5db8e1SLois Curfman McInnes See the nonlinear solvers chapter of the users manual for details. 20065afa06eSLois Curfman McInnes 20165afa06eSLois Curfman McInnes The user should call MatDestroy() when finished with the matrix-free 20265afa06eSLois Curfman McInnes matrix context. 20365afa06eSLois Curfman McInnes 2043d5db8e1SLois Curfman McInnes Options Database Keys: 2053d5db8e1SLois Curfman McInnes $ -snes_mf_err <error_rel> 2063d5db8e1SLois Curfman McInnes $ -snes_mf_unim <umin> 2073d5db8e1SLois Curfman McInnes 20865afa06eSLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix 20965afa06eSLois Curfman McInnes 2103d5db8e1SLois Curfman McInnes .seealso: MatDestroy(), SNESSetMatrixFreeParameters() 2115392566eSBarry Smith @*/ 2125392566eSBarry Smith int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J) 2135392566eSBarry Smith { 2145392566eSBarry Smith MPI_Comm comm; 2155392566eSBarry Smith MFCtx_Private *mfctx; 216eb9086c3SLois Curfman McInnes int n, nloc, ierr, flg; 2173e966953SLois Curfman McInnes char p[64]; 2185392566eSBarry Smith 2190452661fSBarry Smith mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx); 220464493b3SBarry Smith PLogObjectMemory(snes,sizeof(MFCtx_Private)); 221b4fd4287SBarry Smith mfctx->sp = 0; 2225392566eSBarry Smith mfctx->snes = snes; 223eb9086c3SLois Curfman McInnes mfctx->error_rel = 1.e-8; /* assumes double precision */ 224639f9d9dSBarry Smith mfctx->umin = 1.e-6; 225eb9086c3SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr); 226eb9086c3SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr); 227639f9d9dSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 2283e966953SLois Curfman McInnes PetscStrcpy(p,"-"); 2293e966953SLois Curfman McInnes if (snes->prefix) PetscStrcat(p,snes->prefix); 230639f9d9dSBarry Smith if (flg) { 2313e966953SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel); 2323e966953SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin); 233639f9d9dSBarry Smith } 2345392566eSBarry Smith ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 235195ee392SLois Curfman McInnes ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 236195ee392SLois Curfman McInnes ierr = VecGetSize(x,&n); CHKERRQ(ierr); 2377ddc982cSLois Curfman McInnes ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr); 238029af93fSBarry Smith ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr); 2391c1c02c0SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr); 2401c1c02c0SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr); 2411c1c02c0SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr); 242b9fa9cd0SBarry Smith PLogObjectParent(*J,mfctx->w); 243d370d78aSBarry Smith PLogObjectParent(snes,*J); 24481e6777dSBarry Smith return 0; 24581e6777dSBarry Smith } 24681e6777dSBarry Smith 2475615d1e5SSatish Balay #undef __FUNC__ 2485615d1e5SSatish Balay #define __FUNC__ "SNESSetMatrixFreeParameters" 249b4fd4287SBarry Smith /*@ 250eb9086c3SLois Curfman McInnes SNESSetMatrixFreeParameters - Sets the parameters for the approximation of 251eb9086c3SLois Curfman McInnes matrix-vector products using finite differences. 252eb9086c3SLois Curfman McInnes 253639f9d9dSBarry Smith $ J(u)*a = [J(u+h*a) - J(u)]/h where 254639f9d9dSBarry Smith $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 255639f9d9dSBarry Smith $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 256639f9d9dSBarry Smith $ 257eb9086c3SLois Curfman McInnes Input Parameters: 258eb9086c3SLois Curfman McInnes . snes - the SNES context 259ccb6204bSLois Curfman McInnes . error_rel - relative error (should be set to the square root of 260ccb6204bSLois Curfman McInnes the relative error in the function evaluations) 261eb9086c3SLois Curfman McInnes . umin - minimum allowable u-value 262eb9086c3SLois Curfman McInnes 263eb9086c3SLois Curfman McInnes .keywords: SNES, matrix-free, parameters 2643d5db8e1SLois Curfman McInnes 2653d5db8e1SLois Curfman McInnes .seealso: SNESDefaultMatrixFreeMatCreate() 266eb9086c3SLois Curfman McInnes @*/ 267eb9086c3SLois Curfman McInnes int SNESSetMatrixFreeParameters(SNES snes,double error,double umin) 268eb9086c3SLois Curfman McInnes { 269eb9086c3SLois Curfman McInnes MFCtx_Private *ctx; 270eb9086c3SLois Curfman McInnes int ierr; 271eb9086c3SLois Curfman McInnes Mat mat; 272eb9086c3SLois Curfman McInnes 273eb9086c3SLois Curfman McInnes ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 274eb9086c3SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 275eb9086c3SLois Curfman McInnes if (ctx) { 276eb9086c3SLois Curfman McInnes if (error != PETSC_DEFAULT) ctx->error_rel = error; 277eb9086c3SLois Curfman McInnes if (umin != PETSC_DEFAULT) ctx->umin = umin; 278eb9086c3SLois Curfman McInnes } 279eb9086c3SLois Curfman McInnes return 0; 280eb9086c3SLois Curfman McInnes } 281eb9086c3SLois Curfman McInnes 2825615d1e5SSatish Balay #undef __FUNC__ 2835615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace" 284eb9086c3SLois Curfman McInnes /*@ 28521a45821SLois Curfman McInnes SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that 286f525115eSLois Curfman McInnes an operator is supposed to have. Since roundoff will create a 287b4fd4287SBarry Smith small component in the null space, if you know the null space 288b4fd4287SBarry Smith you may have it automatically removed. 289b4fd4287SBarry Smith 290b4fd4287SBarry Smith Input Parameters: 29121a45821SLois Curfman McInnes . J - the matrix-free matrix context 29221a45821SLois Curfman McInnes . has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants 293b4fd4287SBarry Smith . n - number of vectors (excluding constant vector) in null space 29421a45821SLois Curfman McInnes . vecs - the vectors that span the null space (excluding the constant vector); 295b4fd4287SBarry Smith . these vectors must be orthonormal 296b4fd4287SBarry Smith 297b4fd4287SBarry Smith .keywords: SNES, matrix-free, null space 298b4fd4287SBarry Smith @*/ 299b4fd4287SBarry Smith int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs) 300b4fd4287SBarry Smith { 301b4fd4287SBarry Smith int ierr; 302b4fd4287SBarry Smith MFCtx_Private *ctx; 303b4fd4287SBarry Smith MPI_Comm comm; 304b4fd4287SBarry Smith 305b4fd4287SBarry Smith PetscObjectGetComm((PetscObject)J,&comm); 306b4fd4287SBarry Smith 307b4fd4287SBarry Smith ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 308b4fd4287SBarry Smith /* no context indicates that it is not the "matrix free" matrix type */ 309b4fd4287SBarry Smith if (!ctx) return 0; 310b4fd4287SBarry Smith ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr); 311b4fd4287SBarry Smith return 0; 312b4fd4287SBarry Smith } 313b4fd4287SBarry Smith 3145334005bSBarry Smith 3155334005bSBarry Smith 3165334005bSBarry Smith 317