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