xref: /petsc/src/snes/interface/snesob.c (revision 87d70a26ceae4e170d98ea04762233dc00f2c8df)
12a4ee8f2SPeter Brune #include <petsc-private/snesimpl.h>
22a4ee8f2SPeter Brune 
32a4ee8f2SPeter Brune #undef __FUNCT__
42a4ee8f2SPeter Brune #define __FUNCT__ "SNESSetObjective"
52a4ee8f2SPeter Brune /*@C
62a4ee8f2SPeter Brune    SNESSetObjective - Sets the objective function minimized by
72a4ee8f2SPeter Brune    the SNES methods.
82a4ee8f2SPeter Brune 
92a4ee8f2SPeter Brune    Logically Collective on SNES
102a4ee8f2SPeter Brune 
112a4ee8f2SPeter Brune    Input Parameters:
122a4ee8f2SPeter Brune +  snes - the SNES context
132a4ee8f2SPeter Brune .  func - objective evaluation routine
142a4ee8f2SPeter Brune -  ctx - [optional] user-defined context for private data for the
152a4ee8f2SPeter Brune          function evaluation routine (may be PETSC_NULL)
162a4ee8f2SPeter Brune 
172a4ee8f2SPeter Brune    Calling sequence of func:
182a4ee8f2SPeter Brune $    func (SNES snes,Vec x,PetscReal *obj,void *ctx);
192a4ee8f2SPeter Brune 
202a4ee8f2SPeter Brune +  snes - the SNES context
212a4ee8f2SPeter Brune .  X - solution
222a4ee8f2SPeter Brune .  F - current function/gradient
232a4ee8f2SPeter Brune .  obj - real to hold the objective value
242a4ee8f2SPeter Brune -  ctx - optional user-defined objective context
252a4ee8f2SPeter Brune 
262a4ee8f2SPeter Brune    Level: beginner
272a4ee8f2SPeter Brune 
282a4ee8f2SPeter Brune .keywords: SNES, nonlinear, set, objective
292a4ee8f2SPeter Brune 
302a4ee8f2SPeter Brune .seealso: SNESGetObjective(), SNESComputeObjective(), SNESSetFunction(), SNESSetJacobian()
312a4ee8f2SPeter Brune @*/
322a4ee8f2SPeter Brune PetscErrorCode  SNESSetObjective(SNES snes,SNESObjective func,void *ctx)
332a4ee8f2SPeter Brune {
342a4ee8f2SPeter Brune   PetscErrorCode ierr;
352a4ee8f2SPeter Brune   DM             dm;
362a4ee8f2SPeter Brune 
372a4ee8f2SPeter Brune   PetscFunctionBegin;
382a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
392a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
402a4ee8f2SPeter Brune   ierr = DMSNESSetObjective(dm,func,ctx);CHKERRQ(ierr);
412a4ee8f2SPeter Brune   PetscFunctionReturn(0);
422a4ee8f2SPeter Brune }
432a4ee8f2SPeter Brune 
442a4ee8f2SPeter Brune #undef __FUNCT__
452a4ee8f2SPeter Brune #define __FUNCT__ "SNESGetObjective"
462a4ee8f2SPeter Brune /*@C
472a4ee8f2SPeter Brune    SNESGetObjective - Returns the objective function.
482a4ee8f2SPeter Brune 
492a4ee8f2SPeter Brune    Not Collective
502a4ee8f2SPeter Brune 
512a4ee8f2SPeter Brune    Input Parameter:
522a4ee8f2SPeter Brune .  snes - the SNES context
532a4ee8f2SPeter Brune 
542a4ee8f2SPeter Brune    Output Parameter:
552a4ee8f2SPeter Brune +  func - the function (or PETSC_NULL)
562a4ee8f2SPeter Brune -  ctx - the function context (or PETSC_NULL)
572a4ee8f2SPeter Brune 
582a4ee8f2SPeter Brune    Level: advanced
592a4ee8f2SPeter Brune 
602a4ee8f2SPeter Brune .keywords: SNES, nonlinear, get, objective
612a4ee8f2SPeter Brune 
622a4ee8f2SPeter Brune .seealso: SNESSetObjective(), SNESGetSolution()
632a4ee8f2SPeter Brune @*/
642a4ee8f2SPeter Brune PetscErrorCode SNESGetObjective(SNES snes,SNESObjective *func,void **ctx)
652a4ee8f2SPeter Brune {
662a4ee8f2SPeter Brune   PetscErrorCode ierr;
672a4ee8f2SPeter Brune   DM             dm;
682a4ee8f2SPeter Brune   PetscFunctionBegin;
692a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
702a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
712a4ee8f2SPeter Brune   ierr = DMSNESGetObjective(dm,func,ctx);CHKERRQ(ierr);
722a4ee8f2SPeter Brune   PetscFunctionReturn(0);
732a4ee8f2SPeter Brune }
742a4ee8f2SPeter Brune 
752a4ee8f2SPeter Brune #undef __FUNCT__
762a4ee8f2SPeter Brune #define __FUNCT__ "SNESComputeObjective"
772a4ee8f2SPeter Brune /*@C
782a4ee8f2SPeter Brune    SNESComputeObjective - Computes the objective.
792a4ee8f2SPeter Brune 
802a4ee8f2SPeter Brune    Collective on SNES
812a4ee8f2SPeter Brune 
822a4ee8f2SPeter Brune    Input Parameter:
832a4ee8f2SPeter Brune +  snes - the SNES context
842a4ee8f2SPeter Brune -  X    - the state vector
852a4ee8f2SPeter Brune 
862a4ee8f2SPeter Brune    Output Parameter:
872a4ee8f2SPeter Brune .  ob   - the objective value
882a4ee8f2SPeter Brune 
892a4ee8f2SPeter Brune    Level: advanced
902a4ee8f2SPeter Brune 
912a4ee8f2SPeter Brune .keywords: SNES, nonlinear, compute, objective
922a4ee8f2SPeter Brune 
932a4ee8f2SPeter Brune .seealso: SNESSetObjective(), SNESGetSolution()
942a4ee8f2SPeter Brune @*/
952a4ee8f2SPeter Brune PetscErrorCode SNESComputeObjective(SNES snes,Vec X,PetscReal *ob)
962a4ee8f2SPeter Brune {
972a4ee8f2SPeter Brune   PetscErrorCode ierr;
982a4ee8f2SPeter Brune   DM             dm;
99942e3340SBarry Smith   DMSNES         sdm;
1002a4ee8f2SPeter Brune   PetscFunctionBegin;
1012a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
1022a4ee8f2SPeter Brune   PetscValidHeaderSpecific(X,VEC_CLASSID,2);
1032a4ee8f2SPeter Brune   PetscValidPointer(ob,3);
1042a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
105942e3340SBarry Smith   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
10622c6f798SBarry Smith   if (sdm->ops->computeobjective) {
10722c6f798SBarry Smith     ierr = (sdm->ops->computeobjective)(snes,X,ob,sdm->objectivectx);CHKERRQ(ierr);
108*87d70a26SPeter Brune   } else {
109*87d70a26SPeter Brune     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetObjective() before SNESComputeObjective().");
110*87d70a26SPeter Brune   }
1112a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1122a4ee8f2SPeter Brune }
11397584545SPeter Brune 
11497584545SPeter Brune 
11597584545SPeter Brune #undef __FUNCT__
11697584545SPeter Brune #define __FUNCT__ "SNESDefaultObjectiveComputeFunctionFD"
11797584545SPeter Brune PetscErrorCode SNESDefaultObjectiveComputeFunctionFD(SNES snes,Vec X,Vec F,void *ctx) {
11897584545SPeter Brune   /* Quadratically interpolates the change in the objective based upon a change in a particular direction */
11997584545SPeter Brune 
12097584545SPeter Brune   Vec            Xh;
12197584545SPeter Brune   PetscErrorCode ierr;
12297584545SPeter Brune   PetscInt       i,N,start,end;
12397584545SPeter Brune   PetscReal      ob,ob1,ob2,ob3,fob,dx,eps=1e-6;
12497584545SPeter Brune   PetscScalar    fv,xv;
12597584545SPeter Brune 
12697584545SPeter Brune   PetscFunctionBegin;
12797584545SPeter Brune   ierr = VecDuplicate(X,&Xh);CHKERRQ(ierr);
12897584545SPeter Brune   ierr = PetscOptionsReal("-snes_fd_function_eps","Tolerance for nonzero entries in fd function","None",eps,&eps,PETSC_NULL);CHKERRQ(ierr);
12997584545SPeter Brune   ierr = VecSet(F,0.);CHKERRQ(ierr);
13097584545SPeter Brune 
13197584545SPeter Brune   ierr = VecNorm(X,NORM_2,&fob);CHKERRQ(ierr);
13297584545SPeter Brune 
13397584545SPeter Brune   ierr = VecGetSize(X,&N);CHKERRQ(ierr);
13497584545SPeter Brune   ierr = VecGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
13597584545SPeter Brune   ierr = SNESComputeObjective(snes,X,&ob);CHKERRQ(ierr);
13697584545SPeter Brune 
13797584545SPeter Brune   if (fob > 0.) {
13897584545SPeter Brune     dx =1e-6*fob;
13997584545SPeter Brune   } else {
14097584545SPeter Brune     dx = 1e-6;
14197584545SPeter Brune   }
14297584545SPeter Brune 
14397584545SPeter Brune   for (i=0; i<N; i++) {
14497584545SPeter Brune     /* compute the 1st value */
14597584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
14697584545SPeter Brune     if (i>= start && i<end) {
14797584545SPeter Brune       xv = dx;
14897584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
14997584545SPeter Brune     }
15097584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
15197584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
15297584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob1);CHKERRQ(ierr);
15397584545SPeter Brune 
15497584545SPeter Brune     /* compute the 2nd value */
15597584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
15697584545SPeter Brune     if (i>= start && i<end) {
15797584545SPeter Brune       xv = 2.*dx;
15897584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
15997584545SPeter Brune     }
16097584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
16197584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
16297584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob2);CHKERRQ(ierr);
16397584545SPeter Brune 
16497584545SPeter Brune     /* compute the 3rd value */
16597584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
16697584545SPeter Brune     if (i>= start && i<end) {
16797584545SPeter Brune       xv = -dx;
16897584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
16997584545SPeter Brune     }
17097584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
17197584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
17297584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob3);CHKERRQ(ierr);
17397584545SPeter Brune 
17497584545SPeter Brune     if (i >= start && i<end) {
17597584545SPeter Brune       /* set this entry to be the gradient of the objective */
17697584545SPeter Brune       fv = (-ob2 + 6.*ob1 - 3.*ob -2.*ob3) / (6.*dx);
17797584545SPeter Brune       if (PetscAbsScalar(fv) > eps) {
17897584545SPeter Brune         ierr = VecSetValues(F,1,&i,&fv,INSERT_VALUES);CHKERRQ(ierr);
17997584545SPeter Brune       } else {
18097584545SPeter Brune         fv = 0.;
18197584545SPeter Brune         ierr = VecSetValues(F,1,&i,&fv,INSERT_VALUES);CHKERRQ(ierr);
18297584545SPeter Brune       }
18397584545SPeter Brune     }
18497584545SPeter Brune   }
18597584545SPeter Brune 
18697584545SPeter Brune   ierr = VecDestroy(&Xh);CHKERRQ(ierr);
18797584545SPeter Brune 
18897584545SPeter Brune   ierr = VecAssemblyBegin(F);CHKERRQ(ierr);
18997584545SPeter Brune   ierr = VecAssemblyEnd(F);CHKERRQ(ierr);
19097584545SPeter Brune   PetscFunctionReturn(0);
19197584545SPeter Brune }
192