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