Actual source code: fdiff.c

petsc-dev 2014-02-02
Report Typos and Errors
  1: #include <petsctao.h>         /*I  "petsctao.h"  I*/
  2: #include <petsc-private/taoimpl.h>
  3: #include <petscsnes.h>

  5: /*
  6:    For finited difference computations of the Hessian, we use PETSc's SNESComputeJacobianDefault
  7: */

 11: static PetscErrorCode Fsnes(SNES snes ,Vec X,Vec G,void*ctx)
 12: {
 14:   Tao            tao = (Tao)ctx;

 18:   ierr=TaoComputeGradient(tao,X,G);
 19:   return(0);
 20: }

 24: /*@C
 25:   TaoDefaultComputeGradient - computes the gradient using finite differences.

 27:   Collective on Tao

 29:   Input Parameters:
 30: + tao - the Tao context
 31: . X - compute gradient at this point
 32: - dummy - not used

 34:   Output Parameters:
 35: . G - Gradient Vector

 37:    Options Database Key:
 38: +  -tao_fd_gradient - Activates TaoDefaultComputeGradient()
 39: -  -tao_fd_delta <delta> - change in x used to calculate finite differences

 41:    Level: advanced

 43:    Note:
 44:    This routine is slow and expensive, and is not currently optimized
 45:    to take advantage of sparsity in the problem.  Although
 46:    TaoAppDefaultComputeGradient is not recommended for general use
 47:    in large-scale applications, It can be useful in checking the
 48:    correctness of a user-provided gradient.  Use the tao method "tao_fd_test"
 49:    to get an indication of whether your gradient is correct.


 52:    Note:
 53:    This finite difference gradient evaluation can be set using the routine TaoSetGradientRoutine() or by using the command line option -tao_fd_gradient

 55: .seealso: TaoSetGradientRoutine()

 57: @*/
 58: PetscErrorCode TaoDefaultComputeGradient(Tao tao,Vec X,Vec G,void *dummy)
 59: {
 60:   PetscReal      *g;
 61:   PetscReal      f, f2;
 63:   PetscInt       low,high,N,i;
 64:   PetscBool      flg;
 65:   PetscReal      h=PETSC_SQRT_MACHINE_EPSILON;

 68:   TaoComputeObjective(tao, X,&f);
 69:   PetscOptionsGetReal(NULL,"-tao_fd_delta",&h,&flg);
 70:   VecGetSize(X,&N);
 71:   VecGetOwnershipRange(X,&low,&high);
 72:   VecGetArray(G,&g);
 73:   for (i=0;i<N;i++) {
 74:       VecSetValue(X,i,h,ADD_VALUES);
 75:       VecAssemblyBegin(X);
 76:       VecAssemblyEnd(X);

 78:       TaoComputeObjective(tao,X,&f2);

 80:       VecSetValue(X,i,-h,ADD_VALUES);
 81:       VecAssemblyBegin(X);
 82:       VecAssemblyEnd(X);

 84:       if (i>=low && i<high) {
 85:           g[i-low]=(f2-f)/h;
 86:       }
 87:   }
 88:   VecRestoreArray(G,&g);
 89:   return(0);
 90: }

 94: /*@C
 95:    TaoDefaultComputeHessian - Computes the Hessian using finite differences.

 97:    Collective on Tao

 99:    Input Parameters:
100: +  tao - the Tao context
101: .  V - compute Hessian at this point
102: -  dummy - not used

104:    Output Parameters:
105: +  H - Hessian matrix (not altered in this routine)
106: .  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
107: -  flag - flag indicating whether the matrix sparsity structure has changed

109:    Options Database Key:
110: +  -tao_fd - Activates TaoDefaultComputeHessian()
111: -  -tao_view_hessian - view the hessian after each evaluation using PETSC_VIEWER_STDOUT_WORLD

113:    Level: advanced

115:    Notes:
116:    This routine is slow and expensive, and is not currently optimized
117:    to take advantage of sparsity in the problem.  Although
118:    TaoDefaultComputeHessian() is not recommended for general use
119:    in large-scale applications, It can be useful in checking the
120:    correctness of a user-provided Hessian.



124: .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessianColor(), SNESComputeJacobianDefault(), TaoSetGradientRoutine(), TaoDefaultComputeGradient()

126: @*/
127: PetscErrorCode TaoDefaultComputeHessian(Tao tao,Vec V,Mat *H,Mat *B,MatStructure *flag,void *dummy){
128:   PetscErrorCode       ierr;
129:   MPI_Comm             comm;
130:   Vec                  G;
131:   SNES                 snes;

135:   VecDuplicate(V,&G);

137:   PetscInfo(tao,"TAO Using finite differences w/o coloring to compute Hessian matrix\n");

139:   TaoComputeGradient(tao,V,G);

141:   PetscObjectGetComm((PetscObject)(*H),&comm);
142:   SNESCreate(comm,&snes);

144:   SNESSetFunction(snes,G,Fsnes,tao);
145:   SNESComputeJacobianDefault(snes,V,H,B,flag,tao);
146:   SNESDestroy(&snes);
147:   VecDestroy(&G);
148:   return(0);
149: }

153: /*@C
154:    TaoDefaultComputeHessianColor - Computes the Hessian using colored finite differences.

156:    Collective on Tao

158:    Input Parameters:
159: +  tao - the Tao context
160: .  V - compute Hessian at this point
161: -  ctx - the PetscColoring object (must be of type MatFDColoring)

163:    Output Parameters:
164: +  H - Hessian matrix (not altered in this routine)
165: .  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
166: -  flag - flag indicating whether the matrix sparsity structure has changed

168:    Level: advanced


171: .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessian(),SNESComputeJacobianDefaultColor(), TaoSetGradientRoutine()

173: @*/
174: PetscErrorCode TaoDefaultComputeHessianColor(Tao tao, Vec V, Mat *H,Mat *B,MatStructure *flag,void *ctx)
175: {
176:   PetscErrorCode      ierr;
177:   MatFDColoring       coloring = (MatFDColoring)ctx;

181:   *flag = SAME_NONZERO_PATTERN;

183:   ierr=PetscInfo(tao,"TAO computing matrix using finite differences Hessian and coloring\n");
184:   MatFDColoringApply(*B,coloring,V,flag,ctx);

186:   if (*H != *B) {
187:       MatAssemblyBegin(*H, MAT_FINAL_ASSEMBLY);
188:       MatAssemblyEnd(*H, MAT_FINAL_ASSEMBLY);
189:   }
190:   return(0);
191: }