xref: /petsc/src/snes/interface/snesob.c (revision 00677de27e955871ed92da7850af0f7802ab35ae)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h>
22a4ee8f2SPeter Brune 
3075cc632SBarry Smith /*MC
4dc4c0fb0SBarry Smith     SNESObjectiveFunction - functional form used to convey an objective function to the nonlinear solver, that will be used instead of the 2-norm of the residual
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:
11dc4c0fb0SBarry Smith +      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 
18420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESSetObjective()`, `SNESGetObjective()`, `SNESJacobianFunction`, `SNESFunction`
19075cc632SBarry Smith M*/
20075cc632SBarry Smith 
21075cc632SBarry Smith /*@C
22dc4c0fb0SBarry Smith   SNESSetObjective - Sets the objective function minimized by some of the `SNES` linesearch methods, used instead of the 2-norm of the residual
23075cc632SBarry Smith 
24c3339decSBarry Smith   Logically Collective
25075cc632SBarry Smith 
26075cc632SBarry Smith   Input Parameters:
27f6dfbefdSBarry Smith + snes - the `SNES` context
28f6dfbefdSBarry Smith . obj  - objective evaluation routine; see `SNESObjectiveFunction` for details
29075cc632SBarry Smith - ctx  - [optional] user-defined context for private data for the
30dc4c0fb0SBarry Smith          function evaluation routine (may be `NULL`)
31075cc632SBarry Smith 
322b26979fSBarry Smith   Level: intermediate
332a4ee8f2SPeter Brune 
34f6dfbefdSBarry Smith   Note:
35f6dfbefdSBarry Smith   Some of the `SNESLineSearch` methods attempt to minimize a given objective provided by this function to determine a step length.
36eb23ba39SBarry Smith 
37420bcc1bSBarry Smith   If not provided then this defaults to the two-norm of the function evaluation (set with `SNESSetFunction()`)
38075cc632SBarry Smith 
39f6dfbefdSBarry Smith   This is not used in the `SNESLINESEARCHCP` line search.
40f6dfbefdSBarry Smith 
41420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESLineSearch()`, `SNESGetObjective()`, `SNESComputeObjective()`, `SNESSetFunction()`, `SNESSetJacobian()`, `SNESObjectiveFunction`
422a4ee8f2SPeter Brune @*/
43d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESSetObjective(SNES snes, PetscErrorCode (*obj)(SNES, Vec, PetscReal *, void *), void *ctx)
44d71ae5a4SJacob Faibussowitsch {
452a4ee8f2SPeter Brune   DM dm;
462a4ee8f2SPeter Brune 
472a4ee8f2SPeter Brune   PetscFunctionBegin;
482a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
499566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
509566063dSJacob Faibussowitsch   PetscCall(DMSNESSetObjective(dm, obj, ctx));
513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
522a4ee8f2SPeter Brune }
532a4ee8f2SPeter Brune 
542a4ee8f2SPeter Brune /*@C
55dc4c0fb0SBarry Smith   SNESGetObjective - Returns the objective function set with `SNESSetObjective()`
562a4ee8f2SPeter Brune 
572a4ee8f2SPeter Brune   Not Collective
582a4ee8f2SPeter Brune 
592a4ee8f2SPeter Brune   Input Parameter:
60f6dfbefdSBarry Smith . snes - the `SNES` context
612a4ee8f2SPeter Brune 
62d8d19677SJose E. Roman   Output Parameters:
63420bcc1bSBarry Smith + obj - objective evaluation routine (or `NULL`); see `SNESObjectiveFunction` for details
64dc4c0fb0SBarry Smith - ctx - the function context (or `NULL`)
652a4ee8f2SPeter Brune 
662a4ee8f2SPeter Brune   Level: advanced
672a4ee8f2SPeter Brune 
68420bcc1bSBarry Smith .seealso: [](ch_snes), `SNES`, `SNESSetObjective()`, `SNESObjectiveFunction`, `SNESGetSolution()`
692a4ee8f2SPeter Brune @*/
70d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESGetObjective(SNES snes, PetscErrorCode (**obj)(SNES, Vec, PetscReal *, void *), void **ctx)
71d71ae5a4SJacob Faibussowitsch {
722a4ee8f2SPeter Brune   DM dm;
735fd66863SKarl Rupp 
742a4ee8f2SPeter Brune   PetscFunctionBegin;
752a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
769566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
779566063dSJacob Faibussowitsch   PetscCall(DMSNESGetObjective(dm, obj, ctx));
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
792a4ee8f2SPeter Brune }
802a4ee8f2SPeter Brune 
812a4ee8f2SPeter Brune /*@C
82f6dfbefdSBarry Smith   SNESComputeObjective - Computes the objective function that has been provided by `SNESSetObjective()`
832a4ee8f2SPeter Brune 
84c3339decSBarry Smith   Collective
852a4ee8f2SPeter Brune 
86d8d19677SJose E. Roman   Input Parameters:
87f6dfbefdSBarry Smith + snes - the `SNES` context
882a4ee8f2SPeter Brune - X    - the state vector
892a4ee8f2SPeter Brune 
902a4ee8f2SPeter Brune   Output Parameter:
912a4ee8f2SPeter Brune . ob - the objective value
922a4ee8f2SPeter Brune 
93f6dfbefdSBarry Smith   Level: developer
942a4ee8f2SPeter Brune 
95*00677de2SStefano Zampini   Notes:
96*00677de2SStefano Zampini   `SNESComputeObjective()` is typically used within line-search routines,
97*00677de2SStefano Zampini   so users would not generally call this routine themselves.
98*00677de2SStefano Zampini 
99*00677de2SStefano Zampini   When solving for $F(x) = b$, this routine computes $objective(x) - x^T b$ where $objective(x)$ is the function provided with `SNESSetObjective()`
100*00677de2SStefano Zampini 
101420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESLineSearch`, `SNES`, `SNESSetObjective()`, `SNESGetSolution()`
1022a4ee8f2SPeter Brune @*/
103d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESComputeObjective(SNES snes, Vec X, PetscReal *ob)
104d71ae5a4SJacob Faibussowitsch {
1052a4ee8f2SPeter Brune   DM     dm;
106942e3340SBarry Smith   DMSNES sdm;
107f23aa3ddSBarry Smith 
1082a4ee8f2SPeter Brune   PetscFunctionBegin;
1092a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1102a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
1114f572ea9SToby Isaac   PetscAssertPointer(ob, 3);
1129566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
1139566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
1140fdf79fbSJacob Faibussowitsch   PetscCheck(sdm->ops->computeobjective, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetObjective() before SNESComputeObjective().");
1159566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(SNES_ObjectiveEval, snes, X, 0, 0));
1169566063dSJacob Faibussowitsch   PetscCall((sdm->ops->computeobjective)(snes, X, ob, sdm->objectivectx));
1179566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(SNES_ObjectiveEval, snes, X, 0, 0));
118*00677de2SStefano Zampini   if (snes->vec_rhs) {
119*00677de2SStefano Zampini     PetscScalar dot;
120*00677de2SStefano Zampini 
121*00677de2SStefano Zampini     PetscCall(VecDot(snes->vec_rhs, X, &dot));
122*00677de2SStefano Zampini     *ob -= PetscRealPart(dot);
123*00677de2SStefano Zampini   }
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1252a4ee8f2SPeter Brune }
12697584545SPeter Brune 
1275c3e6ab7SPeter Brune /*@C
128f6dfbefdSBarry Smith   SNESObjectiveComputeFunctionDefaultFD - Computes the gradient of a user provided objective function
1295c3e6ab7SPeter Brune 
130c3339decSBarry Smith   Collective
1315c3e6ab7SPeter Brune 
132d8d19677SJose E. Roman   Input Parameters:
133f6dfbefdSBarry Smith + snes - the `SNES` context
1345c3e6ab7SPeter Brune . X    - the state vector
1355c3e6ab7SPeter Brune - ctx  - the (ignored) function context
1365c3e6ab7SPeter Brune 
1375c3e6ab7SPeter Brune   Output Parameter:
1385c3e6ab7SPeter Brune . F - the function value
1395c3e6ab7SPeter Brune 
140f6dfbefdSBarry Smith   Options Database Keys:
141f6dfbefdSBarry Smith + -snes_fd_function_eps - Tolerance for including non-zero entries into the gradient, default is 1.e-6
142f6dfbefdSBarry Smith - -snes_fd_function     - Computes function from user provided objective function (set with `SNESSetObjective()`) with finite difference
1435c3e6ab7SPeter Brune 
144dc4c0fb0SBarry Smith   Level: advanced
145dc4c0fb0SBarry Smith 
1465c3e6ab7SPeter Brune   Notes:
147f6dfbefdSBarry 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
148f6dfbefdSBarry Smith   `SNESObjectiveComputeFunctionDefaultFD()` is similar in character to `SNESComputeJacobianDefault()`.
1495c3e6ab7SPeter Brune   Therefore, it should be used for debugging purposes only.  Using it in conjunction with
150f6dfbefdSBarry Smith   `SNESComputeJacobianDefault()` is excessively costly and produces a Jacobian that is quite
151f6dfbefdSBarry Smith   noisy.  This is often necessary, but should be done with care, even when debugging
1525c3e6ab7SPeter Brune   small problems.
1535c3e6ab7SPeter Brune 
154420bcc1bSBarry Smith   This uses quadratic interpolation of the objective to form each value in the function.
1555c3e6ab7SPeter Brune 
156420bcc1bSBarry Smith .seealso: [](ch_snes), `SNESSetObjective()`, `SNESSetFunction()`, `SNESComputeObjective()`, `SNESComputeJacobianDefault()`
1575c3e6ab7SPeter Brune @*/
158d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESObjectiveComputeFunctionDefaultFD(SNES snes, Vec X, Vec F, void *ctx)
159d71ae5a4SJacob Faibussowitsch {
16097584545SPeter Brune   Vec         Xh;
16197584545SPeter Brune   PetscInt    i, N, start, end;
16297584545SPeter Brune   PetscReal   ob, ob1, ob2, ob3, fob, dx, eps = 1e-6;
16397584545SPeter Brune   PetscScalar fv, xv;
16497584545SPeter Brune 
16597584545SPeter Brune   PetscFunctionBegin;
1669566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &Xh));
167d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->prefix, "Differencing parameters", "SNES");
1689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_fd_function_eps", "Tolerance for nonzero entries in fd function", "None", eps, &eps, NULL));
169d0609cedSBarry Smith   PetscOptionsEnd();
1709566063dSJacob Faibussowitsch   PetscCall(VecSet(F, 0.));
17197584545SPeter Brune 
1729566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, &fob));
17397584545SPeter Brune 
1749566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &N));
1759566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &start, &end));
1769566063dSJacob Faibussowitsch   PetscCall(SNESComputeObjective(snes, X, &ob));
17797584545SPeter Brune 
178f5af7f23SKarl Rupp   if (fob > 0.) dx = 1e-6 * fob;
179f5af7f23SKarl Rupp   else dx = 1e-6;
18097584545SPeter Brune 
18197584545SPeter Brune   for (i = 0; i < N; i++) {
18297584545SPeter Brune     /* compute the 1st value */
1839566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Xh));
18497584545SPeter Brune     if (i >= start && i < end) {
18597584545SPeter Brune       xv = dx;
1869566063dSJacob Faibussowitsch       PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES));
18797584545SPeter Brune     }
1889566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(Xh));
1899566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(Xh));
1909566063dSJacob Faibussowitsch     PetscCall(SNESComputeObjective(snes, Xh, &ob1));
19197584545SPeter Brune 
19297584545SPeter Brune     /* compute the 2nd value */
1939566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Xh));
19497584545SPeter Brune     if (i >= start && i < end) {
19597584545SPeter Brune       xv = 2. * dx;
1969566063dSJacob Faibussowitsch       PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES));
19797584545SPeter Brune     }
1989566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(Xh));
1999566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(Xh));
2009566063dSJacob Faibussowitsch     PetscCall(SNESComputeObjective(snes, Xh, &ob2));
20197584545SPeter Brune 
20297584545SPeter Brune     /* compute the 3rd value */
2039566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Xh));
20497584545SPeter Brune     if (i >= start && i < end) {
20597584545SPeter Brune       xv = -dx;
2069566063dSJacob Faibussowitsch       PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES));
20797584545SPeter Brune     }
2089566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(Xh));
2099566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(Xh));
2109566063dSJacob Faibussowitsch     PetscCall(SNESComputeObjective(snes, Xh, &ob3));
21197584545SPeter Brune 
21297584545SPeter Brune     if (i >= start && i < end) {
21397584545SPeter Brune       /* set this entry to be the gradient of the objective */
21497584545SPeter Brune       fv = (-ob2 + 6. * ob1 - 3. * ob - 2. * ob3) / (6. * dx);
21597584545SPeter Brune       if (PetscAbsScalar(fv) > eps) {
2169566063dSJacob Faibussowitsch         PetscCall(VecSetValues(F, 1, &i, &fv, INSERT_VALUES));
21797584545SPeter Brune       } else {
21897584545SPeter Brune         fv = 0.;
2199566063dSJacob Faibussowitsch         PetscCall(VecSetValues(F, 1, &i, &fv, INSERT_VALUES));
22097584545SPeter Brune       }
22197584545SPeter Brune     }
22297584545SPeter Brune   }
2239566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Xh));
2249566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(F));
2259566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(F));
2263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22797584545SPeter Brune }
228