v0.15.0
Loading...
Searching...
No Matches
MoFEM::TaoCtx Struct Reference

Interface for TAO solvers. More...

#include "src/petsc/TaoCtx.hpp"

Collaboration diagram for MoFEM::TaoCtx:
[legend]

Public Types

using PairNameFEMethodPtr = MoFEM::PairNameFEMethodPtr
 
using FEMethodsSequence = MoFEM::FEMethodsSequence
 
using BasicMethodsSequence = MoFEM::BasicMethodsSequence
 

Public Member Functions

 TaoCtx (Interface &m_field, const std::string &problem_name)
 
virtual ~TaoCtx ()=default
 
FEMethodsSequencegetHessianLoops ()
 
BasicMethodsSequencegetPreProcHessian ()
 
BasicMethodsSequencegetPostProcHessian ()
 
FEMethodsSequencegetObjectiveLoops ()
 
BasicMethodsSequencegetPreProcObjective ()
 
BasicMethodsSequencegetPostProcObjective ()
 
FEMethodsSequencegetGradientLoops ()
 
BasicMethodsSequencegetPreProcGradient ()
 
BasicMethodsSequencegetPostProcGradient ()
 
MoFEMErrorCode copyLoops (const TaoCtx &tao_ctx)
 Copy sequences from another TaoCtx.
 
MoFEMErrorCode clearLoops ()
 Clear loops.
 

Public Attributes

MoFEM::InterfacemField
 database Interface
 
moab::Interface & moab
 moab Interface
 
std::string problemName
 problem name
 
bool vErify
 If true verify vector.
 
FEMethodsSequence loopsHessian
 
FEMethodsSequence loopsObjective
 
FEMethodsSequence loopsGradient
 
BasicMethodsSequence preHessian
 Sequence of methods run before Hessian is assembled.
 
BasicMethodsSequence postHessian
 Sequence of methods run after Hessian is assembled.
 
BasicMethodsSequence preObjective
 Sequence of methods run before objective is evaluated.
 
BasicMethodsSequence postObjective
 Sequence of methods run after objective is evaluated.
 
BasicMethodsSequence preGradient
 Sequence of methods run before gradient is assembled.
 
BasicMethodsSequence postGradient
 Sequence of methods run after gradient is assembled.
 

Private Attributes

boost::movelib::unique_ptr< boolvecAssembleSwitch
 
boost::movelib::unique_ptr< boolmatAssembleSwitch
 
PetscLogEvent MOFEM_EVENT_SnesRhs
 Log events to assemble residual.
 
PetscLogEvent MOFEM_EVENT_SnesMat
 Log events to assemble tangent matrix.
 
PetscLogEvent MOFEM_EVENT_TaoHessian
 Log event for TAO Hessian.
 
PetscLogEvent MOFEM_EVENT_TaoGradient
 Log event for TAO Gradient.
 
PetscLogEvent MOFEM_EVENT_TaoObjective
 Log event for TAO Objective.
 

Friends

PetscErrorCode TaoSetObjective (Tao tao, Vec x, PetscReal *f, void *ctx)
 Sets the objective function value for a TAO optimization context.
 
PetscErrorCode TaoSetGradient (Tao tao, Vec x, Vec g, void *ctx)
 Sets the gradient vector for a TAO optimization context.
 
PetscErrorCode TaoSetObjectiveAndGradient (Tao tao, Vec x, PetscReal *f, Vec g, void *ctx)
 Sets the objective function value and gradient for a TAO optimization solver.
 
PetscErrorCode TaoSetHessian (Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
 Sets the Hessian matrix for a TAO optimization context.
 

Detailed Description

Interface for TAO solvers.

Definition at line 14 of file TaoCtx.hpp.

Member Typedef Documentation

◆ BasicMethodsSequence

◆ FEMethodsSequence

◆ PairNameFEMethodPtr

Constructor & Destructor Documentation

◆ TaoCtx()

MoFEM::TaoCtx::TaoCtx ( Interface & m_field,
const std::string & problem_name )
inline

Definition at line 50 of file TaoCtx.hpp.

51 : mField(m_field), moab(m_field.get_moab()), problemName(problem_name),
52 vErify(false) {
53 PetscLogEventRegister("LoopTAOHessian", 0, &MOFEM_EVENT_TaoHessian);
54 PetscLogEventRegister("LoopTAOObjective", 0, &MOFEM_EVENT_TaoObjective);
55 PetscLogEventRegister("LoopTAOGradient", 0, &MOFEM_EVENT_TaoGradient);
56 if (!LogManager::checkIfChannelExist("TAO_WORLD")) {
57 auto core_log = logging::core::get();
58 core_log->add_sink(
60 LogManager::setLog("TAO_WORLD");
61 MOFEM_LOG_TAG("TAO_WORLD", "TAO");
62 }
63 }
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
static bool checkIfChannelExist(const std::string channel)
Check if channel exist.
MoFEM::Interface & mField
database Interface
Definition TaoCtx.hpp:16
moab::Interface & moab
moab Interface
Definition TaoCtx.hpp:17
PetscLogEvent MOFEM_EVENT_TaoGradient
Log event for TAO Gradient.
Definition TaoCtx.hpp:110
std::string problemName
problem name
Definition TaoCtx.hpp:18
PetscLogEvent MOFEM_EVENT_TaoHessian
Log event for TAO Hessian.
Definition TaoCtx.hpp:109
PetscLogEvent MOFEM_EVENT_TaoObjective
Log event for TAO Objective.
Definition TaoCtx.hpp:111
bool vErify
If true verify vector.
Definition TaoCtx.hpp:19

◆ ~TaoCtx()

virtual MoFEM::TaoCtx::~TaoCtx ( )
virtualdefault

Member Function Documentation

◆ clearLoops()

MoFEMErrorCode MoFEM::TaoCtx::clearLoops ( )

Clear loops.

Returns
MoFEMErrorCode

Definition at line 23 of file TaoCtx.cpp.

23 {
25 loopsHessian.clear();
26 loopsObjective.clear();
27 loopsGradient.clear();
28 preHessian.clear();
29 postHessian.clear();
30 preObjective.clear();
31 postObjective.clear();
32 preGradient.clear();
33 postGradient.clear();
35}
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
BasicMethodsSequence postObjective
Sequence of methods run after objective is evaluated.
Definition TaoCtx.hpp:43
BasicMethodsSequence postGradient
Sequence of methods run after gradient is assembled.
Definition TaoCtx.hpp:48
BasicMethodsSequence preObjective
Sequence of methods run before objective is evaluated.
Definition TaoCtx.hpp:41
FEMethodsSequence loopsObjective
Definition TaoCtx.hpp:29
BasicMethodsSequence preGradient
Sequence of methods run before gradient is assembled.
Definition TaoCtx.hpp:46
BasicMethodsSequence postHessian
Sequence of methods run after Hessian is assembled.
Definition TaoCtx.hpp:38
BasicMethodsSequence preHessian
Sequence of methods run before Hessian is assembled.
Definition TaoCtx.hpp:36
FEMethodsSequence loopsGradient
Definition TaoCtx.hpp:32
FEMethodsSequence loopsHessian
Definition TaoCtx.hpp:26

◆ copyLoops()

MoFEMErrorCode MoFEM::TaoCtx::copyLoops ( const TaoCtx & tao_ctx)

Copy sequences from another TaoCtx.

Parameters
tao_ctxTaoCtx from which sequences are copied
Returns
error code

Definition at line 9 of file TaoCtx.cpp.

9 {
11 loopsHessian = tao_ctx.loopsHessian;
12 loopsObjective = tao_ctx.loopsObjective;
13 loopsGradient = tao_ctx.loopsGradient;
14 preHessian = tao_ctx.preHessian;
15 postHessian = tao_ctx.postHessian;
16 preObjective = tao_ctx.preObjective;
17 postObjective = tao_ctx.postObjective;
18 preGradient = tao_ctx.preGradient;
19 postGradient = tao_ctx.postGradient;
21}

◆ getGradientLoops()

FEMethodsSequence & MoFEM::TaoCtx::getGradientLoops ( )
inline

Definition at line 78 of file TaoCtx.hpp.

78{ return loopsGradient; }

◆ getHessianLoops()

FEMethodsSequence & MoFEM::TaoCtx::getHessianLoops ( )
inline

Definition at line 68 of file TaoCtx.hpp.

68{ return loopsHessian; }

◆ getObjectiveLoops()

FEMethodsSequence & MoFEM::TaoCtx::getObjectiveLoops ( )
inline

Definition at line 73 of file TaoCtx.hpp.

73{ return loopsObjective; }

◆ getPostProcGradient()

BasicMethodsSequence & MoFEM::TaoCtx::getPostProcGradient ( )
inline

Definition at line 80 of file TaoCtx.hpp.

80{ return postGradient; }

◆ getPostProcHessian()

BasicMethodsSequence & MoFEM::TaoCtx::getPostProcHessian ( )
inline

Definition at line 70 of file TaoCtx.hpp.

70{ return postHessian; }

◆ getPostProcObjective()

BasicMethodsSequence & MoFEM::TaoCtx::getPostProcObjective ( )
inline

Definition at line 75 of file TaoCtx.hpp.

75{ return postObjective; }

◆ getPreProcGradient()

BasicMethodsSequence & MoFEM::TaoCtx::getPreProcGradient ( )
inline

Definition at line 79 of file TaoCtx.hpp.

79{ return preGradient; }

◆ getPreProcHessian()

BasicMethodsSequence & MoFEM::TaoCtx::getPreProcHessian ( )
inline

Definition at line 69 of file TaoCtx.hpp.

69{ return preHessian; }

◆ getPreProcObjective()

BasicMethodsSequence & MoFEM::TaoCtx::getPreProcObjective ( )
inline

Definition at line 74 of file TaoCtx.hpp.

74{ return preObjective; }

Friends And Related Symbol Documentation

◆ TaoSetGradient

PetscErrorCode TaoSetGradient ( Tao tao,
Vec x,
Vec g,
void * ctx )
friend

Sets the gradient vector for a TAO optimization context.

Parameters
taoThe TAO optimization solver context.
xThe input vector at which to evaluate the gradient.
gPointer to the gradient vector.
ctxUser-defined context for the gradient function.
Returns
PetscErrorCode Error code (0 on success).

Definition at line 92 of file TaoCtx.cpp.

92 {
93 TaoCtx *tao_ctx = (TaoCtx *)ctx;
95 PetscLogEventBegin(tao_ctx->MOFEM_EVENT_TaoGradient, 0, 0, 0, 0);
96 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
97 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
98 CHKERR tao_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
99 tao_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
100
101 auto zero_ghost_vec = [](Vec g) {
103 Vec l;
104 CHKERR VecGhostGetLocalForm(g, &l);
105 double *a;
106 CHKERR VecGetArray(l, &a);
107 int s;
108 CHKERR VecGetLocalSize(l, &s);
109 for (int i = 0; i != s; ++i)
110 a[i] = 0;
111 CHKERR VecRestoreArray(l, &a);
112 CHKERR VecGhostRestoreLocalForm(g, &l);
114 };
115 CHKERR zero_ghost_vec(f);
116
117 tao_ctx->vecAssembleSwitch = boost::movelib::make_unique<bool>(true);
118 auto cache_ptr = boost::make_shared<CacheTuple>();
119 CHKERR tao_ctx->mField.cache_problem_entities(tao_ctx->problemName,
120 cache_ptr);
121
122 auto set = [&](auto &fe) {
123 fe.tao = tao;
124 fe.tao_x = x;
125 fe.tao_f = f;
127 fe.ksp_ctx = KspMethod::CTX_SETFUNCTION;
128 fe.tao_ctx = TaoMethod::CTX_TAO_GRADIENT;
130 fe.cacheWeakPtr = cache_ptr;
131 };
132
133 auto unset = [&](auto &fe) {
134 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
135 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
136 fe.tao_ctx = TaoMethod::CTX_TAO_GRADIENT;
137 fe.data_ctx = PetscData::CtxSetNone;
138 };
139
140 for (auto &bit : tao_ctx->preGradient) {
141 bit->vecAssembleSwitch = boost::move(tao_ctx->vecAssembleSwitch);
142 set(*bit);
143 CHKERR tao_ctx->mField.problem_basic_method_preProcess(tao_ctx->problemName,
144 *bit);
145 unset(*bit);
146 tao_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
147 }
148
149 for (auto &lit : tao_ctx->loopsGradient) {
150 lit.second->vecAssembleSwitch = boost::move(tao_ctx->vecAssembleSwitch);
151 set(*lit.second);
152 CHKERR tao_ctx->mField.loop_finite_elements(tao_ctx->problemName, lit.first,
153 *lit.second, nullptr,
154 MF_EXIST, cache_ptr);
155 unset(*lit.second);
156 tao_ctx->vecAssembleSwitch = boost::move(lit.second->vecAssembleSwitch);
157 }
158
159 for (auto &bit : tao_ctx->postGradient) {
160 bit->vecAssembleSwitch = boost::move(tao_ctx->vecAssembleSwitch);
161 set(*bit);
162 CHKERR tao_ctx->mField.problem_basic_method_postProcess(
163 tao_ctx->problemName, *bit);
164 unset(*bit);
165 tao_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
166 }
167
168 if (*(tao_ctx->vecAssembleSwitch)) {
169 CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
170 CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
171 CHKERR VecAssemblyBegin(f);
172 CHKERR VecAssemblyEnd(f);
173 }
174 PetscLogEventEnd(tao_ctx->MOFEM_EVENT_TaoGradient, 0, 0, 0, 0);
176}
constexpr double a
@ COL
@ MF_EXIST
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
auto bit
set bit
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'l', 3 > l
constexpr double g
static constexpr Switches CtxSetX
static constexpr Switches CtxSetNone
static constexpr Switches CtxSetF
TaoCtx(Interface &m_field, const std::string &problem_name)
Definition TaoCtx.hpp:50

◆ TaoSetHessian

PetscErrorCode TaoSetHessian ( Tao tao,
Vec x,
Mat H,
Mat Hpre,
void * ctx )
friend

Sets the Hessian matrix for a TAO optimization context.

Parameters
taoThe TAO optimization solver context.
xThe input vector at which to evaluate the Hessian.
HPointer to the Hessian matrix.
HprePreconditioner matrix (can be NULL).
ctxUser-defined context for the Hessian function.
Returns
PetscErrorCode Error code (0 on success).

Definition at line 186 of file TaoCtx.cpp.

186 {
187 TaoCtx *tao_ctx = (TaoCtx *)ctx;
189 PetscLogEventBegin(tao_ctx->MOFEM_EVENT_TaoHessian, 0, 0, 0, 0);
190
191 tao_ctx->matAssembleSwitch = boost::movelib::make_unique<bool>(true);
192 auto cache_ptr = boost::make_shared<CacheTuple>();
193 CHKERR tao_ctx->mField.cache_problem_entities(tao_ctx->problemName,
194 cache_ptr);
195
196 auto set = [&](auto &fe) {
197 fe.tao = tao;
198 fe.tao_A = H;
199 fe.tao_B = Hpre;
200 fe.ksp_ctx = KspMethod::CTX_OPERATORS;
202 fe.tao_ctx = TaoMethod::CTX_TAO_HESSIAN;
204 fe.cacheWeakPtr = cache_ptr;
205 };
206
207 auto unset = [&](auto &fe) {
208 fe.snes_ctx = SnesMethod::CTX_SNESNONE;
209 fe.ksp_ctx = KspMethod::CTX_KSPNONE;
210 fe.tao_ctx = TaoMethod::CTX_TAO_NONE;
211 fe.data_ctx = PetscData::CtxSetNone;
212 };
213
214
215 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
216 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
217 CHKERR tao_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
218 tao_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
219
220 for (auto &bit : tao_ctx->preHessian) {
221 bit->matAssembleSwitch = boost::move(tao_ctx->matAssembleSwitch);
222 set(*bit);
223 CHKERR tao_ctx->mField.problem_basic_method_preProcess(
224 tao_ctx->problemName, *bit);
225 unset(*bit);
226 tao_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
227 }
228
229 for (auto &lit : tao_ctx->loopsHessian) {
230 lit.second->matAssembleSwitch = boost::move(tao_ctx->matAssembleSwitch);
231 set(*lit.second);
232 CHKERR tao_ctx->mField.loop_finite_elements(
233 tao_ctx->problemName, lit.first, *(lit.second), nullptr, MF_EXIST,
234 cache_ptr);
235 unset(*lit.second);
236 tao_ctx->matAssembleSwitch = boost::move(lit.second->matAssembleSwitch);
237 }
238
239 for (auto &bit : tao_ctx->postHessian) {
240 bit->matAssembleSwitch = boost::move(tao_ctx->matAssembleSwitch);
241 set(*bit);
242 CHKERR tao_ctx->mField.problem_basic_method_postProcess(
243 tao_ctx->problemName, *bit);
244 unset(*bit);
245 tao_ctx->matAssembleSwitch = boost::move(bit->matAssembleSwitch);
246 }
247
248 if (*tao_ctx->matAssembleSwitch) {
249 CHKERR MatAssemblyBegin(Hpre, MAT_FINAL_ASSEMBLY);
250 CHKERR MatAssemblyEnd(Hpre, MAT_FINAL_ASSEMBLY);
251 }
252 PetscLogEventEnd(tao_ctx->MOFEM_EVENT_TaoHessian, 0, 0, 0, 0);
254}
double H
Hardening.
Definition plastic.cpp:128
static constexpr Switches CtxSetA
static constexpr Switches CtxSetB

◆ TaoSetObjective

PetscErrorCode TaoSetObjective ( Tao tao,
Vec x,
PetscReal * f,
void * ctx )
friend

Sets the objective function value for a TAO optimization context.

Parameters
taoThe TAO optimization solver context.
xThe input vector at which to evaluate the objective.
fPointer to store the computed objective value.
ctxUser-defined context for the objective function.
Returns
PetscErrorCode Error code (0 on success).

Definition at line 37 of file TaoCtx.cpp.

37 {
38 TaoCtx *tao_ctx = (TaoCtx *)ctx;
40 PetscLogEventBegin(tao_ctx->MOFEM_EVENT_TaoObjective, 0, 0, 0, 0);
41 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
42 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
43 CHKERR tao_ctx->mField.getInterface<VecManager>()->setLocalGhostVector(
44 tao_ctx->problemName, COL, x, INSERT_VALUES, SCATTER_REVERSE);
45
46 auto cache_ptr = boost::make_shared<CacheTuple>();
47 CHKERR tao_ctx->mField.cache_problem_entities(tao_ctx->problemName,
48 cache_ptr);
49
50 auto set = [&](auto &fe) {
51 fe.tao = tao;
52 fe.tao_x = x;
54 fe.data_ctx = PetscData::CtxSetX;
55 fe.cacheWeakPtr = cache_ptr;
56 };
57
58 auto unset = [&](auto &fe) {
59 fe.tao_ctx = TaoMethod::CTX_TAO_NONE;
60 fe.data_ctx = PetscData::CtxSetNone;
61 };
62
63 for (auto &bit : tao_ctx->preObjective) {
64 bit->vecAssembleSwitch = boost::move(tao_ctx->vecAssembleSwitch);
65 set(*bit);
66 CHKERR tao_ctx->mField.problem_basic_method_preProcess(tao_ctx->problemName,
67 *bit);
68 unset(*bit);
69 }
70
71 for (auto &lit : tao_ctx->loopsObjective) {
72 set(*lit.second);
73 CHKERR tao_ctx->mField.loop_finite_elements(tao_ctx->problemName, lit.first,
74 *lit.second, nullptr, MF_EXIST,
75 cache_ptr);
76 unset(*lit.second);
77 }
78
79 for (auto &bit : tao_ctx->postObjective) {
80 bit->vecAssembleSwitch = boost::move(tao_ctx->vecAssembleSwitch);
81 set(*bit);
82 CHKERR tao_ctx->mField.problem_basic_method_postProcess(
83 tao_ctx->problemName, *bit);
84 unset(*bit);
85 tao_ctx->vecAssembleSwitch = boost::move(bit->vecAssembleSwitch);
86 }
87
88 PetscLogEventEnd(tao_ctx->MOFEM_EVENT_TaoObjective, 0, 0, 0, 0);
90}

◆ TaoSetObjectiveAndGradient

PetscErrorCode TaoSetObjectiveAndGradient ( Tao tao,
Vec x,
PetscReal * f,
Vec g,
void * ctx )
friend

Sets the objective function value and gradient for a TAO optimization solver.

Parameters
taoThe TAO solver context.
xInput vector at which to evaluate the objective and gradient.
fPointer to store the computed objective function value.
gVector to store the computed gradient.
ctxUser-defined context for objective and gradient evaluation.
Returns
PetscErrorCode Error code (0 on success).
Note
This function is typically called within a user-defined routine to provide both the objective value and its gradient for optimization algorithms.
See also
https://petsc.org/release/docs/manualpages/Tao/TaoSetObjectiveAndGradient.html

Definition at line 178 of file TaoCtx.cpp.

179 {
181 CHKERR TaoSetObjective(tao, x, f, ctx);
182 CHKERR TaoSetGradient(tao, x, g, ctx);
184}
friend PetscErrorCode TaoSetObjective(Tao tao, Vec x, PetscReal *f, void *ctx)
Sets the objective function value for a TAO optimization context.
Definition TaoCtx.cpp:37
friend PetscErrorCode TaoSetGradient(Tao tao, Vec x, Vec g, void *ctx)
Sets the gradient vector for a TAO optimization context.
Definition TaoCtx.cpp:92

Member Data Documentation

◆ loopsGradient

FEMethodsSequence MoFEM::TaoCtx::loopsGradient

Sequence of finite element methods for assembling the gradient (first derivative) vector in TAO optimization

Definition at line 32 of file TaoCtx.hpp.

◆ loopsHessian

FEMethodsSequence MoFEM::TaoCtx::loopsHessian

Sequence of finite element methods for assembling the Hessian (second derivative) matrix in TAO optimization

Definition at line 26 of file TaoCtx.hpp.

◆ loopsObjective

FEMethodsSequence MoFEM::TaoCtx::loopsObjective

Sequence of finite element methods for evaluating the objective function in TAO optimization

Definition at line 29 of file TaoCtx.hpp.

◆ matAssembleSwitch

boost::movelib::unique_ptr<bool> MoFEM::TaoCtx::matAssembleSwitch
private

Definition at line 106 of file TaoCtx.hpp.

◆ mField

MoFEM::Interface& MoFEM::TaoCtx::mField

database Interface

Definition at line 16 of file TaoCtx.hpp.

◆ moab

moab::Interface& MoFEM::TaoCtx::moab

moab Interface

Definition at line 17 of file TaoCtx.hpp.

◆ MOFEM_EVENT_SnesMat

PetscLogEvent MoFEM::TaoCtx::MOFEM_EVENT_SnesMat
private

Log events to assemble tangent matrix.

Definition at line 108 of file TaoCtx.hpp.

◆ MOFEM_EVENT_SnesRhs

PetscLogEvent MoFEM::TaoCtx::MOFEM_EVENT_SnesRhs
private

Log events to assemble residual.

Definition at line 107 of file TaoCtx.hpp.

◆ MOFEM_EVENT_TaoGradient

PetscLogEvent MoFEM::TaoCtx::MOFEM_EVENT_TaoGradient
private

Log event for TAO Gradient.

Definition at line 110 of file TaoCtx.hpp.

◆ MOFEM_EVENT_TaoHessian

PetscLogEvent MoFEM::TaoCtx::MOFEM_EVENT_TaoHessian
private

Log event for TAO Hessian.

Definition at line 109 of file TaoCtx.hpp.

◆ MOFEM_EVENT_TaoObjective

PetscLogEvent MoFEM::TaoCtx::MOFEM_EVENT_TaoObjective
private

Log event for TAO Objective.

Definition at line 111 of file TaoCtx.hpp.

◆ postGradient

BasicMethodsSequence MoFEM::TaoCtx::postGradient

Sequence of methods run after gradient is assembled.

Definition at line 48 of file TaoCtx.hpp.

◆ postHessian

BasicMethodsSequence MoFEM::TaoCtx::postHessian

Sequence of methods run after Hessian is assembled.

Definition at line 38 of file TaoCtx.hpp.

◆ postObjective

BasicMethodsSequence MoFEM::TaoCtx::postObjective

Sequence of methods run after objective is evaluated.

Definition at line 43 of file TaoCtx.hpp.

◆ preGradient

BasicMethodsSequence MoFEM::TaoCtx::preGradient

Sequence of methods run before gradient is assembled.

Definition at line 46 of file TaoCtx.hpp.

◆ preHessian

BasicMethodsSequence MoFEM::TaoCtx::preHessian

Sequence of methods run before Hessian is assembled.

Definition at line 36 of file TaoCtx.hpp.

◆ preObjective

BasicMethodsSequence MoFEM::TaoCtx::preObjective

Sequence of methods run before objective is evaluated.

Definition at line 41 of file TaoCtx.hpp.

◆ problemName

std::string MoFEM::TaoCtx::problemName

problem name

Definition at line 18 of file TaoCtx.hpp.

◆ vecAssembleSwitch

boost::movelib::unique_ptr<bool> MoFEM::TaoCtx::vecAssembleSwitch
private

Definition at line 105 of file TaoCtx.hpp.

◆ vErify

bool MoFEM::TaoCtx::vErify

If true verify vector.

Definition at line 19 of file TaoCtx.hpp.


The documentation for this struct was generated from the following files: