xref: /petsc/src/tao/interface/fdiff.c (revision ffad99011bdf8bdff5e8540ef3c49b4fd8d6e6bb)
1ba92ff59SBarry Smith #include <petsctao.h>         /*I  "petsctao.h"  I*/
2ba92ff59SBarry Smith #include <petsc-private/taoimpl.h>
3aaa7dc30SBarry Smith #include <petscsnes.h>
4a7e14dcfSSatish Balay 
5a7e14dcfSSatish Balay /*
6a7e14dcfSSatish Balay    For finited difference computations of the Hessian, we use PETSc's SNESComputeJacobianDefault
7a7e14dcfSSatish Balay */
8a7e14dcfSSatish Balay 
9a7e14dcfSSatish Balay #undef __FUNCT__
10a7e14dcfSSatish Balay #define __FUNCT__ "Fsnes"
115ca45b2bSBarry Smith static PetscErrorCode Fsnes(SNES snes ,Vec X,Vec G,void*ctx)
125ca45b2bSBarry Smith {
13a7e14dcfSSatish Balay   PetscErrorCode ierr;
14441846f8SBarry Smith   Tao            tao = (Tao)ctx;
155ca45b2bSBarry Smith 
16a7e14dcfSSatish Balay   PetscFunctionBegin;
17441846f8SBarry Smith   PetscValidHeaderSpecific(ctx,TAO_CLASSID,4);
18a7e14dcfSSatish Balay   ierr=TaoComputeGradient(tao,X,G);CHKERRQ(ierr);
19a7e14dcfSSatish Balay   PetscFunctionReturn(0);
20a7e14dcfSSatish Balay }
21a7e14dcfSSatish Balay 
22a7e14dcfSSatish Balay #undef __FUNCT__
23a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultComputeGradient"
24a7e14dcfSSatish Balay /*@C
25a7e14dcfSSatish Balay   TaoDefaultComputeGradient - computes the gradient using finite differences.
26a7e14dcfSSatish Balay 
27441846f8SBarry Smith   Collective on Tao
28a7e14dcfSSatish Balay 
29a7e14dcfSSatish Balay   Input Parameters:
30441846f8SBarry Smith + tao - the Tao context
31a7e14dcfSSatish Balay . X - compute gradient at this point
32a7e14dcfSSatish Balay - dummy - not used
33a7e14dcfSSatish Balay 
34a7e14dcfSSatish Balay   Output Parameters:
35a7e14dcfSSatish Balay . G - Gradient Vector
36a7e14dcfSSatish Balay 
37a7e14dcfSSatish Balay    Options Database Key:
38a7e14dcfSSatish Balay +  -tao_fd_gradient - Activates TaoDefaultComputeGradient()
39a7e14dcfSSatish Balay -  -tao_fd_delta <delta> - change in x used to calculate finite differences
40a7e14dcfSSatish Balay 
41a7e14dcfSSatish Balay    Level: advanced
42a7e14dcfSSatish Balay 
43a7e14dcfSSatish Balay    Note:
44a7e14dcfSSatish Balay    This routine is slow and expensive, and is not currently optimized
45a7e14dcfSSatish Balay    to take advantage of sparsity in the problem.  Although
46a7e14dcfSSatish Balay    TaoAppDefaultComputeGradient is not recommended for general use
47a7e14dcfSSatish Balay    in large-scale applications, It can be useful in checking the
4858417fe7SBarry Smith    correctness of a user-provided gradient.  Use the tao method TAOTEST
49a7e14dcfSSatish Balay    to get an indication of whether your gradient is correct.
50a7e14dcfSSatish Balay 
51a7e14dcfSSatish Balay 
52a7e14dcfSSatish Balay    Note:
53a7e14dcfSSatish Balay    This finite difference gradient evaluation can be set using the routine TaoSetGradientRoutine() or by using the command line option -tao_fd_gradient
54a7e14dcfSSatish Balay 
55a7e14dcfSSatish Balay .seealso: TaoSetGradientRoutine()
56a7e14dcfSSatish Balay 
57a7e14dcfSSatish Balay @*/
58441846f8SBarry Smith PetscErrorCode TaoDefaultComputeGradient(Tao tao,Vec X,Vec G,void *dummy)
59a7e14dcfSSatish Balay {
60a7e14dcfSSatish Balay   PetscReal      *g;
61a7e14dcfSSatish Balay   PetscReal      f, f2;
62a7e14dcfSSatish Balay   PetscErrorCode ierr;
63a7e14dcfSSatish Balay   PetscInt       low,high,N,i;
64a7e14dcfSSatish Balay   PetscBool      flg;
65a7e14dcfSSatish Balay   PetscReal      h=PETSC_SQRT_MACHINE_EPSILON;
665ca45b2bSBarry Smith 
67a7e14dcfSSatish Balay   PetscFunctionBegin;
68a7e14dcfSSatish Balay   ierr = TaoComputeObjective(tao, X,&f);CHKERRQ(ierr);
696c23d075SBarry Smith   ierr = PetscOptionsGetReal(NULL,"-tao_fd_delta",&h,&flg);CHKERRQ(ierr);
70a7e14dcfSSatish Balay   ierr = VecGetSize(X,&N);CHKERRQ(ierr);
71a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(X,&low,&high);CHKERRQ(ierr);
72a7e14dcfSSatish Balay   ierr = VecGetArray(G,&g);CHKERRQ(ierr);
73a7e14dcfSSatish Balay   for (i=0;i<N;i++) {
745e081366SBarry Smith     if (i>=low && i<high) {
755e081366SBarry Smith       PetscScalar *xx;
765e081366SBarry Smith       ierr = VecGetArray(X,&xx);CHKERRQ(ierr);
775e081366SBarry Smith       xx[i-low] += h;
785e081366SBarry Smith       ierr = VecRestoreArray(X,&xx);CHKERRQ(ierr);
795e081366SBarry Smith     }
80a7e14dcfSSatish Balay 
81a7e14dcfSSatish Balay     ierr = TaoComputeObjective(tao,X,&f2);CHKERRQ(ierr);
82a7e14dcfSSatish Balay 
835e081366SBarry Smith     if (i>=low && i<high) {
845e081366SBarry Smith       PetscScalar *xx;
855e081366SBarry Smith       ierr = VecGetArray(X,&xx);CHKERRQ(ierr);
865e081366SBarry Smith       xx[i-low] -= h;
875e081366SBarry Smith       ierr = VecRestoreArray(X,&xx);CHKERRQ(ierr);
885e081366SBarry Smith     }
89a7e14dcfSSatish Balay     if (i>=low && i<high) {
90a7e14dcfSSatish Balay       g[i-low]=(f2-f)/h;
91a7e14dcfSSatish Balay     }
92a7e14dcfSSatish Balay   }
93a7e14dcfSSatish Balay   ierr = VecRestoreArray(G,&g);CHKERRQ(ierr);
94a7e14dcfSSatish Balay   PetscFunctionReturn(0);
95a7e14dcfSSatish Balay }
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay #undef __FUNCT__
98a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultComputeHessian"
99a7e14dcfSSatish Balay /*@C
100a7e14dcfSSatish Balay    TaoDefaultComputeHessian - Computes the Hessian using finite differences.
101a7e14dcfSSatish Balay 
102441846f8SBarry Smith    Collective on Tao
103a7e14dcfSSatish Balay 
104a7e14dcfSSatish Balay    Input Parameters:
105441846f8SBarry Smith +  tao - the Tao context
106a7e14dcfSSatish Balay .  V - compute Hessian at this point
107a7e14dcfSSatish Balay -  dummy - not used
108a7e14dcfSSatish Balay 
109a7e14dcfSSatish Balay    Output Parameters:
110a7e14dcfSSatish Balay +  H - Hessian matrix (not altered in this routine)
111a7e14dcfSSatish Balay .  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
112a7e14dcfSSatish Balay -  flag - flag indicating whether the matrix sparsity structure has changed
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay    Options Database Key:
115a7e14dcfSSatish Balay +  -tao_fd - Activates TaoDefaultComputeHessian()
116a7e14dcfSSatish Balay -  -tao_view_hessian - view the hessian after each evaluation using PETSC_VIEWER_STDOUT_WORLD
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay    Level: advanced
119a7e14dcfSSatish Balay 
120a7e14dcfSSatish Balay    Notes:
121a7e14dcfSSatish Balay    This routine is slow and expensive, and is not currently optimized
122a7e14dcfSSatish Balay    to take advantage of sparsity in the problem.  Although
123a7e14dcfSSatish Balay    TaoDefaultComputeHessian() is not recommended for general use
124a7e14dcfSSatish Balay    in large-scale applications, It can be useful in checking the
125a7e14dcfSSatish Balay    correctness of a user-provided Hessian.
126a7e14dcfSSatish Balay 
127a7e14dcfSSatish Balay 
128a7e14dcfSSatish Balay 
129a7e14dcfSSatish Balay .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessianColor(), SNESComputeJacobianDefault(), TaoSetGradientRoutine(), TaoDefaultComputeGradient()
130a7e14dcfSSatish Balay 
131a7e14dcfSSatish Balay @*/
132*ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessian(Tao tao,Vec V,Mat H,Mat B,void *dummy)
133d1e9a80fSBarry Smith {
134a7e14dcfSSatish Balay   PetscErrorCode       ierr;
135a7e14dcfSSatish Balay   MPI_Comm             comm;
136a7e14dcfSSatish Balay   Vec                  G;
137a7e14dcfSSatish Balay   SNES                 snes;
138a7e14dcfSSatish Balay 
139a7e14dcfSSatish Balay   PetscFunctionBegin;
140a7e14dcfSSatish Balay   PetscValidHeaderSpecific(V,VEC_CLASSID,2);
141a7e14dcfSSatish Balay   ierr = VecDuplicate(V,&G);CHKERRQ(ierr);
142a7e14dcfSSatish Balay 
143a7e14dcfSSatish Balay   ierr = PetscInfo(tao,"TAO Using finite differences w/o coloring to compute Hessian matrix\n");CHKERRQ(ierr);
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay   ierr = TaoComputeGradient(tao,V,G);CHKERRQ(ierr);
146a7e14dcfSSatish Balay 
14794ab13aaSBarry Smith   ierr = PetscObjectGetComm((PetscObject)H,&comm);CHKERRQ(ierr);
148a7e14dcfSSatish Balay   ierr = SNESCreate(comm,&snes);CHKERRQ(ierr);
149a7e14dcfSSatish Balay 
150a7e14dcfSSatish Balay   ierr = SNESSetFunction(snes,G,Fsnes,tao);CHKERRQ(ierr);
151d1e9a80fSBarry Smith   ierr = SNESComputeJacobianDefault(snes,V,H,B,tao);CHKERRQ(ierr);
152a7e14dcfSSatish Balay   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
153a7e14dcfSSatish Balay   ierr = VecDestroy(&G);CHKERRQ(ierr);
154a7e14dcfSSatish Balay   PetscFunctionReturn(0);
155a7e14dcfSSatish Balay }
156a7e14dcfSSatish Balay 
157a7e14dcfSSatish Balay #undef __FUNCT__
158a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultComputeHessianColor"
159a7e14dcfSSatish Balay /*@C
160a7e14dcfSSatish Balay    TaoDefaultComputeHessianColor - Computes the Hessian using colored finite differences.
161a7e14dcfSSatish Balay 
162441846f8SBarry Smith    Collective on Tao
163a7e14dcfSSatish Balay 
164a7e14dcfSSatish Balay    Input Parameters:
165441846f8SBarry Smith +  tao - the Tao context
166a7e14dcfSSatish Balay .  V - compute Hessian at this point
167a7e14dcfSSatish Balay -  ctx - the PetscColoring object (must be of type MatFDColoring)
168a7e14dcfSSatish Balay 
169a7e14dcfSSatish Balay    Output Parameters:
170a7e14dcfSSatish Balay +  H - Hessian matrix (not altered in this routine)
171a7e14dcfSSatish Balay .  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
172a7e14dcfSSatish Balay -  flag - flag indicating whether the matrix sparsity structure has changed
173a7e14dcfSSatish Balay 
174a7e14dcfSSatish Balay    Level: advanced
175a7e14dcfSSatish Balay 
176a7e14dcfSSatish Balay 
177a7e14dcfSSatish Balay .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessian(),SNESComputeJacobianDefaultColor(), TaoSetGradientRoutine()
178a7e14dcfSSatish Balay 
179a7e14dcfSSatish Balay @*/
180*ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessianColor(Tao tao, Vec V, Mat H,Mat B,void *ctx)
1815ca45b2bSBarry Smith {
182a7e14dcfSSatish Balay   PetscErrorCode      ierr;
183a7e14dcfSSatish Balay   MatFDColoring       coloring = (MatFDColoring)ctx;
184a7e14dcfSSatish Balay 
185a7e14dcfSSatish Balay   PetscFunctionBegin;
186a7e14dcfSSatish Balay   PetscValidHeaderSpecific(ctx,MAT_FDCOLORING_CLASSID,6);
187a7e14dcfSSatish Balay   ierr=PetscInfo(tao,"TAO computing matrix using finite differences Hessian and coloring\n");CHKERRQ(ierr);
188d1e9a80fSBarry Smith   ierr = MatFDColoringApply(B,coloring,V,ctx);CHKERRQ(ierr);
18994ab13aaSBarry Smith   if (H != B) {
19094ab13aaSBarry Smith     ierr = MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19194ab13aaSBarry Smith     ierr = MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
192a7e14dcfSSatish Balay   }
193a7e14dcfSSatish Balay   PetscFunctionReturn(0);
194a7e14dcfSSatish Balay }
195a7e14dcfSSatish Balay 
196a7e14dcfSSatish Balay 
197