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