xref: /petsc/src/tao/interface/fdiff.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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 */
Fsnes(SNES snes,Vec X,Vec G,PetscCtx ctx)9*2a8381b2SBarry Smith static PetscErrorCode Fsnes(SNES snes, Vec X, Vec G, PetscCtx ctx)
10d71ae5a4SJacob Faibussowitsch {
11441846f8SBarry Smith   Tao tao = (Tao)ctx;
125ca45b2bSBarry Smith 
13a7e14dcfSSatish Balay   PetscFunctionBegin;
14f4c1ad5cSStefano Zampini   PetscValidHeaderSpecific(tao, TAO_CLASSID, 4);
159566063dSJacob Faibussowitsch   PetscCall(TaoComputeGradient(tao, X, G));
163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17a7e14dcfSSatish Balay }
18a7e14dcfSSatish Balay 
19a7e14dcfSSatish Balay /*@C
20a7e14dcfSSatish Balay   TaoDefaultComputeGradient - computes the gradient using finite differences.
21a7e14dcfSSatish Balay 
22c3339decSBarry Smith   Collective
23a7e14dcfSSatish Balay 
24a7e14dcfSSatish Balay   Input Parameters:
25441846f8SBarry Smith + tao   - the Tao context
26e056e8ceSJacob Faibussowitsch . Xin   - 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:
3965ba42b6SBarry Smith   This routine is slow and expensive, and is not optimized
40a7e14dcfSSatish Balay   to take advantage of sparsity in the problem.  Although
4165ba42b6SBarry Smith   not recommended for general use
4265ba42b6SBarry Smith   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.
4565ba42b6SBarry Smith   This finite difference gradient evaluation can be set using the routine `TaoSetGradient()` or by using the command line option -tao_fd_gradient
46a7e14dcfSSatish Balay 
4765ba42b6SBarry Smith .seealso: `Tao`, `TaoSetGradient()`
48a7e14dcfSSatish Balay @*/
TaoDefaultComputeGradient(Tao tao,Vec Xin,Vec G,void * dummy)49d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoDefaultComputeGradient(Tao tao, Vec Xin, Vec G, void *dummy)
50d71ae5a4SJacob Faibussowitsch {
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;
599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(((PetscObject)tao)->options, ((PetscObject)tao)->prefix, "-tao_fd_delta", &h, &flg));
609566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(Xin, &X));
619566063dSJacob Faibussowitsch   PetscCall(VecCopy(Xin, X));
629566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &N));
639566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(X, &low, &high));
649566063dSJacob Faibussowitsch   PetscCall(VecSetOption(X, VEC_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(G, &g));
66a7e14dcfSSatish Balay   for (i = 0; i < N; i++) {
679566063dSJacob Faibussowitsch     PetscCall(VecSetValue(X, i, -h, ADD_VALUES));
689566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(X));
699566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(X));
709566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(tao, X, &f));
719566063dSJacob Faibussowitsch     PetscCall(VecSetValue(X, i, 2.0 * h, ADD_VALUES));
729566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(X));
739566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(X));
749566063dSJacob Faibussowitsch     PetscCall(TaoComputeObjective(tao, X, &f2));
759566063dSJacob Faibussowitsch     PetscCall(VecSetValue(X, i, -h, ADD_VALUES));
769566063dSJacob Faibussowitsch     PetscCall(VecAssemblyBegin(X));
779566063dSJacob Faibussowitsch     PetscCall(VecAssemblyEnd(X));
78ad540459SPierre Jolivet     if (i >= low && i < high) g[i - low] = (f2 - f) / (2.0 * h);
79a7e14dcfSSatish Balay   }
809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(G, &g));
819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
83a7e14dcfSSatish Balay }
84a7e14dcfSSatish Balay 
85a7e14dcfSSatish Balay /*@C
86a7e14dcfSSatish Balay   TaoDefaultComputeHessian - Computes the Hessian using finite differences.
87a7e14dcfSSatish Balay 
88c3339decSBarry Smith   Collective
89a7e14dcfSSatish Balay 
90a7e14dcfSSatish Balay   Input Parameters:
91441846f8SBarry Smith + tao   - the Tao context
92a7e14dcfSSatish Balay . V     - compute Hessian at this point
93a7e14dcfSSatish Balay - dummy - not used
94a7e14dcfSSatish Balay 
95a7e14dcfSSatish Balay   Output Parameters:
96a7e14dcfSSatish Balay + H - Hessian matrix (not altered in this routine)
9756744491SBarry Smith - B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
98a7e14dcfSSatish Balay 
99a7e14dcfSSatish Balay   Options Database Key:
100f4c1ad5cSStefano Zampini . -tao_fd_hessian - activates TaoDefaultComputeHessian()
101a7e14dcfSSatish Balay 
102a7e14dcfSSatish Balay   Level: advanced
103a7e14dcfSSatish Balay 
104a7e14dcfSSatish Balay   Notes:
10565ba42b6SBarry Smith   This routine is slow and expensive, and is not optimized
106a7e14dcfSSatish Balay   to take advantage of sparsity in the problem.  Although
10765ba42b6SBarry Smith   it is not recommended for general use
108a7e14dcfSSatish Balay   in large-scale applications, It can be useful in checking the
109a7e14dcfSSatish Balay   correctness of a user-provided Hessian.
110a7e14dcfSSatish Balay 
11165ba42b6SBarry Smith .seealso: `Tao`, `TaoSetHessian()`, `TaoDefaultComputeHessianColor()`, `SNESComputeJacobianDefault()`, `TaoSetGradient()`, `TaoDefaultComputeGradient()`
112a7e14dcfSSatish Balay @*/
TaoDefaultComputeHessian(Tao tao,Vec V,Mat H,Mat B,void * dummy)113d71ae5a4SJacob Faibussowitsch PetscErrorCode TaoDefaultComputeHessian(Tao tao, Vec V, Mat H, Mat B, void *dummy)
114d71ae5a4SJacob Faibussowitsch {
115a7e14dcfSSatish Balay   SNES snes;
116f4c1ad5cSStefano Zampini   DM   dm;
117a7e14dcfSSatish Balay 
118a7e14dcfSSatish Balay   PetscFunctionBegin;
1199566063dSJacob Faibussowitsch   PetscCall(PetscInfo(tao, "TAO Using finite differences w/o coloring to compute Hessian matrix\n"));
1209566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PetscObjectComm((PetscObject)H), &snes));
1219566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes, NULL, Fsnes, tao));
1229566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &dm));
1239566063dSJacob Faibussowitsch   PetscCall(DMShellSetGlobalVector(dm, V));
1249566063dSJacob Faibussowitsch   PetscCall(SNESSetUp(snes));
125f4c1ad5cSStefano Zampini   if (H) {
126f4c1ad5cSStefano Zampini     PetscInt n, N;
127f4c1ad5cSStefano Zampini 
1289566063dSJacob Faibussowitsch     PetscCall(VecGetSize(V, &N));
1299566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(V, &n));
1309566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(H, n, n, N, N));
1319566063dSJacob Faibussowitsch     PetscCall(MatSetUp(H));
132f4c1ad5cSStefano Zampini   }
133f4c1ad5cSStefano Zampini   if (B && B != H) {
134f4c1ad5cSStefano Zampini     PetscInt n, N;
135f4c1ad5cSStefano Zampini 
1369566063dSJacob Faibussowitsch     PetscCall(VecGetSize(V, &N));
1379566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(V, &n));
1389566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(B, n, n, N, N));
1399566063dSJacob Faibussowitsch     PetscCall(MatSetUp(B));
140f4c1ad5cSStefano Zampini   }
1419566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobianDefault(snes, V, H, B, NULL));
1429566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
1433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
144a7e14dcfSSatish Balay }
145a7e14dcfSSatish Balay 
146a7e14dcfSSatish Balay /*@C
147a7e14dcfSSatish Balay   TaoDefaultComputeHessianColor - Computes the Hessian using colored finite differences.
148a7e14dcfSSatish Balay 
149c3339decSBarry Smith   Collective
150a7e14dcfSSatish Balay 
151a7e14dcfSSatish Balay   Input Parameters:
152441846f8SBarry Smith + tao - the Tao context
153a7e14dcfSSatish Balay . V   - compute Hessian at this point
15465ba42b6SBarry Smith - ctx - the color object of type `MatFDColoring`
155a7e14dcfSSatish Balay 
156a7e14dcfSSatish Balay   Output Parameters:
157a7e14dcfSSatish Balay + H - Hessian matrix (not altered in this routine)
15856744491SBarry Smith - B - newly computed Hessian matrix to use with preconditioner (generally the same as H)
159a7e14dcfSSatish Balay 
160a7e14dcfSSatish Balay   Level: advanced
161a7e14dcfSSatish Balay 
16265ba42b6SBarry Smith .seealso: `Tao`, `MatColoring`, `TaoSetHessian()`, `TaoDefaultComputeHessian()`, `SNESComputeJacobianDefaultColor()`, `TaoSetGradient()`
163a7e14dcfSSatish Balay @*/
TaoDefaultComputeHessianColor(Tao tao,Vec V,Mat H,Mat B,PetscCtx ctx)164*2a8381b2SBarry Smith PetscErrorCode TaoDefaultComputeHessianColor(Tao tao, Vec V, Mat H, Mat B, PetscCtx ctx)
165d71ae5a4SJacob Faibussowitsch {
166a7e14dcfSSatish Balay   MatFDColoring coloring = (MatFDColoring)ctx;
167a7e14dcfSSatish Balay 
168a7e14dcfSSatish Balay   PetscFunctionBegin;
169064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(coloring, MAT_FDCOLORING_CLASSID, 5);
1709566063dSJacob Faibussowitsch   PetscCall(PetscInfo(tao, "TAO computing matrix using finite differences Hessian and coloring\n"));
1719566063dSJacob Faibussowitsch   PetscCall(MatFDColoringApply(B, coloring, V, ctx));
17294ab13aaSBarry Smith   if (H != B) {
1739566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
1749566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
175a7e14dcfSSatish Balay   }
1763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
177a7e14dcfSSatish Balay }
178a7e14dcfSSatish Balay 
TaoDefaultComputeHessianMFFD(Tao tao,Vec X,Mat H,Mat B,PetscCtx ctx)179*2a8381b2SBarry Smith PetscErrorCode TaoDefaultComputeHessianMFFD(Tao tao, Vec X, Mat H, Mat B, PetscCtx ctx)
180d71ae5a4SJacob Faibussowitsch {
181f4c1ad5cSStefano Zampini   PetscInt  n, N;
182c5e9d94cSAlp Dener   PetscBool assembled;
183a7e14dcfSSatish Balay 
184f4c1ad5cSStefano Zampini   PetscFunctionBegin;
1853c859ba3SBarry Smith   PetscCheck(!B || B == H, PetscObjectComm((PetscObject)tao), PETSC_ERR_SUP, "Preconditioning Hessian matrix");
1869566063dSJacob Faibussowitsch   PetscCall(MatAssembled(H, &assembled));
187c5e9d94cSAlp Dener   if (!assembled) {
1889566063dSJacob Faibussowitsch     PetscCall(VecGetSize(X, &N));
1899566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(X, &n));
1909566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(H, n, n, N, N));
1919566063dSJacob Faibussowitsch     PetscCall(MatSetType(H, MATMFFD));
1929566063dSJacob Faibussowitsch     PetscCall(MatSetUp(H));
1939566063dSJacob Faibussowitsch     PetscCall(MatMFFDSetFunction(H, (PetscErrorCode (*)(void *, Vec, Vec))TaoComputeGradient, tao));
194c5e9d94cSAlp Dener   }
1959566063dSJacob Faibussowitsch   PetscCall(MatMFFDSetBase(H, X, NULL));
1969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
1979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
1983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
199f4c1ad5cSStefano Zampini }
200