15d83a8b1SBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
22a4ee8f2SPeter Brune
3075cc632SBarry Smith /*@C
49bcc50f1SBarry Smith SNESSetObjective - Sets the objective function minimized by some of the `SNES` linesearch methods, used instead of the 2-norm of the residual in the line search
5075cc632SBarry Smith
6c3339decSBarry Smith Logically Collective
7075cc632SBarry Smith
8075cc632SBarry Smith Input Parameters:
9f6dfbefdSBarry Smith + snes - the `SNES` context
10709e13c7SJeremy L Thompson . obj - objective evaluation routine; see `SNESObjectiveFn` for the calling sequence
11709e13c7SJeremy L Thompson - ctx - [optional] user-defined context for private data for the objective evaluation routine (may be `NULL`)
12075cc632SBarry Smith
132b26979fSBarry Smith Level: intermediate
142a4ee8f2SPeter Brune
15f6dfbefdSBarry Smith Note:
16f6dfbefdSBarry Smith Some of the `SNESLineSearch` methods attempt to minimize a given objective provided by this function to determine a step length.
17eb23ba39SBarry Smith
18420bcc1bSBarry Smith If not provided then this defaults to the two-norm of the function evaluation (set with `SNESSetFunction()`)
19075cc632SBarry Smith
20f6dfbefdSBarry Smith This is not used in the `SNESLINESEARCHCP` line search.
21f6dfbefdSBarry Smith
229bcc50f1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESGetObjective()`, `SNESComputeObjective()`, `SNESSetFunction()`, `SNESSetJacobian()`,
2376c63389SBarry Smith `SNESObjectiveFn`, `SNESSetObjectiveDomainError()`
242a4ee8f2SPeter Brune @*/
SNESSetObjective(SNES snes,SNESObjectiveFn * obj,PetscCtx ctx)25*2a8381b2SBarry Smith PetscErrorCode SNESSetObjective(SNES snes, SNESObjectiveFn *obj, PetscCtx ctx)
26d71ae5a4SJacob Faibussowitsch {
272a4ee8f2SPeter Brune DM dm;
282a4ee8f2SPeter Brune
292a4ee8f2SPeter Brune PetscFunctionBegin;
302a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
319566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm));
329566063dSJacob Faibussowitsch PetscCall(DMSNESSetObjective(dm, obj, ctx));
333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
342a4ee8f2SPeter Brune }
352a4ee8f2SPeter Brune
362a4ee8f2SPeter Brune /*@C
37dc4c0fb0SBarry Smith SNESGetObjective - Returns the objective function set with `SNESSetObjective()`
382a4ee8f2SPeter Brune
392a4ee8f2SPeter Brune Not Collective
402a4ee8f2SPeter Brune
412a4ee8f2SPeter Brune Input Parameter:
42f6dfbefdSBarry Smith . snes - the `SNES` context
432a4ee8f2SPeter Brune
44d8d19677SJose E. Roman Output Parameters:
45709e13c7SJeremy L Thompson + obj - objective evaluation routine (or `NULL`); see `SNESObjectiveFn` for the calling sequence
46dc4c0fb0SBarry Smith - ctx - the function context (or `NULL`)
472a4ee8f2SPeter Brune
482a4ee8f2SPeter Brune Level: advanced
492a4ee8f2SPeter Brune
508434afd1SBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetObjective()`, `SNESGetSolution()`, `SNESObjectiveFn`
512a4ee8f2SPeter Brune @*/
SNESGetObjective(SNES snes,SNESObjectiveFn ** obj,PetscCtxRt ctx)52*2a8381b2SBarry Smith PetscErrorCode SNESGetObjective(SNES snes, SNESObjectiveFn **obj, PetscCtxRt ctx)
53d71ae5a4SJacob Faibussowitsch {
542a4ee8f2SPeter Brune DM dm;
555fd66863SKarl Rupp
562a4ee8f2SPeter Brune PetscFunctionBegin;
572a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
589566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm));
599566063dSJacob Faibussowitsch PetscCall(DMSNESGetObjective(dm, obj, ctx));
603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
612a4ee8f2SPeter Brune }
622a4ee8f2SPeter Brune
635d83a8b1SBarry Smith /*@
64f6dfbefdSBarry Smith SNESComputeObjective - Computes the objective function that has been provided by `SNESSetObjective()`
652a4ee8f2SPeter Brune
66c3339decSBarry Smith Collective
672a4ee8f2SPeter Brune
68d8d19677SJose E. Roman Input Parameters:
69f6dfbefdSBarry Smith + snes - the `SNES` context
702a4ee8f2SPeter Brune - X - the state vector
712a4ee8f2SPeter Brune
722a4ee8f2SPeter Brune Output Parameter:
732a4ee8f2SPeter Brune . ob - the objective value
742a4ee8f2SPeter Brune
75f6dfbefdSBarry Smith Level: developer
762a4ee8f2SPeter Brune
7700677de2SStefano Zampini Notes:
7800677de2SStefano Zampini `SNESComputeObjective()` is typically used within line-search routines,
7900677de2SStefano Zampini so users would not generally call this routine themselves.
8000677de2SStefano Zampini
8100677de2SStefano Zampini When solving for $F(x) = b$, this routine computes $objective(x) - x^T b$ where $objective(x)$ is the function provided with `SNESSetObjective()`
8200677de2SStefano Zampini
83420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESLineSearch`, `SNES`, `SNESSetObjective()`, `SNESGetSolution()`
842a4ee8f2SPeter Brune @*/
SNESComputeObjective(SNES snes,Vec X,PetscReal * ob)85d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESComputeObjective(SNES snes, Vec X, PetscReal *ob)
86d71ae5a4SJacob Faibussowitsch {
872a4ee8f2SPeter Brune DM dm;
88942e3340SBarry Smith DMSNES sdm;
89f23aa3ddSBarry Smith
902a4ee8f2SPeter Brune PetscFunctionBegin;
912a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
922a4ee8f2SPeter Brune PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
934f572ea9SToby Isaac PetscAssertPointer(ob, 3);
949566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm));
959566063dSJacob Faibussowitsch PetscCall(DMGetDMSNES(dm, &sdm));
960fdf79fbSJacob Faibussowitsch PetscCheck(sdm->ops->computeobjective, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetObjective() before SNESComputeObjective().");
979566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(SNES_ObjectiveEval, snes, X, 0, 0));
98f4f49eeaSPierre Jolivet PetscCall(sdm->ops->computeobjective(snes, X, ob, sdm->objectivectx));
999566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(SNES_ObjectiveEval, snes, X, 0, 0));
10000677de2SStefano Zampini if (snes->vec_rhs) {
10100677de2SStefano Zampini PetscScalar dot;
10200677de2SStefano Zampini
10300677de2SStefano Zampini PetscCall(VecDot(snes->vec_rhs, X, &dot));
10400677de2SStefano Zampini *ob -= PetscRealPart(dot);
10500677de2SStefano Zampini }
1063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1072a4ee8f2SPeter Brune }
10897584545SPeter Brune
1095c3e6ab7SPeter Brune /*@C
110f6dfbefdSBarry Smith SNESObjectiveComputeFunctionDefaultFD - Computes the gradient of a user provided objective function
1115c3e6ab7SPeter Brune
112c3339decSBarry Smith Collective
1135c3e6ab7SPeter Brune
114d8d19677SJose E. Roman Input Parameters:
115f6dfbefdSBarry Smith + snes - the `SNES` context
1165c3e6ab7SPeter Brune . X - the state vector
1175c3e6ab7SPeter Brune - ctx - the (ignored) function context
1185c3e6ab7SPeter Brune
1195c3e6ab7SPeter Brune Output Parameter:
1205c3e6ab7SPeter Brune . F - the function value
1215c3e6ab7SPeter Brune
122f6dfbefdSBarry Smith Options Database Keys:
123f6dfbefdSBarry Smith + -snes_fd_function_eps - Tolerance for including non-zero entries into the gradient, default is 1.e-6
124f6dfbefdSBarry Smith - -snes_fd_function - Computes function from user provided objective function (set with `SNESSetObjective()`) with finite difference
1255c3e6ab7SPeter Brune
126dc4c0fb0SBarry Smith Level: advanced
127dc4c0fb0SBarry Smith
1285c3e6ab7SPeter Brune Notes:
129f6dfbefdSBarry 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
130f6dfbefdSBarry Smith `SNESObjectiveComputeFunctionDefaultFD()` is similar in character to `SNESComputeJacobianDefault()`.
1315c3e6ab7SPeter Brune Therefore, it should be used for debugging purposes only. Using it in conjunction with
132f6dfbefdSBarry Smith `SNESComputeJacobianDefault()` is excessively costly and produces a Jacobian that is quite
133f6dfbefdSBarry Smith noisy. This is often necessary, but should be done with care, even when debugging
1345c3e6ab7SPeter Brune small problems.
1355c3e6ab7SPeter Brune
136420bcc1bSBarry Smith This uses quadratic interpolation of the objective to form each value in the function.
1375c3e6ab7SPeter Brune
1388434afd1SBarry Smith .seealso: [](ch_snes), `SNESSetObjective()`, `SNESSetFunction()`, `SNESComputeObjective()`, `SNESComputeJacobianDefault()`, `SNESObjectiveFn`
1395c3e6ab7SPeter Brune @*/
SNESObjectiveComputeFunctionDefaultFD(SNES snes,Vec X,Vec F,PetscCtx ctx)140*2a8381b2SBarry Smith PetscErrorCode SNESObjectiveComputeFunctionDefaultFD(SNES snes, Vec X, Vec F, PetscCtx ctx)
141d71ae5a4SJacob Faibussowitsch {
14297584545SPeter Brune Vec Xh;
14397584545SPeter Brune PetscInt i, N, start, end;
14497584545SPeter Brune PetscReal ob, ob1, ob2, ob3, fob, dx, eps = 1e-6;
14597584545SPeter Brune PetscScalar fv, xv;
14697584545SPeter Brune
14797584545SPeter Brune PetscFunctionBegin;
1489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &Xh));
149d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->prefix, "Differencing parameters", "SNES");
1509566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-snes_fd_function_eps", "Tolerance for nonzero entries in fd function", "None", eps, &eps, NULL));
151d0609cedSBarry Smith PetscOptionsEnd();
1529566063dSJacob Faibussowitsch PetscCall(VecSet(F, 0.));
15397584545SPeter Brune
1549566063dSJacob Faibussowitsch PetscCall(VecNorm(X, NORM_2, &fob));
15597584545SPeter Brune
1569566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &N));
1579566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(X, &start, &end));
1589566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, X, &ob));
15997584545SPeter Brune
160f5af7f23SKarl Rupp if (fob > 0.) dx = 1e-6 * fob;
161f5af7f23SKarl Rupp else dx = 1e-6;
16297584545SPeter Brune
16397584545SPeter Brune for (i = 0; i < N; i++) {
16497584545SPeter Brune /* compute the 1st value */
1659566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xh));
16697584545SPeter Brune if (i >= start && i < end) {
16797584545SPeter Brune xv = dx;
1689566063dSJacob Faibussowitsch PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES));
16997584545SPeter Brune }
1709566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Xh));
1719566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Xh));
1729566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, Xh, &ob1));
17397584545SPeter Brune
17497584545SPeter Brune /* compute the 2nd value */
1759566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xh));
17697584545SPeter Brune if (i >= start && i < end) {
17797584545SPeter Brune xv = 2. * dx;
1789566063dSJacob Faibussowitsch PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES));
17997584545SPeter Brune }
1809566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Xh));
1819566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Xh));
1829566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, Xh, &ob2));
18397584545SPeter Brune
18497584545SPeter Brune /* compute the 3rd value */
1859566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xh));
18697584545SPeter Brune if (i >= start && i < end) {
18797584545SPeter Brune xv = -dx;
1889566063dSJacob Faibussowitsch PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES));
18997584545SPeter Brune }
1909566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(Xh));
1919566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(Xh));
1929566063dSJacob Faibussowitsch PetscCall(SNESComputeObjective(snes, Xh, &ob3));
19397584545SPeter Brune
19497584545SPeter Brune if (i >= start && i < end) {
19597584545SPeter Brune /* set this entry to be the gradient of the objective */
19697584545SPeter Brune fv = (-ob2 + 6. * ob1 - 3. * ob - 2. * ob3) / (6. * dx);
19797584545SPeter Brune if (PetscAbsScalar(fv) > eps) {
1989566063dSJacob Faibussowitsch PetscCall(VecSetValues(F, 1, &i, &fv, INSERT_VALUES));
19997584545SPeter Brune } else {
20097584545SPeter Brune fv = 0.;
2019566063dSJacob Faibussowitsch PetscCall(VecSetValues(F, 1, &i, &fv, INSERT_VALUES));
20297584545SPeter Brune }
20397584545SPeter Brune }
20497584545SPeter Brune }
2059566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Xh));
2069566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(F));
2079566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(F));
2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
20997584545SPeter Brune }
210