1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*2f41ae55SBarry Smith static char vcid[] = "$Id: snesmfj.c,v 1.69 1998/10/26 03:26:40 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 */ 15*2f41ae55SBarry Smith Scalar currenth; /* last differencing parameter used */ 16*2f41ae55SBarry Smith Scalar *historyh; /* history of h */ 1779a81275SBarry Smith int ncurrenth,maxcurrenth; 1839e2f89bSBarry Smith } MFCtx_Private; 1939e2f89bSBarry Smith 205615d1e5SSatish Balay #undef __FUNC__ 21d4bb536fSBarry Smith #define __FUNC__ "SNESMatrixFreeDestroy_Private" 22e1311b90SBarry Smith int SNESMatrixFreeDestroy_Private(Mat mat) 23b9fa9cd0SBarry Smith { 24b9fa9cd0SBarry Smith int ierr; 25fae171e0SBarry Smith MFCtx_Private *ctx; 26fae171e0SBarry Smith 273a40ed3dSBarry Smith PetscFunctionBegin; 280a25c783SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 29b9fa9cd0SBarry Smith ierr = VecDestroy(ctx->w); CHKERRQ(ierr); 30b4fd4287SBarry Smith if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);} 31fae171e0SBarry Smith PetscFree(ctx); 323a40ed3dSBarry Smith PetscFunctionReturn(0); 33b9fa9cd0SBarry Smith } 3450361f65SLois Curfman McInnes 355615d1e5SSatish Balay #undef __FUNC__ 36d4bb536fSBarry Smith #define __FUNC__ "SNESMatrixFreeView_Private" 3739e2f89bSBarry Smith /* 381d1bbde2SLois Curfman McInnes SNESMatrixFreeView_Private - Views matrix-free parameters. 398f6e3e37SBarry Smith 4039e2f89bSBarry Smith */ 41eb9086c3SLois Curfman McInnes int SNESMatrixFreeView_Private(Mat J,Viewer viewer) 42eb9086c3SLois Curfman McInnes { 43eb9086c3SLois Curfman McInnes int ierr; 44eb9086c3SLois Curfman McInnes MFCtx_Private *ctx; 45eb9086c3SLois Curfman McInnes MPI_Comm comm; 46eb9086c3SLois Curfman McInnes FILE *fd; 47eb9086c3SLois Curfman McInnes ViewerType vtype; 48eb9086c3SLois Curfman McInnes 493a40ed3dSBarry Smith PetscFunctionBegin; 50eb9086c3SLois Curfman McInnes PetscObjectGetComm((PetscObject)J,&comm); 51eb9086c3SLois Curfman McInnes ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 52eb9086c3SLois Curfman McInnes ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 53eb9086c3SLois Curfman McInnes ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 54eb9086c3SLois Curfman McInnes if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 55eb9086c3SLois Curfman McInnes PetscFPrintf(comm,fd," SNES matrix-free approximation:\n"); 56eb9086c3SLois Curfman McInnes PetscFPrintf(comm,fd," err=%g (relative error in function evaluation)\n",ctx->error_rel); 57eb9086c3SLois Curfman McInnes PetscFPrintf(comm,fd," umin=%g (minimum iterate parameter)\n",ctx->umin); 585cd90555SBarry Smith } else { 595cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported for this object"); 60eb9086c3SLois Curfman McInnes } 613a40ed3dSBarry Smith PetscFunctionReturn(0); 62eb9086c3SLois Curfman McInnes } 63eb9086c3SLois Curfman McInnes 64e9cf3c21SBarry Smith extern int VecDot_Seq(Vec,Vec,Scalar *); 65c481317fSBarry Smith extern int VecNorm_Seq(Vec,NormType,double *); 66c481317fSBarry Smith 675615d1e5SSatish Balay #undef __FUNC__ 685615d1e5SSatish Balay #define __FUNC__ "SNESMatrixFreeMult_Private" 69eb9086c3SLois Curfman McInnes /* 70eb9086c3SLois Curfman McInnes SNESMatrixFreeMult_Private - Default matrix-free form for Jacobian-vector 71eb9086c3SLois Curfman McInnes product, y = F'(u)*a: 72eb9086c3SLois Curfman McInnes y = ( F(u + ha) - F(u) ) /h, 73eb9086c3SLois Curfman McInnes where F = nonlinear function, as set by SNESSetFunction() 74eb9086c3SLois Curfman McInnes u = current iterate 75eb9086c3SLois Curfman McInnes h = difference interval 76eb9086c3SLois Curfman McInnes */ 77eb9086c3SLois Curfman McInnes int SNESMatrixFreeMult_Private(Mat mat,Vec a,Vec y) 7839e2f89bSBarry Smith { 79fae171e0SBarry Smith MFCtx_Private *ctx; 80fae171e0SBarry Smith SNES snes; 81c31ba22aSSatish Balay double ovalues[3],norm, sum, umin; 825334005bSBarry Smith Scalar h, dot, mone = -1.0; 83fae171e0SBarry Smith Vec w,U,F; 8483e56afdSLois Curfman McInnes int ierr, (*eval_fct)(SNES,Vec,Vec); 85c481317fSBarry Smith MPI_Comm comm; 86c31ba22aSSatish Balay #if !defined(USE_PETSC_COMPLEX) 87c31ba22aSSatish Balay double values[3]; 88c31ba22aSSatish Balay #endif 8939e2f89bSBarry Smith 903a40ed3dSBarry Smith PetscFunctionBegin; 9156cd22aeSBarry Smith PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0); 9256cd22aeSBarry Smith 93c481317fSBarry Smith PetscObjectGetComm((PetscObject)mat,&comm); 946e38a044SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 95fae171e0SBarry Smith snes = ctx->snes; 96fae171e0SBarry Smith w = ctx->w; 97eb9086c3SLois Curfman McInnes umin = ctx->umin; 98fae171e0SBarry Smith 9957c37382SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 10057c37382SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 10157c37382SLois Curfman McInnes storage. We may eventually modify event logging to associate events 10257c37382SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 10357c37382SLois Curfman McInnes 10478b31e54SBarry Smith ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 10583e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 10683e56afdSLois Curfman McInnes eval_fct = SNESComputeFunction; 10778b31e54SBarry Smith ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 10883e56afdSLois Curfman McInnes } 10983e56afdSLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 11083e56afdSLois Curfman McInnes eval_fct = SNESComputeGradient; 11183e56afdSLois Curfman McInnes ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr); 11283e56afdSLois Curfman McInnes } 113a8c6a408SBarry Smith else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class"); 11450361f65SLois Curfman McInnes 115eb9086c3SLois Curfman McInnes /* Determine a "good" step size, h */ 116c481317fSBarry Smith 117c481317fSBarry Smith /* 118eb9086c3SLois Curfman McInnes ierr = VecDot(U,a,&dot); CHKERRQ(ierr); 119eb9086c3SLois Curfman McInnes ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr); 120eb9086c3SLois Curfman McInnes ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr); 121c481317fSBarry Smith */ 122c481317fSBarry Smith 123c481317fSBarry Smith /* 12486f5173aSBarry Smith Call the Seq Vector routines and then do a single reduction 12586f5173aSBarry Smith to reduce the number of communications required 126c481317fSBarry Smith */ 127c481317fSBarry Smith 1283a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 129c481317fSBarry Smith PLogEventBegin(VEC_Dot,U,a,0,0); 130c481317fSBarry Smith ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr); 131c481317fSBarry Smith PLogEventEnd(VEC_Dot,U,a,0,0); 132c481317fSBarry Smith PLogEventBegin(VEC_Norm,a,0,0,0); 133c481317fSBarry Smith ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 134c481317fSBarry Smith ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 135c481317fSBarry Smith ovalues[2] = ovalues[2]*ovalues[2]; 136005c665bSBarry Smith PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 137ca161407SBarry Smith ierr = MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );CHKERRQ(ierr); 138005c665bSBarry Smith PLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm); 139c481317fSBarry Smith dot = values[0]; sum = values[1]; norm = sqrt(values[2]); 140c481317fSBarry Smith PLogEventEnd(VEC_Norm,a,0,0,0); 14186f5173aSBarry Smith #else 14286f5173aSBarry Smith { 14386f5173aSBarry Smith Scalar cvalues[3],covalues[3]; 14486f5173aSBarry Smith 14586f5173aSBarry Smith PLogEventBegin(VEC_Dot,U,a,0,0); 14686f5173aSBarry Smith ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr); 14786f5173aSBarry Smith PLogEventEnd(VEC_Dot,U,a,0,0); 14886f5173aSBarry Smith PLogEventBegin(VEC_Norm,a,0,0,0); 14986f5173aSBarry Smith ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 15086f5173aSBarry Smith ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 15186f5173aSBarry Smith covalues[1] = ovalues[1]; 15286f5173aSBarry Smith covalues[2] = ovalues[2]*ovalues[2]; 153005c665bSBarry Smith PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 154ca161407SBarry Smith ierr = MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 155005c665bSBarry Smith PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 15686f5173aSBarry Smith dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2])); 15786f5173aSBarry Smith PLogEventEnd(VEC_Norm,a,0,0,0); 15886f5173aSBarry Smith } 15986f5173aSBarry Smith #endif 160c481317fSBarry Smith 161eb9086c3SLois Curfman McInnes 162eb9086c3SLois Curfman McInnes /* Safeguard for step sizes too small */ 163edd2f0e1SBarry Smith if (sum == 0.0) {dot = 1.0; norm = 1.0;} 1643a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 165e20fef11SSatish Balay else if (PetscAbsScalar(dot) < umin*sum && PetscReal(dot) >= 0.0) dot = umin*sum; 166e20fef11SSatish Balay else if (PetscAbsScalar(dot) < 0.0 && PetscReal(dot) > -umin*sum) dot = -umin*sum; 16719a167f6SBarry Smith #else 168eb9086c3SLois Curfman McInnes else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 169eb9086c3SLois Curfman McInnes else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 17019a167f6SBarry Smith #endif 171eb9086c3SLois Curfman McInnes h = ctx->error_rel*dot/(norm*norm); 17239e2f89bSBarry Smith 1738f6e3e37SBarry Smith /* keep a record of the current differencing parameter h */ 1748f6e3e37SBarry Smith ctx->currenth = h; 175*2f41ae55SBarry Smith #if defined(USE_PETSC_COMPLEX) 176*2f41ae55SBarry Smith PLogInfo(mat,"Current differencing parameter: %g + %g i\n",PetscReal(h),PetscImaginary(h)); 177*2f41ae55SBarry Smith #else 1788f6e3e37SBarry Smith PLogInfo(mat,"Current differencing parameter: %g\n",h); 179*2f41ae55SBarry Smith #endif 18079a81275SBarry Smith if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 18179a81275SBarry Smith ctx->historyh[ctx->ncurrenth++] = h; 18279a81275SBarry Smith } 1838f6e3e37SBarry Smith 184eb9086c3SLois Curfman McInnes /* Evaluate function at F(u + ha) */ 185eb9086c3SLois Curfman McInnes ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr); 18683e56afdSLois Curfman McInnes ierr = eval_fct(snes,w,y); CHKERRQ(ierr); 1870a25c783SBarry Smith 188195ee392SLois Curfman McInnes ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 1895334005bSBarry Smith h = 1.0/h; 190195ee392SLois Curfman McInnes ierr = VecScale(&h,y); CHKERRQ(ierr); 191b4fd4287SBarry Smith if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);} 19257c37382SLois Curfman McInnes 193eb9086c3SLois Curfman McInnes PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0); 1943a40ed3dSBarry Smith PetscFunctionReturn(0); 19539e2f89bSBarry Smith } 19683e56afdSLois Curfman McInnes 1975615d1e5SSatish Balay #undef __FUNC__ 198*2f41ae55SBarry Smith #define __FUNC__ "SNESDefaultMatrixFreeCreate" 1994b828684SBarry Smith /*@C 200*2f41ae55SBarry Smith SNESDefaultMatrixFreeCreate - Creates a matrix-free matrix 20150361f65SLois Curfman McInnes context for use with a SNES solver. This matrix can be used as 20250361f65SLois Curfman McInnes the Jacobian argument for the routine SNESSetJacobian(). 2035392566eSBarry Smith 204c7afd0dbSLois Curfman McInnes Collective on SNES and Vec 205c7afd0dbSLois Curfman McInnes 2065392566eSBarry Smith Input Parameters: 207c7afd0dbSLois Curfman McInnes + snes - the SNES context 208c7afd0dbSLois Curfman McInnes - x - vector where SNES solution is to be stored. 2095392566eSBarry Smith 210eb9086c3SLois Curfman McInnes Output Parameter: 2115392566eSBarry Smith . J - the matrix-free matrix 2125392566eSBarry Smith 21365afa06eSLois Curfman McInnes Notes: 21465afa06eSLois Curfman McInnes The matrix-free matrix context merely contains the function pointers 21565afa06eSLois Curfman McInnes and work space for performing finite difference approximations of 2163d5db8e1SLois Curfman McInnes Jacobian-vector products, J(u)*a, via 2173d5db8e1SLois Curfman McInnes 218c7afd0dbSLois Curfman McInnes .vb 219c7afd0dbSLois Curfman McInnes J(u)*a = [J(u+h*a) - J(u)]/h where 220c7afd0dbSLois Curfman McInnes h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 221c7afd0dbSLois Curfman McInnes = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 222c7afd0dbSLois Curfman McInnes where 223c7afd0dbSLois Curfman McInnes error_rel = square root of relative error in function evaluation 224c7afd0dbSLois Curfman McInnes umin = minimum iterate parameter 225c7afd0dbSLois Curfman McInnes .ve 2263d5db8e1SLois Curfman McInnes 227*2f41ae55SBarry Smith The user can set these parameters via SNESDefaultMatrixFreeSetParameters(). 2283d5db8e1SLois Curfman McInnes See the nonlinear solvers chapter of the users manual for details. 22965afa06eSLois Curfman McInnes 23065afa06eSLois Curfman McInnes The user should call MatDestroy() when finished with the matrix-free 23165afa06eSLois Curfman McInnes matrix context. 23265afa06eSLois Curfman McInnes 2333d5db8e1SLois Curfman McInnes Options Database Keys: 234c7afd0dbSLois Curfman McInnes + -snes_mf_err <error_rel> - Sets error_rel 235c7afd0dbSLois Curfman McInnes - -snes_mf_unim <umin> - Sets umin 2363d5db8e1SLois Curfman McInnes 23765afa06eSLois Curfman McInnes .keywords: SNES, default, matrix-free, create, matrix 23865afa06eSLois Curfman McInnes 239*2f41ae55SBarry Smith .seealso: MatDestroy(), SNESDefaultMatrixFreeSetParameters() 2405392566eSBarry Smith @*/ 241*2f41ae55SBarry Smith int SNESDefaultMatrixFreeCreate(SNES snes,Vec x, Mat *J) 2425392566eSBarry Smith { 2435392566eSBarry Smith MPI_Comm comm; 2445392566eSBarry Smith MFCtx_Private *mfctx; 245eb9086c3SLois Curfman McInnes int n, nloc, ierr, flg; 2463e966953SLois Curfman McInnes char p[64]; 2475392566eSBarry Smith 2483a40ed3dSBarry Smith PetscFunctionBegin; 2490452661fSBarry Smith mfctx = (MFCtx_Private *) PetscMalloc(sizeof(MFCtx_Private)); CHKPTRQ(mfctx); 250464493b3SBarry Smith PLogObjectMemory(snes,sizeof(MFCtx_Private)); 251b4fd4287SBarry Smith mfctx->sp = 0; 2525392566eSBarry Smith mfctx->snes = snes; 253eb9086c3SLois Curfman McInnes mfctx->error_rel = 1.e-8; /* assumes double precision */ 254639f9d9dSBarry Smith mfctx->umin = 1.e-6; 25579a81275SBarry Smith mfctx->currenth = 0.0; 25679a81275SBarry Smith mfctx->historyh = PETSC_NULL; 25779a81275SBarry Smith mfctx->ncurrenth = 0; 25879a81275SBarry Smith mfctx->maxcurrenth = 0; 25979a81275SBarry Smith 260eb9086c3SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr); 261eb9086c3SLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr); 262639f9d9dSBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 2633e966953SLois Curfman McInnes PetscStrcpy(p,"-"); 2643e966953SLois Curfman McInnes if (snes->prefix) PetscStrcat(p,snes->prefix); 265639f9d9dSBarry Smith if (flg) { 26676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel); 26776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin); 268639f9d9dSBarry Smith } 2695392566eSBarry Smith ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 270195ee392SLois Curfman McInnes ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 271195ee392SLois Curfman McInnes ierr = VecGetSize(x,&n); CHKERRQ(ierr); 2727ddc982cSLois Curfman McInnes ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr); 273029af93fSBarry Smith ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr); 2741c1c02c0SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)SNESMatrixFreeMult_Private);CHKERRQ(ierr); 2751c1c02c0SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)SNESMatrixFreeDestroy_Private);CHKERRQ(ierr); 2761c1c02c0SLois Curfman McInnes ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)SNESMatrixFreeView_Private); CHKERRQ(ierr); 277b9fa9cd0SBarry Smith PLogObjectParent(*J,mfctx->w); 278d370d78aSBarry Smith PLogObjectParent(snes,*J); 2793a40ed3dSBarry Smith PetscFunctionReturn(0); 28081e6777dSBarry Smith } 28181e6777dSBarry Smith 2825615d1e5SSatish Balay #undef __FUNC__ 283*2f41ae55SBarry Smith #define __FUNC__ "SNESDefaultMatrixFreeGetH" 2848f6e3e37SBarry Smith /*@ 285*2f41ae55SBarry Smith SNESDefaultMatrixFreeGetH - Gets the last h that was used as the differencing 2868f6e3e37SBarry Smith parameter. 2878f6e3e37SBarry Smith 2888f6e3e37SBarry Smith Not Collective 2898f6e3e37SBarry Smith 2908f6e3e37SBarry Smith Input Parameters: 291*2f41ae55SBarry Smith . mat - the matrix obtained with SNESDefaultMatrixFreeCreate() 2928f6e3e37SBarry Smith 2938f6e3e37SBarry Smith Output Paramter: 2948f6e3e37SBarry Smith . h - the differencing step size 2958f6e3e37SBarry Smith 2968f6e3e37SBarry Smith .keywords: SNES, matrix-free, parameters 2978f6e3e37SBarry Smith 298*2f41ae55SBarry Smith .seealso: SNESDefaultMatrixFreeCreate() 2998f6e3e37SBarry Smith @*/ 300*2f41ae55SBarry Smith int SNESDefaultMatrixFreeGetH(Mat mat,Scalar *h) 3018f6e3e37SBarry Smith { 3028f6e3e37SBarry Smith MFCtx_Private *ctx; 3038f6e3e37SBarry Smith int ierr; 3048f6e3e37SBarry Smith 3058f6e3e37SBarry Smith PetscFunctionBegin; 3068f6e3e37SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 3078f6e3e37SBarry Smith if (ctx) { 3088f6e3e37SBarry Smith *h = ctx->currenth; 3098f6e3e37SBarry Smith } 3108f6e3e37SBarry Smith PetscFunctionReturn(0); 3118f6e3e37SBarry Smith } 3128f6e3e37SBarry Smith 3138f6e3e37SBarry Smith #undef __FUNC__ 314*2f41ae55SBarry Smith #define __FUNC__ "SNESDefaultMatrixFreeKSPMonitor" 31579a81275SBarry Smith /* 316*2f41ae55SBarry Smith SNESDefaultMatrixFreeKSPMonitor - A KSP monitor for use with the default PETSc 31779a81275SBarry Smith SNES matrix free routines. Prints the h differencing parameter used at each 31879a81275SBarry Smith timestep. 31979a81275SBarry Smith 32079a81275SBarry Smith */ 321*2f41ae55SBarry Smith int SNESDefaultMatrixFreeKSPMonitor(KSP ksp,int n,double rnorm,void *dummy) 32279a81275SBarry Smith { 32379a81275SBarry Smith PC pc; 32479a81275SBarry Smith MFCtx_Private *ctx; 32579a81275SBarry Smith int ierr; 32679a81275SBarry Smith Mat mat; 32779a81275SBarry Smith MPI_Comm comm; 32879a81275SBarry Smith PetscTruth nonzeroinitialguess; 32979a81275SBarry Smith 33079a81275SBarry Smith PetscFunctionBegin; 33179a81275SBarry Smith ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr); 33279a81275SBarry Smith ierr = KSPGetPC(ksp,&pc); CHKERRQ(ierr); 33379a81275SBarry Smith ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr); 33479a81275SBarry Smith ierr = PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 33579a81275SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 33679a81275SBarry Smith if (n > 0 || nonzeroinitialguess) { 337*2f41ae55SBarry Smith #if defined(USE_PETSC_COMPLEX) 338*2f41ae55SBarry Smith PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm, 339*2f41ae55SBarry Smith PetscReal(ctx->currenth),PetscImaginary(ctx->currenth)); 340*2f41ae55SBarry Smith #else 34179a81275SBarry Smith PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth); 342*2f41ae55SBarry Smith #endif 34379a81275SBarry Smith } else { 34479a81275SBarry Smith PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm); 34579a81275SBarry Smith } 34679a81275SBarry Smith PetscFunctionReturn(0); 34779a81275SBarry Smith } 34879a81275SBarry Smith 34979a81275SBarry Smith #undef __FUNC__ 350*2f41ae55SBarry Smith #define __FUNC__ "SNESDefaultMatrixFreeSetParameters" 351b4fd4287SBarry Smith /*@ 352*2f41ae55SBarry Smith SNESDefaultMatrixFreeSetParameters - Sets the parameters for the approximation of 353eb9086c3SLois Curfman McInnes matrix-vector products using finite differences. 354eb9086c3SLois Curfman McInnes 355*2f41ae55SBarry Smith Collective on Mat 356c7afd0dbSLois Curfman McInnes 357eb9086c3SLois Curfman McInnes Input Parameters: 358*2f41ae55SBarry Smith + mat - the matrix free matrix created via SNESDefaultMatrixFreeCreate() 359ccb6204bSLois Curfman McInnes . error_rel - relative error (should be set to the square root of 360ccb6204bSLois Curfman McInnes the relative error in the function evaluations) 361c7afd0dbSLois Curfman McInnes - umin - minimum allowable u-value 362eb9086c3SLois Curfman McInnes 363c7afd0dbSLois Curfman McInnes Notes: 364c7afd0dbSLois Curfman McInnes The default matrix-free matrix-vector product routine computes 365c7afd0dbSLois Curfman McInnes .vb 366c7afd0dbSLois Curfman McInnes J(u)*a = [J(u+h*a) - J(u)]/h where 367c7afd0dbSLois Curfman McInnes h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 368c7afd0dbSLois Curfman McInnes = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 369c7afd0dbSLois Curfman McInnes .ve 370c7afd0dbSLois Curfman McInnes 371c7afd0dbSLois Curfman McInnes Options Database Keys: 372c7afd0dbSLois Curfman McInnes + -snes_mf_err <error_rel> - Sets error_rel 373c7afd0dbSLois Curfman McInnes - -snes_mf_unim <umin> - Sets umin 374fee21e36SBarry Smith 375eb9086c3SLois Curfman McInnes .keywords: SNES, matrix-free, parameters 3763d5db8e1SLois Curfman McInnes 377*2f41ae55SBarry Smith .seealso: SNESDefaultMatrixFreeCreate() 378eb9086c3SLois Curfman McInnes @*/ 379*2f41ae55SBarry Smith int SNESDefaultMatrixFreeSetParameters(Mat mat,double error,double umin) 380eb9086c3SLois Curfman McInnes { 381eb9086c3SLois Curfman McInnes MFCtx_Private *ctx; 382eb9086c3SLois Curfman McInnes int ierr; 383eb9086c3SLois Curfman McInnes 3843a40ed3dSBarry Smith PetscFunctionBegin; 385eb9086c3SLois Curfman McInnes ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 386eb9086c3SLois Curfman McInnes if (ctx) { 387eb9086c3SLois Curfman McInnes if (error != PETSC_DEFAULT) ctx->error_rel = error; 388eb9086c3SLois Curfman McInnes if (umin != PETSC_DEFAULT) ctx->umin = umin; 389eb9086c3SLois Curfman McInnes } 3903a40ed3dSBarry Smith PetscFunctionReturn(0); 391eb9086c3SLois Curfman McInnes } 392eb9086c3SLois Curfman McInnes 3935615d1e5SSatish Balay #undef __FUNC__ 394*2f41ae55SBarry Smith #define __FUNC__ "SNESDefaultMatrixFreeAddNullSpace" 395eb9086c3SLois Curfman McInnes /*@ 396*2f41ae55SBarry Smith SNESDefaultMatrixFreeAddNullSpace - Provides a null space that 397f525115eSLois Curfman McInnes an operator is supposed to have. Since roundoff will create a 398b4fd4287SBarry Smith small component in the null space, if you know the null space 399b4fd4287SBarry Smith you may have it automatically removed. 400b4fd4287SBarry Smith 401c7afd0dbSLois Curfman McInnes Collective on Mat 402c7afd0dbSLois Curfman McInnes 403b4fd4287SBarry Smith Input Parameters: 404c7afd0dbSLois Curfman McInnes + J - the matrix-free matrix context 40521a45821SLois Curfman McInnes . has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants 406b4fd4287SBarry Smith . n - number of vectors (excluding constant vector) in null space 407c7afd0dbSLois Curfman McInnes - vecs - the vectors that span the null space (excluding the constant vector); 408c7afd0dbSLois Curfman McInnes these vectors must be orthonormal 409fee21e36SBarry Smith 410b4fd4287SBarry Smith .keywords: SNES, matrix-free, null space 411b4fd4287SBarry Smith @*/ 412*2f41ae55SBarry Smith int SNESDefaultMatrixFreeAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs) 413b4fd4287SBarry Smith { 414b4fd4287SBarry Smith int ierr; 415b4fd4287SBarry Smith MFCtx_Private *ctx; 416b4fd4287SBarry Smith MPI_Comm comm; 417b4fd4287SBarry Smith 4183a40ed3dSBarry Smith PetscFunctionBegin; 419b4fd4287SBarry Smith PetscObjectGetComm((PetscObject)J,&comm); 420b4fd4287SBarry Smith 421b4fd4287SBarry Smith ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 422b4fd4287SBarry Smith /* no context indicates that it is not the "matrix free" matrix type */ 4233a40ed3dSBarry Smith if (!ctx) PetscFunctionReturn(0); 424b4fd4287SBarry Smith ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr); 4253a40ed3dSBarry Smith PetscFunctionReturn(0); 426b4fd4287SBarry Smith } 427b4fd4287SBarry Smith 4285334005bSBarry Smith 4295334005bSBarry Smith 4305334005bSBarry Smith 431