xref: /petsc/src/snes/interface/snesob.c (revision 534a8f05a7a8aff70dd8cfd53d9cd834400a8dbf)
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 .      F - current function/gradient
142a4ee8f2SPeter Brune .      obj - real to hold the objective value
152a4ee8f2SPeter Brune -      ctx - optional user-defined objective context
162a4ee8f2SPeter Brune 
17878cb397SSatish Balay    Level: advanced
18878cb397SSatish Balay 
19f8b49ee9SBarry Smith .seealso:   SNESSetFunction(), SNESGetFunction(), SNESSetObjective(), SNESGetObjective(), SNESJacobianFunction, SNESFunction
20075cc632SBarry Smith M*/
21075cc632SBarry Smith 
22075cc632SBarry Smith 
23075cc632SBarry Smith /*@C
24eb23ba39SBarry Smith    SNESSetObjective - Sets the objective function minimized by some of the SNES linesearch methods.
25075cc632SBarry Smith 
26075cc632SBarry Smith    Logically Collective on SNES
27075cc632SBarry Smith 
28075cc632SBarry Smith    Input Parameters:
29075cc632SBarry Smith +  snes - the SNES context
30f8b49ee9SBarry Smith .  obj - objective evaluation routine; see SNESObjectiveFunction for details
31075cc632SBarry Smith -  ctx - [optional] user-defined context for private data for the
320298fd71SBarry Smith          function evaluation routine (may be NULL)
33075cc632SBarry Smith 
342b26979fSBarry Smith    Level: intermediate
352a4ee8f2SPeter Brune 
36eb23ba39SBarry Smith    Note: This is not used in the SNESLINESEARCHCP line search.
37eb23ba39SBarry Smith 
38eb23ba39SBarry Smith          If not provided then this defaults to the two norm of the function evaluation (set with SNESSetFunction())
39075cc632SBarry Smith 
40075cc632SBarry Smith .seealso: SNESGetObjective(), SNESComputeObjective(), SNESSetFunction(), SNESSetJacobian(), SNESObjectiveFunction
412a4ee8f2SPeter Brune @*/
42f8b49ee9SBarry Smith PetscErrorCode  SNESSetObjective(SNES snes,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx)
432a4ee8f2SPeter Brune {
442a4ee8f2SPeter Brune   PetscErrorCode ierr;
452a4ee8f2SPeter Brune   DM             dm;
462a4ee8f2SPeter Brune 
472a4ee8f2SPeter Brune   PetscFunctionBegin;
482a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
492a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
50f8b49ee9SBarry Smith   ierr = DMSNESSetObjective(dm,obj,ctx);CHKERRQ(ierr);
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:
602a4ee8f2SPeter Brune .  snes - the SNES context
612a4ee8f2SPeter Brune 
622a4ee8f2SPeter Brune    Output Parameter:
63f8b49ee9SBarry 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 
682a4ee8f2SPeter Brune .seealso: SNESSetObjective(), SNESGetSolution()
692a4ee8f2SPeter Brune @*/
70f8b49ee9SBarry Smith PetscErrorCode SNESGetObjective(SNES snes,PetscErrorCode (**obj)(SNES,Vec,PetscReal*,void*),void **ctx)
712a4ee8f2SPeter Brune {
722a4ee8f2SPeter Brune   PetscErrorCode ierr;
732a4ee8f2SPeter Brune   DM             dm;
745fd66863SKarl Rupp 
752a4ee8f2SPeter Brune   PetscFunctionBegin;
762a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
772a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
78f8b49ee9SBarry Smith   ierr = DMSNESGetObjective(dm,obj,ctx);CHKERRQ(ierr);
792a4ee8f2SPeter Brune   PetscFunctionReturn(0);
802a4ee8f2SPeter Brune }
812a4ee8f2SPeter Brune 
822a4ee8f2SPeter Brune /*@C
832a4ee8f2SPeter Brune    SNESComputeObjective - Computes the objective.
842a4ee8f2SPeter Brune 
852a4ee8f2SPeter Brune    Collective on SNES
862a4ee8f2SPeter Brune 
872a4ee8f2SPeter Brune    Input Parameter:
882a4ee8f2SPeter Brune +  snes - the SNES context
892a4ee8f2SPeter Brune -  X    - the state vector
902a4ee8f2SPeter Brune 
912a4ee8f2SPeter Brune    Output Parameter:
922a4ee8f2SPeter Brune .  ob   - the objective value
932a4ee8f2SPeter Brune 
942a4ee8f2SPeter Brune    Level: advanced
952a4ee8f2SPeter Brune 
962a4ee8f2SPeter Brune .seealso: SNESSetObjective(), SNESGetSolution()
972a4ee8f2SPeter Brune @*/
982a4ee8f2SPeter Brune PetscErrorCode SNESComputeObjective(SNES snes,Vec X,PetscReal *ob)
992a4ee8f2SPeter Brune {
1002a4ee8f2SPeter Brune   PetscErrorCode ierr;
1012a4ee8f2SPeter Brune   DM             dm;
102942e3340SBarry Smith   DMSNES         sdm;
103f23aa3ddSBarry Smith 
1042a4ee8f2SPeter Brune   PetscFunctionBegin;
1052a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1062a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
107*534a8f05SLisandro Dalcin   PetscValidRealPointer(ob,3);
1082a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
109942e3340SBarry Smith   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
11022c6f798SBarry Smith   if (sdm->ops->computeobjective) {
11194db00ebSBarry Smith     ierr = PetscLogEventBegin(SNES_ObjectiveEval,snes,X,0,0);CHKERRQ(ierr);
11222c6f798SBarry Smith     ierr = (sdm->ops->computeobjective)(snes,X,ob,sdm->objectivectx);CHKERRQ(ierr);
11394db00ebSBarry Smith     ierr = PetscLogEventEnd(SNES_ObjectiveEval,snes,X,0,0);CHKERRQ(ierr);
114ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetObjective() before SNESComputeObjective().");
1152a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1162a4ee8f2SPeter Brune }
11797584545SPeter Brune 
11897584545SPeter Brune 
1195c3e6ab7SPeter Brune /*@C
1205c3e6ab7SPeter Brune    SNESObjectiveComputeFunctionDefaultFD - Computes the gradient of a user provided objective
1215c3e6ab7SPeter Brune 
1225c3e6ab7SPeter Brune    Collective on SNES
1235c3e6ab7SPeter Brune 
1245c3e6ab7SPeter Brune    Input Parameter:
1255c3e6ab7SPeter Brune +  snes - the SNES context
1265c3e6ab7SPeter Brune .  X    - the state vector
1275c3e6ab7SPeter Brune -  ctx  - the (ignored) function context
1285c3e6ab7SPeter Brune 
1295c3e6ab7SPeter Brune    Output Parameter:
1305c3e6ab7SPeter Brune .  F   - the function value
1315c3e6ab7SPeter Brune 
1325c3e6ab7SPeter Brune    Options Database Key:
1335c3e6ab7SPeter Brune +  -snes_fd_function_eps - The differencing parameter
1345c3e6ab7SPeter Brune -  -snes_fd_function - Compute function from user provided objective with finite difference
1355c3e6ab7SPeter Brune 
1365c3e6ab7SPeter Brune    Notes:
1375c3e6ab7SPeter Brune    SNESObjectiveComputeFunctionDefaultFD is similar in character to SNESComputeJacobianDefault.
1385c3e6ab7SPeter Brune    Therefore, it should be used for debugging purposes only.  Using it in conjunction with
1395c3e6ab7SPeter Brune    SNESComputeJacobianDefault is excessively costly and produces a Jacobian that is quite
1405c3e6ab7SPeter Brune    noisy.  This is often necessary, but should be done with a grain of salt, even when debugging
1415c3e6ab7SPeter Brune    small problems.
1425c3e6ab7SPeter Brune 
1435c3e6ab7SPeter Brune    Note that this uses quadratic interpolation of the objective to form each value in the function.
1445c3e6ab7SPeter Brune 
1455c3e6ab7SPeter Brune    Level: advanced
1465c3e6ab7SPeter Brune 
1475c3e6ab7SPeter Brune .seealso: SNESSetFunction(), SNESComputeObjective(), SNESComputeJacobianDefault()
1485c3e6ab7SPeter Brune @*/
1498d359177SBarry Smith PetscErrorCode SNESObjectiveComputeFunctionDefaultFD(SNES snes,Vec X,Vec F,void *ctx)
1500adebc6cSBarry Smith {
15197584545SPeter Brune   Vec            Xh;
15297584545SPeter Brune   PetscErrorCode ierr;
15397584545SPeter Brune   PetscInt       i,N,start,end;
15497584545SPeter Brune   PetscReal      ob,ob1,ob2,ob3,fob,dx,eps=1e-6;
15597584545SPeter Brune   PetscScalar    fv,xv;
15697584545SPeter Brune 
15797584545SPeter Brune   PetscFunctionBegin;
15897584545SPeter Brune   ierr = VecDuplicate(X,&Xh);CHKERRQ(ierr);
159302440fdSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"Differencing parameters","SNES");CHKERRQ(ierr);
1600298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_fd_function_eps","Tolerance for nonzero entries in fd function","None",eps,&eps,NULL);CHKERRQ(ierr);
161302440fdSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
16297584545SPeter Brune   ierr = VecSet(F,0.);CHKERRQ(ierr);
16397584545SPeter Brune 
16497584545SPeter Brune   ierr = VecNorm(X,NORM_2,&fob);CHKERRQ(ierr);
16597584545SPeter Brune 
16697584545SPeter Brune   ierr = VecGetSize(X,&N);CHKERRQ(ierr);
16797584545SPeter Brune   ierr = VecGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
16897584545SPeter Brune   ierr = SNESComputeObjective(snes,X,&ob);CHKERRQ(ierr);
16997584545SPeter Brune 
170f5af7f23SKarl Rupp   if (fob > 0.) dx =1e-6*fob;
171f5af7f23SKarl Rupp   else dx = 1e-6;
17297584545SPeter Brune 
17397584545SPeter Brune   for (i=0; i<N; i++) {
17497584545SPeter Brune     /* compute the 1st value */
17597584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
17697584545SPeter Brune     if (i>= start && i<end) {
17797584545SPeter Brune       xv   = dx;
17897584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
17997584545SPeter Brune     }
18097584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
18197584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
18297584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob1);CHKERRQ(ierr);
18397584545SPeter Brune 
18497584545SPeter Brune     /* compute the 2nd value */
18597584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
18697584545SPeter Brune     if (i>= start && i<end) {
18797584545SPeter Brune       xv   = 2.*dx;
18897584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
18997584545SPeter Brune     }
19097584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
19197584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
19297584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob2);CHKERRQ(ierr);
19397584545SPeter Brune 
19497584545SPeter Brune     /* compute the 3rd value */
19597584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
19697584545SPeter Brune     if (i>= start && i<end) {
19797584545SPeter Brune       xv   = -dx;
19897584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
19997584545SPeter Brune     }
20097584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
20197584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
20297584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob3);CHKERRQ(ierr);
20397584545SPeter Brune 
20497584545SPeter Brune     if (i >= start && i<end) {
20597584545SPeter Brune       /* set this entry to be the gradient of the objective */
20697584545SPeter Brune       fv = (-ob2 + 6.*ob1 - 3.*ob -2.*ob3) / (6.*dx);
20797584545SPeter Brune       if (PetscAbsScalar(fv) > eps) {
20897584545SPeter Brune         ierr = VecSetValues(F,1,&i,&fv,INSERT_VALUES);CHKERRQ(ierr);
20997584545SPeter Brune       } else {
21097584545SPeter Brune         fv   = 0.;
21197584545SPeter Brune         ierr = VecSetValues(F,1,&i,&fv,INSERT_VALUES);CHKERRQ(ierr);
21297584545SPeter Brune       }
21397584545SPeter Brune     }
21497584545SPeter Brune   }
21597584545SPeter Brune 
21697584545SPeter Brune   ierr = VecDestroy(&Xh);CHKERRQ(ierr);
21797584545SPeter Brune 
21897584545SPeter Brune   ierr = VecAssemblyBegin(F);CHKERRQ(ierr);
21997584545SPeter Brune   ierr = VecAssemblyEnd(F);CHKERRQ(ierr);
22097584545SPeter Brune   PetscFunctionReturn(0);
22197584545SPeter Brune }
222