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 Optional functions: 17e884886fSBarry Smith ------------------- 18e884886fSBarry Smith MatMFFDView_ - prints information about the parameters being used. 19e884886fSBarry Smith This is called when SNESView() or -snes_view is used. 20e884886fSBarry Smith 21e884886fSBarry Smith MatMFFDSetFromOptions_ - checks the options database for options that 22e884886fSBarry Smith apply to this method. 23e884886fSBarry Smith 24e884886fSBarry Smith MatMFFDDestroy_ - frees any space allocated by the routines above 25e884886fSBarry Smith 26e884886fSBarry Smith */ 27e884886fSBarry Smith 28e884886fSBarry Smith /* 29e884886fSBarry Smith This include file defines the data structure MatMFFD that 30e884886fSBarry Smith includes information about the computation of h. It is shared by 31e884886fSBarry Smith all implementations that people provide 32e884886fSBarry Smith */ 33af0996ceSBarry Smith #include <petsc/private/matimpl.h> 34c6db04a5SJed Brown #include <../src/mat/impls/mffd/mffdimpl.h> /*I "petscmat.h" I*/ 35e884886fSBarry Smith 36e884886fSBarry Smith /* 37e884886fSBarry Smith The method has one parameter that is used to 38e884886fSBarry Smith "cutoff" very small values. This is stored in a data structure 39e884886fSBarry Smith that is only visible to this file. If your method has no parameters 40e884886fSBarry Smith it can omit this, if it has several simply reorganize the data structure. 41e884886fSBarry Smith The data structure is "hung-off" the MatMFFD data structure in 42e884886fSBarry Smith the void *hctx; field. 43e884886fSBarry Smith */ 44e884886fSBarry Smith typedef struct { 45e884886fSBarry Smith PetscReal umin; /* minimum allowable u'a value relative to |u|_1 */ 46e884886fSBarry Smith } MatMFFD_DS; 47e884886fSBarry Smith 48d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDCompute_DS(MatMFFD ctx, Vec U, Vec a, PetscScalar *h, PetscBool *zeroa) 49d71ae5a4SJacob Faibussowitsch { 50e884886fSBarry Smith MatMFFD_DS *hctx = (MatMFFD_DS *)ctx->hctx; 51e884886fSBarry Smith PetscReal nrm, sum, umin = hctx->umin; 52e884886fSBarry Smith PetscScalar dot; 53e884886fSBarry Smith 54e884886fSBarry Smith PetscFunctionBegin; 55e884886fSBarry Smith if (!(ctx->count % ctx->recomputeperiod)) { 56e884886fSBarry Smith /* 57e884886fSBarry Smith This algorithm requires 2 norms and 1 inner product. Rather than 58e884886fSBarry Smith use directly the VecNorm() and VecDot() routines (and thus have 59e884886fSBarry Smith three separate collective operations, we use the VecxxxBegin/End() routines 60e884886fSBarry Smith */ 619566063dSJacob Faibussowitsch PetscCall(VecDotBegin(U, a, &dot)); 629566063dSJacob Faibussowitsch PetscCall(VecNormBegin(a, NORM_1, &sum)); 639566063dSJacob Faibussowitsch PetscCall(VecNormBegin(a, NORM_2, &nrm)); 649566063dSJacob Faibussowitsch PetscCall(VecDotEnd(U, a, &dot)); 659566063dSJacob Faibussowitsch PetscCall(VecNormEnd(a, NORM_1, &sum)); 669566063dSJacob Faibussowitsch PetscCall(VecNormEnd(a, NORM_2, &nrm)); 67e884886fSBarry Smith 68e884886fSBarry Smith if (nrm == 0.0) { 69e884886fSBarry Smith *zeroa = PETSC_TRUE; 703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 71e884886fSBarry Smith } 72e884886fSBarry Smith *zeroa = PETSC_FALSE; 73e884886fSBarry Smith 74e884886fSBarry Smith /* 75e884886fSBarry Smith Safeguard for step sizes that are "too small" 76e884886fSBarry Smith */ 77e884886fSBarry Smith if (PetscAbsScalar(dot) < umin * sum && PetscRealPart(dot) >= 0.0) dot = umin * sum; 78e884886fSBarry Smith else if (PetscAbsScalar(dot) < 0.0 && PetscRealPart(dot) > -umin * sum) dot = -umin * sum; 79e884886fSBarry Smith *h = ctx->error_rel * dot / (nrm * nrm); 80aed4548fSBarry Smith PetscCheck(!PetscIsInfOrNanScalar(*h), 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); 81e884886fSBarry Smith } else { 82e884886fSBarry Smith *h = ctx->currenth; 83e884886fSBarry Smith } 84e884886fSBarry Smith ctx->count++; 853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 86e884886fSBarry Smith } 87e884886fSBarry Smith 88d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDView_DS(MatMFFD ctx, PetscViewer viewer) 89d71ae5a4SJacob Faibussowitsch { 90e884886fSBarry Smith MatMFFD_DS *hctx = (MatMFFD_DS *)ctx->hctx; 91ace3abfcSBarry Smith PetscBool iascii; 92e884886fSBarry Smith 93e884886fSBarry Smith PetscFunctionBegin; 94e884886fSBarry Smith /* 95e884886fSBarry Smith Currently this only handles the ascii file viewers, others 96e884886fSBarry Smith could be added, but for this type of object other viewers 97e884886fSBarry Smith make less sense 98e884886fSBarry Smith */ 999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 10048a46eb9SPierre Jolivet if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, " umin=%g (minimum iterate parameter)\n", (double)hctx->umin)); 1013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 102e884886fSBarry Smith } 103e884886fSBarry Smith 104d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDSetFromOptions_DS(MatMFFD ctx, PetscOptionItems *PetscOptionsObject) 105d71ae5a4SJacob Faibussowitsch { 106e884886fSBarry Smith MatMFFD_DS *hctx = (MatMFFD_DS *)ctx->hctx; 107e884886fSBarry Smith 108e884886fSBarry Smith PetscFunctionBegin; 109*01c1178eSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Finite difference matrix-free parameters"); 1109566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_mffd_umin", "umin", "MatMFFDDSSetUmin", hctx->umin, &hctx->umin, NULL)); 111d0609cedSBarry Smith PetscOptionsHeadEnd(); 1123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 113e884886fSBarry Smith } 114e884886fSBarry Smith 115d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMFFDDestroy_DS(MatMFFD ctx) 116d71ae5a4SJacob Faibussowitsch { 117e884886fSBarry Smith PetscFunctionBegin; 1189566063dSJacob Faibussowitsch PetscCall(PetscFree(ctx->hctx)); 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 120e884886fSBarry Smith } 121e884886fSBarry Smith 122e884886fSBarry Smith /* 123e884886fSBarry Smith The following two routines use the PetscObjectCompose() and PetscObjectQuery() 124e884886fSBarry Smith mechanism to allow the user to change the Umin parameter used in this method. 125e884886fSBarry Smith */ 126d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDDSSetUmin_DS(Mat mat, PetscReal umin) 127d71ae5a4SJacob Faibussowitsch { 128440c47e6SSatish Balay MatMFFD ctx = NULL; 129e884886fSBarry Smith MatMFFD_DS *hctx; 130e884886fSBarry Smith 131e884886fSBarry Smith PetscFunctionBegin; 1329566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat, &ctx)); 13328b400f6SJacob Faibussowitsch PetscCheck(ctx, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "MatMFFDDSSetUmin() attached to non-shell matrix"); 134e884886fSBarry Smith hctx = (MatMFFD_DS *)ctx->hctx; 135e884886fSBarry Smith hctx->umin = umin; 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 137e884886fSBarry Smith } 138e884886fSBarry Smith 139e884886fSBarry Smith /*@ 140e884886fSBarry Smith MatMFFDDSSetUmin - Sets the "umin" parameter used by the 141e884886fSBarry Smith PETSc routine for computing the differencing parameter, h, which is used 14211a5261eSBarry Smith for matrix-free Jacobian-vector products for a `MATMFFD` matrix. 143e884886fSBarry Smith 144e884886fSBarry Smith Input Parameters: 14511a5261eSBarry Smith + A - the `MATMFFD` matrix 146e884886fSBarry Smith - umin - the parameter 147e884886fSBarry Smith 148e884886fSBarry Smith Level: advanced 149e884886fSBarry Smith 15011a5261eSBarry Smith Note: 15111a5261eSBarry Smith See the manual page for `MatCreateSNESMF()` for a complete description of the 152e884886fSBarry Smith algorithm used to compute h. 153e884886fSBarry Smith 15411a5261eSBarry Smith .seealso: `MATMFFD`, `MatMFFDSetFunctionError()`, `MatCreateSNESMF()` 155e884886fSBarry Smith @*/ 156d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMFFDDSSetUmin(Mat A, PetscReal umin) 157d71ae5a4SJacob Faibussowitsch { 158e884886fSBarry Smith PetscFunctionBegin; 1590700a824SBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 160cac4c232SBarry Smith PetscTryMethod(A, "MatMFFDDSSetUmin_C", (Mat, PetscReal), (A, umin)); 1613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 162e884886fSBarry Smith } 163e884886fSBarry Smith 164e884886fSBarry Smith /*MC 16511a5261eSBarry Smith MATMFFD_DS - algorithm for compute the "h" used in the finite difference matrix-free matrix vector product, `MatMult()`. 166e884886fSBarry Smith 167e884886fSBarry Smith Options Database Keys: 16811a5261eSBarry Smith . -mat_mffd_umin <umin> - see `MatMFFDDSSetUmin()` 169e884886fSBarry Smith 170e884886fSBarry Smith Level: intermediate 171e884886fSBarry Smith 17295452b02SPatrick Sanan Notes: 17395452b02SPatrick Sanan Requires 2 norms and 1 inner product, but they are computed together 17411a5261eSBarry Smith so only one parallel collective operation is needed. See `MATMFFD_WP` for a method 17511a5261eSBarry Smith (with `KSPGMRES`) that requires NO collective operations. 176e884886fSBarry Smith 177e884886fSBarry Smith Formula used: 178e884886fSBarry Smith F'(u)*a = [F(u+h*a) - F(u)]/h where 179e884886fSBarry Smith h = error_rel*u'a/||a||^2 if |u'a| > umin*||a||_{1} 180e884886fSBarry Smith = error_rel*umin*sign(u'a)*||a||_{1}/||a||^2 otherwise 181e884886fSBarry Smith where 182e884886fSBarry Smith error_rel = square root of relative error in function evaluation 183e884886fSBarry Smith umin = minimum iterate parameter 184e884886fSBarry Smith 18511a5261eSBarry Smith References: 18611a5261eSBarry Smith . * - Dennis and Schnabel, "Numerical Methods for Unconstrained Optimization and Nonlinear Equations" 187e884886fSBarry Smith 18811a5261eSBarry Smith .seealso: `MATMFFD`, `MATMFFD_WP`, `MatCreateMFFD()`, `MatCreateSNESMF()`, `MATMFFD_WP`, `MatMFFDDSSetUmin()` 189e884886fSBarry Smith M*/ 190d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode MatCreateMFFD_DS(MatMFFD ctx) 191d71ae5a4SJacob Faibussowitsch { 192e884886fSBarry Smith MatMFFD_DS *hctx; 193e884886fSBarry Smith 194e884886fSBarry Smith PetscFunctionBegin; 195e884886fSBarry Smith /* allocate my own private data structure */ 1964dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&hctx)); 197e884886fSBarry Smith ctx->hctx = (void *)hctx; 198e884886fSBarry Smith /* set a default for my parameter */ 199e884886fSBarry Smith hctx->umin = 1.e-6; 200e884886fSBarry Smith 201e884886fSBarry Smith /* set the functions I am providing */ 202e884886fSBarry Smith ctx->ops->compute = MatMFFDCompute_DS; 203e884886fSBarry Smith ctx->ops->destroy = MatMFFDDestroy_DS; 204e884886fSBarry Smith ctx->ops->view = MatMFFDView_DS; 205e884886fSBarry Smith ctx->ops->setfromoptions = MatMFFDSetFromOptions_DS; 206e884886fSBarry Smith 2079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)ctx->mat, "MatMFFDDSSetUmin_C", MatMFFDDSSetUmin_DS)); 2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 209e884886fSBarry Smith } 210