xref: /petsc/src/tao/interface/fdiff.c (revision 9566063d113dddea24716c546802770db7481bc0)
1ba92ff59SBarry Smith #include <petsctao.h>         /*I  "petsctao.h"  I*/
2af0996ceSBarry Smith #include <petsc/private/taoimpl.h>
3aaa7dc30SBarry Smith #include <petscsnes.h>
4f4c1ad5cSStefano Zampini #include <petscdmshell.h>
5a7e14dcfSSatish Balay 
6a7e14dcfSSatish Balay /*
7a7e14dcfSSatish Balay    For finited difference computations of the Hessian, we use PETSc's SNESComputeJacobianDefault
8a7e14dcfSSatish Balay */
95ca45b2bSBarry Smith static PetscErrorCode Fsnes(SNES snes,Vec X,Vec G,void* ctx)
105ca45b2bSBarry Smith {
11441846f8SBarry Smith   Tao            tao = (Tao)ctx;
125ca45b2bSBarry Smith 
13a7e14dcfSSatish Balay   PetscFunctionBegin;
14f4c1ad5cSStefano Zampini   PetscValidHeaderSpecific(tao,TAO_CLASSID,4);
15*9566063dSJacob Faibussowitsch   PetscCall(TaoComputeGradient(tao,X,G));
16a7e14dcfSSatish Balay   PetscFunctionReturn(0);
17a7e14dcfSSatish Balay }
18a7e14dcfSSatish Balay 
19a7e14dcfSSatish Balay /*@C
20a7e14dcfSSatish Balay   TaoDefaultComputeGradient - computes the gradient using finite differences.
21a7e14dcfSSatish Balay 
22441846f8SBarry Smith   Collective on Tao
23a7e14dcfSSatish Balay 
24a7e14dcfSSatish Balay   Input Parameters:
25441846f8SBarry Smith + tao   - the Tao context
26a7e14dcfSSatish Balay . X     - compute gradient at this point
27a7e14dcfSSatish Balay - dummy - not used
28a7e14dcfSSatish Balay 
29f899ff85SJose E. Roman   Output Parameter:
30a7e14dcfSSatish Balay . G - Gradient Vector
31a7e14dcfSSatish Balay 
32a7e14dcfSSatish Balay   Options Database Key:
33f4c1ad5cSStefano Zampini + -tao_fd_gradient      - activates TaoDefaultComputeGradient()
34f4c1ad5cSStefano Zampini - -tao_fd_delta <delta> - change in X used to calculate finite differences
35a7e14dcfSSatish Balay 
36a7e14dcfSSatish Balay   Level: advanced
37a7e14dcfSSatish Balay 
38f4c1ad5cSStefano Zampini   Notes:
39a7e14dcfSSatish Balay   This routine is slow and expensive, and is not currently optimized
40a7e14dcfSSatish Balay   to take advantage of sparsity in the problem.  Although
41f4c1ad5cSStefano Zampini   TaoDefaultComputeGradient is not recommended for general use
42a7e14dcfSSatish Balay   in large-scale applications, It can be useful in checking the
4358417fe7SBarry Smith   correctness of a user-provided gradient.  Use the tao method TAOTEST
44a7e14dcfSSatish Balay   to get an indication of whether your gradient is correct.
45a82e8c82SStefano Zampini   This finite difference gradient evaluation can be set using the routine TaoSetGradient() or by using the command line option -tao_fd_gradient
46a7e14dcfSSatish Balay 
47a82e8c82SStefano Zampini .seealso: TaoSetGradient()
48a7e14dcfSSatish Balay @*/
49f4c1ad5cSStefano Zampini PetscErrorCode TaoDefaultComputeGradient(Tao tao,Vec Xin,Vec G,void *dummy)
50a7e14dcfSSatish Balay {
51f4c1ad5cSStefano Zampini   Vec            X;
52f4c1ad5cSStefano Zampini   PetscScalar    *g;
53a7e14dcfSSatish Balay   PetscReal      f, f2;
54a7e14dcfSSatish Balay   PetscInt       low,high,N,i;
55a7e14dcfSSatish Balay   PetscBool      flg;
56d66999acSBarry Smith   PetscReal      h=.5*PETSC_SQRT_MACHINE_EPSILON;
575ca45b2bSBarry Smith 
58a7e14dcfSSatish Balay   PetscFunctionBegin;
59*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(((PetscObject)tao)->options,((PetscObject)tao)->prefix,"-tao_fd_delta",&h,&flg));
60*9566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(Xin,&X));
61*9566063dSJacob Faibussowitsch   PetscCall(VecCopy(Xin,X));
62*9566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X,&N));
63*9566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X,&low,&high));
64*9566063dSJacob Faibussowitsch   PetscCall(VecSetOption(X,VEC_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
65*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G,&g));
66a7e14dcfSSatish Balay   for (i=0;i<N;i++) {
67*9566063dSJacob Faibussowitsch     PetscCall(VecSetValue(X,i,-h,ADD_VALUES));
68*9566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(X));
69*9566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(X));
70*9566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(tao,X,&f));
71*9566063dSJacob Faibussowitsch     PetscCall(VecSetValue(X,i,2.0*h,ADD_VALUES));
72*9566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(X));
73*9566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(X));
74*9566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(tao,X,&f2));
75*9566063dSJacob Faibussowitsch     PetscCall(VecSetValue(X,i,-h,ADD_VALUES));
76*9566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(X));
77*9566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(X));
78a7e14dcfSSatish Balay     if (i>=low && i<high) {
79d66999acSBarry Smith       g[i-low]=(f2-f)/(2.0*h);
80a7e14dcfSSatish Balay     }
81a7e14dcfSSatish Balay   }
82*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G,&g));
83*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
84a7e14dcfSSatish Balay   PetscFunctionReturn(0);
85a7e14dcfSSatish Balay }
86a7e14dcfSSatish Balay 
87a7e14dcfSSatish Balay /*@C
88a7e14dcfSSatish Balay    TaoDefaultComputeHessian - Computes the Hessian using finite differences.
89a7e14dcfSSatish Balay 
90441846f8SBarry Smith    Collective on Tao
91a7e14dcfSSatish Balay 
92a7e14dcfSSatish Balay    Input Parameters:
93441846f8SBarry Smith +  tao   - the Tao context
94a7e14dcfSSatish Balay .  V     - compute Hessian at this point
95a7e14dcfSSatish Balay -  dummy - not used
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay    Output Parameters:
98a7e14dcfSSatish Balay +  H - Hessian matrix (not altered in this routine)
9956744491SBarry Smith -  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
100a7e14dcfSSatish Balay 
101a7e14dcfSSatish Balay    Options Database Key:
102f4c1ad5cSStefano Zampini .  -tao_fd_hessian - activates TaoDefaultComputeHessian()
103a7e14dcfSSatish Balay 
104a7e14dcfSSatish Balay    Level: advanced
105a7e14dcfSSatish Balay 
106a7e14dcfSSatish Balay    Notes:
107a7e14dcfSSatish Balay    This routine is slow and expensive, and is not currently optimized
108a7e14dcfSSatish Balay    to take advantage of sparsity in the problem.  Although
109a7e14dcfSSatish Balay    TaoDefaultComputeHessian() is not recommended for general use
110a7e14dcfSSatish Balay    in large-scale applications, It can be useful in checking the
111a7e14dcfSSatish Balay    correctness of a user-provided Hessian.
112a7e14dcfSSatish Balay 
113a82e8c82SStefano Zampini .seealso: TaoSetHessian(), TaoDefaultComputeHessianColor(), SNESComputeJacobianDefault(), TaoSetGradient(), TaoDefaultComputeGradient()
114a7e14dcfSSatish Balay @*/
115ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessian(Tao tao,Vec V,Mat H,Mat B,void *dummy)
116d1e9a80fSBarry Smith {
117a7e14dcfSSatish Balay   SNES           snes;
118f4c1ad5cSStefano Zampini   DM             dm;
119a7e14dcfSSatish Balay 
120a7e14dcfSSatish Balay   PetscFunctionBegin;
121*9566063dSJacob Faibussowitsch   PetscCall(PetscInfo(tao,"TAO Using finite differences w/o coloring to compute Hessian matrix\n"));
122*9566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PetscObjectComm((PetscObject)H),&snes));
123*9566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes,NULL,Fsnes,tao));
124*9566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&dm));
125*9566063dSJacob Faibussowitsch   PetscCall(DMShellSetGlobalVector(dm,V));
126*9566063dSJacob Faibussowitsch   PetscCall(SNESSetUp(snes));
127f4c1ad5cSStefano Zampini   if (H) {
128f4c1ad5cSStefano Zampini     PetscInt n,N;
129f4c1ad5cSStefano Zampini 
130*9566063dSJacob Faibussowitsch     PetscCall(VecGetSize(V,&N));
131*9566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(V,&n));
132*9566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(H,n,n,N,N));
133*9566063dSJacob Faibussowitsch     PetscCall(MatSetUp(H));
134f4c1ad5cSStefano Zampini   }
135f4c1ad5cSStefano Zampini   if (B && B != H) {
136f4c1ad5cSStefano Zampini     PetscInt n,N;
137f4c1ad5cSStefano Zampini 
138*9566063dSJacob Faibussowitsch     PetscCall(VecGetSize(V,&N));
139*9566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(V,&n));
140*9566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B,n,n,N,N));
141*9566063dSJacob Faibussowitsch     PetscCall(MatSetUp(B));
142f4c1ad5cSStefano Zampini   }
143*9566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobianDefault(snes,V,H,B,NULL));
144*9566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
145a7e14dcfSSatish Balay   PetscFunctionReturn(0);
146a7e14dcfSSatish Balay }
147a7e14dcfSSatish Balay 
148a7e14dcfSSatish Balay /*@C
149a7e14dcfSSatish Balay    TaoDefaultComputeHessianColor - Computes the Hessian using colored finite differences.
150a7e14dcfSSatish Balay 
151441846f8SBarry Smith    Collective on Tao
152a7e14dcfSSatish Balay 
153a7e14dcfSSatish Balay    Input Parameters:
154441846f8SBarry Smith +  tao - the Tao context
155a7e14dcfSSatish Balay .  V   - compute Hessian at this point
156a7e14dcfSSatish Balay -  ctx - the PetscColoring object (must be of type MatFDColoring)
157a7e14dcfSSatish Balay 
158a7e14dcfSSatish Balay    Output Parameters:
159a7e14dcfSSatish Balay +  H - Hessian matrix (not altered in this routine)
16056744491SBarry Smith -  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
161a7e14dcfSSatish Balay 
162a7e14dcfSSatish Balay    Level: advanced
163a7e14dcfSSatish Balay 
164a82e8c82SStefano Zampini .seealso: TaoSetHessian(), TaoDefaultComputeHessian(),SNESComputeJacobianDefaultColor(), TaoSetGradient()
165a7e14dcfSSatish Balay @*/
166ffad9901SBarry Smith PetscErrorCode TaoDefaultComputeHessianColor(Tao tao,Vec V,Mat H,Mat B,void *ctx)
1675ca45b2bSBarry Smith {
168a7e14dcfSSatish Balay   MatFDColoring       coloring = (MatFDColoring)ctx;
169a7e14dcfSSatish Balay 
170a7e14dcfSSatish Balay   PetscFunctionBegin;
171064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(coloring,MAT_FDCOLORING_CLASSID,5);
172*9566063dSJacob Faibussowitsch   PetscCall(PetscInfo(tao,"TAO computing matrix using finite differences Hessian and coloring\n"));
173*9566063dSJacob Faibussowitsch   PetscCall(MatFDColoringApply(B,coloring,V,ctx));
17494ab13aaSBarry Smith   if (H != B) {
175*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
176*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
177a7e14dcfSSatish Balay   }
178a7e14dcfSSatish Balay   PetscFunctionReturn(0);
179a7e14dcfSSatish Balay }
180a7e14dcfSSatish Balay 
181f4c1ad5cSStefano Zampini PetscErrorCode TaoDefaultComputeHessianMFFD(Tao tao,Vec X,Mat H,Mat B,void *ctx)
182f4c1ad5cSStefano Zampini {
183f4c1ad5cSStefano Zampini   PetscInt       n,N;
184c5e9d94cSAlp Dener   PetscBool      assembled;
185a7e14dcfSSatish Balay 
186f4c1ad5cSStefano Zampini   PetscFunctionBegin;
1873c859ba3SBarry Smith   PetscCheck(!B || B == H,PetscObjectComm((PetscObject)tao),PETSC_ERR_SUP,"Preconditioning Hessian matrix");
188*9566063dSJacob Faibussowitsch   PetscCall(MatAssembled(H, &assembled));
189c5e9d94cSAlp Dener   if (!assembled) {
190*9566063dSJacob Faibussowitsch     PetscCall(VecGetSize(X,&N));
191*9566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(X,&n));
192*9566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(H,n,n,N,N));
193*9566063dSJacob Faibussowitsch     PetscCall(MatSetType(H,MATMFFD));
194*9566063dSJacob Faibussowitsch     PetscCall(MatSetUp(H));
195*9566063dSJacob Faibussowitsch     PetscCall(MatMFFDSetFunction(H,(PetscErrorCode (*)(void*,Vec,Vec))TaoComputeGradient,tao));
196c5e9d94cSAlp Dener   }
197*9566063dSJacob Faibussowitsch   PetscCall(MatMFFDSetBase(H,X,NULL));
198*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY));
199*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY));
200f4c1ad5cSStefano Zampini   PetscFunctionReturn(0);
201f4c1ad5cSStefano Zampini }
202