xref: /petsc/src/tao/interface/fdiff.c (revision 46bdf8c8e7980dd1dca8bd7d868b63fbdf2185ad)
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 {
60*46bdf8c8SLisandro Dalcin   PetscScalar    *x,*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) {
75*46bdf8c8SLisandro Dalcin       ierr = VecGetArray(X,&x);CHKERRQ(ierr);
76*46bdf8c8SLisandro Dalcin       x[i-low] += h;
77*46bdf8c8SLisandro Dalcin       ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
785e081366SBarry Smith     }
79a7e14dcfSSatish Balay 
80a7e14dcfSSatish Balay     ierr = TaoComputeObjective(tao,X,&f2);CHKERRQ(ierr);
81a7e14dcfSSatish Balay 
825e081366SBarry Smith     if (i>=low && i<high) {
83*46bdf8c8SLisandro Dalcin       ierr = VecGetArray(X,&x);CHKERRQ(ierr);
84*46bdf8c8SLisandro Dalcin       x[i-low] -= h;
85*46bdf8c8SLisandro Dalcin       ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
865e081366SBarry Smith     }
87a7e14dcfSSatish Balay     if (i>=low && i<high) {
88a7e14dcfSSatish Balay       g[i-low]=(f2-f)/h;
89a7e14dcfSSatish Balay     }
90a7e14dcfSSatish Balay   }
91a7e14dcfSSatish Balay   ierr = VecRestoreArray(G,&g);CHKERRQ(ierr);
92a7e14dcfSSatish Balay   PetscFunctionReturn(0);
93a7e14dcfSSatish Balay }
94a7e14dcfSSatish Balay 
95a7e14dcfSSatish Balay #undef __FUNCT__
96a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultComputeHessian"
97a7e14dcfSSatish Balay /*@C
98a7e14dcfSSatish Balay    TaoDefaultComputeHessian - Computes the Hessian using finite differences.
99a7e14dcfSSatish Balay 
100441846f8SBarry Smith    Collective on Tao
101a7e14dcfSSatish Balay 
102a7e14dcfSSatish Balay    Input Parameters:
103441846f8SBarry Smith +  tao - the Tao context
104a7e14dcfSSatish Balay .  V - compute Hessian at this point
105a7e14dcfSSatish Balay -  dummy - not used
106a7e14dcfSSatish Balay 
107a7e14dcfSSatish Balay    Output Parameters:
108a7e14dcfSSatish Balay +  H - Hessian matrix (not altered in this routine)
109a7e14dcfSSatish Balay .  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
110a7e14dcfSSatish Balay -  flag - flag indicating whether the matrix sparsity structure has changed
111a7e14dcfSSatish Balay 
112a7e14dcfSSatish Balay    Options Database Key:
113a7e14dcfSSatish Balay +  -tao_fd - Activates TaoDefaultComputeHessian()
114a7e14dcfSSatish Balay -  -tao_view_hessian - view the hessian after each evaluation using PETSC_VIEWER_STDOUT_WORLD
115a7e14dcfSSatish Balay 
116a7e14dcfSSatish Balay    Level: advanced
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay    Notes:
119a7e14dcfSSatish Balay    This routine is slow and expensive, and is not currently optimized
120a7e14dcfSSatish Balay    to take advantage of sparsity in the problem.  Although
121a7e14dcfSSatish Balay    TaoDefaultComputeHessian() is not recommended for general use
122a7e14dcfSSatish Balay    in large-scale applications, It can be useful in checking the
123a7e14dcfSSatish Balay    correctness of a user-provided Hessian.
124a7e14dcfSSatish Balay 
125a7e14dcfSSatish Balay 
126a7e14dcfSSatish Balay 
127a7e14dcfSSatish Balay .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessianColor(), SNESComputeJacobianDefault(), TaoSetGradientRoutine(), TaoDefaultComputeGradient()
128a7e14dcfSSatish Balay 
129a7e14dcfSSatish Balay @*/
130ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessian(Tao tao,Vec V,Mat H,Mat B,void *dummy)
131d1e9a80fSBarry Smith {
132a7e14dcfSSatish Balay   PetscErrorCode       ierr;
133a7e14dcfSSatish Balay   MPI_Comm             comm;
134a7e14dcfSSatish Balay   Vec                  G;
135a7e14dcfSSatish Balay   SNES                 snes;
136a7e14dcfSSatish Balay 
137a7e14dcfSSatish Balay   PetscFunctionBegin;
138a7e14dcfSSatish Balay   PetscValidHeaderSpecific(V,VEC_CLASSID,2);
139a7e14dcfSSatish Balay   ierr = VecDuplicate(V,&G);CHKERRQ(ierr);
140a7e14dcfSSatish Balay 
141a7e14dcfSSatish Balay   ierr = PetscInfo(tao,"TAO Using finite differences w/o coloring to compute Hessian matrix\n");CHKERRQ(ierr);
142a7e14dcfSSatish Balay 
143a7e14dcfSSatish Balay   ierr = TaoComputeGradient(tao,V,G);CHKERRQ(ierr);
144a7e14dcfSSatish Balay 
14594ab13aaSBarry Smith   ierr = PetscObjectGetComm((PetscObject)H,&comm);CHKERRQ(ierr);
146a7e14dcfSSatish Balay   ierr = SNESCreate(comm,&snes);CHKERRQ(ierr);
147a7e14dcfSSatish Balay 
148a7e14dcfSSatish Balay   ierr = SNESSetFunction(snes,G,Fsnes,tao);CHKERRQ(ierr);
149d1e9a80fSBarry Smith   ierr = SNESComputeJacobianDefault(snes,V,H,B,tao);CHKERRQ(ierr);
150a7e14dcfSSatish Balay   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
151a7e14dcfSSatish Balay   ierr = VecDestroy(&G);CHKERRQ(ierr);
152a7e14dcfSSatish Balay   PetscFunctionReturn(0);
153a7e14dcfSSatish Balay }
154a7e14dcfSSatish Balay 
155a7e14dcfSSatish Balay #undef __FUNCT__
156a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultComputeHessianColor"
157a7e14dcfSSatish Balay /*@C
158a7e14dcfSSatish Balay    TaoDefaultComputeHessianColor - Computes the Hessian using colored finite differences.
159a7e14dcfSSatish Balay 
160441846f8SBarry Smith    Collective on Tao
161a7e14dcfSSatish Balay 
162a7e14dcfSSatish Balay    Input Parameters:
163441846f8SBarry Smith +  tao - the Tao context
164a7e14dcfSSatish Balay .  V - compute Hessian at this point
165a7e14dcfSSatish Balay -  ctx - the PetscColoring object (must be of type MatFDColoring)
166a7e14dcfSSatish Balay 
167a7e14dcfSSatish Balay    Output Parameters:
168a7e14dcfSSatish Balay +  H - Hessian matrix (not altered in this routine)
169a7e14dcfSSatish Balay .  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
170a7e14dcfSSatish Balay -  flag - flag indicating whether the matrix sparsity structure has changed
171a7e14dcfSSatish Balay 
172a7e14dcfSSatish Balay    Level: advanced
173a7e14dcfSSatish Balay 
174a7e14dcfSSatish Balay 
175a7e14dcfSSatish Balay .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessian(),SNESComputeJacobianDefaultColor(), TaoSetGradientRoutine()
176a7e14dcfSSatish Balay 
177a7e14dcfSSatish Balay @*/
178ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessianColor(Tao tao, Vec V, Mat H,Mat B,void *ctx)
1795ca45b2bSBarry Smith {
180a7e14dcfSSatish Balay   PetscErrorCode      ierr;
181a7e14dcfSSatish Balay   MatFDColoring       coloring = (MatFDColoring)ctx;
182a7e14dcfSSatish Balay 
183a7e14dcfSSatish Balay   PetscFunctionBegin;
184a7e14dcfSSatish Balay   PetscValidHeaderSpecific(ctx,MAT_FDCOLORING_CLASSID,6);
185a7e14dcfSSatish Balay   ierr=PetscInfo(tao,"TAO computing matrix using finite differences Hessian and coloring\n");CHKERRQ(ierr);
186d1e9a80fSBarry Smith   ierr = MatFDColoringApply(B,coloring,V,ctx);CHKERRQ(ierr);
18794ab13aaSBarry Smith   if (H != B) {
18894ab13aaSBarry Smith     ierr = MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
18994ab13aaSBarry Smith     ierr = MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
190a7e14dcfSSatish Balay   }
191a7e14dcfSSatish Balay   PetscFunctionReturn(0);
192a7e14dcfSSatish Balay }
193a7e14dcfSSatish Balay 
194a7e14dcfSSatish Balay 
195