1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*a4d4d686SBarry Smith static char vcid[] = "$Id: snesmfj.c,v 1.71 1998/10/31 23:58: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" 781e6777dSBarry Smith 8*a4d4d686SBarry Smith typedef struct _p_MatSNESFDMF_HCtx *MatSNESFDMF_HCtx; 9*a4d4d686SBarry Smith typedef struct _p_MatSNESFDMFCtx *MatSNESFDMFCtx; 10*a4d4d686SBarry Smith 11*a4d4d686SBarry Smith typedef struct { 12*a4d4d686SBarry Smith int (*compute)(MatSNESFDMFCtx,Vec,Vec,Scalar *); 13*a4d4d686SBarry Smith int (*view)(MatSNESFDMFCtx); 14*a4d4d686SBarry Smith int (*destroy)(MatSNESFDMFCtx); 15*a4d4d686SBarry Smith } MFHOps; 16*a4d4d686SBarry Smith 17*a4d4d686SBarry Smith struct _p_MatSNESFDMF_HCtx { 18*a4d4d686SBarry Smith MFHOps *ops; 19*a4d4d686SBarry Smith void *data; 20*a4d4d686SBarry Smith }; 21*a4d4d686SBarry Smith 22*a4d4d686SBarry Smith struct _p_MatSNESFDMFCtx { /* context for default matrix-free SNES */ 23eb9086c3SLois Curfman McInnes SNES snes; /* SNES context */ 24eb9086c3SLois Curfman McInnes Vec w; /* work vector */ 25eb9086c3SLois Curfman McInnes PCNullSpace sp; /* null space context */ 26eb9086c3SLois Curfman McInnes double error_rel; /* square root of relative error in computing function */ 27639f9d9dSBarry Smith double umin; /* minimum allowable u'a value relative to |u|_1 */ 282f41ae55SBarry Smith Scalar currenth; /* last differencing parameter used */ 292f41ae55SBarry Smith Scalar *historyh; /* history of h */ 3079a81275SBarry Smith int ncurrenth,maxcurrenth; 31*a4d4d686SBarry Smith MatSNESFDMF_HCtx hctx; 32*a4d4d686SBarry Smith }; 3339e2f89bSBarry Smith 345615d1e5SSatish Balay #undef __FUNC__ 35*a4d4d686SBarry Smith #define __FUNC__ "MatSNESFDMFSetType" 36*a4d4d686SBarry Smith int MatSNESFDMFDestroy_Private(Mat mat) 37b9fa9cd0SBarry Smith { 38b9fa9cd0SBarry Smith int ierr; 39*a4d4d686SBarry Smith MatSNESFDMFCtx ctx; 40*a4d4d686SBarry Smith 41*a4d4d686SBarry Smith PetscFunctionBegin; 42*a4d4d686SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 43*a4d4d686SBarry Smith 44*a4d4d686SBarry Smith 45*a4d4d686SBarry Smith #undef __FUNC__ 46*a4d4d686SBarry Smith #define __FUNC__ "MatSNESFDMFDestroy_Private" 47*a4d4d686SBarry Smith int MatSNESFDMFDestroy_Private(Mat mat) 48*a4d4d686SBarry Smith { 49*a4d4d686SBarry Smith int ierr; 50*a4d4d686SBarry Smith MatSNESFDMFCtx ctx; 51fae171e0SBarry Smith 523a40ed3dSBarry Smith PetscFunctionBegin; 530a25c783SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx);CHKERRQ(ierr); 54b9fa9cd0SBarry Smith ierr = VecDestroy(ctx->w); CHKERRQ(ierr); 55*a4d4d686SBarry Smith if (ctx->hctx->ops->destroy) {ierr = (*ctx->hctx->ops->destroy)(ctx);CHKERRQ(ierr);} 56b4fd4287SBarry Smith if (ctx->sp) {ierr = PCNullSpaceDestroy(ctx->sp);CHKERRQ(ierr);} 57fae171e0SBarry Smith PetscFree(ctx); 583a40ed3dSBarry Smith PetscFunctionReturn(0); 59b9fa9cd0SBarry Smith } 6050361f65SLois Curfman McInnes 615615d1e5SSatish Balay #undef __FUNC__ 62*a4d4d686SBarry Smith #define __FUNC__ "MatSNESFDMFView_Private" 6339e2f89bSBarry Smith /* 64*a4d4d686SBarry Smith MatSNESFDMFView_Private - Views matrix-free parameters. 658f6e3e37SBarry Smith 6639e2f89bSBarry Smith */ 67*a4d4d686SBarry Smith int MatSNESFDMFView_Private(Mat J,Viewer viewer) 68eb9086c3SLois Curfman McInnes { 69eb9086c3SLois Curfman McInnes int ierr; 70*a4d4d686SBarry Smith MatSNESFDMFCtx ctx; 71eb9086c3SLois Curfman McInnes MPI_Comm comm; 72eb9086c3SLois Curfman McInnes FILE *fd; 73eb9086c3SLois Curfman McInnes ViewerType vtype; 74eb9086c3SLois Curfman McInnes 753a40ed3dSBarry Smith PetscFunctionBegin; 76eb9086c3SLois Curfman McInnes PetscObjectGetComm((PetscObject)J,&comm); 77eb9086c3SLois Curfman McInnes ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 78eb9086c3SLois Curfman McInnes ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 79eb9086c3SLois Curfman McInnes ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 80eb9086c3SLois Curfman McInnes if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 81eb9086c3SLois Curfman McInnes PetscFPrintf(comm,fd," SNES matrix-free approximation:\n"); 82eb9086c3SLois Curfman McInnes PetscFPrintf(comm,fd," err=%g (relative error in function evaluation)\n",ctx->error_rel); 83eb9086c3SLois Curfman McInnes PetscFPrintf(comm,fd," umin=%g (minimum iterate parameter)\n",ctx->umin); 845cd90555SBarry Smith } else { 855cd90555SBarry Smith SETERRQ(1,1,"Viewer type not supported for this object"); 86eb9086c3SLois Curfman McInnes } 873a40ed3dSBarry Smith PetscFunctionReturn(0); 88eb9086c3SLois Curfman McInnes } 89eb9086c3SLois Curfman McInnes 90c481317fSBarry Smith 915615d1e5SSatish Balay #undef __FUNC__ 92*a4d4d686SBarry Smith #define __FUNC__ "MatSNESFDMFMult_Private" 93eb9086c3SLois Curfman McInnes /* 94*a4d4d686SBarry Smith MatSNESFDMFMult_Private - Default matrix-free form for Jacobian-vector 95eb9086c3SLois Curfman McInnes product, y = F'(u)*a: 96*a4d4d686SBarry Smith 97eb9086c3SLois Curfman McInnes y = ( F(u + ha) - F(u) ) /h, 98eb9086c3SLois Curfman McInnes where F = nonlinear function, as set by SNESSetFunction() 99eb9086c3SLois Curfman McInnes u = current iterate 100eb9086c3SLois Curfman McInnes h = difference interval 101eb9086c3SLois Curfman McInnes */ 102*a4d4d686SBarry Smith int MatSNESFDMFMult_Private(Mat mat,Vec a,Vec y) 10339e2f89bSBarry Smith { 104*a4d4d686SBarry Smith MatSNESFDMFCtx ctx; 105fae171e0SBarry Smith SNES snes; 106*a4d4d686SBarry Smith Scalar h, mone = -1.0; 107fae171e0SBarry Smith Vec w,U,F; 10883e56afdSLois Curfman McInnes int ierr, (*eval_fct)(SNES,Vec,Vec); 10939e2f89bSBarry Smith 1103a40ed3dSBarry Smith PetscFunctionBegin; 11156cd22aeSBarry Smith PLogEventBegin(MAT_MatrixFreeMult,a,y,0,0); 11256cd22aeSBarry Smith 1136e38a044SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 114fae171e0SBarry Smith snes = ctx->snes; 115fae171e0SBarry Smith w = ctx->w; 116fae171e0SBarry Smith 11757c37382SLois Curfman McInnes /* We log matrix-free matrix-vector products separately, so that we can 11857c37382SLois Curfman McInnes separate the performance monitoring from the cases that use conventional 11957c37382SLois Curfman McInnes storage. We may eventually modify event logging to associate events 12057c37382SLois Curfman McInnes with particular objects, hence alleviating the more general problem. */ 12157c37382SLois Curfman McInnes 12278b31e54SBarry Smith ierr = SNESGetSolution(snes,&U); CHKERRQ(ierr); 12383e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 12483e56afdSLois Curfman McInnes eval_fct = SNESComputeFunction; 12578b31e54SBarry Smith ierr = SNESGetFunction(snes,&F); CHKERRQ(ierr); 126*a4d4d686SBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 12783e56afdSLois Curfman McInnes eval_fct = SNESComputeGradient; 12883e56afdSLois Curfman McInnes ierr = SNESGetGradient(snes,&F); CHKERRQ(ierr); 129*a4d4d686SBarry Smith } else SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Invalid method class"); 13050361f65SLois Curfman McInnes 131*a4d4d686SBarry Smith ierr = (*ctx->hctx->ops->compute)(ctx,U,a,&h); CHKERRQ(ierr); 132*a4d4d686SBarry Smith 133*a4d4d686SBarry Smith /* keep a record of the current differencing parameter h */ 134*a4d4d686SBarry Smith ctx->currenth = h; 135*a4d4d686SBarry Smith #if defined(USE_PETSC_COMPLEX) 136*a4d4d686SBarry Smith PLogInfo(mat,"Current differencing parameter: %g + %g i\n",PetscReal(h),PetscImaginary(h)); 137*a4d4d686SBarry Smith #else 138*a4d4d686SBarry Smith PLogInfo(mat,"Current differencing parameter: %g\n",h); 139*a4d4d686SBarry Smith #endif 140*a4d4d686SBarry Smith if (ctx->historyh && ctx->ncurrenth < ctx->maxcurrenth) { 141*a4d4d686SBarry Smith ctx->historyh[ctx->ncurrenth++] = h; 142*a4d4d686SBarry Smith } 143*a4d4d686SBarry Smith 144*a4d4d686SBarry Smith /* Evaluate function at F(u + ha) */ 145*a4d4d686SBarry Smith ierr = VecWAXPY(&h,a,U,w); CHKERRQ(ierr); 146*a4d4d686SBarry Smith ierr = eval_fct(snes,w,y); CHKERRQ(ierr); 147*a4d4d686SBarry Smith 148*a4d4d686SBarry Smith ierr = VecAXPY(&mone,F,y); CHKERRQ(ierr); 149*a4d4d686SBarry Smith h = 1.0/h; 150*a4d4d686SBarry Smith ierr = VecScale(&h,y); CHKERRQ(ierr); 151*a4d4d686SBarry Smith if (ctx->sp) {ierr = PCNullSpaceRemove(ctx->sp,y); CHKERRQ(ierr);} 152*a4d4d686SBarry Smith 153*a4d4d686SBarry Smith PLogEventEnd(MAT_MatrixFreeMult,a,y,0,0); 154*a4d4d686SBarry Smith PetscFunctionReturn(0); 155*a4d4d686SBarry Smith } 156*a4d4d686SBarry Smith 157*a4d4d686SBarry Smith extern int MatFDMFDefaultComputeH(SNESFDMFCtx,Vec,Vec,Scalar *); 158*a4d4d686SBarry Smith 159*a4d4d686SBarry Smith #undef __FUNC__ 160*a4d4d686SBarry Smith #define __FUNC__ "MatCreateSNESFDMF" 161*a4d4d686SBarry Smith /*@C 162*a4d4d686SBarry Smith MatCreateSNESFDMF - Creates a matrix-free matrix 163*a4d4d686SBarry Smith context for use with a SNES solver. This matrix can be used as 164*a4d4d686SBarry Smith the Jacobian argument for the routine SNESSetJacobian(). 165*a4d4d686SBarry Smith 166*a4d4d686SBarry Smith Collective on SNES and Vec 167*a4d4d686SBarry Smith 168*a4d4d686SBarry Smith Input Parameters: 169*a4d4d686SBarry Smith + snes - the SNES context 170*a4d4d686SBarry Smith - x - vector where SNES solution is to be stored. 171*a4d4d686SBarry Smith 172*a4d4d686SBarry Smith Output Parameter: 173*a4d4d686SBarry Smith . J - the matrix-free matrix 174*a4d4d686SBarry Smith 175*a4d4d686SBarry Smith Notes: 176*a4d4d686SBarry Smith The matrix-free matrix context merely contains the function pointers 177*a4d4d686SBarry Smith and work space for performing finite difference approximations of 178*a4d4d686SBarry Smith Jacobian-vector products, J(u)*a, via 179*a4d4d686SBarry Smith 180*a4d4d686SBarry Smith .vb 181*a4d4d686SBarry Smith J(u)*a = [J(u+h*a) - J(u)]/h where 182*a4d4d686SBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 183*a4d4d686SBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 184*a4d4d686SBarry Smith where 185*a4d4d686SBarry Smith error_rel = square root of relative error in function evaluation 186*a4d4d686SBarry Smith umin = minimum iterate parameter 187*a4d4d686SBarry Smith .ve 188*a4d4d686SBarry Smith 189*a4d4d686SBarry Smith The user can set these parameters via MatSNESFDMFSetParameters(). 190*a4d4d686SBarry Smith See the nonlinear solvers chapter of the users manual for details. 191*a4d4d686SBarry Smith 192*a4d4d686SBarry Smith The user should call MatDestroy() when finished with the matrix-free 193*a4d4d686SBarry Smith matrix context. 194*a4d4d686SBarry Smith 195*a4d4d686SBarry Smith Options Database Keys: 196*a4d4d686SBarry Smith + -snes_mf_err <error_rel> - Sets error_rel 197*a4d4d686SBarry Smith . -snes_mf_unim <umin> - Sets umin 198*a4d4d686SBarry Smith - -snes_mf_ksp_monitor - KSP monitor routine that prints differencing h 199*a4d4d686SBarry Smith 200*a4d4d686SBarry Smith .keywords: SNES, default, matrix-free, create, matrix 201*a4d4d686SBarry Smith 202*a4d4d686SBarry Smith .seealso: MatDestroy(), MatSNESFDMFSetParameters(), 203*a4d4d686SBarry Smith MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(), 204*a4d4d686SBarry Smith MatSNESFDMFGetH(),MatSNESFDMFKSPMonitor() 205*a4d4d686SBarry Smith 206*a4d4d686SBarry Smith @*/ 207*a4d4d686SBarry Smith int MatCreateSNESFDMF(SNES snes,Vec x, Mat *J) 208*a4d4d686SBarry Smith { 209*a4d4d686SBarry Smith MPI_Comm comm; 210*a4d4d686SBarry Smith MatSNESFDMFCtx mfctx; 211*a4d4d686SBarry Smith int n, nloc, ierr, flg; 212*a4d4d686SBarry Smith char p[64]; 213*a4d4d686SBarry Smith 214*a4d4d686SBarry Smith PetscFunctionBegin; 215*a4d4d686SBarry Smith mfctx = (MatSNESFDMFCtx) PetscMalloc(sizeof(struct _p_MatSNESFDMFCtx)); CHKPTRQ(mfctx); 216*a4d4d686SBarry Smith PLogObjectMemory(snes,sizeof(MatSNESFDMFCtx)); 217*a4d4d686SBarry Smith mfctx->sp = 0; 218*a4d4d686SBarry Smith mfctx->snes = snes; 219*a4d4d686SBarry Smith mfctx->error_rel = 1.e-8; /* assumes double precision */ 220*a4d4d686SBarry Smith mfctx->umin = 1.e-6; 221*a4d4d686SBarry Smith mfctx->currenth = 0.0; 222*a4d4d686SBarry Smith mfctx->historyh = PETSC_NULL; 223*a4d4d686SBarry Smith mfctx->ncurrenth = 0; 224*a4d4d686SBarry Smith mfctx->maxcurrenth = 0; 225*a4d4d686SBarry Smith mfctx->hctx = 0; 226*a4d4d686SBarry Smith mfctx->hctx->ops->compute = MatFDMFDefaultComputeH; 227*a4d4d686SBarry Smith mfctx->hctx->ops->destroy = 0; 228*a4d4d686SBarry Smith mfctx->hctx->ops->view = 0; 229*a4d4d686SBarry Smith 230*a4d4d686SBarry Smith ierr = OptionsGetDouble(snes->prefix,"-snes_mf_err",&mfctx->error_rel,&flg); CHKERRQ(ierr); 231*a4d4d686SBarry Smith ierr = OptionsGetDouble(snes->prefix,"-snes_mf_umin",&mfctx->umin,&flg); CHKERRQ(ierr); 232*a4d4d686SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help",&flg); CHKERRQ(ierr); 233*a4d4d686SBarry Smith PetscStrcpy(p,"-"); 234*a4d4d686SBarry Smith if (snes->prefix) PetscStrcat(p,snes->prefix); 235*a4d4d686SBarry Smith if (flg) { 236*a4d4d686SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf_err <err>: set sqrt rel error in function (default %g)\n",p,mfctx->error_rel); 237*a4d4d686SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf_umin <umin> see users manual (default %g)\n",p,mfctx->umin); 238*a4d4d686SBarry Smith } 239*a4d4d686SBarry Smith ierr = VecDuplicate(x,&mfctx->w); CHKERRQ(ierr); 240*a4d4d686SBarry Smith ierr = PetscObjectGetComm((PetscObject)x,&comm); CHKERRQ(ierr); 241*a4d4d686SBarry Smith ierr = VecGetSize(x,&n); CHKERRQ(ierr); 242*a4d4d686SBarry Smith ierr = VecGetLocalSize(x,&nloc); CHKERRQ(ierr); 243*a4d4d686SBarry Smith ierr = MatCreateShell(comm,nloc,n,n,n,mfctx,J); CHKERRQ(ierr); 244*a4d4d686SBarry Smith ierr = MatShellSetOperation(*J,MATOP_MULT,(void*)MatSNESFDMFMult_Private);CHKERRQ(ierr); 245*a4d4d686SBarry Smith ierr = MatShellSetOperation(*J,MATOP_DESTROY,(void *)MatSNESFDMFDestroy_Private);CHKERRQ(ierr); 246*a4d4d686SBarry Smith ierr = MatShellSetOperation(*J,MATOP_VIEW,(void *)MatSNESFDMFView_Private); CHKERRQ(ierr); 247*a4d4d686SBarry Smith PLogObjectParent(*J,mfctx->w); 248*a4d4d686SBarry Smith PLogObjectParent(snes,*J); 249*a4d4d686SBarry Smith PetscFunctionReturn(0); 250*a4d4d686SBarry Smith } 251*a4d4d686SBarry Smith 252*a4d4d686SBarry Smith #undef __FUNC__ 253*a4d4d686SBarry Smith #define __FUNC__ "MatSNESFDMFGetH" 254*a4d4d686SBarry Smith /*@ 255*a4d4d686SBarry Smith MatSNESFDMFGetH - Gets the last h that was used as the differencing 256*a4d4d686SBarry Smith parameter. 257*a4d4d686SBarry Smith 258*a4d4d686SBarry Smith Not Collective 259*a4d4d686SBarry Smith 260*a4d4d686SBarry Smith Input Parameters: 261*a4d4d686SBarry Smith . mat - the matrix obtained with MatCreateSNESFDMF() 262*a4d4d686SBarry Smith 263*a4d4d686SBarry Smith Output Paramter: 264*a4d4d686SBarry Smith . h - the differencing step size 265*a4d4d686SBarry Smith 266*a4d4d686SBarry Smith .keywords: SNES, matrix-free, parameters 267*a4d4d686SBarry Smith 268*a4d4d686SBarry Smith .seealso: MatCreateSNESFDMF(),MatSNESFDMFSetHHistory(), 269*a4d4d686SBarry Smith MatSNESFDMFResetHHistory(),MatSNESFDMFKSPMonitor() 270*a4d4d686SBarry Smith @*/ 271*a4d4d686SBarry Smith int MatSNESFDMFGetH(Mat mat,Scalar *h) 272*a4d4d686SBarry Smith { 273*a4d4d686SBarry Smith MatSNESFDMFCtx ctx; 274*a4d4d686SBarry Smith int ierr; 275*a4d4d686SBarry Smith 276*a4d4d686SBarry Smith PetscFunctionBegin; 277*a4d4d686SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 278*a4d4d686SBarry Smith if (ctx) { 279*a4d4d686SBarry Smith *h = ctx->currenth; 280*a4d4d686SBarry Smith } 281*a4d4d686SBarry Smith PetscFunctionReturn(0); 282*a4d4d686SBarry Smith } 283*a4d4d686SBarry Smith 284*a4d4d686SBarry Smith #undef __FUNC__ 285*a4d4d686SBarry Smith #define __FUNC__ "MatSNESFDMFKSPMonitor" 286*a4d4d686SBarry Smith /* 287*a4d4d686SBarry Smith MatSNESFDMFKSPMonitor - A KSP monitor for use with the default PETSc 288*a4d4d686SBarry Smith SNES matrix free routines. Prints the h differencing parameter used at each 289*a4d4d686SBarry Smith timestep. 290*a4d4d686SBarry Smith 291*a4d4d686SBarry Smith */ 292*a4d4d686SBarry Smith int MatSNESFDMFKSPMonitor(KSP ksp,int n,double rnorm,void *dummy) 293*a4d4d686SBarry Smith { 294*a4d4d686SBarry Smith PC pc; 295*a4d4d686SBarry Smith MatSNESFDMFCtx ctx; 296*a4d4d686SBarry Smith int ierr; 297*a4d4d686SBarry Smith Mat mat; 298*a4d4d686SBarry Smith MPI_Comm comm; 299*a4d4d686SBarry Smith PetscTruth nonzeroinitialguess; 300*a4d4d686SBarry Smith 301*a4d4d686SBarry Smith PetscFunctionBegin; 302*a4d4d686SBarry Smith ierr = PetscObjectGetComm((PetscObject)ksp,&comm);CHKERRQ(ierr); 303*a4d4d686SBarry Smith ierr = KSPGetPC(ksp,&pc); CHKERRQ(ierr); 304*a4d4d686SBarry Smith ierr = KSPGetInitialGuessNonzero(ksp,&nonzeroinitialguess);CHKERRQ(ierr); 305*a4d4d686SBarry Smith ierr = PCGetOperators(pc,&mat,PETSC_NULL,PETSC_NULL); CHKERRQ(ierr); 306*a4d4d686SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 307*a4d4d686SBarry Smith if (n > 0 || nonzeroinitialguess) { 308*a4d4d686SBarry Smith #if defined(USE_PETSC_COMPLEX) 309*a4d4d686SBarry Smith PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g + %g i\n",n,rnorm, 310*a4d4d686SBarry Smith PetscReal(ctx->currenth),PetscImaginary(ctx->currenth)); 311*a4d4d686SBarry Smith #else 312*a4d4d686SBarry Smith PetscPrintf(comm,"%d KSP Residual norm %14.12e h %g \n",n,rnorm,ctx->currenth); 313*a4d4d686SBarry Smith #endif 314*a4d4d686SBarry Smith } else { 315*a4d4d686SBarry Smith PetscPrintf(comm,"%d KSP Residual norm %14.12e\n",n,rnorm); 316*a4d4d686SBarry Smith } 317*a4d4d686SBarry Smith PetscFunctionReturn(0); 318*a4d4d686SBarry Smith } 319*a4d4d686SBarry Smith 320*a4d4d686SBarry Smith #undef __FUNC__ 321*a4d4d686SBarry Smith #define __FUNC__ "MatSNESFDMFSetParameters" 322*a4d4d686SBarry Smith /*@ 323*a4d4d686SBarry Smith MatSNESFDMFSetParameters - Sets the parameters for the approximation of 324*a4d4d686SBarry Smith matrix-vector products using finite differences. 325*a4d4d686SBarry Smith 326*a4d4d686SBarry Smith Collective on Mat 327*a4d4d686SBarry Smith 328*a4d4d686SBarry Smith Input Parameters: 329*a4d4d686SBarry Smith + mat - the matrix free matrix created via MatCreateSNESFDMF() 330*a4d4d686SBarry Smith . error_rel - relative error (should be set to the square root of 331*a4d4d686SBarry Smith the relative error in the function evaluations) 332*a4d4d686SBarry Smith - umin - minimum allowable u-value 333*a4d4d686SBarry Smith 334*a4d4d686SBarry Smith Notes: 335*a4d4d686SBarry Smith The default matrix-free matrix-vector product routine computes 336*a4d4d686SBarry Smith .vb 337*a4d4d686SBarry Smith J(u)*a = [J(u+h*a) - J(u)]/h where 338*a4d4d686SBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 339*a4d4d686SBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 else 340*a4d4d686SBarry Smith .ve 341*a4d4d686SBarry Smith 342*a4d4d686SBarry Smith Options Database Keys: 343*a4d4d686SBarry Smith + -snes_mf_err <error_rel> - Sets error_rel 344*a4d4d686SBarry Smith - -snes_mf_unim <umin> - Sets umin 345*a4d4d686SBarry Smith 346*a4d4d686SBarry Smith .keywords: SNES, matrix-free, parameters 347*a4d4d686SBarry Smith 348*a4d4d686SBarry Smith .seealso: MatCreateSNESFDMF(),MatSNESFDMFGetH(), 349*a4d4d686SBarry Smith MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(), 350*a4d4d686SBarry Smith MatSNESFDMFKSPMonitor() 351*a4d4d686SBarry Smith @*/ 352*a4d4d686SBarry Smith int MatSNESFDMFSetParameters(Mat mat,double error,double umin) 353*a4d4d686SBarry Smith { 354*a4d4d686SBarry Smith MatSNESFDMFCtx ctx; 355*a4d4d686SBarry Smith int ierr; 356*a4d4d686SBarry Smith 357*a4d4d686SBarry Smith PetscFunctionBegin; 358*a4d4d686SBarry Smith ierr = MatShellGetContext(mat,(void **)&ctx); CHKERRQ(ierr); 359*a4d4d686SBarry Smith if (ctx) { 360*a4d4d686SBarry Smith if (error != PETSC_DEFAULT) ctx->error_rel = error; 361*a4d4d686SBarry Smith if (umin != PETSC_DEFAULT) ctx->umin = umin; 362*a4d4d686SBarry Smith } 363*a4d4d686SBarry Smith PetscFunctionReturn(0); 364*a4d4d686SBarry Smith } 365*a4d4d686SBarry Smith 366*a4d4d686SBarry Smith #undef __FUNC__ 367*a4d4d686SBarry Smith #define __FUNC__ "MatSNESFDMFAddNullSpace" 368*a4d4d686SBarry Smith /*@ 369*a4d4d686SBarry Smith MatSNESFDMFAddNullSpace - Provides a null space that 370*a4d4d686SBarry Smith an operator is supposed to have. Since roundoff will create a 371*a4d4d686SBarry Smith small component in the null space, if you know the null space 372*a4d4d686SBarry Smith you may have it automatically removed. 373*a4d4d686SBarry Smith 374*a4d4d686SBarry Smith Collective on Mat 375*a4d4d686SBarry Smith 376*a4d4d686SBarry Smith Input Parameters: 377*a4d4d686SBarry Smith + J - the matrix-free matrix context 378*a4d4d686SBarry Smith . has_cnst - PETSC_TRUE or PETSC_FALSE, indicating if null space has constants 379*a4d4d686SBarry Smith . n - number of vectors (excluding constant vector) in null space 380*a4d4d686SBarry Smith - vecs - the vectors that span the null space (excluding the constant vector); 381*a4d4d686SBarry Smith these vectors must be orthonormal 382*a4d4d686SBarry Smith 383*a4d4d686SBarry Smith .keywords: SNES, matrix-free, null space 384*a4d4d686SBarry Smith 385*a4d4d686SBarry Smith .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(), 386*a4d4d686SBarry Smith MatSNESFDMFSetHHistory(), MatSNESFDMFResetHHistory(), 387*a4d4d686SBarry Smith MatSNESFDMFKSPMonitor(), MatSNESFDMFSetParameters() 388*a4d4d686SBarry Smith 389*a4d4d686SBarry Smith @*/ 390*a4d4d686SBarry Smith int MatSNESFDMFAddNullSpace(Mat J,int has_cnst,int n,Vec *vecs) 391*a4d4d686SBarry Smith { 392*a4d4d686SBarry Smith int ierr; 393*a4d4d686SBarry Smith MatSNESFDMFCtx ctx; 394*a4d4d686SBarry Smith MPI_Comm comm; 395*a4d4d686SBarry Smith 396*a4d4d686SBarry Smith PetscFunctionBegin; 397*a4d4d686SBarry Smith PetscObjectGetComm((PetscObject)J,&comm); 398*a4d4d686SBarry Smith 399*a4d4d686SBarry Smith ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 400*a4d4d686SBarry Smith /* no context indicates that it is not the "matrix free" matrix type */ 401*a4d4d686SBarry Smith if (!ctx) PetscFunctionReturn(0); 402*a4d4d686SBarry Smith ierr = PCNullSpaceCreate(comm,has_cnst,n,vecs,&ctx->sp); CHKERRQ(ierr); 403*a4d4d686SBarry Smith PetscFunctionReturn(0); 404*a4d4d686SBarry Smith } 405*a4d4d686SBarry Smith 406*a4d4d686SBarry Smith #undef __FUNC__ 407*a4d4d686SBarry Smith #define __FUNC__ "MatSNESFDMFSetHHistory" 408*a4d4d686SBarry Smith /*@ 409*a4d4d686SBarry Smith MatSNESFDMFSetHHistory - Sets an array to collect a history 410*a4d4d686SBarry Smith of the differencing values h computed for the matrix free product 411*a4d4d686SBarry Smith 412*a4d4d686SBarry Smith Collective on Mat 413*a4d4d686SBarry Smith 414*a4d4d686SBarry Smith Input Parameters: 415*a4d4d686SBarry Smith + J - the matrix-free matrix context 416*a4d4d686SBarry Smith . histroy - space to hold the h history 417*a4d4d686SBarry Smith - nhistory - number of entries in history, if more h are generated than 418*a4d4d686SBarry Smith nhistory the later ones are discarded 419*a4d4d686SBarry Smith 420*a4d4d686SBarry Smith Notes: 421*a4d4d686SBarry Smith Use MatSNESFDMFResetHHistory() to reset the history counter 422*a4d4d686SBarry Smith and collect a new batch of h. 423*a4d4d686SBarry Smith 424*a4d4d686SBarry Smith .keywords: SNES, matrix-free, h history, differencing history 425*a4d4d686SBarry Smith 426*a4d4d686SBarry Smith .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(), 427*a4d4d686SBarry Smith MatSNESFDMFResetHHistory(), 428*a4d4d686SBarry Smith MatSNESFDMFKSPMonitor(), MatSNESFDMFSetParameters() 429*a4d4d686SBarry Smith 430*a4d4d686SBarry Smith @*/ 431*a4d4d686SBarry Smith int MatSNESFDMFSetHHistory(Mat J,Scalar *history,int nhistory) 432*a4d4d686SBarry Smith { 433*a4d4d686SBarry Smith int ierr; 434*a4d4d686SBarry Smith MatSNESFDMFCtx ctx; 435*a4d4d686SBarry Smith 436*a4d4d686SBarry Smith PetscFunctionBegin; 437*a4d4d686SBarry Smith 438*a4d4d686SBarry Smith ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 439*a4d4d686SBarry Smith /* no context indicates that it is not the "matrix free" matrix type */ 440*a4d4d686SBarry Smith if (!ctx) PetscFunctionReturn(0); 441*a4d4d686SBarry Smith ctx->historyh = history; 442*a4d4d686SBarry Smith ctx->maxcurrenth = nhistory; 443*a4d4d686SBarry Smith ctx->currenth = 0; 444*a4d4d686SBarry Smith 445*a4d4d686SBarry Smith PetscFunctionReturn(0); 446*a4d4d686SBarry Smith } 447*a4d4d686SBarry Smith 448*a4d4d686SBarry Smith #undef __FUNC__ 449*a4d4d686SBarry Smith #define __FUNC__ "MatSNESFDMFResetHHistory" 450*a4d4d686SBarry Smith /*@ 451*a4d4d686SBarry Smith MatSNESFDMFResetHHistory - Resets the counter to zero to begin 452*a4d4d686SBarry Smith collecting a new set of differencing histories. 453*a4d4d686SBarry Smith 454*a4d4d686SBarry Smith Collective on Mat 455*a4d4d686SBarry Smith 456*a4d4d686SBarry Smith Input Parameters: 457*a4d4d686SBarry Smith . J - the matrix-free matrix context 458*a4d4d686SBarry Smith 459*a4d4d686SBarry Smith Notes: 460*a4d4d686SBarry Smith Use MatSNESFDMFSetHHistory() to create the original history counter 461*a4d4d686SBarry Smith 462*a4d4d686SBarry Smith .keywords: SNES, matrix-free, h history, differencing history 463*a4d4d686SBarry Smith 464*a4d4d686SBarry Smith .seealso: MatSNESFDMFGetH(), MatCreateSNESFDMF(), 465*a4d4d686SBarry Smith MatSNESFDMFSetHHistory(), 466*a4d4d686SBarry Smith MatSNESFDMFKSPMonitor(), MatSNESFDMFSetParameters() 467*a4d4d686SBarry Smith 468*a4d4d686SBarry Smith @*/ 469*a4d4d686SBarry Smith int MatSNESFDMFResetHHistory(Mat J) 470*a4d4d686SBarry Smith { 471*a4d4d686SBarry Smith int ierr; 472*a4d4d686SBarry Smith MatSNESFDMFCtx ctx; 473*a4d4d686SBarry Smith 474*a4d4d686SBarry Smith PetscFunctionBegin; 475*a4d4d686SBarry Smith 476*a4d4d686SBarry Smith ierr = MatShellGetContext(J,(void **)&ctx); CHKERRQ(ierr); 477*a4d4d686SBarry Smith /* no context indicates that it is not the "matrix free" matrix type */ 478*a4d4d686SBarry Smith if (!ctx) PetscFunctionReturn(0); 479*a4d4d686SBarry Smith ctx->currenth = 0; 480*a4d4d686SBarry Smith 481*a4d4d686SBarry Smith PetscFunctionReturn(0); 482*a4d4d686SBarry Smith } 483*a4d4d686SBarry Smith 484*a4d4d686SBarry Smith /* ---------------------------------------------------------------------------- */ 485*a4d4d686SBarry Smith extern int VecDot_Seq(Vec,Vec,Scalar *); 486*a4d4d686SBarry Smith extern int VecNorm_Seq(Vec,NormType,double *); 487*a4d4d686SBarry Smith 488*a4d4d686SBarry Smith typedef struct { 489*a4d4d686SBarry Smith double error_rel; /* square root of relative error in computing function */ 490*a4d4d686SBarry Smith double umin; /* minimum allowable u'a value relative to |u|_1 */ 491*a4d4d686SBarry Smith } MatSNESFDMFDefaultComputeH_struct; 492*a4d4d686SBarry Smith 493*a4d4d686SBarry Smith #undef __FUNC__ 494*a4d4d686SBarry Smith #define __FUNC__ "MatSNESFDMFDefaultComputeH" 495*a4d4d686SBarry Smith /* 496*a4d4d686SBarry Smith MatSNESFDMFDefaultComputeH - Standard PETSc code for 497*a4d4d686SBarry Smith computing h with matrix-free finite differences. 498*a4d4d686SBarry Smith 499*a4d4d686SBarry Smith */ 500*a4d4d686SBarry Smith int MatFDMFDefaultComputeH(MatSNESFDMFCtx ctx,Vec U,Vec a,Scalar *h) 501*a4d4d686SBarry Smith { 502*a4d4d686SBarry Smith MPI_Comm comm; 503*a4d4d686SBarry Smith double ovalues[3],norm, sum, umin; 504*a4d4d686SBarry Smith Scalar dot; 505*a4d4d686SBarry Smith int ierr; 506*a4d4d686SBarry Smith #if !defined(USE_PETSC_COMPLEX) 507*a4d4d686SBarry Smith double values[3]; 508*a4d4d686SBarry Smith #endif 509*a4d4d686SBarry Smith 510*a4d4d686SBarry Smith umin = ctx->umin; 511*a4d4d686SBarry Smith 512c481317fSBarry Smith 513c481317fSBarry Smith /* 514eb9086c3SLois Curfman McInnes ierr = VecDot(U,a,&dot); CHKERRQ(ierr); 515eb9086c3SLois Curfman McInnes ierr = VecNorm(a,NORM_1,&sum); CHKERRQ(ierr); 516eb9086c3SLois Curfman McInnes ierr = VecNorm(a,NORM_2,&norm); CHKERRQ(ierr); 517c481317fSBarry Smith */ 518c481317fSBarry Smith 519c481317fSBarry Smith /* 52086f5173aSBarry Smith Call the Seq Vector routines and then do a single reduction 52186f5173aSBarry Smith to reduce the number of communications required 522c481317fSBarry Smith */ 523c481317fSBarry Smith 5243a40ed3dSBarry Smith #if !defined(USE_PETSC_COMPLEX) 525c481317fSBarry Smith PLogEventBegin(VEC_Dot,U,a,0,0); 526c481317fSBarry Smith ierr = VecDot_Seq(U,a,ovalues); CHKERRQ(ierr); 527c481317fSBarry Smith PLogEventEnd(VEC_Dot,U,a,0,0); 528c481317fSBarry Smith PLogEventBegin(VEC_Norm,a,0,0,0); 529c481317fSBarry Smith ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 530c481317fSBarry Smith ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 531c481317fSBarry Smith ovalues[2] = ovalues[2]*ovalues[2]; 532005c665bSBarry Smith PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 533ca161407SBarry Smith ierr = MPI_Allreduce(ovalues,values,3,MPI_DOUBLE,MPI_SUM,comm );CHKERRQ(ierr); 534005c665bSBarry Smith PLogEventBarrierEnd(VEC_NormBarrier,0,0,0,0,comm); 535c481317fSBarry Smith dot = values[0]; sum = values[1]; norm = sqrt(values[2]); 536c481317fSBarry Smith PLogEventEnd(VEC_Norm,a,0,0,0); 53786f5173aSBarry Smith #else 53886f5173aSBarry Smith { 53986f5173aSBarry Smith Scalar cvalues[3],covalues[3]; 54086f5173aSBarry Smith 54186f5173aSBarry Smith PLogEventBegin(VEC_Dot,U,a,0,0); 54286f5173aSBarry Smith ierr = VecDot_Seq(U,a,covalues); CHKERRQ(ierr); 54386f5173aSBarry Smith PLogEventEnd(VEC_Dot,U,a,0,0); 54486f5173aSBarry Smith PLogEventBegin(VEC_Norm,a,0,0,0); 54586f5173aSBarry Smith ierr = VecNorm_Seq(a,NORM_1,ovalues+1); CHKERRQ(ierr); 54686f5173aSBarry Smith ierr = VecNorm_Seq(a,NORM_2,ovalues+2); CHKERRQ(ierr); 54786f5173aSBarry Smith covalues[1] = ovalues[1]; 54886f5173aSBarry Smith covalues[2] = ovalues[2]*ovalues[2]; 549005c665bSBarry Smith PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 550ca161407SBarry Smith ierr = MPI_Allreduce(covalues,cvalues,6,MPI_DOUBLE,MPI_SUM,comm);CHKERRQ(ierr); 551005c665bSBarry Smith PLogEventBarrierBegin(VEC_NormBarrier,0,0,0,0,comm); 55286f5173aSBarry Smith dot = cvalues[0]; sum = PetscReal(cvalues[1]); norm = sqrt(PetscReal(cvalues[2])); 55386f5173aSBarry Smith PLogEventEnd(VEC_Norm,a,0,0,0); 55486f5173aSBarry Smith } 55586f5173aSBarry Smith #endif 556c481317fSBarry Smith 557eb9086c3SLois Curfman McInnes /* Safeguard for step sizes too small */ 558edd2f0e1SBarry Smith if (sum == 0.0) {dot = 1.0; norm = 1.0;} 5593a40ed3dSBarry Smith #if defined(USE_PETSC_COMPLEX) 560e20fef11SSatish Balay else if (PetscAbsScalar(dot) < umin*sum && PetscReal(dot) >= 0.0) dot = umin*sum; 561e20fef11SSatish Balay else if (PetscAbsScalar(dot) < 0.0 && PetscReal(dot) > -umin*sum) dot = -umin*sum; 56219a167f6SBarry Smith #else 563eb9086c3SLois Curfman McInnes else if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 564eb9086c3SLois Curfman McInnes else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 56519a167f6SBarry Smith #endif 566*a4d4d686SBarry Smith *h = ctx->error_rel*dot/(norm*norm); 5673a40ed3dSBarry Smith PetscFunctionReturn(0); 56839e2f89bSBarry Smith } 569