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: 7*aaa7dc30SBarry 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 #undef __FUNCT__ 24075cc632SBarry Smith #define __FUNCT__ "SNESSetObjective" 25075cc632SBarry Smith /*@C 26075cc632SBarry Smith SNESSetObjective - Sets the objective function minimized by the SNES methods. 27075cc632SBarry Smith 28075cc632SBarry Smith Logically Collective on SNES 29075cc632SBarry Smith 30075cc632SBarry Smith Input Parameters: 31075cc632SBarry Smith + snes - the SNES context 32f8b49ee9SBarry Smith . obj - objective evaluation routine; see SNESObjectiveFunction for details 33075cc632SBarry Smith - ctx - [optional] user-defined context for private data for the 340298fd71SBarry Smith function evaluation routine (may be NULL) 35075cc632SBarry Smith 362a4ee8f2SPeter Brune Level: beginner 372a4ee8f2SPeter Brune 38075cc632SBarry Smith Note: If not provided then this defaults to the two norm of the function evaluation (set with SNESSetFunction()) 39075cc632SBarry Smith 402a4ee8f2SPeter Brune .keywords: SNES, nonlinear, set, objective 412a4ee8f2SPeter Brune 42075cc632SBarry Smith .seealso: SNESGetObjective(), SNESComputeObjective(), SNESSetFunction(), SNESSetJacobian(), SNESObjectiveFunction 432a4ee8f2SPeter Brune @*/ 44f8b49ee9SBarry Smith PetscErrorCode SNESSetObjective(SNES snes,PetscErrorCode (*obj)(SNES,Vec,PetscReal*,void*),void *ctx) 452a4ee8f2SPeter Brune { 462a4ee8f2SPeter Brune PetscErrorCode ierr; 472a4ee8f2SPeter Brune DM dm; 482a4ee8f2SPeter Brune 492a4ee8f2SPeter Brune PetscFunctionBegin; 502a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 512a4ee8f2SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 52f8b49ee9SBarry Smith ierr = DMSNESSetObjective(dm,obj,ctx);CHKERRQ(ierr); 532a4ee8f2SPeter Brune PetscFunctionReturn(0); 542a4ee8f2SPeter Brune } 552a4ee8f2SPeter Brune 562a4ee8f2SPeter Brune #undef __FUNCT__ 572a4ee8f2SPeter Brune #define __FUNCT__ "SNESGetObjective" 582a4ee8f2SPeter Brune /*@C 592a4ee8f2SPeter Brune SNESGetObjective - Returns the objective function. 602a4ee8f2SPeter Brune 612a4ee8f2SPeter Brune Not Collective 622a4ee8f2SPeter Brune 632a4ee8f2SPeter Brune Input Parameter: 642a4ee8f2SPeter Brune . snes - the SNES context 652a4ee8f2SPeter Brune 662a4ee8f2SPeter Brune Output Parameter: 67f8b49ee9SBarry Smith + obj - objective evaluation routine (or NULL); see SNESObjectFunction for details 680298fd71SBarry Smith - ctx - the function context (or NULL) 692a4ee8f2SPeter Brune 702a4ee8f2SPeter Brune Level: advanced 712a4ee8f2SPeter Brune 722a4ee8f2SPeter Brune .keywords: SNES, nonlinear, get, objective 732a4ee8f2SPeter Brune 742a4ee8f2SPeter Brune .seealso: SNESSetObjective(), SNESGetSolution() 752a4ee8f2SPeter Brune @*/ 76f8b49ee9SBarry Smith PetscErrorCode SNESGetObjective(SNES snes,PetscErrorCode (**obj)(SNES,Vec,PetscReal*,void*),void **ctx) 772a4ee8f2SPeter Brune { 782a4ee8f2SPeter Brune PetscErrorCode ierr; 792a4ee8f2SPeter Brune DM dm; 805fd66863SKarl Rupp 812a4ee8f2SPeter Brune PetscFunctionBegin; 822a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 832a4ee8f2SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 84f8b49ee9SBarry Smith ierr = DMSNESGetObjective(dm,obj,ctx);CHKERRQ(ierr); 852a4ee8f2SPeter Brune PetscFunctionReturn(0); 862a4ee8f2SPeter Brune } 872a4ee8f2SPeter Brune 882a4ee8f2SPeter Brune #undef __FUNCT__ 892a4ee8f2SPeter Brune #define __FUNCT__ "SNESComputeObjective" 902a4ee8f2SPeter Brune /*@C 912a4ee8f2SPeter Brune SNESComputeObjective - Computes the objective. 922a4ee8f2SPeter Brune 932a4ee8f2SPeter Brune Collective on SNES 942a4ee8f2SPeter Brune 952a4ee8f2SPeter Brune Input Parameter: 962a4ee8f2SPeter Brune + snes - the SNES context 972a4ee8f2SPeter Brune - X - the state vector 982a4ee8f2SPeter Brune 992a4ee8f2SPeter Brune Output Parameter: 1002a4ee8f2SPeter Brune . ob - the objective value 1012a4ee8f2SPeter Brune 1022a4ee8f2SPeter Brune Level: advanced 1032a4ee8f2SPeter Brune 1042a4ee8f2SPeter Brune .keywords: SNES, nonlinear, compute, objective 1052a4ee8f2SPeter Brune 1062a4ee8f2SPeter Brune .seealso: SNESSetObjective(), SNESGetSolution() 1072a4ee8f2SPeter Brune @*/ 1082a4ee8f2SPeter Brune PetscErrorCode SNESComputeObjective(SNES snes,Vec X,PetscReal *ob) 1092a4ee8f2SPeter Brune { 1102a4ee8f2SPeter Brune PetscErrorCode ierr; 1112a4ee8f2SPeter Brune DM dm; 112942e3340SBarry Smith DMSNES sdm; 113f23aa3ddSBarry Smith 1142a4ee8f2SPeter Brune PetscFunctionBegin; 1152a4ee8f2SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 1162a4ee8f2SPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 1172a4ee8f2SPeter Brune PetscValidPointer(ob,3); 1182a4ee8f2SPeter Brune ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 119942e3340SBarry Smith ierr = DMGetDMSNES(dm,&sdm);CHKERRQ(ierr); 12022c6f798SBarry Smith if (sdm->ops->computeobjective) { 12122c6f798SBarry Smith ierr = (sdm->ops->computeobjective)(snes,X,ob,sdm->objectivectx);CHKERRQ(ierr); 122ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetObjective() before SNESComputeObjective()."); 1232a4ee8f2SPeter Brune PetscFunctionReturn(0); 1242a4ee8f2SPeter Brune } 12597584545SPeter Brune 12697584545SPeter Brune 12797584545SPeter Brune #undef __FUNCT__ 1288d359177SBarry Smith #define __FUNCT__ "SNESObjectiveComputeFunctionDefaultFD" 1295c3e6ab7SPeter Brune /*@C 1305c3e6ab7SPeter Brune SNESObjectiveComputeFunctionDefaultFD - Computes the gradient of a user provided objective 1315c3e6ab7SPeter Brune 1325c3e6ab7SPeter Brune Collective on SNES 1335c3e6ab7SPeter Brune 1345c3e6ab7SPeter Brune Input Parameter: 1355c3e6ab7SPeter Brune + snes - the SNES context 1365c3e6ab7SPeter Brune . X - the state vector 1375c3e6ab7SPeter Brune - ctx - the (ignored) function context 1385c3e6ab7SPeter Brune 1395c3e6ab7SPeter Brune Output Parameter: 1405c3e6ab7SPeter Brune . F - the function value 1415c3e6ab7SPeter Brune 1425c3e6ab7SPeter Brune Options Database Key: 1435c3e6ab7SPeter Brune + -snes_fd_function_eps - The differencing parameter 1445c3e6ab7SPeter Brune - -snes_fd_function - Compute function from user provided objective with finite difference 1455c3e6ab7SPeter Brune 1465c3e6ab7SPeter Brune Notes: 1475c3e6ab7SPeter Brune SNESObjectiveComputeFunctionDefaultFD is similar in character to SNESComputeJacobianDefault. 1485c3e6ab7SPeter Brune Therefore, it should be used for debugging purposes only. Using it in conjunction with 1495c3e6ab7SPeter Brune SNESComputeJacobianDefault is excessively costly and produces a Jacobian that is quite 1505c3e6ab7SPeter Brune noisy. This is often necessary, but should be done with a grain of salt, even when debugging 1515c3e6ab7SPeter Brune small problems. 1525c3e6ab7SPeter Brune 1535c3e6ab7SPeter Brune Note that this uses quadratic interpolation of the objective to form each value in the function. 1545c3e6ab7SPeter Brune 1555c3e6ab7SPeter Brune Level: advanced 1565c3e6ab7SPeter Brune 1575c3e6ab7SPeter Brune .keywords: SNES, objective, debugging, finite differences, function 1585c3e6ab7SPeter Brune 1595c3e6ab7SPeter Brune .seealso: SNESSetFunction(), SNESComputeObjective(), SNESComputeJacobianDefault() 1605c3e6ab7SPeter Brune @*/ 1618d359177SBarry Smith PetscErrorCode SNESObjectiveComputeFunctionDefaultFD(SNES snes,Vec X,Vec F,void *ctx) 1620adebc6cSBarry Smith { 16397584545SPeter Brune Vec Xh; 16497584545SPeter Brune PetscErrorCode ierr; 16597584545SPeter Brune PetscInt i,N,start,end; 16697584545SPeter Brune PetscReal ob,ob1,ob2,ob3,fob,dx,eps=1e-6; 16797584545SPeter Brune PetscScalar fv,xv; 16897584545SPeter Brune 16997584545SPeter Brune PetscFunctionBegin; 17097584545SPeter Brune ierr = VecDuplicate(X,&Xh);CHKERRQ(ierr); 1710298fd71SBarry Smith ierr = PetscOptionsReal("-snes_fd_function_eps","Tolerance for nonzero entries in fd function","None",eps,&eps,NULL);CHKERRQ(ierr); 17297584545SPeter Brune ierr = VecSet(F,0.);CHKERRQ(ierr); 17397584545SPeter Brune 17497584545SPeter Brune ierr = VecNorm(X,NORM_2,&fob);CHKERRQ(ierr); 17597584545SPeter Brune 17697584545SPeter Brune ierr = VecGetSize(X,&N);CHKERRQ(ierr); 17797584545SPeter Brune ierr = VecGetOwnershipRange(X,&start,&end);CHKERRQ(ierr); 17897584545SPeter Brune ierr = SNESComputeObjective(snes,X,&ob);CHKERRQ(ierr); 17997584545SPeter Brune 180f5af7f23SKarl Rupp if (fob > 0.) dx =1e-6*fob; 181f5af7f23SKarl Rupp else dx = 1e-6; 18297584545SPeter Brune 18397584545SPeter Brune for (i=0; i<N; i++) { 18497584545SPeter Brune /* compute the 1st value */ 18597584545SPeter Brune ierr = VecCopy(X,Xh);CHKERRQ(ierr); 18697584545SPeter Brune if (i>= start && i<end) { 18797584545SPeter Brune xv = 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,&ob1);CHKERRQ(ierr); 19397584545SPeter Brune 19497584545SPeter Brune /* compute the 2nd value */ 19597584545SPeter Brune ierr = VecCopy(X,Xh);CHKERRQ(ierr); 19697584545SPeter Brune if (i>= start && i<end) { 19797584545SPeter Brune xv = 2.*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,&ob2);CHKERRQ(ierr); 20397584545SPeter Brune 20497584545SPeter Brune /* compute the 3rd value */ 20597584545SPeter Brune ierr = VecCopy(X,Xh);CHKERRQ(ierr); 20697584545SPeter Brune if (i>= start && i<end) { 20797584545SPeter Brune xv = -dx; 20897584545SPeter Brune ierr = VecSetValues(Xh,1,&i,&xv,ADD_VALUES);CHKERRQ(ierr); 20997584545SPeter Brune } 21097584545SPeter Brune ierr = VecAssemblyBegin(Xh);CHKERRQ(ierr); 21197584545SPeter Brune ierr = VecAssemblyEnd(Xh);CHKERRQ(ierr); 21297584545SPeter Brune ierr = SNESComputeObjective(snes,Xh,&ob3);CHKERRQ(ierr); 21397584545SPeter Brune 21497584545SPeter Brune if (i >= start && i<end) { 21597584545SPeter Brune /* set this entry to be the gradient of the objective */ 21697584545SPeter Brune fv = (-ob2 + 6.*ob1 - 3.*ob -2.*ob3) / (6.*dx); 21797584545SPeter Brune if (PetscAbsScalar(fv) > eps) { 21897584545SPeter Brune ierr = VecSetValues(F,1,&i,&fv,INSERT_VALUES);CHKERRQ(ierr); 21997584545SPeter Brune } else { 22097584545SPeter Brune fv = 0.; 22197584545SPeter Brune ierr = VecSetValues(F,1,&i,&fv,INSERT_VALUES);CHKERRQ(ierr); 22297584545SPeter Brune } 22397584545SPeter Brune } 22497584545SPeter Brune } 22597584545SPeter Brune 22697584545SPeter Brune ierr = VecDestroy(&Xh);CHKERRQ(ierr); 22797584545SPeter Brune 22897584545SPeter Brune ierr = VecAssemblyBegin(F);CHKERRQ(ierr); 22997584545SPeter Brune ierr = VecAssemblyEnd(F);CHKERRQ(ierr); 23097584545SPeter Brune PetscFunctionReturn(0); 23197584545SPeter Brune } 232