v0.15.0
Loading...
Searching...
No Matches
TaoCtx.hpp
Go to the documentation of this file.
1/** \file TaoCtx.hpp
2 * \brief Context for PETSc TAO, i.e. optimization solver
3 * \author Anonymous author(s) committing under MIT license
4 */
5
6#ifndef __TAOTX_HPP__
7#define __TAOTX_HPP__
8
9namespace MoFEM {
10
11/** \brief Interface for TAO solvers
12 * \ingroup mofem_petsc_solvers
13 */
14struct TaoCtx {
15
16 MoFEM::Interface &mField; ///< database Interface
17 moab::Interface &moab; ///< moab Interface
18 std::string problemName; ///< problem name
19 bool vErify; ///< If true verify vector
20
24
26 loopsHessian; ///< Sequence of finite element methods for assembling the
27 ///< Hessian (second derivative) matrix in TAO optimization
29 loopsObjective; ///< Sequence of finite element methods for evaluating the
30 ///< objective function in TAO optimization
32 loopsGradient; ///< Sequence of finite element methods for assembling the
33 ///< gradient (first derivative) vector in TAO optimization
34
36 preHessian; ///< Sequence of methods run before Hessian is assembled
38 postHessian; ///< Sequence of methods run after Hessian is assembled
39
41 preObjective; ///< Sequence of methods run before objective is evaluated
43 postObjective; ///< Sequence of methods run after objective is evaluated
44
46 preGradient; ///< Sequence of methods run before gradient is assembled
48 postGradient; ///< Sequence of methods run after gradient is assembled
49
50 TaoCtx(Interface &m_field, const std::string &problem_name)
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 }
64
65 virtual ~TaoCtx() = default;
66
67 // TAO Hessian (second derivative) matrix assembly loop and pre/post sequences
71
72 // TAO Objective function evaluation loop and pre/post sequences
76
77 // TAO Gradient (first derivative) vector assembly loop and pre/post sequences
81
82 /**
83 * \brief Copy sequences from another TaoCtx
84 * @param tao_ctx TaoCtx from which sequences are copied
85 * @return error code
86 */
87 MoFEMErrorCode copyLoops(const TaoCtx &tao_ctx);
88
89 /**
90 * @brief Clear loops
91 *
92 * @return MoFEMErrorCode
93 */
95
96 friend PetscErrorCode TaoSetObjective(Tao tao, Vec x, PetscReal *f,
97 void *ctx);
98 friend PetscErrorCode TaoSetGradient(Tao tao, Vec x, Vec g, void *ctx);
99 friend PetscErrorCode TaoSetObjectiveAndGradient(Tao tao, Vec x, PetscReal *f,
100 Vec g, void *ctx);
101 friend PetscErrorCode TaoSetHessian(Tao tao, Vec x, Mat H, Mat Hpre,
102 void *ctx);
103
104private:
105 boost::movelib::unique_ptr<bool> vecAssembleSwitch;
106 boost::movelib::unique_ptr<bool> matAssembleSwitch;
107 PetscLogEvent MOFEM_EVENT_SnesRhs; ///< Log events to assemble residual
108 PetscLogEvent MOFEM_EVENT_SnesMat; ///< Log events to assemble tangent matrix
109 PetscLogEvent MOFEM_EVENT_TaoHessian; ///< Log event for TAO Hessian
110 PetscLogEvent MOFEM_EVENT_TaoGradient; ///< Log event for TAO Gradient
111 PetscLogEvent MOFEM_EVENT_TaoObjective; ///< Log event for TAO Objective
112};
113
114/**
115 * @brief Sets the objective function value for a TAO optimization context.
116 *
117 * @param tao The TAO optimization solver context.
118 * @param x The input vector at which to evaluate the objective.
119 * @param f Pointer to store the computed objective value.
120 * @param ctx User-defined context for the objective function.
121 *
122 * @return PetscErrorCode Error code (0 on success).
123 */
124PetscErrorCode TaoSetObjective(Tao tao, Vec x, PetscReal *f, void *ctx);
125
126/**
127 * @brief Sets the gradient vector for a TAO optimization context.
128 *
129 * @param tao The TAO optimization solver context.
130 * @param x The input vector at which to evaluate the gradient.
131 * @param g Pointer to the gradient vector.
132 * @param ctx User-defined context for the gradient function.
133 *
134 * @return PetscErrorCode Error code (0 on success).
135 */
136PetscErrorCode TaoSetGradient(Tao tao, Vec x, Vec g, void *ctx);
137
138/**
139 * @brief Sets the objective function value and gradient for a TAO optimization
140 * solver.
141 *
142 * @param tao The TAO solver context.
143 * @param x Input vector at which to evaluate the objective and gradient.
144 * @param f Pointer to store the computed objective function value.
145 * @param g Vector to store the computed gradient.
146 * @param ctx User-defined context for objective and gradient evaluation.
147 *
148 * @return PetscErrorCode Error code (0 on success).
149 *
150 * @note This function is typically called within a user-defined routine to
151 * provide both the objective value and its gradient for optimization
152 * algorithms.
153 *
154 * @see
155 * https://petsc.org/release/docs/manualpages/Tao/TaoSetObjectiveAndGradient.html
156 */
157PetscErrorCode TaoSetObjectiveAndGradient(Tao tao, Vec x, PetscReal *f, Vec g,
158 void *ctx);
159
160/**
161 * @brief Sets the Hessian matrix for a TAO optimization context.
162 *
163 * @param tao The TAO optimization solver context.
164 * @param x The input vector at which to evaluate the Hessian.
165 * @param H Pointer to the Hessian matrix.
166 * @param Hpre Preconditioner matrix (can be NULL).
167 * @param ctx User-defined context for the Hessian function.
168 *
169 * @return PetscErrorCode Error code (0 on success).
170 */
171PetscErrorCode TaoSetHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx);
172
173} // namespace MoFEM
174
175#endif // __TAOTX_HPP__
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
PetscErrorCode TaoSetObjective(Tao tao, Vec x, PetscReal *f, void *ctx)
Sets the objective function value for a TAO optimization context.
Definition TaoCtx.cpp:37
std::deque< BasicMethodPtr > BasicMethodsSequence
Definition AuxPETSc.hpp:57
std::deque< PairNameFEMethodPtr > FEMethodsSequence
Definition AuxPETSc.hpp:56
PetscErrorCode TaoSetObjectiveAndGradient(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx)
Sets the objective function value and gradient for a TAO optimization solver.
Definition TaoCtx.cpp:178
PetscErrorCode TaoSetGradient(Tao tao, Vec x, Vec f, void *ctx)
Sets the gradient vector for a TAO optimization context.
Definition TaoCtx.cpp:92
PetscErrorCode TaoSetHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
Sets the Hessian matrix for a TAO optimization context.
Definition TaoCtx.cpp:186
double H
Hardening.
Definition plastic.cpp:128
constexpr double g
Deprecated interface functions.
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.
Interface for TAO solvers.
Definition TaoCtx.hpp:14
TaoCtx(Interface &m_field, const std::string &problem_name)
Definition TaoCtx.hpp:50
MoFEMErrorCode copyLoops(const TaoCtx &tao_ctx)
Copy sequences from another TaoCtx.
Definition TaoCtx.cpp:9
boost::movelib::unique_ptr< bool > vecAssembleSwitch
Definition TaoCtx.hpp:105
MoFEM::FEMethodsSequence FEMethodsSequence
Definition TaoCtx.hpp:22
MoFEM::Interface & mField
database Interface
Definition TaoCtx.hpp:16
moab::Interface & moab
moab Interface
Definition TaoCtx.hpp:17
virtual ~TaoCtx()=default
FEMethodsSequence & getGradientLoops()
Definition TaoCtx.hpp:78
BasicMethodsSequence postObjective
Sequence of methods run after objective is evaluated.
Definition TaoCtx.hpp:43
BasicMethodsSequence & getPreProcGradient()
Definition TaoCtx.hpp:79
PetscLogEvent MOFEM_EVENT_TaoGradient
Log event for TAO Gradient.
Definition TaoCtx.hpp:110
PetscLogEvent MOFEM_EVENT_SnesMat
Log events to assemble tangent matrix.
Definition TaoCtx.hpp:108
std::string problemName
problem name
Definition TaoCtx.hpp:18
BasicMethodsSequence & getPostProcGradient()
Definition TaoCtx.hpp:80
BasicMethodsSequence & getPostProcObjective()
Definition TaoCtx.hpp:75
BasicMethodsSequence postGradient
Sequence of methods run after gradient is assembled.
Definition TaoCtx.hpp:48
PetscLogEvent MOFEM_EVENT_TaoHessian
Log event for TAO Hessian.
Definition TaoCtx.hpp:109
boost::movelib::unique_ptr< bool > matAssembleSwitch
Definition TaoCtx.hpp:106
BasicMethodsSequence preObjective
Sequence of methods run before objective is evaluated.
Definition TaoCtx.hpp:41
friend PetscErrorCode TaoSetHessian(Tao tao, Vec x, Mat H, Mat Hpre, void *ctx)
Sets the Hessian matrix for a TAO optimization context.
Definition TaoCtx.cpp:186
FEMethodsSequence loopsObjective
Definition TaoCtx.hpp:29
BasicMethodsSequence preGradient
Sequence of methods run before gradient is assembled.
Definition TaoCtx.hpp:46
BasicMethodsSequence & getPostProcHessian()
Definition TaoCtx.hpp:70
PetscLogEvent MOFEM_EVENT_SnesRhs
Log events to assemble residual.
Definition TaoCtx.hpp:107
PetscLogEvent MOFEM_EVENT_TaoObjective
Log event for TAO Objective.
Definition TaoCtx.hpp:111
FEMethodsSequence & getObjectiveLoops()
Definition TaoCtx.hpp:73
friend PetscErrorCode TaoSetObjectiveAndGradient(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx)
Sets the objective function value and gradient for a TAO optimization solver.
Definition TaoCtx.cpp:178
FEMethodsSequence & getHessianLoops()
Definition TaoCtx.hpp:68
MoFEMErrorCode clearLoops()
Clear loops.
Definition TaoCtx.cpp:23
bool vErify
If true verify vector.
Definition TaoCtx.hpp:19
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
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
friend PetscErrorCode TaoSetGradient(Tao tao, Vec x, Vec g, void *ctx)
Sets the gradient vector for a TAO optimization context.
Definition TaoCtx.cpp:92
FEMethodsSequence loopsGradient
Definition TaoCtx.hpp:32
BasicMethodsSequence & getPreProcObjective()
Definition TaoCtx.hpp:74
BasicMethodsSequence & getPreProcHessian()
Definition TaoCtx.hpp:69
MoFEM::BasicMethodsSequence BasicMethodsSequence
Definition TaoCtx.hpp:23
FEMethodsSequence loopsHessian
Definition TaoCtx.hpp:26