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