1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> 22a4ee8f2SPeter Brune 3075cc632SBarry Smith /*MC 4*f6dfbefdSBarry Smith SNESObjectiveFunction - functional form used to convey an objective function to the nonlinear solver 52a4ee8f2SPeter Brune 6075cc632SBarry Smith Synopsis: 7aaa7dc30SBarry Smith #include <petscsnes.h> 8075cc632SBarry Smith SNESObjectiveFunction(SNES snes,Vec x,PetscReal *obj,void *ctx); 92a4ee8f2SPeter Brune 102a4ee8f2SPeter Brune Input Parameters: 112a4ee8f2SPeter Brune + snes - the SNES context 122a4ee8f2SPeter Brune . X - solution 132a4ee8f2SPeter Brune . obj - real to hold the objective value 142a4ee8f2SPeter Brune - ctx - optional user-defined objective context 152a4ee8f2SPeter Brune 16878cb397SSatish Balay Level: advanced 17878cb397SSatish Balay 18*f6dfbefdSBarry Smith .seealso: `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESSetObjective()`, `SNESGetObjective()`, `SNESJacobianFunction`, `SNESFunction` 19075cc632SBarry Smith M*/ 20075cc632SBarry Smith 21075cc632SBarry Smith /*@C 22*f6dfbefdSBarry Smith SNESSetObjective - Sets the objective function minimized by some of the `SNES` linesearch methods. 23075cc632SBarry Smith 24*f6dfbefdSBarry Smith Logically Collective on snes 25075cc632SBarry Smith 26075cc632SBarry Smith Input Parameters: 27*f6dfbefdSBarry Smith + snes - the `SNES` context 28*f6dfbefdSBarry Smith . obj - objective evaluation routine; see `SNESObjectiveFunction` for details 29075cc632SBarry Smith - ctx - [optional] user-defined context for private data for the 300298fd71SBarry Smith function evaluation routine (may be NULL) 31075cc632SBarry Smith 322b26979fSBarry Smith Level: intermediate 332a4ee8f2SPeter Brune 34*f6dfbefdSBarry Smith Note: 35*f6dfbefdSBarry Smith Some of the `SNESLineSearch` methods attempt to minimize a given objective provided by this function to determine a step length. 36eb23ba39SBarry Smith 37*f6dfbefdSBarry Smith If not provided then this defaults to the two norm of the function evaluation (set with `SNESSetFunction()`) 38075cc632SBarry Smith 39*f6dfbefdSBarry Smith This is not used in the `SNESLINESEARCHCP` line search. 40*f6dfbefdSBarry Smith 41*f6dfbefdSBarry Smith .seealso: `SNES`, `SNESLineSearch()`, `SNESGetObjective()`, `SNESComputeObjective()`, `SNESSetFunction()`, `SNESSetJacobian()`, `SNESObjectiveFunction` 422a4ee8f2SPeter Brune @*/ 439371c9d4SSatish Balay PetscErrorCode SNESSetObjective(SNES snes, PetscErrorCode (*obj)(SNES, Vec, PetscReal *, void *), void *ctx) { 442a4ee8f2SPeter Brune DM dm; 452a4ee8f2SPeter Brune 462a4ee8f2SPeter Brune PetscFunctionBegin; 472a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 489566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 499566063dSJacob Faibussowitsch PetscCall(DMSNESSetObjective(dm, obj, ctx)); 502a4ee8f2SPeter Brune PetscFunctionReturn(0); 512a4ee8f2SPeter Brune } 522a4ee8f2SPeter Brune 532a4ee8f2SPeter Brune /*@C 542a4ee8f2SPeter Brune SNESGetObjective - Returns the objective function. 552a4ee8f2SPeter Brune 562a4ee8f2SPeter Brune Not Collective 572a4ee8f2SPeter Brune 582a4ee8f2SPeter Brune Input Parameter: 59*f6dfbefdSBarry Smith . snes - the `SNES` context 602a4ee8f2SPeter Brune 61d8d19677SJose E. Roman Output Parameters: 62*f6dfbefdSBarry Smith + obj - objective evaluation routine (or NULL); see `SNESObjectFunction` for details 630298fd71SBarry Smith - ctx - the function context (or NULL) 642a4ee8f2SPeter Brune 652a4ee8f2SPeter Brune Level: advanced 662a4ee8f2SPeter Brune 67*f6dfbefdSBarry Smith .seealso: `SNES`, `SNESSetObjective()`, `SNESGetSolution()` 682a4ee8f2SPeter Brune @*/ 699371c9d4SSatish Balay PetscErrorCode SNESGetObjective(SNES snes, PetscErrorCode (**obj)(SNES, Vec, PetscReal *, void *), void **ctx) { 702a4ee8f2SPeter Brune DM dm; 715fd66863SKarl Rupp 722a4ee8f2SPeter Brune PetscFunctionBegin; 732a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 749566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 759566063dSJacob Faibussowitsch PetscCall(DMSNESGetObjective(dm, obj, ctx)); 762a4ee8f2SPeter Brune PetscFunctionReturn(0); 772a4ee8f2SPeter Brune } 782a4ee8f2SPeter Brune 792a4ee8f2SPeter Brune /*@C 80*f6dfbefdSBarry Smith SNESComputeObjective - Computes the objective function that has been provided by `SNESSetObjective()` 812a4ee8f2SPeter Brune 82*f6dfbefdSBarry Smith Collective on snes 832a4ee8f2SPeter Brune 84d8d19677SJose E. Roman Input Parameters: 85*f6dfbefdSBarry Smith + snes - the `SNES` context 862a4ee8f2SPeter Brune - X - the state vector 872a4ee8f2SPeter Brune 882a4ee8f2SPeter Brune Output Parameter: 892a4ee8f2SPeter Brune . ob - the objective value 902a4ee8f2SPeter Brune 91*f6dfbefdSBarry Smith Level: developer 922a4ee8f2SPeter Brune 93*f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNES`, `SNESSetObjective()`, `SNESGetSolution()` 942a4ee8f2SPeter Brune @*/ 959371c9d4SSatish Balay PetscErrorCode SNESComputeObjective(SNES snes, Vec X, PetscReal *ob) { 962a4ee8f2SPeter Brune DM dm; 97942e3340SBarry Smith DMSNES sdm; 98f23aa3ddSBarry Smith 992a4ee8f2SPeter Brune PetscFunctionBegin; 1002a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 1012a4ee8f2SPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 102534a8f05SLisandro Dalcin PetscValidRealPointer(ob, 3); 1039566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 1049566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm, &sdm)); 10522c6f798SBarry Smith if (sdm->ops->computeobjective) { 1069566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_ObjectiveEval, snes, X, 0, 0)); 1079566063dSJacob Faibussowitsch PetscCall((sdm->ops->computeobjective)(snes, X, ob, sdm->objectivectx)); 1089566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_ObjectiveEval, snes, X, 0, 0)); 109ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetObjective() before SNESComputeObjective()."); 1102a4ee8f2SPeter Brune PetscFunctionReturn(0); 1112a4ee8f2SPeter Brune } 11297584545SPeter Brune 1135c3e6ab7SPeter Brune /*@C 114*f6dfbefdSBarry Smith SNESObjectiveComputeFunctionDefaultFD - Computes the gradient of a user provided objective function 1155c3e6ab7SPeter Brune 116*f6dfbefdSBarry Smith Collective on snes 1175c3e6ab7SPeter Brune 118d8d19677SJose E. Roman Input Parameters: 119*f6dfbefdSBarry Smith + snes - the `SNES` context 1205c3e6ab7SPeter Brune . X - the state vector 1215c3e6ab7SPeter Brune - ctx - the (ignored) function context 1225c3e6ab7SPeter Brune 1235c3e6ab7SPeter Brune Output Parameter: 1245c3e6ab7SPeter Brune . F - the function value 1255c3e6ab7SPeter Brune 126*f6dfbefdSBarry Smith Options Database Keys: 127*f6dfbefdSBarry Smith + -snes_fd_function_eps - Tolerance for including non-zero entries into the gradient, default is 1.e-6 128*f6dfbefdSBarry Smith - -snes_fd_function - Computes function from user provided objective function (set with `SNESSetObjective()`) with finite difference 1295c3e6ab7SPeter Brune 1305c3e6ab7SPeter Brune Notes: 131*f6dfbefdSBarry Smith This function can be used with `SNESSetFunction()` to have the nonlinear function solved for with `SNES` defined by the gradient of an objective function 132*f6dfbefdSBarry Smith `SNESObjectiveComputeFunctionDefaultFD()` is similar in character to `SNESComputeJacobianDefault()`. 1335c3e6ab7SPeter Brune Therefore, it should be used for debugging purposes only. Using it in conjunction with 134*f6dfbefdSBarry Smith `SNESComputeJacobianDefault()` is excessively costly and produces a Jacobian that is quite 135*f6dfbefdSBarry Smith noisy. This is often necessary, but should be done with care, even when debugging 1365c3e6ab7SPeter Brune small problems. 1375c3e6ab7SPeter Brune 1385c3e6ab7SPeter Brune Note that this uses quadratic interpolation of the objective to form each value in the function. 1395c3e6ab7SPeter Brune 1405c3e6ab7SPeter Brune Level: advanced 1415c3e6ab7SPeter Brune 142*f6dfbefdSBarry Smith .seealso: `SNESSetObjective()`, `SNESSetFunction()`, `SNESComputeObjective()`, `SNESComputeJacobianDefault()` 1435c3e6ab7SPeter Brune @*/ 1449371c9d4SSatish Balay PetscErrorCode SNESObjectiveComputeFunctionDefaultFD(SNES snes, Vec X, Vec F, void *ctx) { 14597584545SPeter Brune Vec Xh; 14697584545SPeter Brune PetscInt i, N, start, end; 14797584545SPeter Brune PetscReal ob, ob1, ob2, ob3, fob, dx, eps = 1e-6; 14897584545SPeter Brune PetscScalar fv, xv; 14997584545SPeter Brune 15097584545SPeter Brune PetscFunctionBegin; 1519566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Xh)); 152d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->prefix, "Differencing parameters", "SNES"); 1539566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_fd_function_eps", "Tolerance for nonzero entries in fd function", "None", eps, &eps, NULL)); 154d0609cedSBarry Smith PetscOptionsEnd(); 1559566063dSJacob Faibussowitsch PetscCall(VecSet(F, 0.)); 15697584545SPeter Brune 1579566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &fob)); 15897584545SPeter Brune 1599566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &N)); 1609566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &start, &end)); 1619566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, X, &ob)); 16297584545SPeter Brune 163f5af7f23SKarl Rupp if (fob > 0.) dx = 1e-6 * fob; 164f5af7f23SKarl Rupp else dx = 1e-6; 16597584545SPeter Brune 16697584545SPeter Brune for (i = 0; i < N; i++) { 16797584545SPeter Brune /* compute the 1st value */ 1689566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xh)); 16997584545SPeter Brune if (i >= start && i < end) { 17097584545SPeter Brune xv = dx; 1719566063dSJacob Faibussowitsch PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES)); 17297584545SPeter Brune } 1739566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Xh)); 1749566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Xh)); 1759566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, Xh, &ob1)); 17697584545SPeter Brune 17797584545SPeter Brune /* compute the 2nd value */ 1789566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xh)); 17997584545SPeter Brune if (i >= start && i < end) { 18097584545SPeter Brune xv = 2. * dx; 1819566063dSJacob Faibussowitsch PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES)); 18297584545SPeter Brune } 1839566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Xh)); 1849566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Xh)); 1859566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, Xh, &ob2)); 18697584545SPeter Brune 18797584545SPeter Brune /* compute the 3rd value */ 1889566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xh)); 18997584545SPeter Brune if (i >= start && i < end) { 19097584545SPeter Brune xv = -dx; 1919566063dSJacob Faibussowitsch PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES)); 19297584545SPeter Brune } 1939566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Xh)); 1949566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Xh)); 1959566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, Xh, &ob3)); 19697584545SPeter Brune 19797584545SPeter Brune if (i >= start && i < end) { 19897584545SPeter Brune /* set this entry to be the gradient of the objective */ 19997584545SPeter Brune fv = (-ob2 + 6. * ob1 - 3. * ob - 2. * ob3) / (6. * dx); 20097584545SPeter Brune if (PetscAbsScalar(fv) > eps) { 2019566063dSJacob Faibussowitsch PetscCall(VecSetValues(F, 1, &i, &fv, INSERT_VALUES)); 20297584545SPeter Brune } else { 20397584545SPeter Brune fv = 0.; 2049566063dSJacob Faibussowitsch PetscCall(VecSetValues(F, 1, &i, &fv, INSERT_VALUES)); 20597584545SPeter Brune } 20697584545SPeter Brune } 20797584545SPeter Brune } 2089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xh)); 20997584545SPeter Brune 2109566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(F)); 2119566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(F)); 21297584545SPeter Brune PetscFunctionReturn(0); 21397584545SPeter Brune } 214