xref: /petsc/src/snes/interface/snesob.c (revision d8d19677bbccf95218448bee62e6b87f4513e133)
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 .      obj - real to hold the objective value
142a4ee8f2SPeter Brune -      ctx - optional user-defined objective context
152a4ee8f2SPeter Brune 
16878cb397SSatish Balay    Level: advanced
17878cb397SSatish Balay 
18f8b49ee9SBarry Smith .seealso:   SNESSetFunction(), SNESGetFunction(), SNESSetObjective(), SNESGetObjective(), SNESJacobianFunction, SNESFunction
19075cc632SBarry Smith M*/
20075cc632SBarry Smith 
21075cc632SBarry Smith /*@C
22eb23ba39SBarry Smith    SNESSetObjective - Sets the objective function minimized by some of the SNES linesearch methods.
23075cc632SBarry Smith 
24075cc632SBarry Smith    Logically Collective on SNES
25075cc632SBarry Smith 
26075cc632SBarry Smith    Input Parameters:
27075cc632SBarry Smith +  snes - the SNES context
28f8b49ee9SBarry 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 
34eb23ba39SBarry Smith    Note: This is not used in the SNESLINESEARCHCP line search.
35eb23ba39SBarry Smith 
36eb23ba39SBarry Smith          If not provided then this defaults to the two norm of the function evaluation (set with SNESSetFunction())
37075cc632SBarry Smith 
38075cc632SBarry Smith .seealso: SNESGetObjective(), SNESComputeObjective(), SNESSetFunction(), SNESSetJacobian(), SNESObjectiveFunction
392a4ee8f2SPeter Brune @*/
40f8b49ee9SBarry Smith PetscErrorCode  SNESSetObjective(SNES snes,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx)
412a4ee8f2SPeter Brune {
422a4ee8f2SPeter Brune   PetscErrorCode ierr;
432a4ee8f2SPeter Brune   DM             dm;
442a4ee8f2SPeter Brune 
452a4ee8f2SPeter Brune   PetscFunctionBegin;
462a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
472a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
48f8b49ee9SBarry Smith   ierr = DMSNESSetObjective(dm,obj,ctx);CHKERRQ(ierr);
492a4ee8f2SPeter Brune   PetscFunctionReturn(0);
502a4ee8f2SPeter Brune }
512a4ee8f2SPeter Brune 
522a4ee8f2SPeter Brune /*@C
532a4ee8f2SPeter Brune    SNESGetObjective - Returns the objective function.
542a4ee8f2SPeter Brune 
552a4ee8f2SPeter Brune    Not Collective
562a4ee8f2SPeter Brune 
572a4ee8f2SPeter Brune    Input Parameter:
582a4ee8f2SPeter Brune .  snes - the SNES context
592a4ee8f2SPeter Brune 
60*d8d19677SJose E. Roman    Output Parameters:
61f8b49ee9SBarry Smith +  obj - objective evaluation routine (or NULL); see SNESObjectFunction for details
620298fd71SBarry Smith -  ctx - the function context (or NULL)
632a4ee8f2SPeter Brune 
642a4ee8f2SPeter Brune    Level: advanced
652a4ee8f2SPeter Brune 
662a4ee8f2SPeter Brune .seealso: SNESSetObjective(), SNESGetSolution()
672a4ee8f2SPeter Brune @*/
68f8b49ee9SBarry Smith PetscErrorCode SNESGetObjective(SNES snes,PetscErrorCode (**obj)(SNES,Vec,PetscReal*,void*),void **ctx)
692a4ee8f2SPeter Brune {
702a4ee8f2SPeter Brune   PetscErrorCode ierr;
712a4ee8f2SPeter Brune   DM             dm;
725fd66863SKarl Rupp 
732a4ee8f2SPeter Brune   PetscFunctionBegin;
742a4ee8f2SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
752a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
76f8b49ee9SBarry Smith   ierr = DMSNESGetObjective(dm,obj,ctx);CHKERRQ(ierr);
772a4ee8f2SPeter Brune   PetscFunctionReturn(0);
782a4ee8f2SPeter Brune }
792a4ee8f2SPeter Brune 
802a4ee8f2SPeter Brune /*@C
812a4ee8f2SPeter Brune    SNESComputeObjective - Computes the objective.
822a4ee8f2SPeter Brune 
832a4ee8f2SPeter Brune    Collective on SNES
842a4ee8f2SPeter Brune 
85*d8d19677SJose E. Roman    Input Parameters:
862a4ee8f2SPeter Brune +  snes - the SNES context
872a4ee8f2SPeter Brune -  X    - the state vector
882a4ee8f2SPeter Brune 
892a4ee8f2SPeter Brune    Output Parameter:
902a4ee8f2SPeter Brune .  ob   - the objective value
912a4ee8f2SPeter Brune 
922a4ee8f2SPeter Brune    Level: advanced
932a4ee8f2SPeter Brune 
942a4ee8f2SPeter Brune .seealso: SNESSetObjective(), SNESGetSolution()
952a4ee8f2SPeter Brune @*/
962a4ee8f2SPeter Brune PetscErrorCode SNESComputeObjective(SNES snes,Vec X,PetscReal *ob)
972a4ee8f2SPeter Brune {
982a4ee8f2SPeter Brune   PetscErrorCode ierr;
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);
1062a4ee8f2SPeter Brune   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
107942e3340SBarry Smith   ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr);
10822c6f798SBarry Smith   if (sdm->ops->computeobjective) {
10994db00ebSBarry Smith     ierr = PetscLogEventBegin(SNES_ObjectiveEval,snes,X,0,0);CHKERRQ(ierr);
11022c6f798SBarry Smith     ierr = (sdm->ops->computeobjective)(snes,X,ob,sdm->objectivectx);CHKERRQ(ierr);
11194db00ebSBarry Smith     ierr = PetscLogEventEnd(SNES_ObjectiveEval,snes,X,0,0);CHKERRQ(ierr);
112ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetObjective() before SNESComputeObjective().");
1132a4ee8f2SPeter Brune   PetscFunctionReturn(0);
1142a4ee8f2SPeter Brune }
11597584545SPeter Brune 
1165c3e6ab7SPeter Brune /*@C
1175c3e6ab7SPeter Brune    SNESObjectiveComputeFunctionDefaultFD - Computes the gradient of a user provided objective
1185c3e6ab7SPeter Brune 
1195c3e6ab7SPeter Brune    Collective on SNES
1205c3e6ab7SPeter Brune 
121*d8d19677SJose E. Roman    Input Parameters:
1225c3e6ab7SPeter Brune +  snes - the SNES context
1235c3e6ab7SPeter Brune .  X    - the state vector
1245c3e6ab7SPeter Brune -  ctx  - the (ignored) function context
1255c3e6ab7SPeter Brune 
1265c3e6ab7SPeter Brune    Output Parameter:
1275c3e6ab7SPeter Brune .  F   - the function value
1285c3e6ab7SPeter Brune 
1295c3e6ab7SPeter Brune    Options Database Key:
1305c3e6ab7SPeter Brune +  -snes_fd_function_eps - The differencing parameter
1315c3e6ab7SPeter Brune -  -snes_fd_function - Compute function from user provided objective with finite difference
1325c3e6ab7SPeter Brune 
1335c3e6ab7SPeter Brune    Notes:
1345c3e6ab7SPeter Brune    SNESObjectiveComputeFunctionDefaultFD is similar in character to SNESComputeJacobianDefault.
1355c3e6ab7SPeter Brune    Therefore, it should be used for debugging purposes only.  Using it in conjunction with
1365c3e6ab7SPeter Brune    SNESComputeJacobianDefault is excessively costly and produces a Jacobian that is quite
1375c3e6ab7SPeter Brune    noisy.  This is often necessary, but should be done with a grain of salt, 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 
1445c3e6ab7SPeter Brune .seealso: SNESSetFunction(), SNESComputeObjective(), SNESComputeJacobianDefault()
1455c3e6ab7SPeter Brune @*/
1468d359177SBarry Smith PetscErrorCode SNESObjectiveComputeFunctionDefaultFD(SNES snes,Vec X,Vec F,void *ctx)
1470adebc6cSBarry Smith {
14897584545SPeter Brune   Vec            Xh;
14997584545SPeter Brune   PetscErrorCode ierr;
15097584545SPeter Brune   PetscInt       i,N,start,end;
15197584545SPeter Brune   PetscReal      ob,ob1,ob2,ob3,fob,dx,eps=1e-6;
15297584545SPeter Brune   PetscScalar    fv,xv;
15397584545SPeter Brune 
15497584545SPeter Brune   PetscFunctionBegin;
15597584545SPeter Brune   ierr = VecDuplicate(X,&Xh);CHKERRQ(ierr);
156302440fdSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,"Differencing parameters","SNES");CHKERRQ(ierr);
1570298fd71SBarry Smith   ierr = PetscOptionsReal("-snes_fd_function_eps","Tolerance for nonzero entries in fd function","None",eps,&eps,NULL);CHKERRQ(ierr);
158302440fdSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
15997584545SPeter Brune   ierr = VecSet(F,0.);CHKERRQ(ierr);
16097584545SPeter Brune 
16197584545SPeter Brune   ierr = VecNorm(X,NORM_2,&fob);CHKERRQ(ierr);
16297584545SPeter Brune 
16397584545SPeter Brune   ierr = VecGetSize(X,&N);CHKERRQ(ierr);
16497584545SPeter Brune   ierr = VecGetOwnershipRange(X,&start,&end);CHKERRQ(ierr);
16597584545SPeter Brune   ierr = SNESComputeObjective(snes,X,&ob);CHKERRQ(ierr);
16697584545SPeter Brune 
167f5af7f23SKarl Rupp   if (fob > 0.) dx =1e-6*fob;
168f5af7f23SKarl Rupp   else dx = 1e-6;
16997584545SPeter Brune 
17097584545SPeter Brune   for (i=0; i<N; i++) {
17197584545SPeter Brune     /* compute the 1st value */
17297584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
17397584545SPeter Brune     if (i>= start && i<end) {
17497584545SPeter Brune       xv   = dx;
17597584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
17697584545SPeter Brune     }
17797584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
17897584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
17997584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob1);CHKERRQ(ierr);
18097584545SPeter Brune 
18197584545SPeter Brune     /* compute the 2nd value */
18297584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
18397584545SPeter Brune     if (i>= start && i<end) {
18497584545SPeter Brune       xv   = 2.*dx;
18597584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
18697584545SPeter Brune     }
18797584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
18897584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
18997584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob2);CHKERRQ(ierr);
19097584545SPeter Brune 
19197584545SPeter Brune     /* compute the 3rd value */
19297584545SPeter Brune     ierr = VecCopy(X,Xh);CHKERRQ(ierr);
19397584545SPeter Brune     if (i>= start && i<end) {
19497584545SPeter Brune       xv   = -dx;
19597584545SPeter Brune       ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr);
19697584545SPeter Brune     }
19797584545SPeter Brune     ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr);
19897584545SPeter Brune     ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr);
19997584545SPeter Brune     ierr = SNESComputeObjective(snes,Xh,&ob3);CHKERRQ(ierr);
20097584545SPeter Brune 
20197584545SPeter Brune     if (i >= start && i<end) {
20297584545SPeter Brune       /* set this entry to be the gradient of the objective */
20397584545SPeter Brune       fv = (-ob2 + 6.*ob1 - 3.*ob -2.*ob3) / (6.*dx);
20497584545SPeter Brune       if (PetscAbsScalar(fv) > eps) {
20597584545SPeter Brune         ierr = VecSetValues(F,1,&i,&fv,INSERT_VALUES);CHKERRQ(ierr);
20697584545SPeter Brune       } else {
20797584545SPeter Brune         fv   = 0.;
20897584545SPeter Brune         ierr = VecSetValues(F,1,&i,&fv,INSERT_VALUES);CHKERRQ(ierr);
20997584545SPeter Brune       }
21097584545SPeter Brune     }
21197584545SPeter Brune   }
21297584545SPeter Brune 
21397584545SPeter Brune   ierr = VecDestroy(&Xh);CHKERRQ(ierr);
21497584545SPeter Brune 
21597584545SPeter Brune   ierr = VecAssemblyBegin(F);CHKERRQ(ierr);
21697584545SPeter Brune   ierr = VecAssemblyEnd(F);CHKERRQ(ierr);
21797584545SPeter Brune   PetscFunctionReturn(0);
21897584545SPeter Brune }
219