xref: /petsc/src/snes/interface/snesob.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
12a4ee8f2SPeter Brune #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:
7075cc632SBarry 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 
17075cc632SBarry Smith .seealso:   SNESSetFunction(), SNESGetFunction(), SNESSetObjectiveFunction(), SNESGetObjectiveFunction()
18075cc632SBarry Smith M*/
19075cc632SBarry Smith 
20075cc632SBarry Smith 
21075cc632SBarry Smith #undef __FUNCT__
22075cc632SBarry Smith #define __FUNCT__ "SNESSetObjective"
23075cc632SBarry Smith /*@C
24075cc632SBarry Smith    SNESSetObjective - Sets the objective function minimized by the SNES methods.
25075cc632SBarry Smith 
26075cc632SBarry Smith    Logically Collective on SNES
27075cc632SBarry Smith 
28075cc632SBarry Smith    Input Parameters:
29075cc632SBarry Smith +  snes - the SNES context
30075cc632SBarry Smith .  SNESObjectiveFunction - objective evaluation routine
31075cc632SBarry Smith -  ctx - [optional] user-defined context for private data for the
320298fd71SBarry Smith          function evaluation routine (may be NULL)
33075cc632SBarry Smith 
342a4ee8f2SPeter Brune    Level: beginner
352a4ee8f2SPeter Brune 
36075cc632SBarry Smith    Note: If not provided then this defaults to the two norm of the function evaluation (set with SNESSetFunction())
37075cc632SBarry Smith 
382a4ee8f2SPeter Brune .keywords: SNES, nonlinear, set, objective
392a4ee8f2SPeter Brune 
40075cc632SBarry Smith .seealso: SNESGetObjective(), SNESComputeObjective(), SNESSetFunction(), SNESSetJacobian(), SNESObjectiveFunction
412a4ee8f2SPeter Brune @*/
42075cc632SBarry Smith PetscErrorCode  SNESSetObjective(SNES snes,PetscErrorCode (*SNESObjectiveFunction)(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);
50075cc632SBarry Smith   ierr = DMSNESSetObjective(dm,SNESObjectiveFunction,ctx);CHKERRQ(ierr);
512a4ee8f2SPeter Brune   PetscFunctionReturn(0);
522a4ee8f2SPeter Brune }
532a4ee8f2SPeter Brune 
542a4ee8f2SPeter Brune #undef __FUNCT__
552a4ee8f2SPeter Brune #define __FUNCT__ "SNESGetObjective"
562a4ee8f2SPeter Brune /*@C
572a4ee8f2SPeter Brune    SNESGetObjective - Returns the objective function.
582a4ee8f2SPeter Brune 
592a4ee8f2SPeter Brune    Not Collective
602a4ee8f2SPeter Brune 
612a4ee8f2SPeter Brune    Input Parameter:
622a4ee8f2SPeter Brune .  snes - the SNES context
632a4ee8f2SPeter Brune 
642a4ee8f2SPeter Brune    Output Parameter:
650298fd71SBarry Smith +  func - the function (or NULL)
660298fd71SBarry Smith -  ctx - the function context (or NULL)
672a4ee8f2SPeter Brune 
682a4ee8f2SPeter Brune    Level: advanced
692a4ee8f2SPeter Brune 
702a4ee8f2SPeter Brune .keywords: SNES, nonlinear, get, objective
712a4ee8f2SPeter Brune 
722a4ee8f2SPeter Brune .seealso: SNESSetObjective(), SNESGetSolution()
732a4ee8f2SPeter Brune @*/
74075cc632SBarry Smith PetscErrorCode SNESGetObjective(SNES snes,PetscErrorCode (**SNESObjectiveFunction)(SNES,Vec,PetscReal*,void*),void **ctx)
752a4ee8f2SPeter Brune {
762a4ee8f2SPeter Brune   PetscErrorCode ierr;
772a4ee8f2SPeter Brune   DM             dm;
785fd66863SKarl Rupp 
792a4ee8f2SPeter Brune   PetscFunctionBegin;
802a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
812a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
82075cc632SBarry Smith   ierr = DMSNESGetObjective(dm,SNESObjectiveFunction,ctx);CHKERRQ(ierr);
832a4ee8f2SPeter Brune   PetscFunctionReturn(0);
842a4ee8f2SPeter Brune }
852a4ee8f2SPeter Brune 
862a4ee8f2SPeter Brune #undef __FUNCT__
872a4ee8f2SPeter Brune #define __FUNCT__ "SNESComputeObjective"
882a4ee8f2SPeter Brune /*@C
892a4ee8f2SPeter Brune    SNESComputeObjective - Computes the objective.
902a4ee8f2SPeter Brune 
912a4ee8f2SPeter Brune    Collective on SNES
922a4ee8f2SPeter Brune 
932a4ee8f2SPeter Brune    Input Parameter:
942a4ee8f2SPeter Brune +  snes - the SNES context
952a4ee8f2SPeter Brune -  X    - the state vector
962a4ee8f2SPeter Brune 
972a4ee8f2SPeter Brune    Output Parameter:
982a4ee8f2SPeter Brune .  ob   - the objective value
992a4ee8f2SPeter Brune 
1002a4ee8f2SPeter Brune    Level: advanced
1012a4ee8f2SPeter Brune 
1022a4ee8f2SPeter Brune .keywords: SNES, nonlinear, compute, objective
1032a4ee8f2SPeter Brune 
1042a4ee8f2SPeter Brune .seealso: SNESSetObjective(), SNESGetSolution()
1052a4ee8f2SPeter Brune @*/
1062a4ee8f2SPeter Brune PetscErrorCode SNESComputeObjective(SNES snes,Vec X,PetscReal *ob)
1072a4ee8f2SPeter Brune {
1082a4ee8f2SPeter Brune   PetscErrorCode ierr;
1092a4ee8f2SPeter Brune   DM             dm;
110942e3340SBarry Smith   DMSNES         sdm;
111f23aa3ddSBarry Smith 
1122a4ee8f2SPeter Brune   PetscFunctionBegin;
1132a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1142a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1152a4ee8f2SPeter Brune   PetscValidPointer(ob,3);
1162a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
117942e3340SBarry Smith   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
11822c6f798SBarry Smith   if (sdm->ops->computeobjective) {
11922c6f798SBarry Smith     ierr = (sdm->ops->computeobjective)(snes,X,ob,sdm->objectivectx);CHKERRQ(ierr);
120*ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetObjective() before SNESComputeObjective().");
1212a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1222a4ee8f2SPeter Brune }
12397584545SPeter Brune 
12497584545SPeter Brune 
12597584545SPeter Brune #undef __FUNCT__
12697584545SPeter Brune #define __FUNCT__ "SNESDefaultObjectiveComputeFunctionFD"
1270adebc6cSBarry Smith PetscErrorCode SNESDefaultObjectiveComputeFunctionFD(SNES snes,Vec X,Vec F,void *ctx)
1280adebc6cSBarry Smith {
12997584545SPeter Brune   /* Quadratically interpolates the change in the objective based upon a change in a particular direction */
13097584545SPeter Brune 
13197584545SPeter Brune   Vec            Xh;
13297584545SPeter Brune   PetscErrorCode ierr;
13397584545SPeter Brune   PetscInt       i,N,start,end;
13497584545SPeter Brune   PetscReal      ob,ob1,ob2,ob3,fob,dx,eps=1e-6;
13597584545SPeter Brune   PetscScalar    fv,xv;
13697584545SPeter Brune 
13797584545SPeter Brune   PetscFunctionBegin;
13897584545SPeter Brune   ierr = VecDuplicate(X,&Xh);CHKERRQ(ierr);
1390298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_fd_function_eps","Tolerance for nonzero entries in fd function","None",eps,&eps,NULL);CHKERRQ(ierr);
14097584545SPeter Brune   ierr = VecSet(F,0.);CHKERRQ(ierr);
14197584545SPeter Brune 
14297584545SPeter Brune   ierr = VecNorm(X,NORM_2,&fob);CHKERRQ(ierr);
14397584545SPeter Brune 
14497584545SPeter Brune   ierr = VecGetSize(X,&N);CHKERRQ(ierr);
14597584545SPeter Brune   ierr = VecGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
14697584545SPeter Brune   ierr = SNESComputeObjective(snes,X,&ob);CHKERRQ(ierr);
14797584545SPeter Brune 
148f5af7f23SKarl Rupp   if (fob > 0.) dx =1e-6*fob;
149f5af7f23SKarl Rupp   else dx = 1e-6;
15097584545SPeter Brune 
15197584545SPeter Brune   for (i=0; i<N; i++) {
15297584545SPeter Brune     /* compute the 1st value */
15397584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
15497584545SPeter Brune     if (i>= start && i<end) {
15597584545SPeter Brune       xv   = dx;
15697584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
15797584545SPeter Brune     }
15897584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
15997584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
16097584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob1);CHKERRQ(ierr);
16197584545SPeter Brune 
16297584545SPeter Brune     /* compute the 2nd value */
16397584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
16497584545SPeter Brune     if (i>= start && i<end) {
16597584545SPeter Brune       xv   = 2.*dx;
16697584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
16797584545SPeter Brune     }
16897584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
16997584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
17097584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob2);CHKERRQ(ierr);
17197584545SPeter Brune 
17297584545SPeter Brune     /* compute the 3rd value */
17397584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
17497584545SPeter Brune     if (i>= start && i<end) {
17597584545SPeter Brune       xv   = -dx;
17697584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
17797584545SPeter Brune     }
17897584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
17997584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
18097584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob3);CHKERRQ(ierr);
18197584545SPeter Brune 
18297584545SPeter Brune     if (i >= start && i<end) {
18397584545SPeter Brune       /* set this entry to be the gradient of the objective */
18497584545SPeter Brune       fv = (-ob2 + 6.*ob1 - 3.*ob -2.*ob3) / (6.*dx);
18597584545SPeter Brune       if (PetscAbsScalar(fv) > eps) {
18697584545SPeter Brune         ierr = VecSetValues(F,1,&i,&fv,INSERT_VALUES);CHKERRQ(ierr);
18797584545SPeter Brune       } else {
18897584545SPeter Brune         fv   = 0.;
18997584545SPeter Brune         ierr = VecSetValues(F,1,&i,&fv,INSERT_VALUES);CHKERRQ(ierr);
19097584545SPeter Brune       }
19197584545SPeter Brune     }
19297584545SPeter Brune   }
19397584545SPeter Brune 
19497584545SPeter Brune   ierr = VecDestroy(&Xh);CHKERRQ(ierr);
19597584545SPeter Brune 
19697584545SPeter Brune   ierr = VecAssemblyBegin(F);CHKERRQ(ierr);
19797584545SPeter Brune   ierr = VecAssemblyEnd(F);CHKERRQ(ierr);
19897584545SPeter Brune   PetscFunctionReturn(0);
19997584545SPeter Brune }
200