1*e884886fSBarry Smith #define PETSCMAT_DLL 2*e884886fSBarry Smith 3*e884886fSBarry Smith /* 4*e884886fSBarry Smith Implements the DS PETSc approach for computing the h 5*e884886fSBarry Smith parameter used with the finite difference based matrix-free 6*e884886fSBarry Smith Jacobian-vector products. 7*e884886fSBarry Smith 8*e884886fSBarry Smith To make your own: clone this file and modify for your needs. 9*e884886fSBarry Smith 10*e884886fSBarry Smith Mandatory functions: 11*e884886fSBarry Smith ------------------- 12*e884886fSBarry Smith MatMFFDCompute_ - for a given point and direction computes h 13*e884886fSBarry Smith 14*e884886fSBarry Smith MatMFFDCreate_ - fills in the MatMFFD data structure 15*e884886fSBarry Smith for this particular implementation 16*e884886fSBarry Smith 17*e884886fSBarry Smith 18*e884886fSBarry Smith Optional functions: 19*e884886fSBarry Smith ------------------- 20*e884886fSBarry Smith MatMFFDView_ - prints information about the parameters being used. 21*e884886fSBarry Smith This is called when SNESView() or -snes_view is used. 22*e884886fSBarry Smith 23*e884886fSBarry Smith MatMFFDSetFromOptions_ - checks the options database for options that 24*e884886fSBarry Smith apply to this method. 25*e884886fSBarry Smith 26*e884886fSBarry Smith MatMFFDDestroy_ - frees any space allocated by the routines above 27*e884886fSBarry Smith 28*e884886fSBarry Smith */ 29*e884886fSBarry Smith 30*e884886fSBarry Smith /* 31*e884886fSBarry Smith This include file defines the data structure MatMFFD that 32*e884886fSBarry Smith includes information about the computation of h. It is shared by 33*e884886fSBarry Smith all implementations that people provide 34*e884886fSBarry Smith */ 35*e884886fSBarry Smith #include "include/private/matimpl.h" 36*e884886fSBarry Smith #include "src/mat/impls/mffd/mffdimpl.h" /*I "petscmat.h" I*/ 37*e884886fSBarry Smith 38*e884886fSBarry Smith /* 39*e884886fSBarry Smith The method has one parameter that is used to 40*e884886fSBarry Smith "cutoff" very small values. This is stored in a data structure 41*e884886fSBarry Smith that is only visible to this file. If your method has no parameters 42*e884886fSBarry Smith it can omit this, if it has several simply reorganize the data structure. 43*e884886fSBarry Smith The data structure is "hung-off" the MatMFFD data structure in 44*e884886fSBarry Smith the void *hctx; field. 45*e884886fSBarry Smith */ 46*e884886fSBarry Smith typedef struct { 47*e884886fSBarry Smith PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */ 48*e884886fSBarry Smith } MatMFFD_DS; 49*e884886fSBarry Smith 50*e884886fSBarry Smith #undef __FUNCT__ 51*e884886fSBarry Smith #define __FUNCT__ "MatMFFDCompute_DS" 52*e884886fSBarry Smith /* 53*e884886fSBarry Smith MatMFFDCompute_DS - Standard PETSc code for computing the 54*e884886fSBarry Smith differencing paramter (h) for use with matrix-free finite differences. 55*e884886fSBarry Smith 56*e884886fSBarry Smith Input Parameters: 57*e884886fSBarry Smith + ctx - the matrix free context 58*e884886fSBarry Smith . U - the location at which you want the Jacobian 59*e884886fSBarry Smith - a - the direction you want the derivative 60*e884886fSBarry Smith 61*e884886fSBarry Smith 62*e884886fSBarry Smith Output Parameter: 63*e884886fSBarry Smith . h - the scale computed 64*e884886fSBarry Smith 65*e884886fSBarry Smith */ 66*e884886fSBarry Smith static PetscErrorCode MatMFFDCompute_DS(MatMFFD ctx,Vec U,Vec a,PetscScalar *h,PetscTruth *zeroa) 67*e884886fSBarry Smith { 68*e884886fSBarry Smith MatMFFD_DS *hctx = (MatMFFD_DS*)ctx->hctx; 69*e884886fSBarry Smith PetscReal nrm,sum,umin = hctx->umin; 70*e884886fSBarry Smith PetscScalar dot; 71*e884886fSBarry Smith PetscErrorCode ierr; 72*e884886fSBarry Smith 73*e884886fSBarry Smith PetscFunctionBegin; 74*e884886fSBarry Smith if (!(ctx->count % ctx->recomputeperiod)) { 75*e884886fSBarry Smith /* 76*e884886fSBarry Smith This algorithm requires 2 norms and 1 inner product. Rather than 77*e884886fSBarry Smith use directly the VecNorm() and VecDot() routines (and thus have 78*e884886fSBarry Smith three separate collective operations, we use the VecxxxBegin/End() routines 79*e884886fSBarry Smith */ 80*e884886fSBarry Smith ierr = VecDotBegin(U,a,&dot);CHKERRQ(ierr); 81*e884886fSBarry Smith ierr = VecNormBegin(a,NORM_1,&sum);CHKERRQ(ierr); 82*e884886fSBarry Smith ierr = VecNormBegin(a,NORM_2,&nrm);CHKERRQ(ierr); 83*e884886fSBarry Smith ierr = VecDotEnd(U,a,&dot);CHKERRQ(ierr); 84*e884886fSBarry Smith ierr = VecNormEnd(a,NORM_1,&sum);CHKERRQ(ierr); 85*e884886fSBarry Smith ierr = VecNormEnd(a,NORM_2,&nrm);CHKERRQ(ierr); 86*e884886fSBarry Smith 87*e884886fSBarry Smith if (nrm == 0.0) { 88*e884886fSBarry Smith *zeroa = PETSC_TRUE; 89*e884886fSBarry Smith PetscFunctionReturn(0); 90*e884886fSBarry Smith } 91*e884886fSBarry Smith *zeroa = PETSC_FALSE; 92*e884886fSBarry Smith 93*e884886fSBarry Smith /* 94*e884886fSBarry Smith Safeguard for step sizes that are "too small" 95*e884886fSBarry Smith */ 96*e884886fSBarry Smith #if defined(PETSC_USE_COMPLEX) 97*e884886fSBarry Smith if (PetscAbsScalar(dot) < umin*sum && PetscRealPart(dot) >= 0.0) dot = umin*sum; 98*e884886fSBarry Smith else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin*sum) dot = -umin*sum; 99*e884886fSBarry Smith #else 100*e884886fSBarry Smith if (dot < umin*sum && dot >= 0.0) dot = umin*sum; 101*e884886fSBarry Smith else if (dot < 0.0 && dot > -umin*sum) dot = -umin*sum; 102*e884886fSBarry Smith #endif 103*e884886fSBarry Smith *h = ctx->error_rel*dot/(nrm*nrm); 104*e884886fSBarry Smith } else { 105*e884886fSBarry Smith *h = ctx->currenth; 106*e884886fSBarry Smith } 107*e884886fSBarry Smith if (*h != *h) SETERRQ3(PETSC_ERR_PLIB,"Differencing parameter is not a number sum = %G dot = %G norm = %G",sum,PetscRealPart(dot),nrm); 108*e884886fSBarry Smith ctx->count++; 109*e884886fSBarry Smith PetscFunctionReturn(0); 110*e884886fSBarry Smith } 111*e884886fSBarry Smith 112*e884886fSBarry Smith #undef __FUNCT__ 113*e884886fSBarry Smith #define __FUNCT__ "MatMFFDView_DS" 114*e884886fSBarry Smith /* 115*e884886fSBarry Smith MatMFFDView_DS - Prints information about this particular 116*e884886fSBarry Smith method for computing h. Note that this does not print the general 117*e884886fSBarry Smith information about the matrix-free method, as such info is printed 118*e884886fSBarry Smith by the calling routine. 119*e884886fSBarry Smith 120*e884886fSBarry Smith Input Parameters: 121*e884886fSBarry Smith + ctx - the matrix free context 122*e884886fSBarry Smith - viewer - the PETSc viewer 123*e884886fSBarry Smith */ 124*e884886fSBarry Smith static PetscErrorCode MatMFFDView_DS(MatMFFD ctx,PetscViewer viewer) 125*e884886fSBarry Smith { 126*e884886fSBarry Smith MatMFFD_DS *hctx = (MatMFFD_DS *)ctx->hctx; 127*e884886fSBarry Smith PetscErrorCode ierr; 128*e884886fSBarry Smith PetscTruth iascii; 129*e884886fSBarry Smith 130*e884886fSBarry Smith PetscFunctionBegin; 131*e884886fSBarry Smith /* 132*e884886fSBarry Smith Currently this only handles the ascii file viewers, others 133*e884886fSBarry Smith could be added, but for this type of object other viewers 134*e884886fSBarry Smith make less sense 135*e884886fSBarry Smith */ 136*e884886fSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 137*e884886fSBarry Smith if (iascii) { 138*e884886fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," umin=%G (minimum iterate parameter)\n",hctx->umin);CHKERRQ(ierr); 139*e884886fSBarry Smith } else { 140*e884886fSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for this SNES matrix free matrix",((PetscObject)viewer)->type_name); 141*e884886fSBarry Smith } 142*e884886fSBarry Smith PetscFunctionReturn(0); 143*e884886fSBarry Smith } 144*e884886fSBarry Smith 145*e884886fSBarry Smith #undef __FUNCT__ 146*e884886fSBarry Smith #define __FUNCT__ "MatMFFDSetFromOptions_DS" 147*e884886fSBarry Smith /* 148*e884886fSBarry Smith MatMFFDSetFromOptions_DS - Looks in the options database for 149*e884886fSBarry Smith any options appropriate for this method. 150*e884886fSBarry Smith 151*e884886fSBarry Smith Input Parameter: 152*e884886fSBarry Smith . ctx - the matrix free context 153*e884886fSBarry Smith 154*e884886fSBarry Smith */ 155*e884886fSBarry Smith static PetscErrorCode MatMFFDSetFromOptions_DS(MatMFFD ctx) 156*e884886fSBarry Smith { 157*e884886fSBarry Smith PetscErrorCode ierr; 158*e884886fSBarry Smith MatMFFD_DS *hctx = (MatMFFD_DS*)ctx->hctx; 159*e884886fSBarry Smith 160*e884886fSBarry Smith PetscFunctionBegin; 161*e884886fSBarry Smith ierr = PetscOptionsHead("Finite difference matrix free parameters");CHKERRQ(ierr); 162*e884886fSBarry Smith ierr = PetscOptionsReal("-mat_mffd_umin","umin","MatMFFDDSSetUmin",hctx->umin,&hctx->umin,0);CHKERRQ(ierr); 163*e884886fSBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 164*e884886fSBarry Smith PetscFunctionReturn(0); 165*e884886fSBarry Smith } 166*e884886fSBarry Smith 167*e884886fSBarry Smith #undef __FUNCT__ 168*e884886fSBarry Smith #define __FUNCT__ "MatMFFDDestroy_DS" 169*e884886fSBarry Smith /* 170*e884886fSBarry Smith MatMFFDDestroy_DS - Frees the space allocated by 171*e884886fSBarry Smith MatMFFDCreate_DS(). 172*e884886fSBarry Smith 173*e884886fSBarry Smith Input Parameter: 174*e884886fSBarry Smith . ctx - the matrix free context 175*e884886fSBarry Smith 176*e884886fSBarry Smith Notes: 177*e884886fSBarry Smith Does not free the ctx, that is handled by the calling routine 178*e884886fSBarry Smith */ 179*e884886fSBarry Smith static PetscErrorCode MatMFFDDestroy_DS(MatMFFD ctx) 180*e884886fSBarry Smith { 181*e884886fSBarry Smith PetscErrorCode ierr; 182*e884886fSBarry Smith 183*e884886fSBarry Smith PetscFunctionBegin; 184*e884886fSBarry Smith ierr = PetscFree(ctx->hctx);CHKERRQ(ierr); 185*e884886fSBarry Smith PetscFunctionReturn(0); 186*e884886fSBarry Smith } 187*e884886fSBarry Smith 188*e884886fSBarry Smith EXTERN_C_BEGIN 189*e884886fSBarry Smith #undef __FUNCT__ 190*e884886fSBarry Smith #define __FUNCT__ "MatMFFDDSSetUmin_Private" 191*e884886fSBarry Smith /* 192*e884886fSBarry Smith The following two routines use the PetscObjectCompose() and PetscObjectQuery() 193*e884886fSBarry Smith mechanism to allow the user to change the Umin parameter used in this method. 194*e884886fSBarry Smith */ 195*e884886fSBarry Smith PetscErrorCode MatMFFDDSSetUmin_Private(Mat mat,PetscReal umin) 196*e884886fSBarry Smith { 197*e884886fSBarry Smith MatMFFD ctx = (MatMFFD)mat->data; 198*e884886fSBarry Smith MatMFFD_DS *hctx; 199*e884886fSBarry Smith 200*e884886fSBarry Smith PetscFunctionBegin; 201*e884886fSBarry Smith if (!ctx) { 202*e884886fSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"MatMFFDDSSetUmin() attached to non-shell matrix"); 203*e884886fSBarry Smith } 204*e884886fSBarry Smith hctx = (MatMFFD_DS*)ctx->hctx; 205*e884886fSBarry Smith hctx->umin = umin; 206*e884886fSBarry Smith PetscFunctionReturn(0); 207*e884886fSBarry Smith } 208*e884886fSBarry Smith EXTERN_C_END 209*e884886fSBarry Smith 210*e884886fSBarry Smith #undef __FUNCT__ 211*e884886fSBarry Smith #define __FUNCT__ "MatMFFDDSSetUmin" 212*e884886fSBarry Smith /*@ 213*e884886fSBarry Smith MatMFFDDSSetUmin - Sets the "umin" parameter used by the 214*e884886fSBarry Smith PETSc routine for computing the differencing parameter, h, which is used 215*e884886fSBarry Smith for matrix-free Jacobian-vector products. 216*e884886fSBarry Smith 217*e884886fSBarry Smith Input Parameters: 218*e884886fSBarry Smith + A - the matrix created with MatCreateSNESMF() 219*e884886fSBarry Smith - umin - the parameter 220*e884886fSBarry Smith 221*e884886fSBarry Smith Level: advanced 222*e884886fSBarry Smith 223*e884886fSBarry Smith Notes: 224*e884886fSBarry Smith See the manual page for MatCreateSNESMF() for a complete description of the 225*e884886fSBarry Smith algorithm used to compute h. 226*e884886fSBarry Smith 227*e884886fSBarry Smith .seealso: MatMFFDSetFunctionError(), MatCreateSNESMF() 228*e884886fSBarry Smith 229*e884886fSBarry Smith @*/ 230*e884886fSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDDSSetUmin(Mat A,PetscReal umin) 231*e884886fSBarry Smith { 232*e884886fSBarry Smith PetscErrorCode ierr,(*f)(Mat,PetscReal); 233*e884886fSBarry Smith 234*e884886fSBarry Smith PetscFunctionBegin; 235*e884886fSBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE,1); 236*e884886fSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)A,"MatMFFDDSSetUmin_C",(void (**)(void))&f);CHKERRQ(ierr); 237*e884886fSBarry Smith if (f) { 238*e884886fSBarry Smith ierr = (*f)(A,umin);CHKERRQ(ierr); 239*e884886fSBarry Smith } 240*e884886fSBarry Smith PetscFunctionReturn(0); 241*e884886fSBarry Smith } 242*e884886fSBarry Smith 243*e884886fSBarry Smith /*MC 244*e884886fSBarry Smith MATMFFD_DS - the code for compute the "h" used in the finite difference 245*e884886fSBarry Smith matrix-free matrix vector product. This code 246*e884886fSBarry Smith implements the strategy in Dennis and Schnabel, "Numerical Methods for Unconstrained 247*e884886fSBarry Smith Optimization and Nonlinear Equations". 248*e884886fSBarry Smith 249*e884886fSBarry Smith Options Database Keys: 250*e884886fSBarry Smith . -mat_mffd_umin <umin> see MatMFFDDSSetUmin() 251*e884886fSBarry Smith 252*e884886fSBarry Smith Level: intermediate 253*e884886fSBarry Smith 254*e884886fSBarry Smith Notes: Requires 2 norms and 1 inner product, but they are computed together 255*e884886fSBarry Smith so only one parallel collective operation is needed. See MATMFFD_WP for a method 256*e884886fSBarry Smith (with GMRES) that requires NO collective operations. 257*e884886fSBarry Smith 258*e884886fSBarry Smith Formula used: 259*e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 260*e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 261*e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 262*e884886fSBarry Smith where 263*e884886fSBarry Smith error_rel = square root of relative error in function evaluation 264*e884886fSBarry Smith umin = minimum iterate parameter 265*e884886fSBarry Smith 266*e884886fSBarry Smith .seealso: MATMFFD, MatCreateMFFD(), MatCreateSNESMF(), MATMFFD_WP, MatMFFDDSSetUmin() 267*e884886fSBarry Smith 268*e884886fSBarry Smith M*/ 269*e884886fSBarry Smith EXTERN_C_BEGIN 270*e884886fSBarry Smith #undef __FUNCT__ 271*e884886fSBarry Smith #define __FUNCT__ "MatMFFDCreate_DS" 272*e884886fSBarry Smith PetscErrorCode PETSCMAT_DLLEXPORT MatMFFDCreate_DS(MatMFFD ctx) 273*e884886fSBarry Smith { 274*e884886fSBarry Smith MatMFFD_DS *hctx; 275*e884886fSBarry Smith PetscErrorCode ierr; 276*e884886fSBarry Smith 277*e884886fSBarry Smith PetscFunctionBegin; 278*e884886fSBarry Smith 279*e884886fSBarry Smith /* allocate my own private data structure */ 280*e884886fSBarry Smith ierr = PetscNew(MatMFFD_DS,&hctx);CHKERRQ(ierr); 281*e884886fSBarry Smith ctx->hctx = (void*)hctx; 282*e884886fSBarry Smith /* set a default for my parameter */ 283*e884886fSBarry Smith hctx->umin = 1.e-6; 284*e884886fSBarry Smith 285*e884886fSBarry Smith /* set the functions I am providing */ 286*e884886fSBarry Smith ctx->ops->compute = MatMFFDCompute_DS; 287*e884886fSBarry Smith ctx->ops->destroy = MatMFFDDestroy_DS; 288*e884886fSBarry Smith ctx->ops->view = MatMFFDView_DS; 289*e884886fSBarry Smith ctx->ops->setfromoptions = MatMFFDSetFromOptions_DS; 290*e884886fSBarry Smith 291*e884886fSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)ctx->mat,"MatMFFDDSSetUmin_C", 292*e884886fSBarry Smith "MatMFFDDSSetUmin_Private", 293*e884886fSBarry Smith MatMFFDDSSetUmin_Private);CHKERRQ(ierr); 294*e884886fSBarry Smith PetscFunctionReturn(0); 295*e884886fSBarry Smith } 296*e884886fSBarry Smith EXTERN_C_END 297*e884886fSBarry Smith 298*e884886fSBarry Smith 299*e884886fSBarry Smith 300*e884886fSBarry Smith 301*e884886fSBarry Smith 302*e884886fSBarry Smith 303*e884886fSBarry Smith 304