xref: /petsc/src/tao/interface/fdiff.c (revision ad540459ab38c4a232710a68077e487eb6627fe2)
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 */
99371c9d4SSatish Balay static PetscErrorCode Fsnes(SNES snes, Vec X, Vec G, void *ctx) {
10441846f8SBarry Smith   Tao tao = (Tao)ctx;
115ca45b2bSBarry Smith 
12a7e14dcfSSatish Balay   PetscFunctionBegin;
13f4c1ad5cSStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 4);
149566063dSJacob Faibussowitsch   PetscCall(TaoComputeGradient(tao, X, G));
15a7e14dcfSSatish Balay   PetscFunctionReturn(0);
16a7e14dcfSSatish Balay }
17a7e14dcfSSatish Balay 
18a7e14dcfSSatish Balay /*@C
19a7e14dcfSSatish Balay   TaoDefaultComputeGradient - computes the gradient using finite differences.
20a7e14dcfSSatish Balay 
2165ba42b6SBarry Smith   Collective on tao
22a7e14dcfSSatish Balay 
23a7e14dcfSSatish Balay   Input Parameters:
24441846f8SBarry Smith + tao   - the Tao context
25a7e14dcfSSatish Balay . X     - compute gradient at this point
26a7e14dcfSSatish Balay - dummy - not used
27a7e14dcfSSatish Balay 
28f899ff85SJose E. Roman   Output Parameter:
29a7e14dcfSSatish Balay . G - Gradient Vector
30a7e14dcfSSatish Balay 
31a7e14dcfSSatish Balay   Options Database Key:
32f4c1ad5cSStefano Zampini + -tao_fd_gradient      - activates TaoDefaultComputeGradient()
33f4c1ad5cSStefano Zampini - -tao_fd_delta <delta> - change in X used to calculate finite differences
34a7e14dcfSSatish Balay 
35a7e14dcfSSatish Balay   Level: advanced
36a7e14dcfSSatish Balay 
37f4c1ad5cSStefano Zampini   Notes:
3865ba42b6SBarry Smith   This routine is slow and expensive, and is not optimized
39a7e14dcfSSatish Balay   to take advantage of sparsity in the problem.  Although
4065ba42b6SBarry Smith   not recommended for general use
4165ba42b6SBarry Smith   in large-scale applications, it can be useful in checking the
4258417fe7SBarry Smith   correctness of a user-provided gradient.  Use the tao method TAOTEST
43a7e14dcfSSatish Balay   to get an indication of whether your gradient is correct.
4465ba42b6SBarry Smith   This finite difference gradient evaluation can be set using the routine `TaoSetGradient()` or by using the command line option -tao_fd_gradient
45a7e14dcfSSatish Balay 
4665ba42b6SBarry Smith .seealso: `Tao`, `TaoSetGradient()`
47a7e14dcfSSatish Balay @*/
489371c9d4SSatish Balay PetscErrorCode TaoDefaultComputeGradient(Tao tao, Vec Xin, Vec G, void *dummy) {
49f4c1ad5cSStefano Zampini   Vec          X;
50f4c1ad5cSStefano Zampini   PetscScalar *g;
51a7e14dcfSSatish Balay   PetscReal    f, f2;
52a7e14dcfSSatish Balay   PetscInt     low, high, N, i;
53a7e14dcfSSatish Balay   PetscBool    flg;
54d66999acSBarry Smith   PetscReal    h = .5 * PETSC_SQRT_MACHINE_EPSILON;
555ca45b2bSBarry Smith 
56a7e14dcfSSatish Balay   PetscFunctionBegin;
579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(((PetscObject)tao)->options, ((PetscObject)tao)->prefix, "-tao_fd_delta", &h, &flg));
589566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(Xin, &X));
599566063dSJacob Faibussowitsch   PetscCall(VecCopy(Xin, X));
609566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &N));
619566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &low, &high));
629566063dSJacob Faibussowitsch   PetscCall(VecSetOption(X, VEC_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
639566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G, &g));
64a7e14dcfSSatish Balay   for (i = 0; i < N; i++) {
659566063dSJacob Faibussowitsch     PetscCall(VecSetValue(X, i, -h, ADD_VALUES));
669566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(X));
679566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(X));
689566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(tao, X, &f));
699566063dSJacob Faibussowitsch     PetscCall(VecSetValue(X, i, 2.0 * h, ADD_VALUES));
709566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(X));
719566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(X));
729566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(tao, X, &f2));
739566063dSJacob Faibussowitsch     PetscCall(VecSetValue(X, i, -h, ADD_VALUES));
749566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(X));
759566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(X));
76*ad540459SPierre Jolivet     if (i >= low && i < high) g[i - low] = (f2 - f) / (2.0 * h);
77a7e14dcfSSatish Balay   }
789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G, &g));
799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
80a7e14dcfSSatish Balay   PetscFunctionReturn(0);
81a7e14dcfSSatish Balay }
82a7e14dcfSSatish Balay 
83a7e14dcfSSatish Balay /*@C
84a7e14dcfSSatish Balay    TaoDefaultComputeHessian - Computes the Hessian using finite differences.
85a7e14dcfSSatish Balay 
8665ba42b6SBarry Smith    Collective on tao
87a7e14dcfSSatish Balay 
88a7e14dcfSSatish Balay    Input Parameters:
89441846f8SBarry Smith +  tao   - the Tao context
90a7e14dcfSSatish Balay .  V     - compute Hessian at this point
91a7e14dcfSSatish Balay -  dummy - not used
92a7e14dcfSSatish Balay 
93a7e14dcfSSatish Balay    Output Parameters:
94a7e14dcfSSatish Balay +  H - Hessian matrix (not altered in this routine)
9556744491SBarry Smith -  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
96a7e14dcfSSatish Balay 
97a7e14dcfSSatish Balay    Options Database Key:
98f4c1ad5cSStefano Zampini .  -tao_fd_hessian - activates TaoDefaultComputeHessian()
99a7e14dcfSSatish Balay 
100a7e14dcfSSatish Balay    Level: advanced
101a7e14dcfSSatish Balay 
102a7e14dcfSSatish Balay    Notes:
10365ba42b6SBarry Smith    This routine is slow and expensive, and is not optimized
104a7e14dcfSSatish Balay    to take advantage of sparsity in the problem.  Although
10565ba42b6SBarry Smith    it is not recommended for general use
106a7e14dcfSSatish Balay    in large-scale applications, It can be useful in checking the
107a7e14dcfSSatish Balay    correctness of a user-provided Hessian.
108a7e14dcfSSatish Balay 
10965ba42b6SBarry Smith .seealso: `Tao`, `TaoSetHessian()`, `TaoDefaultComputeHessianColor()`, `SNESComputeJacobianDefault()`, `TaoSetGradient()`, `TaoDefaultComputeGradient()`
110a7e14dcfSSatish Balay @*/
1119371c9d4SSatish Balay PetscErrorCode TaoDefaultComputeHessian(Tao tao, Vec V, Mat H, Mat B, void *dummy) {
112a7e14dcfSSatish Balay   SNES snes;
113f4c1ad5cSStefano Zampini   DM   dm;
114a7e14dcfSSatish Balay 
115a7e14dcfSSatish Balay   PetscFunctionBegin;
1169566063dSJacob Faibussowitsch   PetscCall(PetscInfo(tao, "TAO Using finite differences w/o coloring to compute Hessian matrix\n"));
1179566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PetscObjectComm((PetscObject)H), &snes));
1189566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, Fsnes, tao));
1199566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
1209566063dSJacob Faibussowitsch   PetscCall(DMShellSetGlobalVector(dm, V));
1219566063dSJacob Faibussowitsch   PetscCall(SNESSetUp(snes));
122f4c1ad5cSStefano Zampini   if (H) {
123f4c1ad5cSStefano Zampini     PetscInt n, N;
124f4c1ad5cSStefano Zampini 
1259566063dSJacob Faibussowitsch     PetscCall(VecGetSize(V, &N));
1269566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(V, &n));
1279566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(H, n, n, N, N));
1289566063dSJacob Faibussowitsch     PetscCall(MatSetUp(H));
129f4c1ad5cSStefano Zampini   }
130f4c1ad5cSStefano Zampini   if (B && B != H) {
131f4c1ad5cSStefano Zampini     PetscInt n, N;
132f4c1ad5cSStefano Zampini 
1339566063dSJacob Faibussowitsch     PetscCall(VecGetSize(V, &N));
1349566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(V, &n));
1359566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, n, n, N, N));
1369566063dSJacob Faibussowitsch     PetscCall(MatSetUp(B));
137f4c1ad5cSStefano Zampini   }
1389566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobianDefault(snes, V, H, B, NULL));
1399566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
140a7e14dcfSSatish Balay   PetscFunctionReturn(0);
141a7e14dcfSSatish Balay }
142a7e14dcfSSatish Balay 
143a7e14dcfSSatish Balay /*@C
144a7e14dcfSSatish Balay    TaoDefaultComputeHessianColor - Computes the Hessian using colored finite differences.
145a7e14dcfSSatish Balay 
146441846f8SBarry Smith    Collective on Tao
147a7e14dcfSSatish Balay 
148a7e14dcfSSatish Balay    Input Parameters:
149441846f8SBarry Smith +  tao - the Tao context
150a7e14dcfSSatish Balay .  V   - compute Hessian at this point
15165ba42b6SBarry Smith -  ctx - the color object of type `MatFDColoring`
152a7e14dcfSSatish Balay 
153a7e14dcfSSatish Balay    Output Parameters:
154a7e14dcfSSatish Balay +  H - Hessian matrix (not altered in this routine)
15556744491SBarry Smith -  B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
156a7e14dcfSSatish Balay 
157a7e14dcfSSatish Balay    Level: advanced
158a7e14dcfSSatish Balay 
15965ba42b6SBarry Smith .seealso: `Tao`, `MatColoring`, `TaoSetHessian()`, `TaoDefaultComputeHessian()`, `SNESComputeJacobianDefaultColor()`, `TaoSetGradient()`
160a7e14dcfSSatish Balay @*/
1619371c9d4SSatish Balay PetscErrorCode TaoDefaultComputeHessianColor(Tao tao, Vec V, Mat H, Mat B, void *ctx) {
162a7e14dcfSSatish Balay   MatFDColoring coloring = (MatFDColoring)ctx;
163a7e14dcfSSatish Balay 
164a7e14dcfSSatish Balay   PetscFunctionBegin;
165064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 5);
1669566063dSJacob Faibussowitsch   PetscCall(PetscInfo(tao, "TAO computing matrix using finite differences Hessian and coloring\n"));
1679566063dSJacob Faibussowitsch   PetscCall(MatFDColoringApply(B, coloring, V, ctx));
16894ab13aaSBarry Smith   if (H != B) {
1699566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
1709566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
171a7e14dcfSSatish Balay   }
172a7e14dcfSSatish Balay   PetscFunctionReturn(0);
173a7e14dcfSSatish Balay }
174a7e14dcfSSatish Balay 
1759371c9d4SSatish Balay PetscErrorCode TaoDefaultComputeHessianMFFD(Tao tao, Vec X, Mat H, Mat B, void *ctx) {
176f4c1ad5cSStefano Zampini   PetscInt  n, N;
177c5e9d94cSAlp Dener   PetscBool assembled;
178a7e14dcfSSatish Balay 
179f4c1ad5cSStefano Zampini   PetscFunctionBegin;
1803c859ba3SBarry Smith   PetscCheck(!B || B == H, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Preconditioning Hessian matrix");
1819566063dSJacob Faibussowitsch   PetscCall(MatAssembled(H, &assembled));
182c5e9d94cSAlp Dener   if (!assembled) {
1839566063dSJacob Faibussowitsch     PetscCall(VecGetSize(X, &N));
1849566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(X, &n));
1859566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(H, n, n, N, N));
1869566063dSJacob Faibussowitsch     PetscCall(MatSetType(H, MATMFFD));
1879566063dSJacob Faibussowitsch     PetscCall(MatSetUp(H));
1889566063dSJacob Faibussowitsch     PetscCall(MatMFFDSetFunction(H, (PetscErrorCode(*)(void *, Vec, Vec))TaoComputeGradient, tao));
189c5e9d94cSAlp Dener   }
1909566063dSJacob Faibussowitsch   PetscCall(MatMFFDSetBase(H, X, NULL));
1919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
1929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
193f4c1ad5cSStefano Zampini   PetscFunctionReturn(0);
194f4c1ad5cSStefano Zampini }
195