xref: /petsc/src/snes/interface/snesob.c (revision 0fdf79fb08699bf9be0aa4d8ba0185e387a216c8)
1af0996ceSBarry Smith #include <petsc/private/snesimpl.h>
22a4ee8f2SPeter Brune 
3075cc632SBarry Smith /*MC
4f6dfbefdSBarry 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 
18f6dfbefdSBarry Smith .seealso: `SNES`, `SNESSetFunction()`, `SNESGetFunction()`, `SNESSetObjective()`, `SNESGetObjective()`, `SNESJacobianFunction`, `SNESFunction`
19075cc632SBarry Smith M*/
20075cc632SBarry Smith 
21075cc632SBarry Smith /*@C
22f6dfbefdSBarry Smith    SNESSetObjective - Sets the objective function minimized by some of the `SNES` linesearch methods.
23075cc632SBarry Smith 
24f6dfbefdSBarry Smith    Logically Collective on snes
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
300298fd71SBarry 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 
37f6dfbefdSBarry 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 
41f6dfbefdSBarry Smith .seealso: `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));
512a4ee8f2SPeter Brune   PetscFunctionReturn(0);
522a4ee8f2SPeter Brune }
532a4ee8f2SPeter Brune 
542a4ee8f2SPeter Brune /*@C
552a4ee8f2SPeter Brune    SNESGetObjective - Returns the objective function.
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:
63f6dfbefdSBarry Smith +  obj - objective evaluation routine (or NULL); see `SNESObjectFunction` for details
640298fd71SBarry Smith -  ctx - the function context (or NULL)
652a4ee8f2SPeter Brune 
662a4ee8f2SPeter Brune    Level: advanced
672a4ee8f2SPeter Brune 
68f6dfbefdSBarry Smith .seealso: `SNES`, `SNESSetObjective()`, `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));
782a4ee8f2SPeter Brune   PetscFunctionReturn(0);
792a4ee8f2SPeter Brune }
802a4ee8f2SPeter Brune 
812a4ee8f2SPeter Brune /*@C
82f6dfbefdSBarry Smith    SNESComputeObjective - Computes the objective function that has been provided by `SNESSetObjective()`
832a4ee8f2SPeter Brune 
84f6dfbefdSBarry Smith    Collective on snes
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 
95f6dfbefdSBarry Smith .seealso: `SNESLineSearch`, `SNES`, `SNESSetObjective()`, `SNESGetSolution()`
962a4ee8f2SPeter Brune @*/
97d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESComputeObjective(SNES snes, Vec X, PetscReal *ob)
98d71ae5a4SJacob Faibussowitsch {
992a4ee8f2SPeter Brune   DM     dm;
100942e3340SBarry Smith   DMSNES sdm;
101f23aa3ddSBarry Smith 
1022a4ee8f2SPeter Brune   PetscFunctionBegin;
1032a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
1042a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X, VEC_CLASSID, 2);
105534a8f05SLisandro Dalcin   PetscValidRealPointer(ob, 3);
1069566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
1079566063dSJacob Faibussowitsch   PetscCall(DMGetDMSNES(dm, &sdm));
108*0fdf79fbSJacob Faibussowitsch   PetscCheck(sdm->ops->computeobjective, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetObjective() before SNESComputeObjective().");
1099566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(SNES_ObjectiveEval, snes, X, 0, 0));
1109566063dSJacob Faibussowitsch   PetscCall((sdm->ops->computeobjective)(snes, X, ob, sdm->objectivectx));
1119566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(SNES_ObjectiveEval, snes, X, 0, 0));
1122a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1132a4ee8f2SPeter Brune }
11497584545SPeter Brune 
1155c3e6ab7SPeter Brune /*@C
116f6dfbefdSBarry Smith    SNESObjectiveComputeFunctionDefaultFD - Computes the gradient of a user provided objective function
1175c3e6ab7SPeter Brune 
118f6dfbefdSBarry Smith    Collective on snes
1195c3e6ab7SPeter Brune 
120d8d19677SJose E. Roman    Input Parameters:
121f6dfbefdSBarry Smith +  snes - the `SNES` context
1225c3e6ab7SPeter Brune .  X    - the state vector
1235c3e6ab7SPeter Brune -  ctx  - the (ignored) function context
1245c3e6ab7SPeter Brune 
1255c3e6ab7SPeter Brune    Output Parameter:
1265c3e6ab7SPeter Brune .  F   - the function value
1275c3e6ab7SPeter Brune 
128f6dfbefdSBarry Smith    Options Database Keys:
129f6dfbefdSBarry Smith +  -snes_fd_function_eps - Tolerance for including non-zero entries into the gradient, default is 1.e-6
130f6dfbefdSBarry Smith -  -snes_fd_function - Computes function from user provided objective function (set with `SNESSetObjective()`) with finite difference
1315c3e6ab7SPeter Brune 
1325c3e6ab7SPeter Brune    Notes:
133f6dfbefdSBarry 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
134f6dfbefdSBarry Smith    `SNESObjectiveComputeFunctionDefaultFD()` is similar in character to `SNESComputeJacobianDefault()`.
1355c3e6ab7SPeter Brune    Therefore, it should be used for debugging purposes only.  Using it in conjunction with
136f6dfbefdSBarry Smith    `SNESComputeJacobianDefault()` is excessively costly and produces a Jacobian that is quite
137f6dfbefdSBarry Smith    noisy.  This is often necessary, but should be done with care, even when debugging
1385c3e6ab7SPeter Brune    small problems.
1395c3e6ab7SPeter Brune 
1405c3e6ab7SPeter Brune    Note that this uses quadratic interpolation of the objective to form each value in the function.
1415c3e6ab7SPeter Brune 
1425c3e6ab7SPeter Brune    Level: advanced
1435c3e6ab7SPeter Brune 
144f6dfbefdSBarry Smith .seealso: `SNESSetObjective()`, `SNESSetFunction()`, `SNESComputeObjective()`, `SNESComputeJacobianDefault()`
1455c3e6ab7SPeter Brune @*/
146d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESObjectiveComputeFunctionDefaultFD(SNES snes, Vec X, Vec F, void *ctx)
147d71ae5a4SJacob Faibussowitsch {
14897584545SPeter Brune   Vec         Xh;
14997584545SPeter Brune   PetscInt    i, N, start, end;
15097584545SPeter Brune   PetscReal   ob, ob1, ob2, ob3, fob, dx, eps = 1e-6;
15197584545SPeter Brune   PetscScalar fv, xv;
15297584545SPeter Brune 
15397584545SPeter Brune   PetscFunctionBegin;
1549566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &Xh));
155d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->prefix, "Differencing parameters", "SNES");
1569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-snes_fd_function_eps", "Tolerance for nonzero entries in fd function", "None", eps, &eps, NULL));
157d0609cedSBarry Smith   PetscOptionsEnd();
1589566063dSJacob Faibussowitsch   PetscCall(VecSet(F, 0.));
15997584545SPeter Brune 
1609566063dSJacob Faibussowitsch   PetscCall(VecNorm(X, NORM_2, &fob));
16197584545SPeter Brune 
1629566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &N));
1639566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &start, &end));
1649566063dSJacob Faibussowitsch   PetscCall(SNESComputeObjective(snes, X, &ob));
16597584545SPeter Brune 
166f5af7f23SKarl Rupp   if (fob > 0.) dx = 1e-6 * fob;
167f5af7f23SKarl Rupp   else dx = 1e-6;
16897584545SPeter Brune 
16997584545SPeter Brune   for (i = 0; i < N; i++) {
17097584545SPeter Brune     /* compute the 1st value */
1719566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Xh));
17297584545SPeter Brune     if (i >= start && i < end) {
17397584545SPeter Brune       xv = dx;
1749566063dSJacob Faibussowitsch       PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES));
17597584545SPeter Brune     }
1769566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(Xh));
1779566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(Xh));
1789566063dSJacob Faibussowitsch     PetscCall(SNESComputeObjective(snes, Xh, &ob1));
17997584545SPeter Brune 
18097584545SPeter Brune     /* compute the 2nd value */
1819566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Xh));
18297584545SPeter Brune     if (i >= start && i < end) {
18397584545SPeter Brune       xv = 2. * dx;
1849566063dSJacob Faibussowitsch       PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES));
18597584545SPeter Brune     }
1869566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(Xh));
1879566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(Xh));
1889566063dSJacob Faibussowitsch     PetscCall(SNESComputeObjective(snes, Xh, &ob2));
18997584545SPeter Brune 
19097584545SPeter Brune     /* compute the 3rd value */
1919566063dSJacob Faibussowitsch     PetscCall(VecCopy(X, Xh));
19297584545SPeter Brune     if (i >= start && i < end) {
19397584545SPeter Brune       xv = -dx;
1949566063dSJacob Faibussowitsch       PetscCall(VecSetValues(Xh, 1, &i, &xv, ADD_VALUES));
19597584545SPeter Brune     }
1969566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(Xh));
1979566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(Xh));
1989566063dSJacob Faibussowitsch     PetscCall(SNESComputeObjective(snes, Xh, &ob3));
19997584545SPeter Brune 
20097584545SPeter Brune     if (i >= start && i < end) {
20197584545SPeter Brune       /* set this entry to be the gradient of the objective */
20297584545SPeter Brune       fv = (-ob2 + 6. * ob1 - 3. * ob - 2. * ob3) / (6. * dx);
20397584545SPeter Brune       if (PetscAbsScalar(fv) > eps) {
2049566063dSJacob Faibussowitsch         PetscCall(VecSetValues(F, 1, &i, &fv, INSERT_VALUES));
20597584545SPeter Brune       } else {
20697584545SPeter Brune         fv = 0.;
2079566063dSJacob Faibussowitsch         PetscCall(VecSetValues(F, 1, &i, &fv, INSERT_VALUES));
20897584545SPeter Brune       }
20997584545SPeter Brune     }
21097584545SPeter Brune   }
2119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Xh));
21297584545SPeter Brune 
2139566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(F));
2149566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(F));
21597584545SPeter Brune   PetscFunctionReturn(0);
21697584545SPeter Brune }
217