v0.15.5
Loading...
Searching...
No Matches
Public Member Functions | Private Member Functions | Private Attributes | List of all members
ShapeOptimization::ObjectiveFunctionDataImpl Struct Reference

Implementation of ObjectiveFunctionData interface using Python integration. More...

Inheritance diagram for ShapeOptimization::ObjectiveFunctionDataImpl:
[legend]
Collaboration diagram for ShapeOptimization::ObjectiveFunctionDataImpl:
[legend]

Public Member Functions

 ObjectiveFunctionDataImpl ()=default
 
virtual ~ObjectiveFunctionDataImpl ()=default
 
MoFEMErrorCode initPython (const std::string py_file)
 Initialize Python interpreter and load objective function script.
 
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)
 Evaluate objective function at current state.
 
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)
 Compute gradient of objective function with respect to stress.
 
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)
 Compute gradient of objective function with respect to strain.
 
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)
 Compute gradient of objective function with respect to displacement.
 
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)
 
MoFEMErrorCode evalBoundaryObjectiveGradientTraction (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > t_ptr, boost::shared_ptr< MatrixDouble > o_ptr)
 Evaluate gradient of objective function w.r.t. traction-like vector field.
 
MoFEMErrorCode evalBoundaryObjectiveGradientU (MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > t_ptr, boost::shared_ptr< MatrixDouble > o_ptr)
 
MoFEMErrorCode numberOfModes (int block_id, int &modes)
 Return number of topology optimization modes for given material block.
 
MoFEMErrorCode blockModes (int block_id, MatrixDouble &coords, std::array< double, 3 > &centroid, std::array< double, 6 > &bbodx, MatrixDouble &o_ptr)
 Define spatial topology modes for design optimization.
 
- Public Member Functions inherited from ShapeOptimization::ObjectiveFunctionData
virtual ~ObjectiveFunctionData ()=default
 

Private Member Functions

MoFEMErrorCode interiorObjectiveFunctionImpl (np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
 Internal implementation for objective function evaluation.
 
MoFEMErrorCode interiorObjectiveGradientStressImpl (np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
 Internal implementation for stress gradient computation.
 
MoFEMErrorCode interiorObjectiveGradientStrainImpl (np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
 Internal implementation for strain gradient computation.
 
MoFEMErrorCode interiorObjectiveGradientUImpl (np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
 Internal implementation for displacement gradient computation.
 
MoFEMErrorCode boundaryObjectiveGradientTractionImpl (np::ndarray coords, np::ndarray u, np::ndarray t, np::ndarray &o)
 
MoFEMErrorCode boundaryObjectiveFunctionImpl (np::ndarray coords, np::ndarray u, np::ndarray t, np::ndarray &o)
 
MoFEMErrorCode boundaryObjectiveGradientUImpl (np::ndarray coords, np::ndarray u, np::ndarray t, np::ndarray &o)
 
MoFEMErrorCode blockModesImpl (int block_id, np::ndarray coords, np::ndarray centroid, np::ndarray bbodx, np::ndarray &o_ptr)
 Internal implementation for topology mode generation.
 
np::ndarray convertToNumPy (std::vector< double > &data, int rows, int nb_gauss_pts)
 Convert std::vector to NumPy array for Python interface.
 
np::ndarray convertToNumPy (double *ptr, int s)
 Convert raw pointer to NumPy array for Python interface.
 
MatrixDouble copyToFull (MatrixDouble &s)
 Convert symmetric tensor storage to full matrix format.
 
void copyToSymmetric (double *ptr, MatrixDouble &s)
 Convert full matrix to symmetric tensor storage format.
 

Private Attributes

bp::object mainNamespace
 Main Python namespace for script execution.
 

Detailed Description

Implementation of ObjectiveFunctionData interface using Python integration.

This class provides a concrete implementation of the ObjectiveFunctionData interface that bridges MoFEM C++ data structures with Python-defined objective functions. It enables flexible definition of optimization objectives through Python scripting while maintaining high-performance computation in the C++ finite element framework.

Key features:

The class handles:

  1. Loading and executing Python objective function scripts
  2. Converting MoFEM data structures to NumPy arrays
  3. Calling Python functions for objective evaluation
  4. Converting Python results back to MoFEM format
  5. Managing Python interpreter state and namespace

Example Python interface functions that must be defined:

Definition at line 51 of file ObjectiveFunctionData.cpp.

Constructor & Destructor Documentation

◆ ObjectiveFunctionDataImpl()

ShapeOptimization::ObjectiveFunctionDataImpl::ObjectiveFunctionDataImpl ( )
default

◆ ~ObjectiveFunctionDataImpl()

virtual ShapeOptimization::ObjectiveFunctionDataImpl::~ObjectiveFunctionDataImpl ( )
virtualdefault

Member Function Documentation

◆ blockModes()

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

Define spatial topology modes for design optimization.

Generate spatial topology modes for design optimization.

Generates basis functions that define how the geometry can be modified during topology optimization. These modes serve as design variables and define the design space for optimization.

Parameters
block_idMaterial block identifier
coordsElement coordinates
centroidBlock centroid coordinates
bbodxBounding box dimensions [xmin,xmax,ymin,ymax,zmin,zmax]
o_ptrOutput mode vectors
Returns
MoFEMErrorCode Success or error code

This method defines the design parameterization for topology optimization by generating spatial basis functions (modes) that describe how the geometry can be modified during optimization. These modes serve as design variables and define the feasible design space for the optimization problem.

Mathematical context: The geometry modification is parameterized as: x_new = x_original + Σ(αᵢ * φᵢ(x)) where αᵢ are design variables and φᵢ(x) are spatial mode functions

Common mode types:

  • Radial basis functions: φ(x) = exp(-||x-c||²/σ²) for localized changes
  • Polynomial modes: φ(x) = xⁿyᵐzᵖ for global shape changes
  • Sinusoidal modes: φ(x) = sin(kx)cos(ly) for periodic patterns
  • Principal component modes: Derived from geometric sensitivity analysis

Process:

  1. Query Python function for number of modes for this material block
  2. Convert coordinate data and geometric information to NumPy format
  3. Call Python function block_modes(block_id, coords, centroid, bbox)
  4. Python returns mode vectors for each coordinate at each mode
  5. Reshape and store modes for use as design variables in optimization

The modes enable efficient design space exploration and gradient-based optimization while maintaining geometric feasibility and smoothness.

Parameters
block_idMaterial block identifier for mode generation
coordsElement coordinates where modes are evaluated
centroidGeometric centroid of the material block [x,y,z]
bbodxBounding box dimensions [xmin,xmax,ymin,ymax,zmin,zmax]
o_ptrOutput matrix: modes × (coordinates × spatial_dimension)
Returns
MoFEMErrorCode Success or error code

Implements ShapeOptimization::ObjectiveFunctionData.

Definition at line 956 of file ObjectiveFunctionData.cpp.

958 {
960 try {
961
962 // Query Python function for number of topology modes for this block
963 int nb_modes;
964 CHKERR numberOfModes(block_id, nb_modes);
965
966 // Convert coordinate matrix to NumPy format for Python processing
967 auto np_coords =
968 convertToNumPy(coords.data(), coords.size1(), coords.size2());
969
970 // Convert geometric information to NumPy arrays
971 auto np_centroid =
972 convertToNumPy(centroid.data(), 3); // Block centroid [x,y,z]
973 auto np_bbodx = convertToNumPy(
974 bbodx.data(), 6); // Bounding box [xmin,xmax,ymin,ymax,zmin,zmax]
975
976 // Prepare output array: [modes × (coordinates * spatial_dimensions)]
977 np::ndarray np_output =
978 np::empty(bp::make_tuple(nb_modes, coords.size1(), coords.size2()),
979 np::dtype::get_builtin<double>());
980
981 // Call Python implementation to generate topology modes
982 CHKERR blockModesImpl(block_id, np_coords, np_centroid, np_bbodx,
983 np_output);
984
985 // Check the shape of returned array
986 if (np_output.get_shape()[0] != nb_modes ||
987 np_output.get_shape()[1] != coords.size1() ||
988 np_output.get_shape()[2] != coords.size2()) {
990 "Wrong shape of Modes from python expected (" +
991 std::to_string(nb_modes) + ", " +
992 std::to_string(coords.size1()) + ", " +
993 std::to_string(coords.size2()) + "), got (" +
994 std::to_string(np_output.get_shape()[0]) + ", " +
995 std::to_string(np_output.get_shape()[1]) + ", " +
996 std::to_string(np_output.get_shape()[2]) + ")");
997 }
998
999 // Reshape output matrix for MoFEM format: [modes × (coordinates * spatial_dimensions)]
1000 o_ptr.resize(nb_modes, coords.size1() * coords.size2(), false);
1001 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
1002 // Copy flattened mode data to output matrix
1003 std::copy(val_ptr, val_ptr + coords.size1() * coords.size2() * nb_modes,
1004 o_ptr.data().begin());
1005
1006 } catch (bp::error_already_set const &) {
1007 // Handle Python errors in mode generation
1008 PyErr_Print();
1010 }
1012}
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
MoFEMErrorCode blockModesImpl(int block_id, np::ndarray coords, np::ndarray centroid, np::ndarray bbodx, np::ndarray &o_ptr)
Internal implementation for topology mode generation.
np::ndarray convertToNumPy(std::vector< double > &data, int rows, int nb_gauss_pts)
Convert std::vector to NumPy array for Python interface.
MoFEMErrorCode numberOfModes(int block_id, int &modes)
Return number of topology optimization modes for given material block.

◆ blockModesImpl()

MoFEMErrorCode ShapeOptimization::ObjectiveFunctionDataImpl::blockModesImpl ( int  block_id,
np::ndarray  coords,
np::ndarray  centroid,
np::ndarray  bbodx,
np::ndarray &  o_ptr 
)
private

Internal implementation for topology mode generation.

Calls Python function to generate spatial basis functions for topology optimization. These modes define the design space and parametrize allowable geometry modifications during optimization.

Parameters
block_idMaterial block identifier
coordsNumPy array of element coordinates
centroidNumPy array of block centroid
bbodxNumPy array of bounding box dimensions
o_ptrOutput NumPy array for mode vectors
Returns
MoFEMErrorCode Success or error code

Definition at line 1243 of file ObjectiveFunctionData.cpp.

1247 {
1249 try {
1250 if (bp::extract<bool>(mainNamespace.attr("__contains__")("block_modes"))) {
1251 o = bp::extract<np::ndarray>(
1252 mainNamespace["block_modes"](block_id, coords, centroid, bbodx));
1253 } else {
1256 "Python function block_modes(block_id,coords,centroid,bbox) is not "
1257 "defined");
1258 }
1259 } catch (bp::error_already_set const &) {
1260 // print all other errors to stderr
1261 PyErr_Print();
1263 }
1265}
bp::object mainNamespace
Main Python namespace for script execution.

◆ boundaryObjectiveFunctionImpl()

MoFEMErrorCode ShapeOptimization::ObjectiveFunctionDataImpl::boundaryObjectiveFunctionImpl ( np::ndarray  coords,
np::ndarray  u,
np::ndarray  t,
np::ndarray &  o 
)
private

Definition at line 1167 of file ObjectiveFunctionData.cpp.

1173 {
1175 try {
1176
1177 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_boundary"))) {
1178 o = bp::extract<np::ndarray>(mainNamespace["f_boundary"](coords, u, t));
1179 } else if (bp::extract<bool>(
1180 mainNamespace.attr("__contains__")("f_boundary_function"))) {
1181 o = bp::extract<np::ndarray>(
1182 mainNamespace["f_boundary_function"](coords, u, t));
1183 } else {
1186 "Python function f_boundary(coords,u,t) is not defined");
1187 }
1188
1189 } catch (bp::error_already_set const &) {
1190 PyErr_Print();
1192 }
1194}
constexpr double t
plate stiffness
Definition plate.cpp:58

◆ boundaryObjectiveGradientTractionImpl()

MoFEMErrorCode ShapeOptimization::ObjectiveFunctionDataImpl::boundaryObjectiveGradientTractionImpl ( np::ndarray  coords,
np::ndarray  u,
np::ndarray  t,
np::ndarray &  o 
)
private

Definition at line 1142 of file ObjectiveFunctionData.cpp.

1148 {
1150 try {
1151
1152 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_boundary_t"))) {
1153 o = bp::extract<np::ndarray>(mainNamespace["f_boundary_t"](coords, u, t));
1154 } else {
1157 "Python function f_boundary_t(coords,u,t) is not defined");
1158 }
1159
1160 } catch (bp::error_already_set const &) {
1161 PyErr_Print();
1163 }
1165}

◆ boundaryObjectiveGradientUImpl()

MoFEMErrorCode ShapeOptimization::ObjectiveFunctionDataImpl::boundaryObjectiveGradientUImpl ( np::ndarray  coords,
np::ndarray  u,
np::ndarray  t,
np::ndarray &  o 
)
private

Definition at line 1196 of file ObjectiveFunctionData.cpp.

1202 {
1204 try {
1205
1206 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_boundary_u"))) {
1207 o = bp::extract<np::ndarray>(mainNamespace["f_boundary_u"](coords, u, t));
1208 } else {
1211 "Python function f_boundary_u(coords,u,t) is not defined");
1212 }
1213
1214 } catch (bp::error_already_set const &) {
1215 PyErr_Print();
1217 }
1219}

◆ convertToNumPy() [1/2]

np::ndarray ShapeOptimization::ObjectiveFunctionDataImpl::convertToNumPy ( double ptr,
int  s 
)
inlineprivate

Convert raw pointer to NumPy array for Python interface.

Low-level conversion for direct memory access to create NumPy arrays. Provides zero-copy conversion when possible for performance.

Parameters
ptrRaw data pointer
sArray size
Returns
np::ndarray NumPy array for Python use

Definition at line 1296 of file ObjectiveFunctionData.cpp.

1297 {
1298 auto dtype = np::dtype::get_builtin<double>();
1299 auto size = bp::make_tuple(s);
1300 auto stride = bp::make_tuple(sizeof(double));
1301 return (np::from_data(ptr, dtype, size, stride, bp::object()));
1302}

◆ convertToNumPy() [2/2]

np::ndarray ShapeOptimization::ObjectiveFunctionDataImpl::convertToNumPy ( std::vector< double > &  data,
int  rows,
int  nb_gauss_pts 
)
inlineprivate

Convert std::vector to NumPy array for Python interface.

Converts a std::vector<double> to a NumPy ndarray.

Efficient conversion from MoFEM data structures to NumPy arrays for seamless Python function calls without data copying.

Parameters
dataSource vector data
rowsNumber of rows in resulting array
nb_gauss_ptsNumber of Gauss points (affects array structure)
Returns
np::ndarray NumPy array for Python use

This function wraps the given vector data into a NumPy array with the specified number of rows and Gauss points. The resulting ndarray shares memory with the input vector, so changes to one will affect the other.

Parameters
dataReference to the vector containing double values to be converted.
rowsNumber of rows in the resulting NumPy array.
nb_gauss_ptsNumber of Gauss points (columns) in the resulting NumPy array.
Returns
np::ndarray NumPy array view of the input data.
Note
  • size specifies the shape of the resulting ndarray as a tuple (rows, nb_gauss_pts).
  • stride specifies the step size in bytes to move to the next element in memory. Here, it is set to sizeof(double), indicating contiguous storage for each element.

Definition at line 1288 of file ObjectiveFunctionData.cpp.

1289 {
1290 auto dtype = np::dtype::get_builtin<double>();
1291 auto size = bp::make_tuple(rows, nb_gauss_pts);
1292 auto stride = bp::make_tuple(nb_gauss_pts * sizeof(double), sizeof(double));
1293 return (np::from_data(data.data(), dtype, size, stride, bp::object()));
1294}

◆ copyToFull()

MatrixDouble ShapeOptimization::ObjectiveFunctionDataImpl::copyToFull ( MatrixDouble &  s)
private

Convert symmetric tensor storage to full matrix format.

Definition at line 418 of file ObjectiveFunctionData.cpp.

418 {
419 MatrixDouble f(9, s.size2());
420 f.clear();
423 auto t_f = getFTensor2FromMat<3, 3>(f);
424 auto t_s = getFTensor2SymmetricFromMat<SPACE_DIM>(s);
425 for (int ii = 0; ii != s.size2(); ++ii) {
426 t_f(i, j) = t_s(i, j);
427 ++t_f;
428 ++t_s;
429 }
430 return f;
431};
#define FTENSOR_INDEX(DIM, I)
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
constexpr int SPACE_DIM
Space dimension of problem (2D or 3D), set at compile time.

◆ copyToSymmetric()

void ShapeOptimization::ObjectiveFunctionDataImpl::copyToSymmetric ( double ptr,
MatrixDouble &  s 
)
private

Convert full matrix to symmetric tensor storage format.

Definition at line 433 of file ObjectiveFunctionData.cpp.

433 {
437
438 ptr + 0 * s.size2(), ptr + 1 * s.size2(),
439 ptr + 2 * s.size2(), ptr + 3 * s.size2(),
440 ptr + 4 * s.size2(), ptr + 5 * s.size2(),
441 ptr + 6 * s.size2(), ptr + 7 * s.size2(),
442 ptr + 8 * s.size2()
443
444 };
445 auto t_s = getFTensor2SymmetricFromMat<SPACE_DIM>(s);
446 for (int ii = 0; ii != s.size2(); ++ii) {
447 t_s(i, j) = (t_f(i, j) || t_f(j, i)) / 2.0;
448 ++t_f;
449 ++t_s;
450 }
451}

◆ evalBoundaryObjectiveFunction()

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

Implements ShapeOptimization::ObjectiveFunctionData.

Definition at line 842 of file ObjectiveFunctionData.cpp.

845 {
847 try {
848 (void)symmetrize;
849
850 auto np_coords =
851 convertToNumPy(coords.data(), coords.size1(), coords.size2());
852 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
853 auto np_t = convertToNumPy(t_ptr->data(), t_ptr->size1(), t_ptr->size2());
854
855 np::ndarray np_output = np::empty(bp::make_tuple(t_ptr->size2()),
856 np::dtype::get_builtin<double>());
857
858 CHKERR boundaryObjectiveFunctionImpl(np_coords, np_u, np_t, np_output);
859
860 // Check the shape of returned array
861 if (np_output.get_nd() != 1 || np_output.get_shape()[0] != t_ptr->size2()) {
864 "Wrong shape of Objective Function from python expected (" +
865 std::to_string(t_ptr->size2()) + "), got (" +
866 std::to_string(np_output.get_shape()[0]) + ")");
867 }
868
869 o_ptr->resize(t_ptr->size2(), false);
870 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
871 std::copy(val_ptr, val_ptr + t_ptr->size2(), o_ptr->data().begin());
872
873 } catch (bp::error_already_set const &) {
874 PyErr_Print();
876 }
878}
MoFEMErrorCode boundaryObjectiveFunctionImpl(np::ndarray coords, np::ndarray u, np::ndarray t, np::ndarray &o)

◆ evalBoundaryObjectiveGradientTraction()

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

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

Implements ShapeOptimization::ObjectiveFunctionData.

Definition at line 800 of file ObjectiveFunctionData.cpp.

803 {
805 try {
806
807 auto np_coords =
808 convertToNumPy(coords.data(), coords.size1(), coords.size2());
809 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
810 auto np_t = convertToNumPy(t_ptr->data(), t_ptr->size1(), t_ptr->size2());
811
812 np::ndarray np_output =
813 np::empty(bp::make_tuple(u_ptr->size1(), u_ptr->size2()),
814 np::dtype::get_builtin<double>());
815
816 CHKERR boundaryObjectiveGradientTractionImpl(np_coords, np_u, np_t, np_output);
817
818 // Check the shape of returned array
819 if (np_output.get_shape()[0] != u_ptr->size1() ||
820 np_output.get_shape()[1] != u_ptr->size2()) {
823 "Wrong shape of Objective Gradient from python expected (" +
824 std::to_string(u_ptr->size1()) + ", " +
825 std::to_string(u_ptr->size2()) + "), got (" +
826 std::to_string(np_output.get_shape()[0]) + ", " +
827 std::to_string(np_output.get_shape()[1]) + ")");
828 }
829
830 o_ptr->resize(u_ptr->size1(), u_ptr->size2(), false);
831 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
832 std::copy(val_ptr, val_ptr + u_ptr->size1() * u_ptr->size2(),
833 o_ptr->data().begin());
834
835 } catch (bp::error_already_set const &) {
836 PyErr_Print();
838 }
840}
MoFEMErrorCode boundaryObjectiveGradientTractionImpl(np::ndarray coords, np::ndarray u, np::ndarray t, np::ndarray &o)

◆ evalBoundaryObjectiveGradientU()

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

Implements ShapeOptimization::ObjectiveFunctionData.

Definition at line 880 of file ObjectiveFunctionData.cpp.

883 {
885 try {
886 auto np_coords =
887 convertToNumPy(coords.data(), coords.size1(), coords.size2());
888 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
889 auto np_t = convertToNumPy(t_ptr->data(), t_ptr->size1(), t_ptr->size2());
890
891 np::ndarray np_output =
892 np::empty(bp::make_tuple(u_ptr->size1(), u_ptr->size2()),
893 np::dtype::get_builtin<double>());
894
895 CHKERR boundaryObjectiveGradientUImpl(np_coords, np_u, np_t, np_output);
896
897 // Check the shape of returned array
898 if (np_output.get_shape()[0] != u_ptr->size1() ||
899 np_output.get_shape()[1] != u_ptr->size2()) {
902 "Wrong shape of Objective Gradient from python expected (" +
903 std::to_string(u_ptr->size1()) + ", " +
904 std::to_string(u_ptr->size2()) + "), got (" +
905 std::to_string(np_output.get_shape()[0]) + ", " +
906 std::to_string(np_output.get_shape()[1]) + ")");
907 }
908
909 o_ptr->resize(u_ptr->size1(), u_ptr->size2(), false);
910 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
911 std::copy(val_ptr, val_ptr + u_ptr->size1() * u_ptr->size2(),
912 o_ptr->data().begin());
913
914 } catch (bp::error_already_set const &) {
915 PyErr_Print();
917 }
919}
MoFEMErrorCode boundaryObjectiveGradientUImpl(np::ndarray coords, np::ndarray u, np::ndarray t, np::ndarray &o)

◆ evalInteriorObjectiveFunction()

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

Evaluate objective function at current state.

Evaluate objective function at current finite element state.

Calls Python-defined objective function with current displacement, stress, and strain fields. Used during optimization to compute the objective value that drives the optimization process.

Parameters
coordsGauss point coordinates
u_ptrDisplacement field values
stress_ptrStress tensor values
strain_ptrStrain tensor values
o_ptrOutput objective function values
Returns
MoFEMErrorCode Success or error code

This method bridges MoFEM finite element data with Python-defined objective functions for topology optimization. It handles the complete data conversion workflow from MoFEM matrices to NumPy arrays, calls the Python objective function, and converts results back to MoFEM format.

Process:

  1. Convert coordinate matrix to NumPy format for Python access
  2. Convert displacement field data to NumPy arrays
  3. Convert symmetric stress/strain tensors to full 3x3 matrix format
  4. Call Python objective function: f(coords, u, stress, strain)
  5. Extract results and copy back to MoFEM vector format

The objective function typically computes scalar quantities like:

  • Compliance: ∫ u^T * f dΩ (minimize structural deformation)
  • Stress constraints: ∫ ||σ - σ_target||² dΩ (control stress distribution)
  • Volume constraints: ∫ ρ dΩ (material usage limitations)
Parameters
coordsGauss point coordinates for current element
u_ptrDisplacement field values at Gauss points
stress_ptrCauchy stress tensor values (symmetric storage)
strain_ptrStrain tensor values (symmetric storage)
o_ptrOutput objective function values at each Gauss point
Returns
MoFEMErrorCode Success or error code

Implements ShapeOptimization::ObjectiveFunctionData.

Definition at line 480 of file ObjectiveFunctionData.cpp.

484 {
486 try {
487
488 // Convert coordinates to NumPy array for Python function
489 auto np_coords =
490 convertToNumPy(coords.data(), coords.size1(), coords.size2());
491 // Convert displacement field to NumPy array
492 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
493
494 // Convert symmetric tensor storage to full matrix format for Python
495 // MoFEM stores symmetric tensors in like Voigt notation, Python expects full 3x3 matrices
496 auto full_stress = symmetrize ? copyToFull(*(stress_ptr)) : *(stress_ptr);
497 auto full_strain = symmetrize ? copyToFull(*(strain_ptr)) : *(strain_ptr);
498
499 // Create NumPy arrays for stress and strain tensors
500 auto np_stress = convertToNumPy(full_stress.data(), full_stress.size1(),
501 full_stress.size2());
502 auto np_strain = convertToNumPy(full_strain.data(), full_strain.size1(),
503 full_strain.size2());
504
505 // Prepare output array for objective function values
506 np::ndarray np_output = np::empty(bp::make_tuple(strain_ptr->size2()),
507 np::dtype::get_builtin<double>());
508
509 // Call Python objective function implementation
510 CHKERR interiorObjectiveFunctionImpl(np_coords, np_u, np_stress, np_strain,
511 np_output);
512
513 //Check the shape of returned array
514 if (np_output.get_nd() != 1 ||
515 np_output.get_shape()[0] != strain_ptr->size2()) {
518 "Wrong shape of Objective Function from python expected (" +
519 std::to_string(strain_ptr->size2()) + "), got (" +
520 std::to_string(np_output.get_shape()[0]) + ")");
521 }
522
523 // Copy Python results back to MoFEM vector format
524 o_ptr->resize(stress_ptr->size2(), false);
525 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
526 std::copy(val_ptr, val_ptr + strain_ptr->size2(), o_ptr->data().begin());
527
528 } catch (bp::error_already_set const &) {
529 // Handle Python errors with detailed error reporting
530 PyErr_Print();
532 }
534}
MoFEMErrorCode interiorObjectiveFunctionImpl(np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
Internal implementation for objective function evaluation.
MatrixDouble copyToFull(MatrixDouble &s)
Convert symmetric tensor storage to full matrix format.

◆ evalInteriorObjectiveGradientStrain()

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

Compute gradient of objective function with respect to strain.

Compute gradient of objective function with respect to strain tensor.

Evaluates ∂f/∂ε where f is objective function and ε is strain tensor. Used when objective function depends directly on strain measures, complementing stress-based gradients in adjoint sensitivity analysis.

Parameters
coordsGauss point coordinates
u_ptrDisplacement field values
stress_ptrStress tensor values
strain_ptrStrain tensor values
o_ptrOutput gradient values
Returns
MoFEMErrorCode Success or error code

This method evaluates ∂f/∂ε, the partial derivative of the objective function with respect to the strain tensor. While many structural objectives depend primarily on stress, strain-based gradients are important for certain optimization formulations and provide additional sensitivity information.

Mathematical context: For strain energy-based objectives: f = ½ε:C:ε The gradient is: ∂f/∂ε = C:ε = σ (stress tensor)

For strain-based constraints or objectives like strain concentration: ∂f/∂ε = ∂/∂ε[∫(ε_vm - ε_target)² dΩ] = 2(ε_vm - ε_target) * ∂ε_vm/∂ε

Process:

  1. Convert field data to NumPy format for Python compatibility
  2. Call Python function f_strain(coords, u, stress, strain)
  3. Python returns gradient matrices ∂f/∂ε for each Gauss point
  4. Convert from full 3x3 format back to symmetric storage
  5. Results used in adjoint sensitivity analysis

This gradient complements stress-based gradients in comprehensive sensitivity analysis for topology optimization problems.

Parameters
coordsGauss point coordinates
u_ptrDisplacement field values
stress_ptrCurrent stress tensor values (symmetric storage)
strain_ptrCurrent strain tensor values (symmetric storage)
o_ptrOutput strain gradients ∂f/∂ε (symmetric storage)
Returns
MoFEMErrorCode Success or error code

Implements ShapeOptimization::ObjectiveFunctionData.

Definition at line 654 of file ObjectiveFunctionData.cpp.

658 {
660 try {
661
662 // Convert coordinates and displacement data to NumPy format
663 auto np_coords =
664 convertToNumPy(coords.data(), coords.size1(), coords.size2());
665 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
666
667 // Convert symmetric tensor data to full 3x3 matrices for Python
668 auto full_stress = symmetrize ? copyToFull(*(stress_ptr)) : *(stress_ptr);
669 auto full_strain = symmetrize ? copyToFull(*(strain_ptr)) : *(strain_ptr);
670
671 auto np_stress = convertToNumPy(full_stress.data(), full_stress.size1(),
672 full_stress.size2());
673 auto np_strain = convertToNumPy(full_strain.data(), full_strain.size1(),
674 full_strain.size2());
675
676 // Prepare output array for strain gradients
677 np::ndarray np_output =
678 np::empty(bp::make_tuple(full_strain.size1(), full_strain.size2()),
679 np::dtype::get_builtin<double>());
680
681 // Call Python implementation for strain gradient computation
682 CHKERR interiorObjectiveGradientStrainImpl(np_coords, np_u, np_stress, np_strain,
683 np_output);
684
685 // Check the shape of returned array
686 if (np_output.get_shape()[0] != full_strain.size1() ||
687 np_output.get_shape()[1] != full_strain.size2()) {
690 "Wrong shape of Objective Gradient from python expected (" +
691 std::to_string(full_strain.size1()) + ", " +
692 std::to_string(full_strain.size2()) + "), got (" +
693 std::to_string(np_output.get_shape()[0]) + ", " +
694 std::to_string(np_output.get_shape()[1]) + ")");
695 }
696
697 o_ptr->resize(strain_ptr->size1(), strain_ptr->size2(), false);
698 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
699 copyToSymmetric(val_ptr, *(o_ptr));
700
701 } catch (bp::error_already_set const &) {
702 PyErr_Print();
704 }
706}
MoFEMErrorCode interiorObjectiveGradientStrainImpl(np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
Internal implementation for strain gradient computation.
void copyToSymmetric(double *ptr, MatrixDouble &s)
Convert full matrix to symmetric tensor storage format.

◆ evalInteriorObjectiveGradientStress()

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

Compute gradient of objective function with respect to stress.

Compute gradient of objective function with respect to stress tensor.

Evaluates ∂f/∂σ where f is objective function and σ is stress tensor. This gradient is used in the adjoint method to compute sensitivities efficiently. Essential for gradient-based topology optimization.

Parameters
coordsGauss point coordinates
u_ptrDisplacement field values
stress_ptrStress tensor values
strain_ptrStrain tensor values
o_ptrOutput gradient values
Returns
MoFEMErrorCode Success or error code

This method evaluates ∂f/∂σ, the partial derivative of the objective function with respect to the Cauchy stress tensor. This gradient is fundamental to the adjoint method for topology optimization, as it provides the driving force for the adjoint equation solution.

Mathematical context: The adjoint method requires ∂f/∂σ to compute sensitivities efficiently. For stress-based objectives like von Mises stress constraints: ∂f/∂σ = ∂/∂σ[∫(σ_vm - σ_target)² dΩ] = 2(σ_vm - σ_target) * ∂σ_vm/∂σ

Process:

  1. Convert all field data to NumPy format for Python processing
  2. Call Python function f_stress(coords, u, stress, strain)
  3. Python returns full 3x3 gradient matrices for each Gauss point
  4. Convert back to symmetric tensor storage used by MoFEM
  5. Store results for use in adjoint equation assembly

The resulting gradients drive the adjoint solution that enables efficient computation of design sensitivities independent of design variable count.

Parameters
coordsGauss point coordinates
u_ptrDisplacement field values
stress_ptrCurrent stress tensor values (symmetric storage)
strain_ptrCurrent strain tensor values (symmetric storage)
o_ptrOutput stress gradients ∂f/∂σ (symmetric storage)
Returns
MoFEMErrorCode Success or error code

Implements ShapeOptimization::ObjectiveFunctionData.

Definition at line 566 of file ObjectiveFunctionData.cpp.

570 {
572 try {
573
574 // Convert coordinates and displacement field to NumPy format
575 auto np_coords =
576 convertToNumPy(coords.data(), coords.size1(), coords.size2());
577 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
578
579 // Convert symmetric tensors to full 3x3 format for Python processing
580 auto full_stress = symmetrize ? copyToFull(*(stress_ptr)) : *(stress_ptr);
581 auto full_strain = symmetrize ? copyToFull(*(strain_ptr)) : *(strain_ptr);
582
583 // Create NumPy arrays for stress and strain tensors
584 auto np_stress = convertToNumPy(full_stress.data(), full_stress.size1(),
585 full_stress.size2());
586 auto np_strain = convertToNumPy(full_strain.data(), full_strain.size1(),
587 full_strain.size2());
588 // Prepare output array for stress gradients (full matrix format)
589 np::ndarray np_output =
590 np::empty(bp::make_tuple(full_strain.size1(), full_strain.size2()),
591 np::dtype::get_builtin<double>());
592
593 // Call Python implementation for stress gradient computation
594 CHKERR interiorObjectiveGradientStressImpl(np_coords, np_u, np_stress, np_strain,
595 np_output);
596
597 // Check the shape of returned array
598 if (np_output.get_shape()[0] != full_strain.size1() ||
599 np_output.get_shape()[1] != full_strain.size2()) {
602 "Wrong shape of Objective Gradient from python expected (" +
603 std::to_string(full_strain.size1()) + ", " +
604 std::to_string(full_strain.size2()) + "), got (" +
605 std::to_string(np_output.get_shape()[0]) + ", " +
606 std::to_string(np_output.get_shape()[1]) + ")");
607 }
608
609 // Prepare output matrix in symmetric format
610 o_ptr->resize(stress_ptr->size1(), stress_ptr->size2(), false);
611 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
612 // Convert full matrix results back to symmetric tensor storage
613 copyToSymmetric(val_ptr, *(o_ptr));
614
615 } catch (bp::error_already_set const &) {
616 PyErr_Print();
618 }
620}
MoFEMErrorCode interiorObjectiveGradientStressImpl(np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
Internal implementation for stress gradient computation.

◆ evalInteriorObjectiveGradientU()

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

Compute gradient of objective function with respect to displacement.

Compute gradient of objective function with respect to displacement field.

Evaluates ∂f/∂u where f is objective function and u is displacement field. This provides direct sensitivity of objective to displacement changes, used as right-hand side in adjoint equation K^T * λ = ∂f/∂u.

Parameters
coordsGauss point coordinates
u_ptrDisplacement field values
stress_ptrStress tensor values
strain_ptrStrain tensor values
o_ptrOutput gradient values
Returns
MoFEMErrorCode Success or error code

This method evaluates ∂f/∂u, the partial derivative of the objective function
with respect to the displacement field. This gradient is crucial for the adjoint method as it forms the right-hand side of the adjoint equation: K^T * λ = ∂f/∂u

Mathematical context: The adjoint method solves: K^T * λ = ∂f/∂u where λ are the adjoint variables (Lagrange multipliers)

For compliance minimization: f = ½u^T * K * u The gradient is: ∂f/∂u = K * u (applied forces)

For displacement-based constraints: f = ||u - u_target||² The gradient is: ∂f/∂u = 2(u - u_target)

Process:

  1. Convert all field data to NumPy format for Python processing
  2. Call Python function f_u(coords, u, stress, strain)
  3. Python returns displacement gradients for each component
  4. Copy results directly (no tensor conversion needed for vectors)
  5. Results drive adjoint equation solution for sensitivity analysis

This gradient is fundamental to adjoint-based topology optimization, enabling efficient sensitivity computation for any number of design variables.

Parameters
coordsGauss point coordinates
u_ptrDisplacement field values
stress_ptrCurrent stress tensor values
strain_ptrCurrent strain tensor values
o_ptrOutput displacement gradients ∂f/∂u
Returns
MoFEMErrorCode Success or error code

Implements ShapeOptimization::ObjectiveFunctionData.

Definition at line 742 of file ObjectiveFunctionData.cpp.

746 {
748 try {
749
750 // Convert coordinates and displacement field to NumPy format
751 auto np_coords =
752 convertToNumPy(coords.data(), coords.size1(), coords.size2());
753 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
754
755 // Convert stress and strain tensors to full matrix format
756 auto full_stress = symmetrize ? copyToFull(*(stress_ptr)) : *(stress_ptr);
757 auto full_strain = symmetrize ? copyToFull(*(strain_ptr)) : *(strain_ptr);
758
759 auto np_stress = convertToNumPy(full_stress.data(), full_stress.size1(),
760 full_stress.size2());
761 auto np_strain = convertToNumPy(full_strain.data(), full_strain.size1(),
762 full_strain.size2());
763
764 // Prepare output array for displacement gradients (same size as displacement field)
765 np::ndarray np_output =
766 np::empty(bp::make_tuple(u_ptr->size1(), u_ptr->size2()),
767 np::dtype::get_builtin<double>());
768
769 // Call Python implementation for displacement gradient computation
770 // Note: This should call interiorObjectiveGradientUImpl, not interiorObjectiveGradientStrainImpl
771 CHKERR interiorObjectiveGradientUImpl(np_coords, np_u, np_stress, np_strain,
772 np_output);
773
774 // Check the shape of returned array
775 if (np_output.get_shape()[0] != u_ptr->size1() ||
776 np_output.get_shape()[1] != u_ptr->size2()) {
779 "Wrong shape of Objective Gradient from python expected (" +
780 std::to_string(u_ptr->size1()) + ", " +
781 std::to_string(u_ptr->size2()) + "), got (" +
782 std::to_string(np_output.get_shape()[0]) + ", " +
783 std::to_string(np_output.get_shape()[1]) + ")");
784 }
785
786 // Copy results directly to output matrix (no tensor conversion needed for vectors)
787 o_ptr->resize(u_ptr->size1(), u_ptr->size2(), false);
788 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
789 std::copy(val_ptr, val_ptr + u_ptr->size1() * u_ptr->size2(),
790 o_ptr->data().begin());
791
792 } catch (bp::error_already_set const &) {
793 // Handle Python errors with detailed reporting
794 PyErr_Print();
796 }
798}
MoFEMErrorCode interiorObjectiveGradientUImpl(np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
Internal implementation for displacement gradient computation.

◆ initPython()

MoFEMErrorCode ShapeOptimization::ObjectiveFunctionDataImpl::initPython ( const std::string  py_file)

Initialize Python interpreter and load objective function script.

This method sets up the Python environment for objective function evaluation by loading a Python script that defines the required optimization functions. It establishes the bridge between MoFEM's C++ finite element computations and user-defined Python objective functions.

The Python script must define the following functions:

  • f(coords, u, stress, strain): Main objective function
  • f_stress(coords, u, stress, strain): Gradient w.r.t. stress ∂f/∂σ
  • f_strain(coords, u, stress, strain): Gradient w.r.t. strain ∂f/∂ε
  • f_u(coords, u, stress, strain): Gradient w.r.t. displacement ∂f/∂u
  • number_of_modes(block_id): Return number of topology modes
  • block_modes(block_id, coords, centroid, bbox): Define topology modes

All functions receive NumPy arrays and must return NumPy arrays of appropriate dimensions for the finite element computation.

Parameters
py_filePath to Python script containing objective function definitions
Returns
MoFEMErrorCode Success or error code
Exceptions
MOFEM_OPERATION_UNSUCCESSFULif Python script has errors

Definition at line 397 of file ObjectiveFunctionData.cpp.

397 {
399 try {
400
401 // Create main Python module and namespace for script execution
402 auto main_module = bp::import("__main__");
403 mainNamespace = main_module.attr("__dict__");
404
405 // Execute the Python script in the main namespace
406 bp::exec_file(py_file.c_str(), mainNamespace, mainNamespace);
407
408 // Python callbacks are resolved lazily with explicit existence checks
409
410 } catch (bp::error_already_set const &) {
411 // Handle Python errors by printing to stderr and throwing MoFEM exception
412 PyErr_Print();
414 }
416}

◆ interiorObjectiveFunctionImpl()

MoFEMErrorCode ShapeOptimization::ObjectiveFunctionDataImpl::interiorObjectiveFunctionImpl ( np::ndarray  coords,
np::ndarray  u,
np::ndarray  stress,
np::ndarray  strain,
np::ndarray &  o 
)
private

Internal implementation for objective function evaluation.

Handles low-level Python function call with NumPy array conversion. Converts MoFEM matrices to NumPy arrays, calls Python function, and handles return value conversion.

Parameters
coordsNumPy array of coordinates
uNumPy array of displacements
stressNumPy array of stress tensors
strainNumPy array of strain tensors
oOutput NumPy array for objective values
Returns
MoFEMErrorCode Success or error code

Definition at line 1014 of file ObjectiveFunctionData.cpp.

1020 {
1022 try {
1023
1024 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f"))) {
1025 // Deprecated: Check for main objective function 'f' first for backward compatibility
1026 o = bp::extract<np::ndarray>(
1027 mainNamespace["f"](coords, u, stress, strain));
1028 } else if (bp::extract<bool>(
1029 mainNamespace.attr("__contains__")("f_interior"))) {
1030 o = bp::extract<np::ndarray>(
1031 mainNamespace["f_interior"](coords, u, stress, strain));
1032 } else {
1035 "Python function f_interior(coords,u,stress,strain) is not defined");
1036 }
1037
1038 } catch (bp::error_already_set const &) {
1039 // print all other errors to stderr
1040 PyErr_Print();
1042 }
1044}

◆ interiorObjectiveGradientStrainImpl()

MoFEMErrorCode ShapeOptimization::ObjectiveFunctionDataImpl::interiorObjectiveGradientStrainImpl ( np::ndarray  coords,
np::ndarray  u,
np::ndarray  stress,
np::ndarray  strain,
np::ndarray &  o 
)
private

Internal implementation for strain gradient computation.

Evaluates ∂f/∂ε through Python interface with NumPy array handling. Provides strain-based sensitivities for comprehensive gradient computation.

Parameters
coordsNumPy array of coordinates
uNumPy array of displacements
stressNumPy array of stress tensors
strainNumPy array of strain tensors
oOutput NumPy array for gradient values
Returns
MoFEMErrorCode Success or error code

Definition at line 1078 of file ObjectiveFunctionData.cpp.

1084 {
1086 try {
1087
1088 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_strain"))) {
1089 // Deprecated: keep legacy name first for backward compatibility
1090 o = bp::extract<np::ndarray>(
1091 mainNamespace["f_strain"](coords, u, stress, strain));
1092 } else if (bp::extract<bool>(
1093 mainNamespace.attr("__contains__")("f_interior_strain"))) {
1094 o = bp::extract<np::ndarray>(
1095 mainNamespace["f_interior_strain"](coords, u, stress, strain));
1096 } else {
1099 "Python function f_interior_strain(coords,u,stress,strain) is not defined");
1100 }
1101
1102 } catch (bp::error_already_set const &) {
1103 // print all other errors to stderr
1104 PyErr_Print();
1106 }
1108}

◆ interiorObjectiveGradientStressImpl()

MoFEMErrorCode ShapeOptimization::ObjectiveFunctionDataImpl::interiorObjectiveGradientStressImpl ( np::ndarray  coords,
np::ndarray  u,
np::ndarray  stress,
np::ndarray  strain,
np::ndarray &  o 
)
private

Internal implementation for stress gradient computation.

Calls Python function to compute ∂f/∂σ with automatic array conversion. Essential for adjoint-based sensitivity analysis in topology optimization.

Parameters
coordsNumPy array of coordinates
uNumPy array of displacements
stressNumPy array of stress tensors
strainNumPy array of strain tensors
oOutput NumPy array for gradient values
Returns
MoFEMErrorCode Success or error code

Definition at line 1046 of file ObjectiveFunctionData.cpp.

1052 {
1054 try {
1055
1056 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_stress"))) {
1057 // Deprecated: keep legacy name first for backward compatibility
1058 o = bp::extract<np::ndarray>(
1059 mainNamespace["f_stress"](coords, u, stress, strain));
1060 } else if (bp::extract<bool>(
1061 mainNamespace.attr("__contains__")("f_interior_stress"))) {
1062 o = bp::extract<np::ndarray>(
1063 mainNamespace["f_interior_stress"](coords, u, stress, strain));
1064 } else {
1067 "Python function f_interior_stress(coords,u,stress,strain) is not defined");
1068 }
1069
1070 } catch (bp::error_already_set const &) {
1071 // print all other errors to stderr
1072 PyErr_Print();
1074 }
1076}

◆ interiorObjectiveGradientUImpl()

MoFEMErrorCode ShapeOptimization::ObjectiveFunctionDataImpl::interiorObjectiveGradientUImpl ( np::ndarray  coords,
np::ndarray  u,
np::ndarray  stress,
np::ndarray  strain,
np::ndarray &  o 
)
private

Internal implementation for displacement gradient computation.

Computes ∂f/∂u through Python interface for adjoint equation right-hand side. This gradient drives the adjoint solution that enables efficient sensitivity computation independent of design variable count.

Parameters
coordsNumPy array of coordinates
uNumPy array of displacements
stressNumPy array of stress tensors
strainNumPy array of strain tensors
oOutput NumPy array for gradient values
Returns
MoFEMErrorCode Success or error code

Definition at line 1110 of file ObjectiveFunctionData.cpp.

1116 {
1118 try {
1119
1120 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_u"))) {
1121 // Deprecated: keep legacy name first for backward compatibility
1122 o = bp::extract<np::ndarray>(
1123 mainNamespace["f_u"](coords, u, stress, strain));
1124 } else if (bp::extract<bool>(
1125 mainNamespace.attr("__contains__")("f_interior_u"))) {
1126 o = bp::extract<np::ndarray>(
1127 mainNamespace["f_interior_u"](coords, u, stress, strain));
1128 } else {
1131 "Python function f_interior_u(coords,u,stress,strain) is not defined");
1132 }
1133
1134 } catch (bp::error_already_set const &) {
1135 // print all other errors to stderr
1136 PyErr_Print();
1138 }
1140}

◆ numberOfModes()

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

Return number of topology optimization modes for given material block.

Implements ShapeOptimization::ObjectiveFunctionData.

Definition at line 1221 of file ObjectiveFunctionData.cpp.

1222 {
1224 try {
1225
1226 if (bp::extract<bool>(
1227 mainNamespace.attr("__contains__")("number_of_modes"))) {
1228 modes = bp::extract<int>(mainNamespace["number_of_modes"](block_id));
1229 } else {
1232 "Python function number_of_modes(block_id) is not defined");
1233 }
1234
1235 } catch (bp::error_already_set const &) {
1236 // print all other errors to stderr
1237 PyErr_Print();
1239 }
1241}

Member Data Documentation

◆ mainNamespace

bp::object ShapeOptimization::ObjectiveFunctionDataImpl::mainNamespace
private

Main Python namespace for script execution.

Definition at line 180 of file ObjectiveFunctionData.cpp.


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