v0.15.5
Loading...
Searching...
No Matches
ObjectiveFunctionData.hpp
Go to the documentation of this file.
1/**
2 * @file ObjectiveFunctionData.hpp
3 * @brief Interface for Python-based objective function evaluation in topology optimization
4 * @details Provides abstract interface and implementation for objective function evaluation
5 * using Python with NumPy integration, enabling flexible objective function
6 * definition for adjoint-based topology optimization
7 */
8
9#ifndef __OBJECTIVE_FUNCTION_DATA_HPP__
10#define __OBJECTIVE_FUNCTION_DATA_HPP__
11
12namespace ShapeOptimization {
13
14/**
15 * @brief Abstract interface for Python-defined objective functions
16 *
17 * This interface allows users to define custom objective functions in Python
18 * that can be called from C++ during the optimization process. The objective
19 * function and its gradients are evaluated at integration points using:
20 * - Coordinates of integration points
21 * - Displacement field values
22 * - Stress and strain tensors
23 *
24 * Key methods:
25 * - evalInteriorObjectiveFunction: Compute f(x,u,σ,ε)
26 * - evalInteriorObjectiveGradientStress: Compute ∂f/∂σ
27 * - evalInteriorObjectiveGradientStrain: Compute ∂f/∂ε
28 * - numberOfModes: Return number of topology modes for optimization
29 * - blockModes: Define spatial topology modes for shape optimization
30 */
32
33 /// Evaluate objective function f(coords, u, stress, strain)
34 virtual MoFEMErrorCode evalInteriorObjectiveFunction(
35 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
36 boost::shared_ptr<MatrixDouble> stress_ptr,
37 boost::shared_ptr<MatrixDouble> strain_ptr,
38 boost::shared_ptr<VectorDouble> o_ptr, bool symmetrize = true) = 0;
39
40 /// Evaluate gradient of objective function w.r.t. stress: ∂f/∂σ
42 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
43 boost::shared_ptr<MatrixDouble> stress_ptr,
44 boost::shared_ptr<MatrixDouble> strain_ptr,
45 boost::shared_ptr<MatrixDouble> o_ptr, bool symmetrize = true) = 0;
46
47 /// Evaluate gradient of objective function w.r.t. strain: ∂f/∂ε
49 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
50 boost::shared_ptr<MatrixDouble> stress_ptr,
51 boost::shared_ptr<MatrixDouble> strain_ptr,
52 boost::shared_ptr<MatrixDouble> o_ptr, bool symmetrize = true) = 0;
53
54 /// Evaluate gradient of objective function w.r.t. displacement: ∂f/∂u
55 virtual MoFEMErrorCode evalInteriorObjectiveGradientU(
56 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
57 boost::shared_ptr<MatrixDouble> stress_ptr,
58 boost::shared_ptr<MatrixDouble> strain_ptr,
59 boost::shared_ptr<MatrixDouble> o_ptr, bool symmetrize = true) = 0;
60
61 virtual MoFEMErrorCode evalBoundaryObjectiveFunction(
62 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
63 boost::shared_ptr<MatrixDouble> t_ptr,
64 boost::shared_ptr<VectorDouble> o_ptr, bool symmetrize = true) = 0;
65
66 /// Evaluate gradient of objective function w.r.t. traction-like vector field
68 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
69 boost::shared_ptr<MatrixDouble> t_ptr,
70 boost::shared_ptr<MatrixDouble> o_ptr) = 0;
71
72 virtual MoFEMErrorCode evalBoundaryObjectiveGradientU(
73 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
74 boost::shared_ptr<MatrixDouble> t_ptr,
75 boost::shared_ptr<MatrixDouble> o_ptr) = 0;
76
77 /// Return number of topology optimization modes for given block
78 virtual MoFEMErrorCode numberOfModes(int block_id, int &modes) = 0;
79
80 /// Define spatial modes for topology optimization
81 virtual MoFEMErrorCode blockModes(int block_id, MatrixDouble &coords,
82 std::array<double, 3> &centroid,
83 std::array<double, 6> &bbox,
84 MatrixDouble &o_ptr) = 0;
85
86 virtual ~ObjectiveFunctionData() = default;
87};
88
89boost::shared_ptr<ObjectiveFunctionData>
90create_python_objective_function(std::string py_file);
91
92} // namespace ShapeOptimization
93
94#endif // __OBJECTIVE_FUNCTION_DATA_HPP__
boost::shared_ptr< ObjectiveFunctionData > create_python_objective_function(std::string py_file)
Factory function to create Python-integrated objective function interface.
Abstract interface for Python-defined objective 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 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 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 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 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 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 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 blockModes(int block_id, MatrixDouble &coords, std::array< double, 3 > &centroid, std::array< double, 6 > &bbox, MatrixDouble &o_ptr)=0
Define spatial modes for topology optimization.