|  | v0.15.0 | 
Abstract interface for Python-defined objective functions. More...
| Public Member Functions | |
| virtual MoFEMErrorCode | evalObjectiveFunction (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)=0 | 
| Evaluate objective function f(coords, u, stress, strain) | |
| virtual MoFEMErrorCode | evalObjectiveGradientStress (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)=0 | 
| Evaluate gradient of objective function w.r.t. stress: ∂f/∂σ | |
| virtual MoFEMErrorCode | evalObjectiveGradientStrain (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)=0 | 
| Evaluate gradient of objective function w.r.t. strain: ∂f/∂ε | |
| 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 124 of file adjoint.cpp.
| 
 | virtualdefault | 
| 
 | pure virtual | 
Define spatial modes for topology optimization.
Implemented in ObjectiveFunctionDataImpl.
| 
 | pure virtual | 
Evaluate objective function f(coords, u, stress, strain)
Implemented in ObjectiveFunctionDataImpl.
| 
 | pure virtual | 
Evaluate gradient of objective function w.r.t. strain: ∂f/∂ε 
 
Implemented in ObjectiveFunctionDataImpl.
| 
 | pure virtual | 
Evaluate gradient of objective function w.r.t. stress: ∂f/∂σ
Implemented in ObjectiveFunctionDataImpl.
| 
 | pure virtual | 
Return number of topology optimization modes for given block.
Implemented in ObjectiveFunctionDataImpl.
