![]() |
v0.15.5 |
Abstract interface for Python-defined objective functions. More...
#include "tutorials/vec-7_shape_optimisation/src/ObjectiveFunctionData.hpp"
Public Member Functions | |
| virtual MoFEMErrorCode | evalInteriorObjectiveFunction (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > stress_ptr, boost::shared_ptr< MatrixDouble > strain_ptr, boost::shared_ptr< VectorDouble > o_ptr, bool symmetrize=true)=0 |
| Evaluate objective function f(coords, u, stress, strain) | |
| virtual MoFEMErrorCode | evalInteriorObjectiveGradientStress (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > stress_ptr, boost::shared_ptr< MatrixDouble > strain_ptr, boost::shared_ptr< MatrixDouble > o_ptr, bool symmetrize=true)=0 |
| Evaluate gradient of objective function w.r.t. stress: ∂f/∂σ | |
| virtual MoFEMErrorCode | evalInteriorObjectiveGradientStrain (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > stress_ptr, boost::shared_ptr< MatrixDouble > strain_ptr, boost::shared_ptr< MatrixDouble > o_ptr, bool symmetrize=true)=0 |
| Evaluate gradient of objective function w.r.t. strain: ∂f/∂ε | |
| virtual MoFEMErrorCode | evalInteriorObjectiveGradientU (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > stress_ptr, boost::shared_ptr< MatrixDouble > strain_ptr, boost::shared_ptr< MatrixDouble > o_ptr, bool symmetrize=true)=0 |
| Evaluate gradient of objective function w.r.t. displacement: ∂f/∂u. | |
| virtual MoFEMErrorCode | evalBoundaryObjectiveFunction (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > t_ptr, boost::shared_ptr< VectorDouble > o_ptr, bool symmetrize=true)=0 |
| virtual MoFEMErrorCode | evalBoundaryObjectiveGradientTraction (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > t_ptr, boost::shared_ptr< MatrixDouble > o_ptr)=0 |
| Evaluate gradient of objective function w.r.t. traction-like vector field. | |
| virtual MoFEMErrorCode | evalBoundaryObjectiveGradientU (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > t_ptr, boost::shared_ptr< MatrixDouble > o_ptr)=0 |
| virtual MoFEMErrorCode | numberOfModes (int block_id, int &modes)=0 |
| Return number of topology optimization modes for given block. | |
| virtual MoFEMErrorCode | blockModes (int block_id, MatrixDouble &coords, std::array< double, 3 > ¢roid, std::array< double, 6 > &bbox, MatrixDouble &o_ptr)=0 |
| Define spatial modes for topology optimization. | |
| virtual | ~ObjectiveFunctionData ()=default |
Abstract interface for Python-defined objective functions.
This interface allows users to define custom objective functions in Python that can be called from C++ during the optimization process. The objective function and its gradients are evaluated at integration points using:
Key methods:
Definition at line 31 of file ObjectiveFunctionData.hpp.
|
virtualdefault |
|
pure virtual |
Define spatial modes for topology optimization.
Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.
|
pure virtual |
Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.
|
pure virtual |
Evaluate gradient of objective function w.r.t. traction-like vector field.
Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.
|
pure virtual |
Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.
|
pure virtual |
Evaluate objective function f(coords, u, stress, strain)
Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.
|
pure virtual |
Evaluate gradient of objective function w.r.t. strain: ∂f/∂ε
Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.
|
pure virtual |
Evaluate gradient of objective function w.r.t. stress: ∂f/∂σ
Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.
|
pure virtual |
Evaluate gradient of objective function w.r.t. displacement: ∂f/∂u.
Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.
|
pure virtual |
Return number of topology optimization modes for given block.
Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.