1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*3a40ed3dSBarry Smith static char vcid[] = "$Id: snesmfj.c,v 1.56 1997/08/22 15:17:50 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__ 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 25*3a40ed3dSBarry 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); 30*3a40ed3dSBarry 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 46*3a40ed3dSBarry 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 } 56*3a40ed3dSBarry 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; 76c481317fSBarry Smith double ovalues[3],values[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; 8139e2f89bSBarry Smith 82*3a40ed3dSBarry Smith PetscFunctionBegin; 8356cd22aeSBarry Smith PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0); 8456cd22aeSBarry Smith 85c481317fSBarry Smith PetscObjectGetComm((PetscObject)mat,&comm); 866e38a044SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 87fae171e0SBarry Smith snes = ctx->snes; 88fae171e0SBarry Smith w = ctx->w; 89eb9086c3SLois Curfman McInnes umin = ctx->umin; 90fae171e0SBarry Smith 9157c37382SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 9257c37382SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 9357c37382SLois Curfman McInnes storage. We may eventually modify event logging to associate events 9457c37382SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 9557c37382SLois Curfman McInnes 9678b31e54SBarry Smith ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 9783e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 9883e56afdSLois Curfman McInnes eval_fct = SNESComputeFunction; 9978b31e54SBarry Smith ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 10083e56afdSLois Curfman McInnes } 10183e56afdSLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 10283e56afdSLois Curfman McInnes eval_fct = SNESComputeGradient; 10383e56afdSLois Curfman McInnes ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr); 10483e56afdSLois Curfman McInnes } 105e3372554SBarry Smith else SETERRQ(1,0,"Invalid method class"); 10650361f65SLois Curfman McInnes 107eb9086c3SLois Curfman McInnes /* Determine a "good" step size, h */ 108c481317fSBarry Smith 109c481317fSBarry Smith /* 110eb9086c3SLois Curfman McInnes ierr = VecDot(U,a,&dot); CHKERRQ(ierr); 111eb9086c3SLois Curfman McInnes ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr); 112eb9086c3SLois Curfman McInnes ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr); 113c481317fSBarry Smith */ 114c481317fSBarry Smith 115c481317fSBarry Smith /* 11686f5173aSBarry Smith Call the Seq Vector routines and then do a single reduction 11786f5173aSBarry Smith to reduce the number of communications required 118c481317fSBarry Smith */ 119c481317fSBarry Smith 120*3a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 121c481317fSBarry Smith PLogEventBegin(VEC_Dot,U,a,0,0); 122c481317fSBarry Smith ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr); 123c481317fSBarry Smith PLogEventEnd(VEC_Dot,U,a,0,0); 124c481317fSBarry Smith PLogEventBegin(VEC_Norm,a,0,0,0); 125c481317fSBarry Smith ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 126c481317fSBarry Smith ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 127c481317fSBarry Smith ovalues[2] = ovalues[2]*ovalues[2]; 128005c665bSBarry Smith PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 129c481317fSBarry Smith MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm ); 130005c665bSBarry Smith PLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm); 131c481317fSBarry Smith dot = values[0]; sum = values[1]; norm = sqrt(values[2]); 132c481317fSBarry Smith PLogEventEnd(VEC_Norm,a,0,0,0); 13386f5173aSBarry Smith #else 13486f5173aSBarry Smith { 13586f5173aSBarry Smith Scalar cvalues[3],covalues[3]; 13686f5173aSBarry Smith 13786f5173aSBarry Smith PLogEventBegin(VEC_Dot,U,a,0,0); 13886f5173aSBarry Smith ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr); 13986f5173aSBarry Smith PLogEventEnd(VEC_Dot,U,a,0,0); 14086f5173aSBarry Smith PLogEventBegin(VEC_Norm,a,0,0,0); 14186f5173aSBarry Smith ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 14286f5173aSBarry Smith ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 14386f5173aSBarry Smith covalues[1] = ovalues[1]; 14486f5173aSBarry Smith covalues[2] = ovalues[2]*ovalues[2]; 145005c665bSBarry Smith PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 14686f5173aSBarry Smith MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm ); 147005c665bSBarry Smith PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 14886f5173aSBarry Smith dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2])); 14986f5173aSBarry Smith PLogEventEnd(VEC_Norm,a,0,0,0); 15086f5173aSBarry Smith } 15186f5173aSBarry Smith #endif 152c481317fSBarry Smith 153eb9086c3SLois Curfman McInnes 154eb9086c3SLois Curfman McInnes /* Safeguard for step sizes too small */ 155edd2f0e1SBarry Smith if (sum == 0.0) {dot = 1.0; norm = 1.0;} 156*3a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 157eb9086c3SLois Curfman McInnes else if (abs(dot) < umin*sum && real(dot) >= 0.0) dot = umin*sum; 158eb9086c3SLois Curfman McInnes else if (abs(dot) < 0.0 && real(dot) > -umin*sum) dot = -umin*sum; 15919a167f6SBarry Smith #else 160eb9086c3SLois Curfman McInnes else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 161eb9086c3SLois Curfman McInnes else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 16219a167f6SBarry Smith #endif 163eb9086c3SLois Curfman McInnes h = ctx->error_rel*dot/(norm*norm); 16439e2f89bSBarry Smith 165eb9086c3SLois Curfman McInnes /* Evaluate function at F(u + ha) */ 166eb9086c3SLois Curfman McInnes ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr); 16783e56afdSLois Curfman McInnes ierr = eval_fct(snes,w,y); CHKERRQ(ierr); 168195ee392SLois Curfman McInnes ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 1695334005bSBarry Smith h = 1.0/h; 170195ee392SLois Curfman McInnes ierr = VecScale(&h,y); CHKERRQ(ierr); 171b4fd4287SBarry Smith if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);} 17257c37382SLois Curfman McInnes 173eb9086c3SLois Curfman McInnes PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0); 174*3a40ed3dSBarry Smith PetscFunctionReturn(0); 17539e2f89bSBarry Smith } 17683e56afdSLois Curfman McInnes 1775615d1e5SSatish Balay #undef __FUNC__ 1785615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatCreate" 1794b828684SBarry Smith /*@C 1805392566eSBarry Smith SNESDefaultMatrixFreeMatCreate - Creates a matrix-free matrix 18150361f65SLois Curfman McInnes context for use with a SNES solver. This matrix can be used as 18250361f65SLois Curfman McInnes the Jacobian argument for the routine SNESSetJacobian(). 1835392566eSBarry Smith 1845392566eSBarry Smith Input Parameters: 185eb9086c3SLois Curfman McInnes . snes - the SNES context 1865392566eSBarry Smith . x - vector where SNES solution is to be stored. 1875392566eSBarry Smith 188eb9086c3SLois Curfman McInnes Output Parameter: 1895392566eSBarry Smith . J - the matrix-free matrix 1905392566eSBarry Smith 19165afa06eSLois Curfman McInnes Notes: 19265afa06eSLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 19365afa06eSLois Curfman McInnes and work space for performing finite difference approximations of 1943d5db8e1SLois Curfman McInnes Jacobian-vector products, J(u)*a, via 1953d5db8e1SLois Curfman McInnes 1963d5db8e1SLois Curfman McInnes $ J(u)*a = [J(u+h*a) - J(u)]/h where 1973d5db8e1SLois Curfman McInnes $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 1983d5db8e1SLois Curfman McInnes $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 1993d5db8e1SLois Curfman McInnes $ where 2003d5db8e1SLois Curfman McInnes $ error_rel = square root of relative error in 2013d5db8e1SLois Curfman McInnes $ function evaluation 2023d5db8e1SLois Curfman McInnes $ umin = minimum iterate parameter 2033d5db8e1SLois Curfman McInnes 2043d5db8e1SLois Curfman McInnes The user can set these parameters via SNESSetMatrixFreeParameters(). 2053d5db8e1SLois Curfman McInnes See the nonlinear solvers chapter of the users manual for details. 20665afa06eSLois Curfman McInnes 20765afa06eSLois Curfman McInnes The user should call MatDestroy() when finished with the matrix-free 20865afa06eSLois Curfman McInnes matrix context. 20965afa06eSLois Curfman McInnes 2103d5db8e1SLois Curfman McInnes Options Database Keys: 2113d5db8e1SLois Curfman McInnes $ -snes_mf_err <error_rel> 2123d5db8e1SLois Curfman McInnes $ -snes_mf_unim <umin> 2133d5db8e1SLois Curfman McInnes 21465afa06eSLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix 21565afa06eSLois Curfman McInnes 2163d5db8e1SLois Curfman McInnes .seealso: MatDestroy(), SNESSetMatrixFreeParameters() 2175392566eSBarry Smith @*/ 2185392566eSBarry Smith int SNESDefaultMatrixFreeMatCreate(SNES snes,Vec x, Mat *J) 2195392566eSBarry Smith { 2205392566eSBarry Smith MPI_Comm comm; 2215392566eSBarry Smith MFCtx_Private *mfctx; 222eb9086c3SLois Curfman McInnes int n, nloc, ierr, flg; 2233e966953SLois Curfman McInnes char p[64]; 2245392566eSBarry Smith 225*3a40ed3dSBarry Smith PetscFunctionBegin; 2260452661fSBarry Smith mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx); 227464493b3SBarry Smith PLogObjectMemory(snes,sizeof(MFCtx_Private)); 228b4fd4287SBarry Smith mfctx->sp = 0; 2295392566eSBarry Smith mfctx->snes = snes; 230eb9086c3SLois Curfman McInnes mfctx->error_rel = 1.e-8; /* assumes double precision */ 231639f9d9dSBarry Smith mfctx->umin = 1.e-6; 232eb9086c3SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr); 233eb9086c3SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr); 234639f9d9dSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 2353e966953SLois Curfman McInnes PetscStrcpy(p,"-"); 2363e966953SLois Curfman McInnes if (snes->prefix) PetscStrcat(p,snes->prefix); 237639f9d9dSBarry Smith if (flg) { 2383e966953SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel); 2393e966953SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin); 240639f9d9dSBarry Smith } 2415392566eSBarry Smith ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 242195ee392SLois Curfman McInnes ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 243195ee392SLois Curfman McInnes ierr = VecGetSize(x,&n); CHKERRQ(ierr); 2447ddc982cSLois Curfman McInnes ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr); 245029af93fSBarry Smith ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr); 2461c1c02c0SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr); 2471c1c02c0SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr); 2481c1c02c0SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr); 249b9fa9cd0SBarry Smith PLogObjectParent(*J,mfctx->w); 250d370d78aSBarry Smith PLogObjectParent(snes,*J); 251*3a40ed3dSBarry Smith PetscFunctionReturn(0); 25281e6777dSBarry Smith } 25381e6777dSBarry Smith 2545615d1e5SSatish Balay #undef __FUNC__ 2555615d1e5SSatish Balay #define __FUNC__ "SNESSetMatrixFreeParameters" 256b4fd4287SBarry Smith /*@ 257eb9086c3SLois Curfman McInnes SNESSetMatrixFreeParameters - Sets the parameters for the approximation of 258eb9086c3SLois Curfman McInnes matrix-vector products using finite differences. 259eb9086c3SLois Curfman McInnes 260639f9d9dSBarry Smith $ J(u)*a = [J(u+h*a) - J(u)]/h where 261639f9d9dSBarry Smith $ h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 262639f9d9dSBarry Smith $ = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 263639f9d9dSBarry Smith $ 264eb9086c3SLois Curfman McInnes Input Parameters: 265eb9086c3SLois Curfman McInnes . snes - the SNES context 266ccb6204bSLois Curfman McInnes . error_rel - relative error (should be set to the square root of 267ccb6204bSLois Curfman McInnes the relative error in the function evaluations) 268eb9086c3SLois Curfman McInnes . umin - minimum allowable u-value 269eb9086c3SLois Curfman McInnes 270eb9086c3SLois Curfman McInnes .keywords: SNES, matrix-free, parameters 2713d5db8e1SLois Curfman McInnes 2723d5db8e1SLois Curfman McInnes .seealso: SNESDefaultMatrixFreeMatCreate() 273eb9086c3SLois Curfman McInnes @*/ 274eb9086c3SLois Curfman McInnes int SNESSetMatrixFreeParameters(SNES snes,double error,double umin) 275eb9086c3SLois Curfman McInnes { 276eb9086c3SLois Curfman McInnes MFCtx_Private *ctx; 277eb9086c3SLois Curfman McInnes int ierr; 278eb9086c3SLois Curfman McInnes Mat mat; 279eb9086c3SLois Curfman McInnes 280*3a40ed3dSBarry Smith PetscFunctionBegin; 281eb9086c3SLois Curfman McInnes ierr = SNESGetJacobian(snes,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 282eb9086c3SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 283eb9086c3SLois Curfman McInnes if (ctx) { 284eb9086c3SLois Curfman McInnes if (error != PETSC_DEFAULT) ctx->error_rel = error; 285eb9086c3SLois Curfman McInnes if (umin != PETSC_DEFAULT) ctx->umin = umin; 286eb9086c3SLois Curfman McInnes } 287*3a40ed3dSBarry Smith PetscFunctionReturn(0); 288eb9086c3SLois Curfman McInnes } 289eb9086c3SLois Curfman McInnes 2905615d1e5SSatish Balay #undef __FUNC__ 2915615d1e5SSatish Balay #define __FUNC__ "SNESDefaultMatrixFreeMatAddNullSpace" 292eb9086c3SLois Curfman McInnes /*@ 29321a45821SLois Curfman McInnes SNESDefaultMatrixFreeMatAddNullSpace - Provides a null space that 294f525115eSLois Curfman McInnes an operator is supposed to have. Since roundoff will create a 295b4fd4287SBarry Smith small component in the null space, if you know the null space 296b4fd4287SBarry Smith you may have it automatically removed. 297b4fd4287SBarry Smith 298b4fd4287SBarry Smith Input Parameters: 29921a45821SLois Curfman McInnes . J - the matrix-free matrix context 30021a45821SLois Curfman McInnes . has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants 301b4fd4287SBarry Smith . n - number of vectors (excluding constant vector) in null space 30221a45821SLois Curfman McInnes . vecs - the vectors that span the null space (excluding the constant vector); 303b4fd4287SBarry Smith . these vectors must be orthonormal 304b4fd4287SBarry Smith 305b4fd4287SBarry Smith .keywords: SNES, matrix-free, null space 306b4fd4287SBarry Smith @*/ 307b4fd4287SBarry Smith int SNESDefaultMatrixFreeMatAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs) 308b4fd4287SBarry Smith { 309b4fd4287SBarry Smith int ierr; 310b4fd4287SBarry Smith MFCtx_Private *ctx; 311b4fd4287SBarry Smith MPI_Comm comm; 312b4fd4287SBarry Smith 313*3a40ed3dSBarry Smith PetscFunctionBegin; 314b4fd4287SBarry Smith PetscObjectGetComm((PetscObject)J,&comm); 315b4fd4287SBarry Smith 316b4fd4287SBarry Smith ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 317b4fd4287SBarry Smith /* no context indicates that it is not the "matrix free" matrix type */ 318*3a40ed3dSBarry Smith if (!ctx) PetscFunctionReturn(0); 319b4fd4287SBarry Smith ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr); 320*3a40ed3dSBarry Smith PetscFunctionReturn(0); 321b4fd4287SBarry Smith } 322b4fd4287SBarry Smith 3235334005bSBarry Smith 3245334005bSBarry Smith 3255334005bSBarry Smith 326