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 876 of file ObjectiveFunctionData.cpp.

878 {
880 try {
881
882 // Query Python function for number of topology modes for this block
883 int nb_modes;
884 CHKERR numberOfModes(block_id, nb_modes);
885
886 // Convert coordinate matrix to NumPy format for Python processing
887 auto np_coords =
888 convertToNumPy(coords.data(), coords.size1(), coords.size2());
889
890 // Convert geometric information to NumPy arrays
891 auto np_centroid =
892 convertToNumPy(centroid.data(), 3); // Block centroid [x,y,z]
893 auto np_bbodx = convertToNumPy(
894 bbodx.data(), 6); // Bounding box [xmin,xmax,ymin,ymax,zmin,zmax]
895
896 // Prepare output array: [coordinates × modes × spatial_dimensions]
897 np::ndarray np_output =
898 np::empty(bp::make_tuple(coords.size1(), nb_modes, 3),
899 np::dtype::get_builtin<double>());
900
901 // Call Python implementation to generate topology modes
902 CHKERR blockModesImpl(block_id, np_coords, np_centroid, np_bbodx,
903 np_output);
904
905 // Reshape output matrix for MoFEM format: [modes × (coordinates * spatial_dimensions)]
906 o_ptr.resize(nb_modes, coords.size1() * coords.size2(), false);
907 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
908 // Copy flattened mode data to output matrix
909 std::copy(val_ptr, val_ptr + coords.size1() * coords.size2() * nb_modes,
910 o_ptr.data().begin());
911
912 } catch (bp::error_already_set const &) {
913 // Handle Python errors in mode generation
914 PyErr_Print();
916 }
918}
#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
#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 1149 of file ObjectiveFunctionData.cpp.

1153 {
1155 try {
1156 if (bp::extract<bool>(mainNamespace.attr("__contains__")("block_modes"))) {
1157 o = bp::extract<np::ndarray>(
1158 mainNamespace["block_modes"](block_id, coords, centroid, bbodx));
1159 } else {
1162 "Python function block_modes(block_id,coords,centroid,bbox) is not "
1163 "defined");
1164 }
1165 } catch (bp::error_already_set const &) {
1166 // print all other errors to stderr
1167 PyErr_Print();
1169 }
1171}
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 1073 of file ObjectiveFunctionData.cpp.

1079 {
1081 try {
1082
1083 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_boundary"))) {
1084 o = bp::extract<np::ndarray>(mainNamespace["f_boundary"](coords, u, t));
1085 } else if (bp::extract<bool>(
1086 mainNamespace.attr("__contains__")("f_boundary_function"))) {
1087 o = bp::extract<np::ndarray>(
1088 mainNamespace["f_boundary_function"](coords, u, t));
1089 } else {
1092 "Python function f_boundary(coords,u,t) is not defined");
1093 }
1094
1095 } catch (bp::error_already_set const &) {
1096 PyErr_Print();
1098 }
1100}
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 1048 of file ObjectiveFunctionData.cpp.

1054 {
1056 try {
1057
1058 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_boundary_t"))) {
1059 o = bp::extract<np::ndarray>(mainNamespace["f_boundary_t"](coords, u, t));
1060 } else {
1063 "Python function f_boundary_t(coords,u,t) is not defined");
1064 }
1065
1066 } catch (bp::error_already_set const &) {
1067 PyErr_Print();
1069 }
1071}

◆ boundaryObjectiveGradientUImpl()

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

Definition at line 1102 of file ObjectiveFunctionData.cpp.

1108 {
1110 try {
1111
1112 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_boundary_u"))) {
1113 o = bp::extract<np::ndarray>(mainNamespace["f_boundary_u"](coords, u, t));
1114 } else {
1117 "Python function f_boundary_u(coords,u,t) is not defined");
1118 }
1119
1120 } catch (bp::error_already_set const &) {
1121 PyErr_Print();
1123 }
1125}

◆ 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 1202 of file ObjectiveFunctionData.cpp.

1203 {
1204 auto dtype = np::dtype::get_builtin<double>();
1205 auto size = bp::make_tuple(s);
1206 auto stride = bp::make_tuple(sizeof(double));
1207 return (np::from_data(ptr, dtype, size, stride, bp::object()));
1208}

◆ 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 1194 of file ObjectiveFunctionData.cpp.

1195 {
1196 auto dtype = np::dtype::get_builtin<double>();
1197 auto size = bp::make_tuple(rows, nb_gauss_pts);
1198 auto stride = bp::make_tuple(nb_gauss_pts * sizeof(double), sizeof(double));
1199 return (np::from_data(data.data(), dtype, size, stride, bp::object()));
1200}

◆ 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 783 of file ObjectiveFunctionData.cpp.

786 {
788 try {
789 (void)symmetrize;
790
791 auto np_coords =
792 convertToNumPy(coords.data(), coords.size1(), coords.size2());
793 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
794 auto np_t = convertToNumPy(t_ptr->data(), t_ptr->size1(), t_ptr->size2());
795
796 np::ndarray np_output = np::empty(bp::make_tuple(t_ptr->size2()),
797 np::dtype::get_builtin<double>());
798
799 CHKERR boundaryObjectiveFunctionImpl(np_coords, np_u, np_t, np_output);
800
801 o_ptr->resize(t_ptr->size2(), false);
802 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
803 std::copy(val_ptr, val_ptr + t_ptr->size2(), o_ptr->data().begin());
804
805 } catch (bp::error_already_set const &) {
806 PyErr_Print();
808 }
810}
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 753 of file ObjectiveFunctionData.cpp.

756 {
758 try {
759
760 auto np_coords =
761 convertToNumPy(coords.data(), coords.size1(), coords.size2());
762 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
763 auto np_t = convertToNumPy(t_ptr->data(), t_ptr->size1(), t_ptr->size2());
764
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 CHKERR boundaryObjectiveGradientTractionImpl(np_coords, np_u, np_t, np_output);
770
771 o_ptr->resize(u_ptr->size1(), u_ptr->size2(), false);
772 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
773 std::copy(val_ptr, val_ptr + u_ptr->size1() * u_ptr->size2(),
774 o_ptr->data().begin());
775
776 } catch (bp::error_already_set const &) {
777 PyErr_Print();
779 }
781}
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 812 of file ObjectiveFunctionData.cpp.

815 {
817 try {
818 auto np_coords =
819 convertToNumPy(coords.data(), coords.size1(), coords.size2());
820 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
821 auto np_t = convertToNumPy(t_ptr->data(), t_ptr->size1(), t_ptr->size2());
822
823 np::ndarray np_output =
824 np::empty(bp::make_tuple(u_ptr->size1(), u_ptr->size2()),
825 np::dtype::get_builtin<double>());
826
827 CHKERR boundaryObjectiveGradientUImpl(np_coords, np_u, np_t, np_output);
828
829 o_ptr->resize(u_ptr->size1(), u_ptr->size2(), false);
830 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
831 std::copy(val_ptr, val_ptr + u_ptr->size1() * u_ptr->size2(),
832 o_ptr->data().begin());
833
834 } catch (bp::error_already_set const &) {
835 PyErr_Print();
837 }
839}
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 // Copy Python results back to MoFEM vector format
514 o_ptr->resize(stress_ptr->size2(), false);
515 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
516 std::copy(val_ptr, val_ptr + strain_ptr->size2(), o_ptr->data().begin());
517
518 } catch (bp::error_already_set const &) {
519 // Handle Python errors with detailed error reporting
520 PyErr_Print();
522 }
524}
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 632 of file ObjectiveFunctionData.cpp.

636 {
638 try {
639
640 // Convert coordinates and displacement data to NumPy format
641 auto np_coords =
642 convertToNumPy(coords.data(), coords.size1(), coords.size2());
643 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
644
645 // Convert symmetric tensor data to full 3x3 matrices for Python
646 auto full_stress = symmetrize ? copyToFull(*(stress_ptr)) : *(stress_ptr);
647 auto full_strain = symmetrize ? copyToFull(*(strain_ptr)) : *(strain_ptr);
648
649 auto np_stress = convertToNumPy(full_stress.data(), full_stress.size1(),
650 full_stress.size2());
651 auto np_strain = convertToNumPy(full_strain.data(), full_strain.size1(),
652 full_strain.size2());
653
654 // Prepare output array for strain gradients
655 np::ndarray np_output =
656 np::empty(bp::make_tuple(full_strain.size1(), full_strain.size2()),
657 np::dtype::get_builtin<double>());
658
659 // Call Python implementation for strain gradient computation
660 CHKERR interiorObjectiveGradientStrainImpl(np_coords, np_u, np_stress, np_strain,
661 np_output);
662 o_ptr->resize(strain_ptr->size1(), strain_ptr->size2(), false);
663 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
664 copyToSymmetric(val_ptr, *(o_ptr));
665
666 } catch (bp::error_already_set const &) {
667 PyErr_Print();
669 }
671}
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 556 of file ObjectiveFunctionData.cpp.

560 {
562 try {
563
564 // Convert coordinates and displacement field to NumPy format
565 auto np_coords =
566 convertToNumPy(coords.data(), coords.size1(), coords.size2());
567 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
568
569 // Convert symmetric tensors to full 3x3 format for Python processing
570 auto full_stress = symmetrize ? copyToFull(*(stress_ptr)) : *(stress_ptr);
571 auto full_strain = symmetrize ? copyToFull(*(strain_ptr)) : *(strain_ptr);
572
573 // Create NumPy arrays for stress and strain tensors
574 auto np_stress = convertToNumPy(full_stress.data(), full_stress.size1(),
575 full_stress.size2());
576 auto np_strain = convertToNumPy(full_strain.data(), full_strain.size1(),
577 full_strain.size2());
578 // Prepare output array for stress gradients (full matrix format)
579 np::ndarray np_output =
580 np::empty(bp::make_tuple(full_strain.size1(), full_strain.size2()),
581 np::dtype::get_builtin<double>());
582
583 // Call Python implementation for stress gradient computation
584 CHKERR interiorObjectiveGradientStressImpl(np_coords, np_u, np_stress, np_strain,
585 np_output);
586
587 // Prepare output matrix in symmetric format
588 o_ptr->resize(stress_ptr->size1(), stress_ptr->size2(), false);
589 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
590 // Convert full matrix results back to symmetric tensor storage
591 copyToSymmetric(val_ptr, *(o_ptr));
592
593 } catch (bp::error_already_set const &) {
594 PyErr_Print();
596 }
598}
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 707 of file ObjectiveFunctionData.cpp.

711 {
713 try {
714
715 // Convert coordinates and displacement field to NumPy format
716 auto np_coords =
717 convertToNumPy(coords.data(), coords.size1(), coords.size2());
718 auto np_u = convertToNumPy(u_ptr->data(), u_ptr->size1(), u_ptr->size2());
719
720 // Convert stress and strain tensors to full matrix format
721 auto full_stress = symmetrize ? copyToFull(*(stress_ptr)) : *(stress_ptr);
722 auto full_strain = symmetrize ? copyToFull(*(strain_ptr)) : *(strain_ptr);
723
724 auto np_stress = convertToNumPy(full_stress.data(), full_stress.size1(),
725 full_stress.size2());
726 auto np_strain = convertToNumPy(full_strain.data(), full_strain.size1(),
727 full_strain.size2());
728
729 // Prepare output array for displacement gradients (same size as displacement field)
730 np::ndarray np_output =
731 np::empty(bp::make_tuple(u_ptr->size1(), u_ptr->size2()),
732 np::dtype::get_builtin<double>());
733
734 // Call Python implementation for displacement gradient computation
735 // Note: This should call interiorObjectiveGradientUImpl, not interiorObjectiveGradientStrainImpl
736 CHKERR interiorObjectiveGradientUImpl(np_coords, np_u, np_stress, np_strain,
737 np_output);
738
739 // Copy results directly to output matrix (no tensor conversion needed for vectors)
740 o_ptr->resize(u_ptr->size1(), u_ptr->size2(), false);
741 double *val_ptr = reinterpret_cast<double *>(np_output.get_data());
742 std::copy(val_ptr, val_ptr + u_ptr->size1() * u_ptr->size2(),
743 o_ptr->data().begin());
744
745 } catch (bp::error_already_set const &) {
746 // Handle Python errors with detailed reporting
747 PyErr_Print();
749 }
751}
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 920 of file ObjectiveFunctionData.cpp.

926 {
928 try {
929
930 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f"))) {
931 // Deprecated: Check for main objective function 'f' first for backward compatibility
932 o = bp::extract<np::ndarray>(
933 mainNamespace["f"](coords, u, stress, strain));
934 } else if (bp::extract<bool>(
935 mainNamespace.attr("__contains__")("f_interior"))) {
936 o = bp::extract<np::ndarray>(
937 mainNamespace["f_interior"](coords, u, stress, strain));
938 } else {
941 "Python function f_interior(coords,u,stress,strain) is not defined");
942 }
943
944 } catch (bp::error_already_set const &) {
945 // print all other errors to stderr
946 PyErr_Print();
948 }
950}

◆ 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 984 of file ObjectiveFunctionData.cpp.

990 {
992 try {
993
994 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_strain"))) {
995 // Deprecated: keep legacy name first for backward compatibility
996 o = bp::extract<np::ndarray>(
997 mainNamespace["f_strain"](coords, u, stress, strain));
998 } else if (bp::extract<bool>(
999 mainNamespace.attr("__contains__")("f_interior_strain"))) {
1000 o = bp::extract<np::ndarray>(
1001 mainNamespace["f_interior_strain"](coords, u, stress, strain));
1002 } else {
1005 "Python function f_interior_strain(coords,u,stress,strain) is not defined");
1006 }
1007
1008 } catch (bp::error_already_set const &) {
1009 // print all other errors to stderr
1010 PyErr_Print();
1012 }
1014}

◆ 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 952 of file ObjectiveFunctionData.cpp.

958 {
960 try {
961
962 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_stress"))) {
963 // Deprecated: keep legacy name first for backward compatibility
964 o = bp::extract<np::ndarray>(
965 mainNamespace["f_stress"](coords, u, stress, strain));
966 } else if (bp::extract<bool>(
967 mainNamespace.attr("__contains__")("f_interior_stress"))) {
968 o = bp::extract<np::ndarray>(
969 mainNamespace["f_interior_stress"](coords, u, stress, strain));
970 } else {
973 "Python function f_interior_stress(coords,u,stress,strain) is not defined");
974 }
975
976 } catch (bp::error_already_set const &) {
977 // print all other errors to stderr
978 PyErr_Print();
980 }
982}

◆ 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 1016 of file ObjectiveFunctionData.cpp.

1022 {
1024 try {
1025
1026 if (bp::extract<bool>(mainNamespace.attr("__contains__")("f_u"))) {
1027 // Deprecated: keep legacy name first for backward compatibility
1028 o = bp::extract<np::ndarray>(
1029 mainNamespace["f_u"](coords, u, stress, strain));
1030 } else if (bp::extract<bool>(
1031 mainNamespace.attr("__contains__")("f_interior_u"))) {
1032 o = bp::extract<np::ndarray>(
1033 mainNamespace["f_interior_u"](coords, u, stress, strain));
1034 } else {
1037 "Python function f_interior_u(coords,u,stress,strain) is not defined");
1038 }
1039
1040 } catch (bp::error_already_set const &) {
1041 // print all other errors to stderr
1042 PyErr_Print();
1044 }
1046}

◆ 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 1127 of file ObjectiveFunctionData.cpp.

1128 {
1130 try {
1131
1132 if (bp::extract<bool>(
1133 mainNamespace.attr("__contains__")("number_of_modes"))) {
1134 modes = bp::extract<int>(mainNamespace["number_of_modes"](block_id));
1135 } else {
1138 "Python function number_of_modes(block_id) is not defined");
1139 }
1140
1141 } catch (bp::error_already_set const &) {
1142 // print all other errors to stderr
1143 PyErr_Print();
1145 }
1147}

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: