xref: /petsc/src/tao/snes/taosnes.c (revision 517f4e3460cd8c7e73c19b8bf533f1efe47f30a7)
1 #include <petsc/private/taoimpl.h> /*I "petsctao.h" I*/
2 
3 typedef struct {
4   SNES snes;
5 } Tao_SNES;
6 
7 static PetscErrorCode TaoSolve_SNES(Tao tao)
8 {
9   Tao_SNES *taosnes = (Tao_SNES *)tao->data;
10 
11   PetscFunctionBegin;
12   PetscCall(SNESSolve(taosnes->snes, NULL, tao->solution));
13   /* TODO REASONS */
14   tao->reason = TAO_CONVERGED_USER;
15   PetscFunctionReturn(PETSC_SUCCESS);
16 }
17 
18 static PetscErrorCode TaoDestroy_SNES(Tao tao)
19 {
20   Tao_SNES *taosnes = (Tao_SNES *)tao->data;
21 
22   PetscFunctionBegin;
23   PetscCall(SNESDestroy(&taosnes->snes));
24   PetscCall(PetscFree(tao->data));
25   PetscFunctionReturn(PETSC_SUCCESS);
26 }
27 
28 static PetscErrorCode TAOSNESObj(SNES snes, Vec X, PetscReal *f, void *ctx)
29 {
30   Tao tao = (Tao)ctx;
31 
32   PetscFunctionBegin;
33   PetscCall(TaoComputeObjective(tao, X, f));
34   PetscFunctionReturn(PETSC_SUCCESS);
35 }
36 
37 static PetscErrorCode TAOSNESFunc(SNES snes, Vec X, Vec F, void *ctx)
38 {
39   Tao tao = (Tao)ctx;
40 
41   PetscFunctionBegin;
42   PetscCall(TaoComputeGradient(tao, X, F));
43   PetscFunctionReturn(PETSC_SUCCESS);
44 }
45 
46 static PetscErrorCode TAOSNESJac(SNES snes, Vec X, Mat A, Mat P, void *ctx)
47 {
48   Tao tao = (Tao)ctx;
49 
50   PetscFunctionBegin;
51   PetscCall(TaoComputeHessian(tao, X, A, P));
52   PetscFunctionReturn(PETSC_SUCCESS);
53 }
54 
55 static PetscErrorCode TaoSetUp_SNES(Tao tao)
56 {
57   Tao_SNES *taosnes = (Tao_SNES *)tao->data;
58   Mat       A, P;
59 
60   PetscFunctionBegin;
61   PetscCall(SNESSetSolution(taosnes->snes, tao->solution));
62   PetscCall(SNESSetObjective(taosnes->snes, TAOSNESObj, tao));
63   PetscCall(SNESSetFunction(taosnes->snes, NULL, TAOSNESFunc, tao));
64   PetscCall(TaoGetHessian(tao, &A, &P, NULL, NULL));
65   if (A) PetscCall(SNESSetJacobian(taosnes->snes, A, P, TAOSNESJac, tao));
66   /* TODO TYPES */
67   PetscCall(SNESSetUp(taosnes->snes));
68   PetscFunctionReturn(PETSC_SUCCESS);
69 }
70 
71 static PetscErrorCode TaoSetFromOptions_SNES(Tao tao, PetscOptionItems *PetscOptionsObject)
72 {
73   Tao_SNES *taosnes = (Tao_SNES *)tao->data;
74 
75   PetscFunctionBegin;
76   PetscCall(SNESSetFromOptions(taosnes->snes));
77   PetscFunctionReturn(PETSC_SUCCESS);
78 }
79 
80 static PetscErrorCode TaoView_SNES(Tao tao, PetscViewer viewer)
81 {
82   Tao_SNES *taosnes = (Tao_SNES *)tao->data;
83 
84   PetscFunctionBegin;
85   PetscCall(SNESView(taosnes->snes, viewer));
86   PetscFunctionReturn(PETSC_SUCCESS);
87 }
88 
89 /*MC
90   TAOSNES - nonlinear solver using SNES
91 
92    Level: advanced
93 
94 .seealso: `TaoCreate()`, `Tao`, `TaoSetType()`, `TaoType`
95 M*/
96 PETSC_EXTERN PetscErrorCode TaoCreate_SNES(Tao tao)
97 {
98   Tao_SNES *taosnes;
99 
100   PetscFunctionBegin;
101   tao->ops->destroy        = TaoDestroy_SNES;
102   tao->ops->setup          = TaoSetUp_SNES;
103   tao->ops->setfromoptions = TaoSetFromOptions_SNES;
104   tao->ops->view           = TaoView_SNES;
105   tao->ops->solve          = TaoSolve_SNES;
106 
107   PetscCall(PetscNew(&taosnes));
108   tao->data = (void *)taosnes;
109   PetscCall(SNESCreate(PetscObjectComm((PetscObject)tao), &taosnes->snes));
110   PetscCall(PetscObjectIncrementTabLevel((PetscObject)taosnes->snes, (PetscObject)tao, 1));
111   PetscFunctionReturn(PETSC_SUCCESS);
112 }
113