v0.15.5
Loading...
Searching...
No Matches
ObjectiveFunctionData.cpp
Go to the documentation of this file.
1/**
2 * @file ObjectiveFunctionData.cpp
3 * @brief Implementation of Python-integrated objective function interface
4 * @details Implements Python-C++ bridge for objective function evaluation
5 * using Boost.Python and NumPy for topology optimization
6 */
7
8#include <boost/python.hpp>
9#include <boost/python/def.hpp>
10#include <boost/python/numpy.hpp>
11namespace bp = boost::python;
12namespace np = boost::python::numpy;
13
14#include <MoFEM.hpp>
15using namespace MoFEM;
17
19
20constexpr int SPACE_DIM =
21 EXECUTABLE_DIMENSION; ///< Space dimension of problem (2D or 3D), set at compile time
22
23/**
24 * @brief Implementation of ObjectiveFunctionData interface using Python integration
25 *
26 * This class provides a concrete implementation of the ObjectiveFunctionData interface
27 * that bridges MoFEM C++ data structures with Python-defined objective functions.
28 * It enables flexible definition of optimization objectives through Python scripting
29 * while maintaining high-performance computation in the C++ finite element framework.
30 *
31 * Key features:
32 * - Python function evaluation with automatic data conversion
33 * - Support for objective function and gradient computations
34 * - Topology optimization mode definition through Python
35 * - Efficient NumPy array interfacing for large data sets
36 * - Automatic memory management between C++ and Python
37 *
38 * The class handles:
39 * 1. Loading and executing Python objective function scripts
40 * 2. Converting MoFEM data structures to NumPy arrays
41 * 3. Calling Python functions for objective evaluation
42 * 4. Converting Python results back to MoFEM format
43 * 5. Managing Python interpreter state and namespace
44 *
45 * Example Python interface functions that must be defined:
46 * - objectiveInteriorFunction(coords, u, stress, strain) -> objective_value
47 * - objectiveInteriorGradientStress(coords, u, stress, strain) -> gradient_wrt_stress
48 * - numberOfModes(block_id) -> number_of_topology_modes
49 * - blockModes(block_id, coords, centroid, bbox) -> mode_vectors
50 */
53 virtual ~ObjectiveFunctionDataImpl() = default;
54
55 /// Initialize Python interpreter and load objective function script
56 MoFEMErrorCode initPython(const std::string py_file);
57
58 /**
59 * @brief Evaluate objective function at current state
60 *
61 * Calls Python-defined objective function with current displacement,
62 * stress, and strain fields. Used during optimization to compute
63 * the objective value that drives the optimization process.
64 *
65 * @param coords Gauss point coordinates
66 * @param u_ptr Displacement field values
67 * @param stress_ptr Stress tensor values
68 * @param strain_ptr Strain tensor values
69 * @param o_ptr Output objective function values
70 * @return MoFEMErrorCode Success or error code
71 */
73 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
74 boost::shared_ptr<MatrixDouble> stress_ptr,
75 boost::shared_ptr<MatrixDouble> strain_ptr,
76 boost::shared_ptr<VectorDouble> o_ptr, bool symmetrize = true);
77
78 /**
79 * @brief Compute gradient of objective function with respect to stress
80 *
81 * Evaluates ∂f/∂σ where f is objective function and σ is stress tensor.
82 * This gradient is used in the adjoint method to compute sensitivities
83 * efficiently. Essential for gradient-based topology optimization.
84 *
85 * @param coords Gauss point coordinates
86 * @param u_ptr Displacement field values
87 * @param stress_ptr Stress tensor values
88 * @param strain_ptr Strain tensor values
89 * @param o_ptr Output gradient values
90 * @return MoFEMErrorCode Success or error code
91 */
93 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
94 boost::shared_ptr<MatrixDouble> stress_ptr,
95 boost::shared_ptr<MatrixDouble> strain_ptr,
96 boost::shared_ptr<MatrixDouble> o_ptr, bool symmetrize = true);
97
98 /**
99 * @brief Compute gradient of objective function with respect to strain
100 *
101 * Evaluates ∂f/∂ε where f is objective function and ε is strain tensor.
102 * Used when objective function depends directly on strain measures,
103 * complementing stress-based gradients in adjoint sensitivity analysis.
104 *
105 * @param coords Gauss point coordinates
106 * @param u_ptr Displacement field values
107 * @param stress_ptr Stress tensor values
108 * @param strain_ptr Strain tensor values
109 * @param o_ptr Output gradient values
110 * @return MoFEMErrorCode Success or error code
111 */
113 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
114 boost::shared_ptr<MatrixDouble> stress_ptr,
115 boost::shared_ptr<MatrixDouble> strain_ptr,
116 boost::shared_ptr<MatrixDouble> o_ptr, bool symmetrize = true);
117
118 /**
119 * @brief Compute gradient of objective function with respect to displacement
120 *
121 * Evaluates ∂f/∂u where f is objective function and u is displacement field.
122 * This provides direct sensitivity of objective to displacement changes,
123 * used as right-hand side in adjoint equation K^T * λ = ∂f/∂u.
124 *
125 * @param coords Gauss point coordinates
126 * @param u_ptr Displacement field values
127 * @param stress_ptr Stress tensor values
128 * @param strain_ptr Strain tensor values
129 * @param o_ptr Output gradient values
130 * @return MoFEMErrorCode Success or error code
131 */
133 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
134 boost::shared_ptr<MatrixDouble> stress_ptr,
135 boost::shared_ptr<MatrixDouble> strain_ptr,
136 boost::shared_ptr<MatrixDouble> o_ptr, bool symmetrize = true);
137
140 boost::shared_ptr<MatrixDouble> u_ptr,
141 boost::shared_ptr<MatrixDouble> t_ptr,
142 boost::shared_ptr<VectorDouble> o_ptr,
143 bool symmetrize = true);
144
147 boost::shared_ptr<MatrixDouble> u_ptr,
148 boost::shared_ptr<MatrixDouble> t_ptr,
149 boost::shared_ptr<MatrixDouble> o_ptr);
150
153 boost::shared_ptr<MatrixDouble> u_ptr,
154 boost::shared_ptr<MatrixDouble> t_ptr,
155 boost::shared_ptr<MatrixDouble> o_ptr);
156
157 /// Return number of topology optimization modes for given material block
158 MoFEMErrorCode numberOfModes(int block_id, int &modes);
159
160 /**
161 * @brief Define spatial topology modes for design optimization
162 *
163 * Generates basis functions that define how the geometry can be modified
164 * during topology optimization. These modes serve as design variables
165 * and define the design space for optimization.
166 *
167 * @param block_id Material block identifier
168 * @param coords Element coordinates
169 * @param centroid Block centroid coordinates
170 * @param bbodx Bounding box dimensions [xmin,xmax,ymin,ymax,zmin,zmax]
171 * @param o_ptr Output mode vectors
172 * @return MoFEMErrorCode Success or error code
173 */
174 MoFEMErrorCode blockModes(int block_id, MatrixDouble &coords,
175 std::array<double, 3> &centroid,
176 std::array<double, 6> &bbodx, MatrixDouble &o_ptr);
177
178private:
179 // Python interpreter objects for objective function evaluation
180 bp::object mainNamespace; ///< Main Python namespace for script execution
181
182 /**
183 * @brief Internal implementation for objective function evaluation
184 *
185 * Handles low-level Python function call with NumPy array conversion.
186 * Converts MoFEM matrices to NumPy arrays, calls Python function,
187 * and handles return value conversion.
188 *
189 * @param coords NumPy array of coordinates
190 * @param u NumPy array of displacements
191 * @param stress NumPy array of stress tensors
192 * @param strain NumPy array of strain tensors
193 * @param o Output NumPy array for objective values
194 * @return MoFEMErrorCode Success or error code
195 */
197
198 np::ndarray coords, np::ndarray u,
199
200 np::ndarray stress, np::ndarray strain, np::ndarray &o
201
202 );
203
204 /**
205 * @brief Internal implementation for stress gradient computation
206 *
207 * Calls Python function to compute ∂f/∂σ with automatic array conversion.
208 * Essential for adjoint-based sensitivity analysis in topology optimization.
209 *
210 * @param coords NumPy array of coordinates
211 * @param u NumPy array of displacements
212 * @param stress NumPy array of stress tensors
213 * @param strain NumPy array of strain tensors
214 * @param o Output NumPy array for gradient values
215 * @return MoFEMErrorCode Success or error code
216 */
218
219 np::ndarray coords, np::ndarray u,
220
221 np::ndarray stress, np::ndarray strain, np::ndarray &o
222
223 );
224
225 /**
226 * @brief Internal implementation for strain gradient computation
227 *
228 * Evaluates ∂f/∂ε through Python interface with NumPy array handling.
229 * Provides strain-based sensitivities for comprehensive gradient computation.
230 *
231 * @param coords NumPy array of coordinates
232 * @param u NumPy array of displacements
233 * @param stress NumPy array of stress tensors
234 * @param strain NumPy array of strain tensors
235 * @param o Output NumPy array for gradient values
236 * @return MoFEMErrorCode Success or error code
237 */
239
240 np::ndarray coords, np::ndarray u,
241
242 np::ndarray stress, np::ndarray strain, np::ndarray &o
243
244 );
245
246 /**
247 * @brief Internal implementation for displacement gradient computation
248 *
249 * Computes ∂f/∂u through Python interface for adjoint equation right-hand side.
250 * This gradient drives the adjoint solution that enables efficient sensitivity
251 * computation independent of design variable count.
252 *
253 * @param coords NumPy array of coordinates
254 * @param u NumPy array of displacements
255 * @param stress NumPy array of stress tensors
256 * @param strain NumPy array of strain tensors
257 * @param o Output NumPy array for gradient values
258 * @return MoFEMErrorCode Success or error code
259 */
261
262 np::ndarray coords, np::ndarray u,
263
264 np::ndarray stress, np::ndarray strain, np::ndarray &o
265
266 );
267
269
270 np::ndarray coords, np::ndarray u,
271
272 np::ndarray t, np::ndarray &o
273
274 );
275
277
278 np::ndarray coords, np::ndarray u,
279
280 np::ndarray t, np::ndarray &o
281
282 );
283
285
286 np::ndarray coords, np::ndarray u,
287
288 np::ndarray t, np::ndarray &o
289
290 );
291
292 /**
293 * @brief Internal implementation for topology mode generation
294 *
295 * Calls Python function to generate spatial basis functions for topology
296 * optimization. These modes define the design space and parametrize
297 * allowable geometry modifications during optimization.
298 *
299 * @param block_id Material block identifier
300 * @param coords NumPy array of element coordinates
301 * @param centroid NumPy array of block centroid
302 * @param bbodx NumPy array of bounding box dimensions
303 * @param o_ptr Output NumPy array for mode vectors
304 * @return MoFEMErrorCode Success or error code
305 */
307
308 int block_id, np::ndarray coords, np::ndarray centroid, np::ndarray bbodx,
309 np::ndarray &o_ptr
310
311 );
312
313 /**
314 * @brief Convert std::vector to NumPy array for Python interface
315 *
316 * Efficient conversion from MoFEM data structures to NumPy arrays
317 * for seamless Python function calls without data copying.
318 *
319 * @param data Source vector data
320 * @param rows Number of rows in resulting array
321 * @param nb_gauss_pts Number of Gauss points (affects array structure)
322 * @return np::ndarray NumPy array for Python use
323 */
324 np::ndarray convertToNumPy(std::vector<double> &data, int rows,
325 int nb_gauss_pts);
326
327 /**
328 * @brief Convert raw pointer to NumPy array for Python interface
329 *
330 * Low-level conversion for direct memory access to create NumPy arrays.
331 * Provides zero-copy conversion when possible for performance.
332 *
333 * @param ptr Raw data pointer
334 * @param s Array size
335 * @return np::ndarray NumPy array for Python use
336 */
337 np::ndarray convertToNumPy(double *ptr, int s);
338
339 /// Convert symmetric tensor storage to full matrix format
341
342 /// Convert full matrix to symmetric tensor storage format
343 void copyToSymmetric(double *ptr, MatrixDouble &s);
344};
345
346/**
347 * @brief Factory function to create Python-integrated objective function interface
348 *
349 * Creates and initializes an ObjectiveFunctionDataImpl instance that bridges
350 * MoFEM finite element computations with Python-defined objective functions.
351 * This enables flexible objective function definition for topology optimization
352 * while maintaining computational efficiency.
353 *
354 * The Python file must define specific functions with correct signatures:
355 * - objectiveInteriorFunction(coords, u, stress, strain) -> objective_value
356 * - objectiveInteriorGradientStress(coords, u, stress, strain) -> gradient_array
357 * - objectiveInteriorGradientStrain(coords, u, stress, strain) -> gradient_array
358 * - objectiveInteriorGradientU(coords, u, stress, strain) -> gradient_array
359 * - numberOfModes(block_id) -> integer
360 * - blockModes(block_id, coords, centroid, bbox) -> mode_array
361 *
362 * @param py_file Path to Python file containing objective function definitions
363 * @return boost::shared_ptr<ObjectiveFunctionData> Configured objective function interface
364 * @throws MoFEM exception if Python initialization fails
365 */
366boost::shared_ptr<ObjectiveFunctionData>
368 auto ptr = boost::make_shared<ObjectiveFunctionDataImpl>();
369 CHK_THROW_MESSAGE(ptr->initPython(py_file), "init python");
370 return ptr;
371}
372
373/**
374 * @brief Initialize Python interpreter and load objective function script
375 *
376 * This method sets up the Python environment for objective function evaluation
377 * by loading a Python script that defines the required optimization functions.
378 * It establishes the bridge between MoFEM's C++ finite element computations
379 * and user-defined Python objective functions.
380 *
381 * The Python script must define the following functions:
382 * - f(coords, u, stress, strain): Main objective function
383 * - f_stress(coords, u, stress, strain): Gradient w.r.t. stress ∂f/∂σ
384 * - f_strain(coords, u, stress, strain): Gradient w.r.t. strain ∂f/∂ε
385 * - f_u(coords, u, stress, strain): Gradient w.r.t. displacement ∂f/∂u
386 * - number_of_modes(block_id): Return number of topology modes
387 * - block_modes(block_id, coords, centroid, bbox): Define topology modes
388 *
389 * All functions receive NumPy arrays and must return NumPy arrays of
390 * appropriate dimensions for the finite element computation.
391 *
392 * @param py_file Path to Python script containing objective function definitions
393 * @return MoFEMErrorCode Success or error code
394 * @throws MOFEM_OPERATION_UNSUCCESSFUL if Python script has errors
395 */
397ObjectiveFunctionDataImpl::initPython(const std::string py_file) {
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}
417
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};
432
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}
452
453/**
454 * @brief Evaluate objective function at current finite element state
455 *
456 * This method bridges MoFEM finite element data with Python-defined objective
457 * functions for topology optimization. It handles the complete data conversion
458 * workflow from MoFEM matrices to NumPy arrays, calls the Python objective
459 * function, and converts results back to MoFEM format.
460 *
461 * Process:
462 * 1. Convert coordinate matrix to NumPy format for Python access
463 * 2. Convert displacement field data to NumPy arrays
464 * 3. Convert symmetric stress/strain tensors to full 3x3 matrix format
465 * 4. Call Python objective function: f(coords, u, stress, strain)
466 * 5. Extract results and copy back to MoFEM vector format
467 *
468 * The objective function typically computes scalar quantities like:
469 * - Compliance: ∫ u^T * f dΩ (minimize structural deformation)
470 * - Stress constraints: ∫ ||σ - σ_target||² dΩ (control stress distribution)
471 * - Volume constraints: ∫ ρ dΩ (material usage limitations)
472 *
473 * @param coords Gauss point coordinates for current element
474 * @param u_ptr Displacement field values at Gauss points
475 * @param stress_ptr Cauchy stress tensor values (symmetric storage)
476 * @param strain_ptr Strain tensor values (symmetric storage)
477 * @param o_ptr Output objective function values at each Gauss point
478 * @return MoFEMErrorCode Success or error code
479 */
481 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
482 boost::shared_ptr<MatrixDouble> stress_ptr,
483 boost::shared_ptr<MatrixDouble> strain_ptr,
484 boost::shared_ptr<VectorDouble> o_ptr, bool symmetrize) {
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}
535
536/**
537 * @brief Compute gradient of objective function with respect to stress tensor
538 *
539 * This method evaluates ∂f/∂σ, the partial derivative of the objective function
540 * with respect to the Cauchy stress tensor. This gradient is fundamental to the
541 * adjoint method for topology optimization, as it provides the driving force
542 * for the adjoint equation solution.
543 *
544 * Mathematical context:
545 * The adjoint method requires ∂f/∂σ to compute sensitivities efficiently.
546 * For stress-based objectives like von Mises stress constraints:
547 * ∂f/∂σ = ∂/∂σ[∫(σ_vm - σ_target)² dΩ] = 2(σ_vm - σ_target) * ∂σ_vm/∂σ
548 *
549 * Process:
550 * 1. Convert all field data to NumPy format for Python processing
551 * 2. Call Python function f_stress(coords, u, stress, strain)
552 * 3. Python returns full 3x3 gradient matrices for each Gauss point
553 * 4. Convert back to symmetric tensor storage used by MoFEM
554 * 5. Store results for use in adjoint equation assembly
555 *
556 * The resulting gradients drive the adjoint solution that enables efficient
557 * computation of design sensitivities independent of design variable count.
558 *
559 * @param coords Gauss point coordinates
560 * @param u_ptr Displacement field values
561 * @param stress_ptr Current stress tensor values (symmetric storage)
562 * @param strain_ptr Current strain tensor values (symmetric storage)
563 * @param o_ptr Output stress gradients ∂f/∂σ (symmetric storage)
564 * @return MoFEMErrorCode Success or error code
565 */
567 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
568 boost::shared_ptr<MatrixDouble> stress_ptr,
569 boost::shared_ptr<MatrixDouble> strain_ptr,
570 boost::shared_ptr<MatrixDouble> o_ptr, bool symmetrize) {
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}
621
622/**
623 * @brief Compute gradient of objective function with respect to strain tensor
624 *
625 * This method evaluates ∂f/∂ε, the partial derivative of the objective function
626 * with respect to the strain tensor. While many structural objectives depend
627 * primarily on stress, strain-based gradients are important for certain
628 * optimization formulations and provide additional sensitivity information.
629 *
630 * Mathematical context:
631 * For strain energy-based objectives: f = ½ε:C:ε
632 * The gradient is: ∂f/∂ε = C:ε = σ (stress tensor)
633 *
634 * For strain-based constraints or objectives like strain concentration:
635 * ∂f/∂ε = ∂/∂ε[∫(ε_vm - ε_target)² dΩ] = 2(ε_vm - ε_target) * ∂ε_vm/∂ε
636 *
637 * Process:
638 * 1. Convert field data to NumPy format for Python compatibility
639 * 2. Call Python function f_strain(coords, u, stress, strain)
640 * 3. Python returns gradient matrices ∂f/∂ε for each Gauss point
641 * 4. Convert from full 3x3 format back to symmetric storage
642 * 5. Results used in adjoint sensitivity analysis
643 *
644 * This gradient complements stress-based gradients in comprehensive
645 * sensitivity analysis for topology optimization problems.
646 *
647 * @param coords Gauss point coordinates
648 * @param u_ptr Displacement field values
649 * @param stress_ptr Current stress tensor values (symmetric storage)
650 * @param strain_ptr Current strain tensor values (symmetric storage)
651 * @param o_ptr Output strain gradients ∂f/∂ε (symmetric storage)
652 * @return MoFEMErrorCode Success or error code
653 */
655 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
656 boost::shared_ptr<MatrixDouble> stress_ptr,
657 boost::shared_ptr<MatrixDouble> strain_ptr,
658 boost::shared_ptr<MatrixDouble> o_ptr, bool symmetrize) {
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}
707
708/**
709 * @brief Compute gradient of objective function with respect to displacement field
710 *
711 * This method evaluates ∂f/∂u, the partial derivative of the objective function
712 * with respect to the displacement field. This gradient is crucial for the adjoint
713 * method as it forms the right-hand side of the adjoint equation: K^T * λ = ∂f/∂u
714 *
715 * Mathematical context:
716 * The adjoint method solves: K^T * λ = ∂f/∂u
717 * where λ are the adjoint variables (Lagrange multipliers)
718 *
719 * For compliance minimization: f = ½u^T * K * u
720 * The gradient is: ∂f/∂u = K * u (applied forces)
721 *
722 * For displacement-based constraints: f = ||u - u_target||²
723 * The gradient is: ∂f/∂u = 2(u - u_target)
724 *
725 * Process:
726 * 1. Convert all field data to NumPy format for Python processing
727 * 2. Call Python function f_u(coords, u, stress, strain)
728 * 3. Python returns displacement gradients for each component
729 * 4. Copy results directly (no tensor conversion needed for vectors)
730 * 5. Results drive adjoint equation solution for sensitivity analysis
731 *
732 * This gradient is fundamental to adjoint-based topology optimization,
733 * enabling efficient sensitivity computation for any number of design variables.
734 *
735 * @param coords Gauss point coordinates
736 * @param u_ptr Displacement field values
737 * @param stress_ptr Current stress tensor values
738 * @param strain_ptr Current strain tensor values
739 * @param o_ptr Output displacement gradients ∂f/∂u
740 * @return MoFEMErrorCode Success or error code
741 */
743 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
744 boost::shared_ptr<MatrixDouble> stress_ptr,
745 boost::shared_ptr<MatrixDouble> strain_ptr,
746 boost::shared_ptr<MatrixDouble> o_ptr, bool symmetrize) {
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}
799
801 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
802 boost::shared_ptr<MatrixDouble> t_ptr,
803 boost::shared_ptr<MatrixDouble> o_ptr) {
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}
841
843 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
844 boost::shared_ptr<MatrixDouble> t_ptr,
845 boost::shared_ptr<VectorDouble> o_ptr, bool symmetrize) {
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}
879
881 MatrixDouble &coords, boost::shared_ptr<MatrixDouble> u_ptr,
882 boost::shared_ptr<MatrixDouble> t_ptr,
883 boost::shared_ptr<MatrixDouble> o_ptr) {
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}
920
921/**
922 * @brief Generate spatial topology modes for design optimization
923 *
924 * This method defines the design parameterization for topology optimization by
925 * generating spatial basis functions (modes) that describe how the geometry
926 * can be modified during optimization. These modes serve as design variables
927 * and define the feasible design space for the optimization problem.
928 *
929 * Mathematical context:
930 * The geometry modification is parameterized as: x_new = x_original + Σ(αᵢ * φᵢ(x))
931 * where αᵢ are design variables and φᵢ(x) are spatial mode functions
932 *
933 * Common mode types:
934 * - Radial basis functions: φ(x) = exp(-||x-c||²/σ²) for localized changes
935 * - Polynomial modes: φ(x) = xⁿyᵐzᵖ for global shape changes
936 * - Sinusoidal modes: φ(x) = sin(kx)cos(ly) for periodic patterns
937 * - Principal component modes: Derived from geometric sensitivity analysis
938 *
939 * Process:
940 * 1. Query Python function for number of modes for this material block
941 * 2. Convert coordinate data and geometric information to NumPy format
942 * 3. Call Python function block_modes(block_id, coords, centroid, bbox)
943 * 4. Python returns mode vectors for each coordinate at each mode
944 * 5. Reshape and store modes for use as design variables in optimization
945 *
946 * The modes enable efficient design space exploration and gradient-based
947 * optimization while maintaining geometric feasibility and smoothness.
948 *
949 * @param block_id Material block identifier for mode generation
950 * @param coords Element coordinates where modes are evaluated
951 * @param centroid Geometric centroid of the material block [x,y,z]
952 * @param bbodx Bounding box dimensions [xmin,xmax,ymin,ymax,zmin,zmax]
953 * @param o_ptr Output matrix: modes × (coordinates × spatial_dimension)
954 * @return MoFEMErrorCode Success or error code
955 */
957 int block_id, MatrixDouble &coords, std::array<double, 3> &centroid,
958 std::array<double, 6> &bbodx, MatrixDouble &o_ptr) {
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}
1013
1015
1016 np::ndarray coords, np::ndarray u,
1017
1018 np::ndarray stress, np::ndarray strain, np::ndarray &o
1019
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}
1045
1047
1048 np::ndarray coords, np::ndarray u,
1049
1050 np::ndarray stress, np::ndarray strain, np::ndarray &o
1051
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}
1077
1079
1080 np::ndarray coords, np::ndarray u,
1081
1082 np::ndarray stress, np::ndarray strain, np::ndarray &o
1083
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}
1109
1111
1112 np::ndarray coords, np::ndarray u,
1113
1114 np::ndarray stress, np::ndarray strain, np::ndarray &o
1115
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}
1141
1143
1144 np::ndarray coords, np::ndarray u,
1145
1146 np::ndarray t, np::ndarray &o
1147
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}
1166
1168
1169 np::ndarray coords, np::ndarray u,
1170
1171 np::ndarray t, np::ndarray &o
1172
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}
1195
1197
1198 np::ndarray coords, np::ndarray u,
1199
1200 np::ndarray t, np::ndarray &o
1201
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}
1220
1222 int &modes) {
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}
1242
1244 np::ndarray coords,
1245 np::ndarray centroid,
1246 np::ndarray bbodx,
1247 np::ndarray &o) {
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}
1266
1267/**
1268 * @brief Converts a std::vector<double> to a NumPy ndarray.
1269 *
1270 * This function wraps the given vector data into a NumPy array with the
1271 * specified number of rows and Gauss points. The resulting ndarray shares
1272 * memory with the input vector, so changes to one will affect the other.
1273 *
1274 * @param data Reference to the vector containing double values to be converted.
1275 * @param rows Number of rows in the resulting NumPy array.
1276 * @param nb_gauss_pts Number of Gauss points (columns) in the resulting NumPy
1277 * array.
1278 * @return np::ndarray NumPy array view of the input data.
1279 *
1280 * @note
1281 * - `size` specifies the shape of the resulting ndarray as a tuple (rows,
1282 * nb_gauss_pts).
1283 * - `stride` specifies the step size in bytes to move to the next element in
1284 * memory. Here, it is set to sizeof(double), indicating contiguous storage for
1285 * each element.
1286 */
1287inline np::ndarray
1288ObjectiveFunctionDataImpl::convertToNumPy(std::vector<double> &data, int rows,
1289 int nb_gauss_pts) {
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}
1295
1296inline np::ndarray ObjectiveFunctionDataImpl::convertToNumPy(double *ptr,
1297 int s) {
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}
1303
1304} // namespace ShapeOptimization
Interface for Python-based objective function evaluation in topology optimization.
#define FTENSOR_INDEX(DIM, I)
#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.
FTensor::Index< 'i', SPACE_DIM > i
FTensor::Index< 'j', 3 > j
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
constexpr int SPACE_DIM
Space dimension of problem (2D or 3D), set at compile time.
boost::shared_ptr< ObjectiveFunctionData > create_python_objective_function(std::string py_file)
Factory function to create Python-integrated objective function interface.
constexpr double t
plate stiffness
Definition plate.cpp:58
Implementation of ObjectiveFunctionData interface using Python integration.
MoFEMErrorCode blockModesImpl(int block_id, np::ndarray coords, np::ndarray centroid, np::ndarray bbodx, np::ndarray &o_ptr)
Internal implementation for topology mode generation.
MoFEMErrorCode boundaryObjectiveFunctionImpl(np::ndarray coords, np::ndarray u, np::ndarray t, np::ndarray &o)
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 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 boundaryObjectiveGradientUImpl(np::ndarray coords, np::ndarray u, np::ndarray t, np::ndarray &o)
bp::object mainNamespace
Main Python namespace for script execution.
MoFEMErrorCode interiorObjectiveGradientStrainImpl(np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
Internal implementation for strain gradient computation.
MoFEMErrorCode initPython(const std::string py_file)
Initialize Python interpreter and load objective function script.
MoFEMErrorCode interiorObjectiveFunctionImpl(np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
Internal implementation for objective function evaluation.
void copyToSymmetric(double *ptr, MatrixDouble &s)
Convert full matrix to symmetric tensor storage format.
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 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.
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)
np::ndarray convertToNumPy(std::vector< double > &data, int rows, int nb_gauss_pts)
Convert std::vector to NumPy array for Python interface.
MatrixDouble copyToFull(MatrixDouble &s)
Convert symmetric tensor storage to full matrix format.
MoFEMErrorCode interiorObjectiveGradientUImpl(np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
Internal implementation for displacement gradient computation.
MoFEMErrorCode numberOfModes(int block_id, int &modes)
Return number of topology optimization modes for given material block.
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 interiorObjectiveGradientStressImpl(np::ndarray coords, np::ndarray u, np::ndarray stress, np::ndarray strain, np::ndarray &o)
Internal implementation for stress gradient computation.
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 evalBoundaryObjectiveGradientU(MatrixDouble &coords, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< MatrixDouble > t_ptr, boost::shared_ptr< MatrixDouble > o_ptr)
MoFEMErrorCode boundaryObjectiveGradientTractionImpl(np::ndarray coords, np::ndarray u, np::ndarray t, np::ndarray &o)
Abstract interface for Python-defined objective functions.
#define EXECUTABLE_DIMENSION
Definition plastic.cpp:13