xref: /petsc/src/tao/interface/fdiff.c (revision c5929fdf3082647d199855a5c1d0286204349b03)
1ba92ff59SBarry Smith #include <petsctao.h>         /*I  "petsctao.h"  I*/
2af0996ceSBarry 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 {
6046bdf8c8SLisandro Dalcin   PetscScalar    *x,*g;
61a7e14dcfSSatish Balay   PetscReal      f, f2;
62a7e14dcfSSatish Balay   PetscErrorCode ierr;
63a7e14dcfSSatish Balay   PetscInt       low,high,N,i;
64a7e14dcfSSatish Balay   PetscBool      flg;
65d66999acSBarry Smith   PetscReal      h=.5*PETSC_SQRT_MACHINE_EPSILON;
665ca45b2bSBarry Smith 
67a7e14dcfSSatish Balay   PetscFunctionBegin;
68*c5929fdfSBarry Smith   ierr = PetscOptionsGetReal(((PetscObject)tao)->options,((PetscObject)tao)->prefix,"-tao_fd_delta",&h,&flg);CHKERRQ(ierr);
69a7e14dcfSSatish Balay   ierr = VecGetSize(X,&N);CHKERRQ(ierr);
70a7e14dcfSSatish Balay   ierr = VecGetOwnershipRange(X,&low,&high);CHKERRQ(ierr);
71a7e14dcfSSatish Balay   ierr = VecGetArray(G,&g);CHKERRQ(ierr);
72a7e14dcfSSatish Balay   for (i=0;i<N;i++) {
735e081366SBarry Smith     if (i>=low && i<high) {
7446bdf8c8SLisandro Dalcin       ierr = VecGetArray(X,&x);CHKERRQ(ierr);
75d66999acSBarry Smith       x[i-low] -= h;
76d66999acSBarry Smith       ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
77d66999acSBarry Smith     }
78d66999acSBarry Smith 
79d66999acSBarry Smith     ierr = TaoComputeObjective(tao, X,&f);CHKERRQ(ierr);
80d66999acSBarry Smith 
81d66999acSBarry Smith     if (i>=low && i<high) {
82d66999acSBarry Smith       ierr = VecGetArray(X,&x);CHKERRQ(ierr);
83d66999acSBarry Smith       x[i-low] += 2*h;
8446bdf8c8SLisandro Dalcin       ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
855e081366SBarry Smith     }
86a7e14dcfSSatish Balay 
87a7e14dcfSSatish Balay     ierr = TaoComputeObjective(tao,X,&f2);CHKERRQ(ierr);
88a7e14dcfSSatish Balay 
895e081366SBarry Smith     if (i>=low && i<high) {
9046bdf8c8SLisandro Dalcin       ierr = VecGetArray(X,&x);CHKERRQ(ierr);
9146bdf8c8SLisandro Dalcin       x[i-low] -= h;
9246bdf8c8SLisandro Dalcin       ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
935e081366SBarry Smith     }
94a7e14dcfSSatish Balay     if (i>=low && i<high) {
95d66999acSBarry Smith       g[i-low]=(f2-f)/(2.0*h);
96a7e14dcfSSatish Balay     }
97a7e14dcfSSatish Balay   }
98a7e14dcfSSatish Balay   ierr = VecRestoreArray(G,&g);CHKERRQ(ierr);
99a7e14dcfSSatish Balay   PetscFunctionReturn(0);
100a7e14dcfSSatish Balay }
101a7e14dcfSSatish Balay 
102a7e14dcfSSatish Balay #undef __FUNCT__
103a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultComputeHessian"
104a7e14dcfSSatish Balay /*@C
105a7e14dcfSSatish Balay    TaoDefaultComputeHessian - Computes the Hessian using finite differences.
106a7e14dcfSSatish Balay 
107441846f8SBarry Smith    Collective on Tao
108a7e14dcfSSatish Balay 
109a7e14dcfSSatish Balay    Input Parameters:
110441846f8SBarry Smith +  tao - the Tao context
111a7e14dcfSSatish Balay .  V - compute Hessian at this point
112a7e14dcfSSatish Balay -  dummy - not used
113a7e14dcfSSatish Balay 
114a7e14dcfSSatish Balay    Output Parameters:
115a7e14dcfSSatish Balay +  H - Hessian matrix (not altered in this routine)
11656744491SBarry Smith -  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay    Options Database Key:
119a7e14dcfSSatish Balay +  -tao_fd - Activates TaoDefaultComputeHessian()
120a7e14dcfSSatish Balay -  -tao_view_hessian - view the hessian after each evaluation using PETSC_VIEWER_STDOUT_WORLD
121a7e14dcfSSatish Balay 
122a7e14dcfSSatish Balay    Level: advanced
123a7e14dcfSSatish Balay 
124a7e14dcfSSatish Balay    Notes:
125a7e14dcfSSatish Balay    This routine is slow and expensive, and is not currently optimized
126a7e14dcfSSatish Balay    to take advantage of sparsity in the problem.  Although
127a7e14dcfSSatish Balay    TaoDefaultComputeHessian() is not recommended for general use
128a7e14dcfSSatish Balay    in large-scale applications, It can be useful in checking the
129a7e14dcfSSatish Balay    correctness of a user-provided Hessian.
130a7e14dcfSSatish Balay 
131a7e14dcfSSatish Balay .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessianColor(), SNESComputeJacobianDefault(), TaoSetGradientRoutine(), TaoDefaultComputeGradient()
132a7e14dcfSSatish Balay 
133a7e14dcfSSatish Balay @*/
134ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessian(Tao tao,Vec V,Mat H,Mat B,void *dummy)
135d1e9a80fSBarry Smith {
136a7e14dcfSSatish Balay   PetscErrorCode       ierr;
137a7e14dcfSSatish Balay   MPI_Comm             comm;
138a7e14dcfSSatish Balay   Vec                  G;
139a7e14dcfSSatish Balay   SNES                 snes;
140a7e14dcfSSatish Balay 
141a7e14dcfSSatish Balay   PetscFunctionBegin;
142a7e14dcfSSatish Balay   PetscValidHeaderSpecific(V,VEC_CLASSID,2);
143a7e14dcfSSatish Balay   ierr = VecDuplicate(V,&G);CHKERRQ(ierr);
144a7e14dcfSSatish Balay 
145a7e14dcfSSatish Balay   ierr = PetscInfo(tao,"TAO Using finite differences w/o coloring to compute Hessian matrix\n");CHKERRQ(ierr);
146a7e14dcfSSatish Balay 
147a7e14dcfSSatish Balay   ierr = TaoComputeGradient(tao,V,G);CHKERRQ(ierr);
148a7e14dcfSSatish Balay 
14994ab13aaSBarry Smith   ierr = PetscObjectGetComm((PetscObject)H,&comm);CHKERRQ(ierr);
150a7e14dcfSSatish Balay   ierr = SNESCreate(comm,&snes);CHKERRQ(ierr);
151a7e14dcfSSatish Balay 
152a7e14dcfSSatish Balay   ierr = SNESSetFunction(snes,G,Fsnes,tao);CHKERRQ(ierr);
153d1e9a80fSBarry Smith   ierr = SNESComputeJacobianDefault(snes,V,H,B,tao);CHKERRQ(ierr);
154a7e14dcfSSatish Balay   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
155a7e14dcfSSatish Balay   ierr = VecDestroy(&G);CHKERRQ(ierr);
156a7e14dcfSSatish Balay   PetscFunctionReturn(0);
157a7e14dcfSSatish Balay }
158a7e14dcfSSatish Balay 
159a7e14dcfSSatish Balay #undef __FUNCT__
160a7e14dcfSSatish Balay #define __FUNCT__ "TaoDefaultComputeHessianColor"
161a7e14dcfSSatish Balay /*@C
162a7e14dcfSSatish Balay    TaoDefaultComputeHessianColor - Computes the Hessian using colored finite differences.
163a7e14dcfSSatish Balay 
164441846f8SBarry Smith    Collective on Tao
165a7e14dcfSSatish Balay 
166a7e14dcfSSatish Balay    Input Parameters:
167441846f8SBarry Smith +  tao - the Tao context
168a7e14dcfSSatish Balay .  V - compute Hessian at this point
169a7e14dcfSSatish Balay -  ctx - the PetscColoring object (must be of type MatFDColoring)
170a7e14dcfSSatish Balay 
171a7e14dcfSSatish Balay    Output Parameters:
172a7e14dcfSSatish Balay +  H - Hessian matrix (not altered in this routine)
17356744491SBarry Smith -  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
174a7e14dcfSSatish Balay 
175a7e14dcfSSatish Balay    Level: advanced
176a7e14dcfSSatish Balay 
177a7e14dcfSSatish Balay 
178a7e14dcfSSatish Balay .seealso: TaoSetHessianRoutine(), TaoDefaultComputeHessian(),SNESComputeJacobianDefaultColor(), TaoSetGradientRoutine()
179a7e14dcfSSatish Balay 
180a7e14dcfSSatish Balay @*/
181ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessianColor(Tao tao, Vec V, Mat H,Mat B,void *ctx)
1825ca45b2bSBarry Smith {
183a7e14dcfSSatish Balay   PetscErrorCode      ierr;
184a7e14dcfSSatish Balay   MatFDColoring       coloring = (MatFDColoring)ctx;
185a7e14dcfSSatish Balay 
186a7e14dcfSSatish Balay   PetscFunctionBegin;
187a7e14dcfSSatish Balay   PetscValidHeaderSpecific(ctx,MAT_FDCOLORING_CLASSID,6);
188a7e14dcfSSatish Balay   ierr=PetscInfo(tao,"TAO computing matrix using finite differences Hessian and coloring\n");CHKERRQ(ierr);
189d1e9a80fSBarry Smith   ierr = MatFDColoringApply(B,coloring,V,ctx);CHKERRQ(ierr);
19094ab13aaSBarry Smith   if (H != B) {
19194ab13aaSBarry Smith     ierr = MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19294ab13aaSBarry Smith     ierr = MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
193a7e14dcfSSatish Balay   }
194a7e14dcfSSatish Balay   PetscFunctionReturn(0);
195a7e14dcfSSatish Balay }
196a7e14dcfSSatish Balay 
197a7e14dcfSSatish Balay 
198