1af0996ceSBarry Smith #include <petsc/private/snesimpl.h> 22a4ee8f2SPeter Brune 3075cc632SBarry Smith /*MC 4075cc632SBarry Smith SNESObjectiveFunction - functional form used to convey the 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 18db781477SPatrick Sanan .seealso: `SNESSetFunction()`, `SNESGetFunction()`, `SNESSetObjective()`, `SNESGetObjective()`, `SNESJacobianFunction`, `SNESFunction` 19075cc632SBarry Smith M*/ 20075cc632SBarry Smith 21075cc632SBarry Smith /*@C 22eb23ba39SBarry Smith SNESSetObjective - Sets the objective function minimized by some of the SNES linesearch methods. 23075cc632SBarry Smith 24075cc632SBarry Smith Logically Collective on SNES 25075cc632SBarry Smith 26075cc632SBarry Smith Input Parameters: 27075cc632SBarry Smith + snes - the SNES context 28f8b49ee9SBarry 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 34eb23ba39SBarry Smith Note: This is not used in the SNESLINESEARCHCP line search. 35eb23ba39SBarry Smith 36eb23ba39SBarry Smith If not provided then this defaults to the two norm of the function evaluation (set with SNESSetFunction()) 37075cc632SBarry Smith 38db781477SPatrick Sanan .seealso: `SNESGetObjective()`, `SNESComputeObjective()`, `SNESSetFunction()`, `SNESSetJacobian()`, `SNESObjectiveFunction` 392a4ee8f2SPeter Brune @*/ 40*9371c9d4SSatish Balay PetscErrorCode SNESSetObjective(SNES snes, PetscErrorCode (*obj)(SNES, Vec, PetscReal *, void *), void *ctx) { 412a4ee8f2SPeter Brune DM dm; 422a4ee8f2SPeter Brune 432a4ee8f2SPeter Brune PetscFunctionBegin; 442a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 459566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 469566063dSJacob Faibussowitsch PetscCall(DMSNESSetObjective(dm, obj, ctx)); 472a4ee8f2SPeter Brune PetscFunctionReturn(0); 482a4ee8f2SPeter Brune } 492a4ee8f2SPeter Brune 502a4ee8f2SPeter Brune /*@C 512a4ee8f2SPeter Brune SNESGetObjective - Returns the objective function. 522a4ee8f2SPeter Brune 532a4ee8f2SPeter Brune Not Collective 542a4ee8f2SPeter Brune 552a4ee8f2SPeter Brune Input Parameter: 562a4ee8f2SPeter Brune . snes - the SNES context 572a4ee8f2SPeter Brune 58d8d19677SJose E. Roman Output Parameters: 59f8b49ee9SBarry Smith + obj - objective evaluation routine (or NULL); see SNESObjectFunction for details 600298fd71SBarry Smith - ctx - the function context (or NULL) 612a4ee8f2SPeter Brune 622a4ee8f2SPeter Brune Level: advanced 632a4ee8f2SPeter Brune 64db781477SPatrick Sanan .seealso: `SNESSetObjective()`, `SNESGetSolution()` 652a4ee8f2SPeter Brune @*/ 66*9371c9d4SSatish Balay PetscErrorCode SNESGetObjective(SNES snes, PetscErrorCode (**obj)(SNES, Vec, PetscReal *, void *), void **ctx) { 672a4ee8f2SPeter Brune DM dm; 685fd66863SKarl Rupp 692a4ee8f2SPeter Brune PetscFunctionBegin; 702a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 719566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 729566063dSJacob Faibussowitsch PetscCall(DMSNESGetObjective(dm, obj, ctx)); 732a4ee8f2SPeter Brune PetscFunctionReturn(0); 742a4ee8f2SPeter Brune } 752a4ee8f2SPeter Brune 762a4ee8f2SPeter Brune /*@C 772a4ee8f2SPeter Brune SNESComputeObjective - Computes the objective. 782a4ee8f2SPeter Brune 792a4ee8f2SPeter Brune Collective on SNES 802a4ee8f2SPeter Brune 81d8d19677SJose E. Roman Input Parameters: 822a4ee8f2SPeter Brune + snes - the SNES context 832a4ee8f2SPeter Brune - X - the state vector 842a4ee8f2SPeter Brune 852a4ee8f2SPeter Brune Output Parameter: 862a4ee8f2SPeter Brune . ob - the objective value 872a4ee8f2SPeter Brune 882a4ee8f2SPeter Brune Level: advanced 892a4ee8f2SPeter Brune 90db781477SPatrick Sanan .seealso: `SNESSetObjective()`, `SNESGetSolution()` 912a4ee8f2SPeter Brune @*/ 92*9371c9d4SSatish Balay PetscErrorCode SNESComputeObjective(SNES snes, Vec X, PetscReal *ob) { 932a4ee8f2SPeter Brune DM dm; 94942e3340SBarry Smith DMSNES sdm; 95f23aa3ddSBarry Smith 962a4ee8f2SPeter Brune PetscFunctionBegin; 972a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1); 982a4ee8f2SPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2); 99534a8f05SLisandro Dalcin PetscValidRealPointer(ob, 3); 1009566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 1019566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm, &sdm)); 10222c6f798SBarry Smith if (sdm->ops->computeobjective) { 1039566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_ObjectiveEval, snes, X, 0, 0)); 1049566063dSJacob Faibussowitsch PetscCall((sdm->ops->computeobjective)(snes, X, ob, sdm->objectivectx)); 1059566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_ObjectiveEval, snes, X, 0, 0)); 106ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetObjective() before SNESComputeObjective()."); 1072a4ee8f2SPeter Brune PetscFunctionReturn(0); 1082a4ee8f2SPeter Brune } 10997584545SPeter Brune 1105c3e6ab7SPeter Brune /*@C 1115c3e6ab7SPeter Brune SNESObjectiveComputeFunctionDefaultFD - Computes the gradient of a user provided objective 1125c3e6ab7SPeter Brune 1135c3e6ab7SPeter Brune Collective on SNES 1145c3e6ab7SPeter Brune 115d8d19677SJose E. Roman Input Parameters: 1165c3e6ab7SPeter Brune + snes - the SNES context 1175c3e6ab7SPeter Brune . X - the state vector 1185c3e6ab7SPeter Brune - ctx - the (ignored) function context 1195c3e6ab7SPeter Brune 1205c3e6ab7SPeter Brune Output Parameter: 1215c3e6ab7SPeter Brune . F - the function value 1225c3e6ab7SPeter Brune 1235c3e6ab7SPeter Brune Options Database Key: 1245c3e6ab7SPeter Brune + -snes_fd_function_eps - The differencing parameter 1255c3e6ab7SPeter Brune - -snes_fd_function - Compute function from user provided objective with finite difference 1265c3e6ab7SPeter Brune 1275c3e6ab7SPeter Brune Notes: 1285c3e6ab7SPeter Brune SNESObjectiveComputeFunctionDefaultFD is similar in character to SNESComputeJacobianDefault. 1295c3e6ab7SPeter Brune Therefore, it should be used for debugging purposes only. Using it in conjunction with 1305c3e6ab7SPeter Brune SNESComputeJacobianDefault is excessively costly and produces a Jacobian that is quite 1315c3e6ab7SPeter Brune noisy. This is often necessary, but should be done with a grain of salt, even when debugging 1325c3e6ab7SPeter Brune small problems. 1335c3e6ab7SPeter Brune 1345c3e6ab7SPeter Brune Note that this uses quadratic interpolation of the objective to form each value in the function. 1355c3e6ab7SPeter Brune 1365c3e6ab7SPeter Brune Level: advanced 1375c3e6ab7SPeter Brune 138db781477SPatrick Sanan .seealso: `SNESSetFunction()`, `SNESComputeObjective()`, `SNESComputeJacobianDefault()` 1395c3e6ab7SPeter Brune @*/ 140*9371c9d4SSatish Balay PetscErrorCode SNESObjectiveComputeFunctionDefaultFD(SNES snes, Vec X, Vec F, void *ctx) { 14197584545SPeter Brune Vec Xh; 14297584545SPeter Brune PetscInt i, N, start, end; 14397584545SPeter Brune PetscReal ob, ob1, ob2, ob3, fob, dx, eps = 1e-6; 14497584545SPeter Brune PetscScalar fv, xv; 14597584545SPeter Brune 14697584545SPeter Brune PetscFunctionBegin; 1479566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Xh)); 148d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->prefix, "Differencing parameters", "SNES"); 1499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_fd_function_eps", "Tolerance for nonzero entries in fd function", "None", eps, &eps, NULL)); 150d0609cedSBarry Smith PetscOptionsEnd(); 1519566063dSJacob Faibussowitsch PetscCall(VecSet(F, 0.)); 15297584545SPeter Brune 1539566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &fob)); 15497584545SPeter Brune 1559566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &N)); 1569566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &start, &end)); 1579566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, X, &ob)); 15897584545SPeter Brune 159f5af7f23SKarl Rupp if (fob > 0.) dx = 1e-6 * fob; 160f5af7f23SKarl Rupp else dx = 1e-6; 16197584545SPeter Brune 16297584545SPeter Brune for (i = 0; i < N; i++) { 16397584545SPeter Brune /* compute the 1st value */ 1649566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xh)); 16597584545SPeter Brune if (i >= start && i < end) { 16697584545SPeter Brune xv = dx; 1679566063dSJacob Faibussowitsch PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES)); 16897584545SPeter Brune } 1699566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Xh)); 1709566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Xh)); 1719566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, Xh, &ob1)); 17297584545SPeter Brune 17397584545SPeter Brune /* compute the 2nd value */ 1749566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xh)); 17597584545SPeter Brune if (i >= start && i < end) { 17697584545SPeter Brune xv = 2. * dx; 1779566063dSJacob Faibussowitsch PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES)); 17897584545SPeter Brune } 1799566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Xh)); 1809566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Xh)); 1819566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, Xh, &ob2)); 18297584545SPeter Brune 18397584545SPeter Brune /* compute the 3rd value */ 1849566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xh)); 18597584545SPeter Brune if (i >= start && i < end) { 18697584545SPeter Brune xv = -dx; 1879566063dSJacob Faibussowitsch PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES)); 18897584545SPeter Brune } 1899566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Xh)); 1909566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Xh)); 1919566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, Xh, &ob3)); 19297584545SPeter Brune 19397584545SPeter Brune if (i >= start && i < end) { 19497584545SPeter Brune /* set this entry to be the gradient of the objective */ 19597584545SPeter Brune fv = (-ob2 + 6. * ob1 - 3. * ob - 2. * ob3) / (6. * dx); 19697584545SPeter Brune if (PetscAbsScalar(fv) > eps) { 1979566063dSJacob Faibussowitsch PetscCall(VecSetValues(F, 1, &i, &fv, INSERT_VALUES)); 19897584545SPeter Brune } else { 19997584545SPeter Brune fv = 0.; 2009566063dSJacob Faibussowitsch PetscCall(VecSetValues(F, 1, &i, &fv, INSERT_VALUES)); 20197584545SPeter Brune } 20297584545SPeter Brune } 20397584545SPeter Brune } 20497584545SPeter Brune 2059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xh)); 20697584545SPeter Brune 2079566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(F)); 2089566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(F)); 20997584545SPeter Brune PetscFunctionReturn(0); 21097584545SPeter Brune } 211