v0.15.5
Loading...
Searching...
No Matches
Public Member Functions | List of all members
ShapeOptimization::ObjectiveFunctionData Struct Referenceabstract

Abstract interface for Python-defined objective functions. More...

#include "tutorials/vec-7_shape_optimisation/src/ObjectiveFunctionData.hpp"

Inheritance diagram for ShapeOptimization::ObjectiveFunctionData:
[legend]

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 > &centroid, std::array< double, 6 > &bbox, MatrixDouble &o_ptr)=0
 Define spatial modes for topology optimization.
 
virtual ~ObjectiveFunctionData ()=default
 

Detailed Description

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.

Constructor & Destructor Documentation

◆ ~ObjectiveFunctionData()

virtual ShapeOptimization::ObjectiveFunctionData::~ObjectiveFunctionData ( )
virtualdefault

Member Function Documentation

◆ blockModes()

virtual MoFEMErrorCode ShapeOptimization::ObjectiveFunctionData::blockModes ( int  block_id,
MatrixDouble &  coords,
std::array< double, 3 > &  centroid,
std::array< double, 6 > &  bbox,
MatrixDouble &  o_ptr 
)
pure virtual

Define spatial modes for topology optimization.

Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.

◆ evalBoundaryObjectiveFunction()

virtual MoFEMErrorCode ShapeOptimization::ObjectiveFunctionData::evalBoundaryObjectiveFunction ( MatrixDouble &  coords,
boost::shared_ptr< MatrixDouble >  u_ptr,
boost::shared_ptr< MatrixDouble >  t_ptr,
boost::shared_ptr< VectorDouble >  o_ptr,
bool  symmetrize = true 
)
pure virtual

◆ evalBoundaryObjectiveGradientTraction()

virtual MoFEMErrorCode ShapeOptimization::ObjectiveFunctionData::evalBoundaryObjectiveGradientTraction ( MatrixDouble &  coords,
boost::shared_ptr< MatrixDouble >  u_ptr,
boost::shared_ptr< MatrixDouble >  t_ptr,
boost::shared_ptr< MatrixDouble >  o_ptr 
)
pure virtual

Evaluate gradient of objective function w.r.t. traction-like vector field.

Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.

◆ evalBoundaryObjectiveGradientU()

virtual MoFEMErrorCode ShapeOptimization::ObjectiveFunctionData::evalBoundaryObjectiveGradientU ( MatrixDouble &  coords,
boost::shared_ptr< MatrixDouble >  u_ptr,
boost::shared_ptr< MatrixDouble >  t_ptr,
boost::shared_ptr< MatrixDouble >  o_ptr 
)
pure virtual

◆ evalInteriorObjectiveFunction()

virtual MoFEMErrorCode ShapeOptimization::ObjectiveFunctionData::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 
)
pure virtual

Evaluate objective function f(coords, u, stress, strain)

Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.

◆ evalInteriorObjectiveGradientStrain()

virtual MoFEMErrorCode ShapeOptimization::ObjectiveFunctionData::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 
)
pure virtual

Evaluate gradient of objective function w.r.t. strain: ∂f/∂ε

Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.

◆ evalInteriorObjectiveGradientStress()

virtual MoFEMErrorCode ShapeOptimization::ObjectiveFunctionData::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 
)
pure virtual

Evaluate gradient of objective function w.r.t. stress: ∂f/∂σ

Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.

◆ evalInteriorObjectiveGradientU()

virtual MoFEMErrorCode ShapeOptimization::ObjectiveFunctionData::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 
)
pure virtual

Evaluate gradient of objective function w.r.t. displacement: ∂f/∂u.

Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.

◆ numberOfModes()

virtual MoFEMErrorCode ShapeOptimization::ObjectiveFunctionData::numberOfModes ( int  block_id,
int &  modes 
)
pure virtual

Return number of topology optimization modes for given block.

Implemented in ShapeOptimization::ObjectiveFunctionDataImpl.


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