v0.15.5
Loading...
Searching...
No Matches
adjoint.cpp
Go to the documentation of this file.
1/**
2 * @file adjoint.cpp
3 * @brief Topology optimization using adjoint method with Python objective functions
4 * @details This tutorial demonstrates:
5 * - Adjoint-based sensitivity analysis for structural topology optimization
6 * - Integration with Python for flexible objective function definition
7 * - Higher-order finite element analysis with geometry fields
8 * - Automatic differentiation using MoFEM's adjoint capabilities
9 * @author Anonymous author(s) committing under MIT license
10 * @example mofem/tutorials/vec-7_shape_optimisation/adjoint.cpp
11 *
12 */
13
14#include <boost/python.hpp>
15#include <boost/python/def.hpp>
16#include <boost/python/numpy.hpp>
17namespace bp = boost::python;
18namespace np = boost::python::numpy;
19
20#include <MoFEM.hpp>
21
22using namespace MoFEM;
23
24//! [Constants and material properties]
25constexpr int BASE_DIM = 1; ///< Dimension of the base functions
26
27//! [Define dimension]
28constexpr int SPACE_DIM =
29 EXECUTABLE_DIMENSION; ///< Space dimension of problem (2D or 3D), set at compile time
30
31//! [Define dimension]
32constexpr AssemblyType A = AssemblyType::PETSC; ///< Use PETSc for matrix/vector assembly
33constexpr IntegrationType I =
34 IntegrationType::GAUSS; ///< Use Gauss quadrature for integration
35
36//! [Material properties for linear elasticity]
37constexpr double young_modulus = 1; ///< Young's modulus E
38constexpr double poisson_ratio = 0.3; ///< Poisson's ratio ν
39constexpr double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio)); ///< Bulk modulus K = E/(3(1-2ν))
40constexpr double shear_modulus_G = young_modulus / (2 * (1 + poisson_ratio)); ///< Shear modulus G = E/(2(1+ν))
41
42PetscBool is_plane_strain = PETSC_FALSE; ///< Flag for plane strain vs plane stress in 2D
43//! [Constants and material properties]
44
45//! [Define finite element types and operators]
46using EntData = EntitiesFieldData::EntData; ///< Entity data for field operations
49using DomainEleOp = DomainEle::UserDataOperator; ///< Domain element operators
50using BoundaryEleOp = BoundaryEle::UserDataOperator; ///< Boundary element operators
51//! [Define finite element types and operators]
52
53//! [Boundary condition types]
54struct DomainBCs {}; ///< Domain boundary conditions marker
55struct BoundaryBCs {}; ///< Boundary conditions marker
56
57//! [Natural boundary condition operators]
64
65template <int DIM> struct PostProcEleByDim;
66
67template <> struct PostProcEleByDim<2> {
71};
72
73template <> struct PostProcEleByDim<3> {
77};
78// Here you can see how the template is being used
82// Forward declaration
83template <int SPACE_DIM, IntegrationType I, typename OpBase>
85
86/// Forward declaration of operator for gradient times symmetric tensor operations
87template <int SPACE_DIM, IntegrationType I, typename OpBase>
89
94
96
97#include <ElasticSpring.hpp>
98#include <FluidLevel.hpp>
99#include <CalculateTraction.hpp>
100#include <NaturalDomainBC.hpp>
101#include <NaturalBoundaryBC.hpp>
102#include <HookeOps.hpp>
103#include <ElasticPostProc.hpp>
104
106using namespace ShapeOptimization;
107
109 const std::string block_name, int dim);
110
111MoFEMErrorCode save_range(moab::Interface &moab, const std::string name,
112 const Range r);
113struct Example {
114
115 Example(MoFEM::Interface &m_field) : mField(m_field) {}
116
117 /// Main driver function for the optimization process
119
120private:
121 MoFEM::Interface &mField; ///< Reference to MoFEM interface
122
123 boost::shared_ptr<MatrixDouble> vectorFieldPtr = nullptr; ///< Field values at evaluation points
124
125 // Problem setup methods
126 MoFEMErrorCode readMesh(); ///< Read mesh from file and setup meshsets
127 MoFEMErrorCode setupProblem(); ///< Setup fields, approximation spaces and DOFs
128 MoFEMErrorCode setupAdJoint(); ///< Setup adjoint fields and finite elements
129 MoFEMErrorCode boundaryCondition(); ///< Apply essential boundary conditions
130 MoFEMErrorCode topologyModes(); ///< Compute topology optimization modes
131 MoFEMErrorCode assembleSystem(); ///< Setup operators in finite element pipeline
132
133 // Analysis methods
134 MoFEMErrorCode solveElastic(); ///< Solve forward elastic problem
135 MoFEMErrorCode postprocessElastic(int iter, SmartPetscObj<Vec> adjoint_vector = nullptr); ///< Post-process and output results
136
137 /// Calculate objective function gradient using adjoint method
138 MoFEMErrorCode calculateGradient(PetscReal *objective_function_value,
139 Vec objective_function_gradient,
140 Vec adjoint_vector);
141
142 // Problem configuration
143 FieldApproximationBase base; ///< Choice of finite element basis functions
144 int fieldOrder = 2; ///< Polynomial order for approximation
145
146 // Solver and data management
147 SmartPetscObj<KSP> kspElastic; ///< Linear solver for elastic problem
148 SmartPetscObj<DM> adjointDM; ///< Data manager for adjoint problem
149 boost::shared_ptr<ObjectiveFunctionData> pythonPtr; ///< Interface to Python objective function
150
151 // Topology optimization data
152 std::vector<SmartPetscObj<Vec>> modeVecs; ///< Topology mode vectors (design variables)
153 std::vector<std::array<double, 3>> modeCentroids; ///< Centroids of optimization blocks
154 std::vector<std::array<double, 6>> modeBBoxes; ///< Bounding boxes of optimization blocks
155 SmartPetscObj<Vec> initialGeometry; ///< Initial geometry field
156};
157
158//! [Run topology optimization problem]
159/**
160 * @brief Main driver for topology optimization using adjoint sensitivity analysis
161 *
162 * This function orchestrates the complete topology optimization workflow:
163 * 1. Initialize Python objective function interface
164 * 2. Setup finite element problems (forward and adjoint)
165 * 3. Compute initial elastic solution
166 * 4. Generate topology optimization modes
167 * 5. Run TAO optimization loop with adjoint-based gradients
168 *
169 * The optimization uses TAO (Toolkit for Advanced Optimization) with L-BFGS
170 * algorithm. At each iteration:
171 * - Update geometry based on current design variables
172 * - Solve forward elastic problem
173 * - Compute objective function and gradients using adjoint method
174 * - Post-process results
175 *
176 * @return MoFEMErrorCode Success or error code
177 */
180
181 // Read objective function from Python file
182 char objective_function_file_name[255] = "objective_function.py";
184 PETSC_NULLPTR, PETSC_NULLPTR, "-objective_function",
185 objective_function_file_name, 255, PETSC_NULLPTR);
186
187 // Verify that the Python objective function file exists
188 auto file_exists = [](std::string myfile) {
189 std::ifstream file(myfile.c_str());
190 if (file) {
191 return true;
192 }
193 return false;
194 };
195 if (!file_exists(objective_function_file_name)) {
196 MOFEM_LOG("WORLD", Sev::error) << "Objective function file NOT found: "
197 << objective_function_file_name;
198 CHK_THROW_MESSAGE(MOFEM_NOT_FOUND, "file NOT found");
199 }
200
201 // Create Python interface for objective function
202 pythonPtr = create_python_objective_function(objective_function_file_name);
203
204 char sensitivity_method_name[32] = "adjoint";
205 CHKERR PetscOptionsGetString(PETSC_NULLPTR, PETSC_NULLPTR,
206 "-sensitivity_method", sensitivity_method_name,
207 sizeof(sensitivity_method_name),
208 PETSC_NULLPTR);
209 std::string sensitivity_method = sensitivity_method_name;
210 std::transform(sensitivity_method.begin(), sensitivity_method.end(),
211 sensitivity_method.begin(),
212 [](unsigned char c) { return std::tolower(c); });
213 if (sensitivity_method == "direct") {
215 } else if (sensitivity_method == "adjoint") {
217 } else {
218 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
219 "Unknown -sensitivity_method. Use 'direct' or 'adjoint'.");
220 }
221 MOFEM_LOG("WORLD", Sev::inform)
222 << "Sensitivity method: " << sensitivity_method;
223
224 // Setup finite element problems
225 CHKERR readMesh(); // Read mesh and meshsets
226 CHKERR setupProblem(); // Setup displacement field and geometry field
227 CHKERR setupAdJoint(); // Setup adjoint field for sensitivity analysis
228 CHKERR boundaryCondition(); // Apply essential boundary conditions
229 CHKERR assembleSystem(); // Setup finite element operators
230
231 // Create linear solver for elastic problem
232 auto pip = mField.getInterface<PipelineManager>();
233 kspElastic = pip->createKSP();
234 CHKERR KSPSetFromOptions(kspElastic);
235
236 // Solve initial elastic problem
238 CHKERR postprocessElastic(-1); // Post-process initial solution
239
240 auto create_vec_modes = [&](auto block_name) {
242 auto mesh_mng = mField.getInterface<MeshsetsManager>();
243 auto bcs = mesh_mng->getCubitMeshsetPtr(
244
245 std::regex((boost::format("%s(.*)") % block_name).str())
246
247 );
248
249 int nb_total_modes = 0;
250 for (auto &bc : bcs) {
251 auto id = bc->getMeshsetId();
252 int nb_modes;
253 CHKERR pythonPtr->numberOfModes(id, nb_modes);
254 nb_total_modes += nb_modes;
255 }
256
257 MOFEM_LOG("WORLD", Sev::inform)
258 << "Total number of modes to apply: " << nb_total_modes;
259
260 modeVecs.resize(nb_total_modes);
261
263 };
264
265 auto get_modes_bounding_boxes = [&](auto block_name) {
267 auto mesh_mng = mField.getInterface<MeshsetsManager>();
268 auto bcs = mesh_mng->getCubitMeshsetPtr(
269
270 std::regex((boost::format("%s(.*)") % block_name).str())
271
272 );
273
274 for (auto &bc : bcs) {
275 auto meshset = bc->getMeshset();
276 Range ents;
277 CHKERR mField.get_moab().get_entities_by_handle(meshset, ents, true);
278 Range verts;
279 CHKERR mField.get_moab().get_connectivity(ents, verts, false);
280 std::vector<double> x(verts.size());
281 std::vector<double> y(verts.size());
282 std::vector<double> z(verts.size());
283 CHKERR mField.get_moab().get_coords(verts, x.data(), y.data(), z.data());
284 std::array<double, 3> centroid = {0, 0, 0};
285 for (int i = 0; i != verts.size(); ++i) {
286 centroid[0] += x[i];
287 centroid[1] += y[i];
288 centroid[2] += z[i];
289 }
290 MPI_Allreduce(MPI_IN_PLACE, centroid.data(), 3, MPI_DOUBLE, MPI_SUM,
291 mField.get_comm());
292 int nb_vertex = verts.size();
293 MPI_Allreduce(MPI_IN_PLACE, &nb_vertex, 1, MPI_INT, MPI_SUM,
294 mField.get_comm());
295 if (nb_vertex) {
296 centroid[0] /= nb_vertex;
297 centroid[1] /= nb_vertex;
298 centroid[2] /= nb_vertex;
299 }
300 std::array<double, 6> bbox = {centroid[0], centroid[1], centroid[2],
301 centroid[0], centroid[1], centroid[2]};
302 for (int i = 0; i != verts.size(); ++i) {
303 bbox[0] = std::min(bbox[0], x[i]);
304 bbox[1] = std::min(bbox[1], y[i]);
305 bbox[2] = std::min(bbox[2], z[i]);
306 bbox[3] = std::max(bbox[3], x[i]);
307 bbox[4] = std::max(bbox[4], y[i]);
308 bbox[5] = std::max(bbox[5], z[i]);
309 }
310 MPI_Allreduce(MPI_IN_PLACE, &bbox[0], 3, MPI_DOUBLE, MPI_MIN,
311 mField.get_comm());
312 MPI_Allreduce(MPI_IN_PLACE, &bbox[3], 3, MPI_DOUBLE, MPI_MAX,
313 mField.get_comm());
314
315 MOFEM_LOG("WORLD", Sev::inform)
316 << "Block: " << bc->getName() << " centroid: " << centroid[0] << " "
317 << centroid[1] << " " << centroid[2];
318 MOFEM_LOG("WORLD", Sev::inform)
319 << "Block: " << bc->getName() << " bbox: " << bbox[0] << " "
320 << bbox[1] << " " << bbox[2] << " " << bbox[3] << " " << bbox[4]
321 << " " << bbox[5];
322
323 modeCentroids.push_back(centroid);
324 modeBBoxes.push_back(bbox);
325 }
326
328 };
329
330 auto eval_objective_and_gradient = [](Tao tao, Vec x, PetscReal *f, Vec g,
331 void *ctx) -> PetscErrorCode {
333 auto *ex_ptr = (Example *)ctx;
334
335 int iter;
336 CHKERR TaoGetIterationNumber(tao, &iter);
337
338 auto set_geometry = [&](Vec x) {
340
341 VecScatter ctx;
342 Vec x_local;
343 CHKERR VecScatterCreateToAll(x, &ctx, &x_local);
344 // scatter as many times as you need
345 CHKERR VecScatterBegin(ctx, x, x_local, INSERT_VALUES, SCATTER_FORWARD);
346 CHKERR VecScatterEnd(ctx, x, x_local, INSERT_VALUES, SCATTER_FORWARD);
347 // destroy scatter context and local vector when no longer needed
348 CHKERR VecScatterDestroy(&ctx);
349
350 auto current_geometry = vectorDuplicate(ex_ptr->initialGeometry);
351 CHKERR VecCopy(ex_ptr->initialGeometry, current_geometry);
352 const double *a;
353 CHKERR VecGetArrayRead(x_local, &a);
354 const double *coeff = a;
355 for (auto &mode_vec : ex_ptr->modeVecs) {
356 MOFEM_LOG("WORLD", Sev::verbose)
357 << "Adding mode with coeff: " << *coeff;
358 CHKERR VecAXPY(current_geometry, (*coeff), mode_vec);
359 ++coeff;
360 }
361 CHKERR VecRestoreArrayRead(x_local, &a);
362 CHKERR VecGhostUpdateBegin(current_geometry, INSERT_VALUES,
363 SCATTER_FORWARD);
364 CHKERR VecGhostUpdateEnd(current_geometry, INSERT_VALUES,
365 SCATTER_FORWARD);
366 CHKERR ex_ptr->mField.getInterface<VecManager>()
367 ->setOtherLocalGhostVector("ADJOINT", "ADJOINT_FIELD", "GEOMETRY",
368 RowColData::ROW, current_geometry,
369 INSERT_VALUES, SCATTER_REVERSE);
370
371 CHKERR VecDestroy(&x_local);
373 };
374
375 CHKERR set_geometry(x);
376 CHKERR KSPReset(ex_ptr->kspElastic);
377 CHKERR ex_ptr->solveElastic();
378 auto simple = ex_ptr->mField.getInterface<Simple>();
379 auto adjoint_vector = createDMVector(simple->getDM());
380 CHKERR ex_ptr->calculateGradient(f, g, adjoint_vector);
381 CHKERR ex_ptr->postprocessElastic(iter, adjoint_vector);
382
384 };
385
386 // Create topology optimization modes and setup TAO solver
387 CHKERR create_vec_modes("OPTIMISE");
388 CHKERR get_modes_bounding_boxes("OPTIMISE");
389
390 /**
391 * Setup TAO (Toolkit for Advanced Optimization) solver for topology optimization
392 *
393 * TAO provides various optimization algorithms. Here we use TAOLMVM (Limited Memory
394 * Variable Metric) which is a quasi-Newton method suitable for unconstrained
395 * optimization with gradient information.
396 *
397 * The optimization variables are coefficients for the topology modes, and
398 * gradients are computed using the adjoint method for efficiency.
399 */
400 auto tao = createTao(mField.get_comm());
401 CHKERR TaoSetType(tao, TAOLMVM);
402
403 auto rank = mField.get_comm_rank();
404 auto g = createVectorMPI(mField.get_comm(), (!rank) ? modeVecs.size() : 0,
405 PETSC_DECIDE);
406 CHKERR TaoSetObjectiveAndGradient(tao, g, eval_objective_and_gradient,
407 (void *)this);
408
409 // Store initial geometry for reference during optimization
411 CHKERR mField.getInterface<VecManager>()->setOtherLocalGhostVector(
412 "ADJOINT", "ADJOINT_FIELD", "GEOMETRY", RowColData::ROW, initialGeometry,
413 INSERT_VALUES, SCATTER_FORWARD);
414
415 // Setup topology optimization modes
416 /**
417 * Generate modes for topology optimization design parameterization
418 *
419 * These modes represent perturbations to the geometry that can be used
420 * as design variables in topology optimization. The modes are defined
421 * through the Python interface and provide spatial basis functions
422 * for shape modifications.
423 */
425
426 // Initialize optimization variables and run solver
427 /**
428 * Start optimization with zero initial guess for design variables
429 *
430 * The TAO solver will iteratively:
431 * 1. Evaluate objective function at current design point
432 * 2. Compute gradients using adjoint sensitivity analysis
433 * 3. Update design variables using L-BFGS algorithm
434 * 4. Check convergence criteria
435 *
436 * Each function evaluation involves solving the forward elastic problem
437 * and the adjoint problem to compute sensitivities efficiently.
438 */
439 auto x0 = vectorDuplicate(g);
440
441 CHKERR VecSet(x0, 0.0);
442 CHKERR TaoSetSolution(tao, x0);
443 CHKERR TaoSetFromOptions(tao);
444 CHKERR TaoSolve(tao);
445
446 // Optimization complete - results available in solution vectors
447 MOFEM_LOG("WORLD", Sev::inform) << "Topology optimization completed";
448
450}
451//! [Run problem]
452
453//! [Read mesh]
454/**
455 * @brief Read mesh from file and setup material/boundary condition meshsets
456 *
457 * This function loads the finite element mesh from file and processes
458 * associated meshsets that define material properties and boundary conditions.
459 * The mesh is typically generated using CUBIT and exported in .h5m format.
460 *
461 * Meshsets are used to group elements/faces by:
462 * - Material properties (for different material blocks)
463 * - Boundary conditions (for applying loads and constraints)
464 * - Optimization regions (for topology optimization)
465 *
466 * @return MoFEMErrorCode Success or error code
467 */
471 CHKERR simple->getOptions(); // Read command line options
472 CHKERR simple->loadFile(); // Load mesh from file
473 // Add meshsets if config file provided
474 CHKERR mField.getInterface<MeshsetsManager>()->setMeshsetFromFile();
476}
477//! [Read mesh]
478
479//! [Set up problem]
480/**
481 * @brief Setup finite element fields, approximation spaces and degrees of freedom
482 *
483 * This function configures the finite element problem by:
484 * 1. Setting up the displacement field "U" with vector approximation
485 * 2. Setting up the geometry field "GEOMETRY" for mesh deformation
486 * 3. Defining polynomial approximation order and basis functions
487 * 4. Creating degrees of freedom on mesh entities
488 *
489 * The displacement field uses H1 vector space for standard elasticity.
490 * The geometry field allows mesh modification during topology optimization.
491 * Different basis functions (Ainsworth-Legendre vs Demkowicz) can be selected.
492 *
493 * @return MoFEMErrorCode Success or error code
494 */
498
499 // Select basis functions for finite element approximation
500 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
501 const char *list_bases[LASBASETOPT] = {"ainsworth", "demkowicz"};
502 PetscInt choice_base_value = AINSWORTH;
503 CHKERR PetscOptionsGetEList(PETSC_NULLPTR, NULL, "-base", list_bases,
504 LASBASETOPT, &choice_base_value, PETSC_NULLPTR);
505
506 switch (choice_base_value) {
507 case AINSWORTH:
509 MOFEM_LOG("WORLD", Sev::inform)
510 << "Set AINSWORTH_LEGENDRE_BASE for displacements";
511 break;
512 case DEMKOWICZ:
514 MOFEM_LOG("WORLD", Sev::inform)
515 << "Set DEMKOWICZ_JACOBI_BASE for displacements";
516 break;
517 default:
518 base = LASTBASE;
519 break;
520 }
521
522 // Add finite element fields
523 /**
524 * Setup displacement field "U" - the primary unknown in elasticity
525 * This field represents displacement vector at each node/DOF
526 */
527 CHKERR simple->addDomainField("U", H1, base, SPACE_DIM);
528 CHKERR simple->addBoundaryField("U", H1, base, SPACE_DIM);
529 CHKERR PetscOptionsGetInt(PETSC_NULLPTR, "", "-order", &fieldOrder,
530 PETSC_NULLPTR);
531
532 /**
533 * Setup geometry field "GEOMETRY" - used for mesh deformation in optimization
534 * This field stores current nodal coordinates and can be modified
535 * during topology optimization to represent design changes
536 */
537 CHKERR simple->addDataField("GEOMETRY", H1, base, SPACE_DIM);
538
539 // Set polynomial approximation order for both fields
540 CHKERR simple->setFieldOrder("U", fieldOrder);
541 CHKERR simple->setFieldOrder("GEOMETRY", fieldOrder);
542 CHKERR simple->setUp();
543
544 // Project higher-order geometry representation onto geometry field
545 /**
546 * For higher-order elements, this projects the exact geometry
547 * onto the geometry field to maintain curved boundaries accurately
548 */
549 auto project_ho_geometry = [&]() {
550 Projection10NodeCoordsOnField ent_method(mField, "GEOMETRY");
551 return mField.loop_dofs("GEOMETRY", ent_method);
552 };
553 CHKERR project_ho_geometry();
554
555 // Check if plane strain assumption should be used
556 CHKERR PetscOptionsGetBool(PETSC_NULLPTR, "", "-plane_strain",
557 &is_plane_strain, PETSC_NULLPTR);
558
560}
561//! [Set up problem]
562
563//! [Setup adjoint]
564/**
565 * @brief Setup adjoint fields and finite elements for sensitivity analysis
566 *
567 * The adjoint method is used to efficiently compute gradients of the objective
568 * function with respect to design variables. This function sets up:
569 *
570 * 1. ADJOINT_FIELD - stores adjoint variables (Lagrange multipliers)
571 * 2. ADJOINT_DM - data manager for adjoint problem
572 * 3. Adjoint finite elements for domain and boundary
573 *
574 * The adjoint equation is: K^T * λ = ∂f/∂u
575 * where λ are adjoint variables, K is stiffness matrix, f is objective
576 *
577 * @return MoFEMErrorCode Success or error code
578 */
582
583 // Create adjoint data manager and field
584 auto create_adjoint_dm = [&]() {
585 auto adjoint_dm = createDM(mField.get_comm(), "DMMOFEM");
586
587 auto add_field = [&]() {
589 CHKERR mField.add_field("ADJOINT_FIELD", H1, base, SPACE_DIM);
591 "ADJOINT_FIELD");
592 for (auto tt = MBEDGE; tt <= moab::CN::TypeDimensionMap[SPACE_DIM].second;
593 ++tt)
594 CHKERR mField.set_field_order(simple->getMeshset(), tt, "ADJOINT_FIELD",
595 fieldOrder);
596 CHKERR mField.set_field_order(simple->getMeshset(), MBVERTEX,
597 "ADJOINT_FIELD", 1);
600 };
601
602 auto add_adjoint_fe_impl = [&]() {
604 CHKERR mField.add_finite_element("ADJOINT_DOMAIN_FE");
606 "ADJOINT_FIELD");
608 "ADJOINT_FIELD");
610 "ADJOINT_FIELD");
612 "GEOMETRY");
614 simple->getMeshset(), SPACE_DIM, "ADJOINT_DOMAIN_FE");
615 CHKERR mField.build_finite_elements("ADJOINT_DOMAIN_FE");
616
617 CHKERR mField.add_finite_element("ADJOINT_BOUNDARY_FE");
619 "ADJOINT_FIELD");
621 "ADJOINT_FIELD");
623 "ADJOINT_FIELD");
625 "GEOMETRY");
626
627 auto block_name = "OPTIMISE";
628 auto mesh_mng = mField.getInterface<MeshsetsManager>();
629 auto bcs = mesh_mng->getCubitMeshsetPtr(
630
631 std::regex((boost::format("%s(.*)") % block_name).str())
632
633 );
634
635 for (auto bc : bcs) {
637 bc->getMeshset(), SPACE_DIM - 1, "ADJOINT_BOUNDARY_FE");
638 }
639
640 CHKERR mField.build_finite_elements("ADJOINT_BOUNDARY_FE");
641
642 CHKERR mField.build_adjacencies(simple->getBitRefLevel(),
643 simple->getBitRefLevelMask());
644
646 };
647
648 auto set_adjoint_dm_imp = [&]() {
650 CHKERR DMMoFEMCreateMoFEM(adjoint_dm, &mField, "ADJOINT",
651 simple->getBitRefLevel(),
652 simple->getBitRefLevelMask());
653 CHKERR DMMoFEMSetDestroyProblem(adjoint_dm, PETSC_TRUE);
654 CHKERR DMSetFromOptions(adjoint_dm);
655 CHKERR DMMoFEMAddElement(adjoint_dm, "ADJOINT_DOMAIN_FE");
656 CHKERR DMMoFEMAddElement(adjoint_dm, "ADJOINT_BOUNDARY_FE");
657 CHKERR DMMoFEMSetSquareProblem(adjoint_dm, PETSC_TRUE);
658 CHKERR DMMoFEMSetIsPartitioned(adjoint_dm, PETSC_TRUE);
659 mField.getInterface<ProblemsManager>()->buildProblemFromFields =
660 PETSC_TRUE;
661 CHKERR DMSetUp(adjoint_dm);
662 mField.getInterface<ProblemsManager>()->buildProblemFromFields =
663 PETSC_FALSE;
665 };
666
667 CHK_THROW_MESSAGE(add_field(), "add adjoint field");
668 CHK_THROW_MESSAGE(add_adjoint_fe_impl(), "add adjoint fe");
669 CHK_THROW_MESSAGE(set_adjoint_dm_imp(), "set adjoint dm");
670
671 return adjoint_dm;
672 };
673
674 adjointDM = create_adjoint_dm();
675
677}
678//
679
680//! [Boundary condition]
684 auto bc_mng = mField.getInterface<BcManager>();
685
686 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_X",
687 "U", 0, 0);
688 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Y",
689 "U", 1, 1);
690 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(), "REMOVE_Z",
691 "U", 2, 2);
692 CHKERR bc_mng->removeBlockDOFsOnEntities(simple->getProblemName(),
693 "REMOVE_ALL", "U", 0, 3);
694 CHKERR bc_mng->pushMarkDOFsOnEntities<DisplacementCubitBcData>(
695 simple->getProblemName(), "U");
696
698}
699//! [Boundary condition]
700
701//! [Adjoint modes]
704
705 auto opt_ents = get_range_from_block(mField, "OPTIMISE", SPACE_DIM - 1);
706 auto subset_dm_bdy = createDM(mField.get_comm(), "DMMOFEM");
707 CHKERR DMMoFEMSetSquareProblem(subset_dm_bdy, PETSC_TRUE);
708 CHKERR DMMoFEMCreateSubDM(subset_dm_bdy, adjointDM, "SUBSET_BDY");
709 CHKERR DMMoFEMAddElement(subset_dm_bdy, "ADJOINT_BOUNDARY_FE");
710 CHKERR DMMoFEMAddSubFieldRow(subset_dm_bdy, "ADJOINT_FIELD",
711 boost::make_shared<Range>(opt_ents));
712 CHKERR DMMoFEMAddSubFieldCol(subset_dm_bdy, "ADJOINT_FIELD",
713 boost::make_shared<Range>(opt_ents));
714 CHKERR DMSetUp(subset_dm_bdy);
715
716 auto subset_dm_domain = createDM(mField.get_comm(), "DMMOFEM");
717 CHKERR DMMoFEMSetSquareProblem(subset_dm_domain, PETSC_TRUE);
718 CHKERR DMMoFEMCreateSubDM(subset_dm_domain, adjointDM, "SUBSET_DOMAIN");
719 CHKERR DMMoFEMAddElement(subset_dm_domain, "ADJOINT_DOMAIN_FE");
720 CHKERR DMMoFEMAddSubFieldRow(subset_dm_domain, "ADJOINT_FIELD");
721 CHKERR DMMoFEMAddSubFieldCol(subset_dm_domain, "ADJOINT_FIELD");
722 CHKERR DMSetUp(subset_dm_domain);
723
724 // remove dofs on boundary of the domain
725 auto remove_dofs = [&]() {
727
728 std::array<Range, 3> remove_dim_ents;
729 remove_dim_ents[0] =
730 get_range_from_block(mField, "OPT_REMOVE_X", SPACE_DIM - 1);
731 remove_dim_ents[1] =
732 get_range_from_block(mField, "OPT_REMOVE_Y", SPACE_DIM - 1);
733 remove_dim_ents[2] =
734 get_range_from_block(mField, "OPT_REMOVE_Z", SPACE_DIM - 1);
735
736 for (int d = 0; d != 3; ++d) {
737 MOFEM_LOG("WORLD", Sev::inform)
738 << "Removing topology modes on block OPT_REMOVE_" << (char)('X' + d)
739 << " with " << remove_dim_ents[d].size() << " entities";
740 }
741
742 Range body_ents;
743 CHKERR mField.get_moab().get_entities_by_dimension(0, SPACE_DIM, body_ents,
744 true);
745 auto skin = moab::Skinner(&mField.get_moab());
746 Range boundary_ents;
747 CHKERR skin.find_skin(0, body_ents, false, boundary_ents);
748 for (int d = 0; d != 3; ++d) {
749 boundary_ents = subtract(boundary_ents, remove_dim_ents[d]);
750 }
751 ParallelComm *pcomm =
752 ParallelComm::get_pcomm(&mField.get_moab(), MYPCOMM_INDEX);
753 CHKERR pcomm->filter_pstatus(boundary_ents,
754 PSTATUS_SHARED | PSTATUS_MULTISHARED,
755 PSTATUS_NOT, -1, nullptr);
756 for (auto d = SPACE_DIM - 2; d >= 0; --d) {
757 if (d >= 0) {
758 Range ents;
759 CHKERR mField.get_moab().get_adjacencies(boundary_ents, d, false, ents,
760 moab::Interface::UNION);
761 boundary_ents.merge(ents);
762 } else {
763 Range verts;
764 CHKERR mField.get_moab().get_connectivity(boundary_ents, verts);
765 boundary_ents.merge(verts);
766 }
767 CHKERR mField.getInterface<CommInterface>()->synchroniseEntities(
768 boundary_ents);
769 }
770 boundary_ents.merge(opt_ents);
771 CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
772 "SUBSET_DOMAIN", "ADJOINT_FIELD", boundary_ents);
773 for (int d = 0; d != 3; ++d) {
774 CHKERR mField.getInterface<ProblemsManager>()->removeDofsOnEntities(
775 "SUBSET_DOMAIN", "ADJOINT_FIELD", remove_dim_ents[d], d, d);
776 }
777
778 // #ifndef NDEBUG
779 if (mField.get_comm_rank() == 0) {
780 CHKERR save_range(mField.get_moab(), "topoMode_boundary_ents.vtk",
781 boundary_ents);
782 }
783 // #endif
784
786 };
787
788 CHKERR remove_dofs();
789
790 auto get_lhs_fe = [&]() {
791 auto fe_lhs = boost::make_shared<BoundaryEle>(mField);
792 fe_lhs->getRuleHook = [](int, int, int p_data) {
793 return 2 * p_data + p_data - 1;
794 };
795 auto &pip = fe_lhs->getOpPtrVector();
797 "GEOMETRY");
800 pip.push_back(new OpMass("ADJOINT_FIELD", "ADJOINT_FIELD",
801 [](double, double, double) { return 1.; }));
802 return fe_lhs;
803 };
804
805 auto get_rhs_fe = [&]() {
806 auto fe_rhs = boost::make_shared<BoundaryEle>(mField);
807 fe_rhs->getRuleHook = [](int, int, int p_data) {
808 return 2 * p_data + p_data - 1;
809 };
810 auto &pip = fe_rhs->getOpPtrVector();
812 "GEOMETRY");
813
814 return fe_rhs;
815 };
816
817 auto block_name = "OPTIMISE";
818 auto mesh_mng = mField.getInterface<MeshsetsManager>();
819 auto bcs = mesh_mng->getCubitMeshsetPtr(
820
821 std::regex((boost::format("%s(.*)") % block_name).str())
822
823 );
824
825 for (auto &v : modeVecs) {
826 v = createDMVector(subset_dm_bdy);
827 }
828
830 struct OpMode : public OP {
831 OpMode(const std::string name,
832 boost::shared_ptr<ObjectiveFunctionData> python_ptr, int id,
833 std::vector<SmartPetscObj<Vec>> mode_vecs,
834 std::vector<std::array<double, 3>> mode_centroids,
835 std::vector<std::array<double, 6>> mode_bboxes, int block_counter,
836 int mode_counter, boost::shared_ptr<Range> range = nullptr)
837 : OP(name, name, OP::OPROW, range), pythonPtr(python_ptr), iD(id),
838 modeVecs(mode_vecs), modeCentroids(mode_centroids),
839 modeBboxes(mode_bboxes), blockCounter(block_counter),
840 modeCounter(mode_counter) {}
841
842 MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
844
845 if (OP::entsPtr) {
846 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
848 }
849
850 auto nb_rows = data.getIndices().size();
851 if (!nb_rows) {
853 }
854 auto nb_base_functions = data.getN().size2();
855
857 CHKERR pythonPtr->blockModes(iD, OP::getCoordsAtGaussPts(),
858 modeCentroids[blockCounter],
859 modeBboxes[blockCounter], blockModes);
860
861 auto nb_integration_pts = getGaussPts().size2();
862 if (blockModes.size2() != 3 * nb_integration_pts) {
863 MOFEM_LOG("WORLD", Sev::error)
864 << "Number of modes does not match number of integration points: "
865 << blockModes.size2() << "!=" << 3 * nb_integration_pts;
866 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "modes/integration points");
867 }
868
869 VectorDouble nf(nb_rows);
870
871 int nb_modes = blockModes.size1();
872 for (auto mode = 0; mode != nb_modes; ++mode) {
873 nf.clear();
874 // get mode
875 auto t_mode = getFTensor1FromPtr<3>(&blockModes(mode, 0));
876 // get element volume
877 const double vol = OP::getMeasure();
878 // get integration weights
879 auto t_w = OP::getFTensor0IntegrationWeight();
880 // get base function gradient on rows
881 auto t_base = data.getFTensor0N();
882 // loop over integration points
883 for (int gg = 0; gg != nb_integration_pts; gg++) {
884
885 // take into account Jacobian
886 const double alpha = t_w * vol;
887 // loop over rows base functions
888 auto t_nf = getFTensor1FromPtr<SPACE_DIM>(nf.data().data());
889 int rr = 0;
890 for (; rr != nb_rows / SPACE_DIM; ++rr) {
891 t_nf(i) += alpha * t_base * t_mode(i);
892 ++t_base;
893 ++t_nf;
894 }
895 for (; rr < nb_base_functions; ++rr)
896 ++t_base;
897 ++t_w; // move to another integration weight
898 ++t_mode; // move to another mode
899 }
900 Vec vec = modeVecs[modeCounter + mode];
901 auto size = data.getIndices().size();
902 auto *indices = data.getIndices().data().data();
903 auto *nf_data = nf.data().data();
904 CHKERR VecSetValues(vec, size, indices, nf_data, ADD_VALUES);
905 }
906
908 }
909
910 private:
911 boost::shared_ptr<ObjectiveFunctionData> pythonPtr;
912 MatrixDouble blockModes;
913 std::vector<std::array<double, 3>> modeCentroids;
914 std::vector<std::array<double, 6>> modeBboxes;
915 int iD;
916 std::vector<SmartPetscObj<Vec>> modeVecs;
917 int blockCounter;
918 int modeCounter;
919 };
920
921 auto solve_bdy = [&]() {
923
924 auto fe_lhs = get_lhs_fe();
925 auto fe_rhs = get_rhs_fe();
926 int block_counter = 0;
927 int mode_counter = 0;
928 for (auto &bc : bcs) {
929 auto id = bc->getMeshsetId();
930 Range ents;
931 CHKERR mField.get_moab().get_entities_by_handle(bc->getMeshset(), ents,
932 true);
933 auto range = boost::make_shared<Range>(ents);
934 auto &pip_rhs = fe_rhs->getOpPtrVector();
935 pip_rhs.push_back(new OpMode("ADJOINT_FIELD", pythonPtr, id, modeVecs,
936 modeCentroids, modeBBoxes, block_counter,
937 mode_counter, range));
938 CHKERR DMoFEMLoopFiniteElements(subset_dm_bdy, "ADJOINT_BOUNDARY_FE",
939 fe_rhs);
940 pip_rhs.pop_back();
941 int nb_modes;
942 CHKERR pythonPtr->numberOfModes(id, nb_modes);
943 ++block_counter;
944 mode_counter += nb_modes;
945 MOFEM_LOG("WORLD", Sev::inform)
946 << "Setting mode block block: " << bc->getName()
947 << " with ID: " << bc->getMeshsetId()
948 << " total modes: " << mode_counter;
949 }
950
951 for (auto &v : modeVecs) {
952 CHKERR VecAssemblyBegin(v);
953 CHKERR VecAssemblyEnd(v);
954 CHKERR VecGhostUpdateBegin(v, ADD_VALUES, SCATTER_REVERSE);
955 CHKERR VecGhostUpdateEnd(v, ADD_VALUES, SCATTER_REVERSE);
956 }
957
958 auto M = createDMMatrix(subset_dm_bdy);
959 fe_lhs->B = M;
960 CHKERR DMoFEMLoopFiniteElements(subset_dm_bdy, "ADJOINT_BOUNDARY_FE",
961 fe_lhs);
962 CHKERR MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY);
963 CHKERR MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY);
964
965 auto solver = createKSP(mField.get_comm());
966 CHKERR KSPSetOperators(solver, M, M);
967 CHKERR KSPSetFromOptions(solver);
968 CHKERR KSPSetUp(solver);
969 auto v = createDMVector(subset_dm_bdy);
970 for (auto &f : modeVecs) {
971 CHKERR KSPSolve(solver, f, v);
972 CHKERR VecSwap(f, v);
973 }
974
975 for (auto &v : modeVecs) {
976 CHKERR VecGhostUpdateBegin(v, INSERT_VALUES, SCATTER_FORWARD);
977 CHKERR VecGhostUpdateEnd(v, INSERT_VALUES, SCATTER_FORWARD);
978 }
979
981 };
982
983 CHKERR solve_bdy();
984
985 auto get_elastic_fe_lhs = [&]() {
986 auto fe = boost::make_shared<DomainEle>(mField);
987 fe->getRuleHook = [](int, int, int p_data) {
988 return 2 * p_data + p_data - 1;
989 };
990 auto &pip = fe->getOpPtrVector();
992 "GEOMETRY");
993 CHKERR HookeOps::opFactoryDomainLhs<SPACE_DIM, A, I, DomainEleOp>(
994 mField, pip, "ADJOINT_FIELD", "MAT_ADJOINT", Sev::noisy);
995 return fe;
996 };
997
998 auto get_elastic_fe_rhs = [&]() {
999 auto fe = boost::make_shared<DomainEle>(mField);
1000 fe->getRuleHook = [](int, int, int p_data) {
1001 return 2 * p_data + p_data - 1;
1002 };
1003 auto &pip = fe->getOpPtrVector();
1005 "GEOMETRY");
1006 CHKERR HookeOps::opFactoryDomainRhs<SPACE_DIM, A, I, DomainEleOp>(
1007 mField, pip, "ADJOINT_FIELD", "MAT_ADJOINT", Sev::noisy);
1008 return fe;
1009 };
1010
1011 auto adjoint_gradient_postprocess = [&](auto mode) {
1013 auto post_proc_mesh = boost::make_shared<moab::Core>();
1014 auto post_proc_begin =
1015 boost::make_shared<PostProcBrokenMeshInMoabBaseBegin>(mField,
1016 post_proc_mesh);
1017 auto post_proc_end = boost::make_shared<PostProcBrokenMeshInMoabBaseEnd>(
1018 mField, post_proc_mesh);
1019
1020 auto geom_vec = boost::make_shared<MatrixDouble>();
1021
1022 auto post_proc_fe =
1023 boost::make_shared<PostProcEleDomain>(mField, post_proc_mesh);
1025 post_proc_fe->getOpPtrVector(), {H1}, "GEOMETRY");
1026 post_proc_fe->getOpPtrVector().push_back(
1027 new OpCalculateVectorFieldValues<SPACE_DIM>("ADJOINT_FIELD", geom_vec,
1028 modeVecs[mode]));
1029
1031
1032 post_proc_fe->getOpPtrVector().push_back(
1033
1034 new OpPPMap(
1035
1036 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
1037
1038 {},
1039
1040 {{"MODE", geom_vec}},
1041
1042 {},
1043
1044 {}
1045
1046 )
1047
1048 );
1049
1051 post_proc_begin->getFEMethod());
1052 CHKERR DMoFEMLoopFiniteElements(adjointDM, "ADJOINT_DOMAIN_FE",
1053 post_proc_fe);
1055 post_proc_begin->getFEMethod());
1056
1057 CHKERR post_proc_end->writeFile("mode_" + std::to_string(mode) + ".h5m");
1058
1060 };
1061
1062 auto solve_domain = [&]() {
1064 auto fe_lhs = get_elastic_fe_lhs();
1065 auto fe_rhs = get_elastic_fe_rhs();
1066 auto v = createDMVector(subset_dm_domain);
1067 auto F = vectorDuplicate(v);
1068 fe_rhs->f = F;
1069
1070 auto M = createDMMatrix(subset_dm_domain);
1071 fe_lhs->B = M;
1072 CHKERR DMoFEMLoopFiniteElements(subset_dm_domain, "ADJOINT_DOMAIN_FE",
1073 fe_lhs);
1074 CHKERR MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY);
1075 CHKERR MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY);
1076
1077 auto solver = createKSP(mField.get_comm());
1078 CHKERR KSPSetOperators(solver, M, M);
1079 CHKERR KSPSetFromOptions(solver);
1080 CHKERR KSPSetUp(solver);
1081
1082 int mode_counter = 0;
1083 for (auto &f : modeVecs) {
1084 CHKERR mField.getInterface<FieldBlas>()->setField(0, "ADJOINT_FIELD");
1085 CHKERR DMoFEMMeshToLocalVector(subset_dm_bdy, f, INSERT_VALUES,
1086 SCATTER_REVERSE);
1087 CHKERR VecZeroEntries(F);
1088 CHKERR DMoFEMLoopFiniteElements(subset_dm_domain, "ADJOINT_DOMAIN_FE",
1089 fe_rhs);
1090 CHKERR VecAssemblyBegin(F);
1091 CHKERR VecAssemblyEnd(F);
1092 CHKERR VecGhostUpdateBegin(F, ADD_VALUES, SCATTER_REVERSE);
1093 CHKERR VecGhostUpdateEnd(F, ADD_VALUES, SCATTER_REVERSE);
1094 CHKERR KSPSolve(solver, F, v);
1095 CHKERR VecGhostUpdateBegin(v, INSERT_VALUES, SCATTER_FORWARD);
1096 CHKERR VecGhostUpdateEnd(v, INSERT_VALUES, SCATTER_FORWARD);
1097 CHKERR DMoFEMMeshToLocalVector(subset_dm_domain, v, INSERT_VALUES,
1098 SCATTER_REVERSE);
1099 auto m = createDMVector(adjointDM);
1100 CHKERR DMoFEMMeshToLocalVector(adjointDM, m, INSERT_VALUES,
1101 SCATTER_FORWARD);
1102 f = m;
1103 ++mode_counter;
1104 }
1106 };
1107
1108 CHKERR solve_domain();
1109
1110 for (int i = 0; i < modeVecs.size(); ++i) {
1111 CHKERR adjoint_gradient_postprocess(i);
1112 }
1113
1115}
1116//! [Adjoint modes]
1117
1118//! [Push operators to pipeline]
1121 auto pip = mField.getInterface<PipelineManager>();
1122
1123 //! [Integration rule]
1124 auto integration_rule = [](int, int, int approx_order) {
1125 return 2 * approx_order + 1;
1126 };
1127
1128 auto integration_rule_bc = [](int, int, int approx_order) {
1129 return 2 * approx_order + 1;
1130 };
1131
1133 CHKERR pip->setDomainLhsIntegrationRule(integration_rule);
1134 CHKERR pip->setBoundaryRhsIntegrationRule(integration_rule_bc);
1135 CHKERR pip->setBoundaryLhsIntegrationRule(integration_rule_bc);
1136 //! [Integration rule]
1137
1139 pip->getOpDomainLhsPipeline(), {H1}, "GEOMETRY");
1141 pip->getOpDomainRhsPipeline(), {H1}, "GEOMETRY");
1143 pip->getOpBoundaryRhsPipeline(), {NOSPACE}, "GEOMETRY");
1145 pip->getOpBoundaryLhsPipeline(), {NOSPACE}, "GEOMETRY");
1146
1147 //! [Push domain stiffness matrix to pipeline]
1148 // Add LHS operator for elasticity (stiffness matrix)
1149 CHKERR HookeOps::opFactoryDomainLhs<SPACE_DIM, A, I, DomainEleOp>(
1150 mField, pip->getOpDomainLhsPipeline(), "U", "MAT_ELASTIC", Sev::noisy);
1151 //! [Push domain stiffness matrix to pipeline]
1152
1153 // Add RHS operator for internal forces
1154 CHKERR HookeOps::opFactoryDomainRhs<SPACE_DIM, A, I, DomainEleOp>(
1155 mField, pip->getOpDomainRhsPipeline(), "U", "MAT_ELASTIC", Sev::noisy);
1156
1157 //! [Push Body forces]
1159 pip->getOpDomainRhsPipeline(), mField, "U", Sev::inform);
1160 //! [Push Body forces]
1161
1162 //! [Push natural boundary conditions]
1163 // Add force boundary condition
1165 pip->getOpBoundaryRhsPipeline(), mField, "U", -1, Sev::inform);
1166 // Add case for mix type of BCs
1168 pip->getOpBoundaryLhsPipeline(), mField, "U", Sev::noisy);
1169 //! [Push natural boundary conditions]
1171}
1172//! [Push operators to pipeline]
1173
1174//! [Solve]
1177 auto simple = mField.getInterface<Simple>();
1178 auto dm = simple->getDM();
1179 auto f = createDMVector(dm);
1180 auto d = vectorDuplicate(f);
1181 CHKERR VecZeroEntries(d);
1182 CHKERR DMoFEMMeshToLocalVector(dm, d, INSERT_VALUES, SCATTER_REVERSE);
1183
1184 auto set_essential_bc = [&]() {
1186 // This is low level pushing finite elements (pipelines) to solver
1187
1188 auto ksp_ctx_ptr = getDMKspCtx(dm);
1189 auto pre_proc_rhs = boost::make_shared<FEMethod>();
1190 auto post_proc_rhs = boost::make_shared<FEMethod>();
1191 auto post_proc_lhs = boost::make_shared<FEMethod>();
1192
1193 auto get_pre_proc_hook = [&]() {
1195 {});
1196 };
1197 pre_proc_rhs->preProcessHook = get_pre_proc_hook();
1198
1199 auto get_post_proc_hook_rhs = [this, post_proc_rhs]() {
1201
1203 post_proc_rhs, 1.)();
1205 };
1206
1207 auto get_post_proc_hook_lhs = [this, post_proc_lhs]() {
1209
1211 post_proc_lhs, 1.)();
1213 };
1214
1215 post_proc_rhs->postProcessHook = get_post_proc_hook_rhs;
1216 post_proc_lhs->postProcessHook = get_post_proc_hook_lhs;
1217
1218 ksp_ctx_ptr->getPreProcComputeRhs().push_front(pre_proc_rhs);
1219 ksp_ctx_ptr->getPostProcComputeRhs().push_back(post_proc_rhs);
1220 ksp_ctx_ptr->getPostProcSetOperators().push_back(post_proc_lhs);
1222 };
1223
1224 auto setup_and_solve = [&](auto solver) {
1226 BOOST_LOG_SCOPED_THREAD_ATTR("Timeline", attrs::timer());
1227 MOFEM_LOG("TIMER", Sev::noisy) << "KSPSetUp";
1228 CHKERR KSPSetUp(solver);
1229 MOFEM_LOG("TIMER", Sev::noisy) << "KSPSetUp <= Done";
1230 MOFEM_LOG("TIMER", Sev::noisy) << "KSPSolve";
1231 CHKERR KSPSolve(solver, f, d);
1232 MOFEM_LOG("TIMER", Sev::noisy) << "KSPSolve <= Done";
1234 };
1235
1236 MOFEM_LOG_CHANNEL("TIMER");
1237 MOFEM_LOG_TAG("TIMER", "timer");
1238
1239 CHKERR set_essential_bc();
1240 CHKERR setup_and_solve(kspElastic);
1241
1242 CHKERR VecGhostUpdateBegin(d, INSERT_VALUES, SCATTER_FORWARD);
1243 CHKERR VecGhostUpdateEnd(d, INSERT_VALUES, SCATTER_FORWARD);
1244 CHKERR DMoFEMMeshToLocalVector(dm, d, INSERT_VALUES, SCATTER_REVERSE);
1245
1246 auto evaluate_field_at_the_point = [&]() {
1248
1249 int coords_dim = 3;
1250 std::array<double, 3> field_eval_coords{0.0, 0.0, 0.0};
1251 PetscBool do_eval_field = PETSC_FALSE;
1252 CHKERR PetscOptionsGetRealArray(NULL, NULL, "-field_eval_coords",
1253 field_eval_coords.data(), &coords_dim,
1254 &do_eval_field);
1255
1256 if (do_eval_field) {
1257
1258 vectorFieldPtr = boost::make_shared<MatrixDouble>();
1259 auto field_eval_data =
1260 mField.getInterface<FieldEvaluatorInterface>()->getData<DomainEle>();
1261
1263 ->buildTree<SPACE_DIM>(field_eval_data, simple->getDomainFEName());
1264
1265 field_eval_data->setEvalPoints(field_eval_coords.data(), 1);
1266 auto no_rule = [](int, int, int) { return -1; };
1267 auto field_eval_fe_ptr = field_eval_data->feMethodPtr.lock();
1268 field_eval_fe_ptr->getRuleHook = no_rule;
1269
1270 field_eval_fe_ptr->getOpPtrVector().push_back(
1272
1274 ->evalFEAtThePoint<SPACE_DIM>(
1275 field_eval_coords.data(), 1e-12, simple->getProblemName(),
1276 simple->getDomainFEName(), field_eval_data,
1278 QUIET);
1279
1280 if (vectorFieldPtr->size1()) {
1281 auto t_disp = getFTensor1FromMat<SPACE_DIM>(*vectorFieldPtr);
1282 if constexpr (SPACE_DIM == 2)
1283 MOFEM_LOG("FieldEvaluator", Sev::inform)
1284 << "U_X: " << t_disp(0) << " U_Y: " << t_disp(1);
1285 else
1286 MOFEM_LOG("FieldEvaluator", Sev::inform)
1287 << "U_X: " << t_disp(0) << " U_Y: " << t_disp(1)
1288 << " U_Z: " << t_disp(2);
1289 }
1290
1292 }
1294 };
1295
1296 CHKERR evaluate_field_at_the_point();
1297
1299}
1300//! [Solve]
1301
1302//! [Postprocess results]
1304 SmartPetscObj<Vec> adjoint_vector) {
1306 auto simple = mField.getInterface<Simple>();
1309 mField, simple->getDM(), simple->getDomainFEName(),
1310 "out_elastic_" + std::to_string(iter) + ".h5m", adjoint_vector,
1311 "ADJOINT", Sev::noisy);
1313}
1314//! [Postprocess results]
1315
1316//! [calculateGradient]
1317
1319
1320
1321
1322template <int SPACE_DIM>
1324 DomainBaseOp> : public DomainBaseOp {
1325
1327
1329 boost::shared_ptr<HookeOps::CommonData> comm_ptr,
1330 boost::shared_ptr<MatrixDouble> jac,
1331 boost::shared_ptr<MatrixDouble> diff_jac,
1332 boost::shared_ptr<VectorDouble> cof_vals)
1333 : OP(field_name, field_name, DomainEleOp::OPROW), commPtr(comm_ptr),
1334 jac(jac), diffJac(diff_jac), cofVals(cof_vals) {}
1335
1336protected:
1337 boost::shared_ptr<HookeOps::CommonData> commPtr;
1338 boost::shared_ptr<MatrixDouble> jac;
1339 boost::shared_ptr<MatrixDouble> diffJac;
1340 boost::shared_ptr<VectorDouble> cofVals;
1342};
1343
1344template <int DIM> inline auto diff_symmetrize(FTensor::Number<DIM>) {
1345
1346 FTensor::Index<'i', DIM> i;
1347 FTensor::Index<'j', DIM> j;
1348 FTensor::Index<'k', DIM> k;
1349 FTensor::Index<'l', DIM> l;
1350
1352
1353 t_diff(i, j, k, l) = 0;
1354 t_diff(0, 0, 0, 0) = 1;
1355 t_diff(1, 1, 1, 1) = 1;
1356
1357 t_diff(1, 0, 1, 0) = 0.5;
1358 t_diff(1, 0, 0, 1) = 0.5;
1359
1360 t_diff(0, 1, 0, 1) = 0.5;
1361 t_diff(0, 1, 1, 0) = 0.5;
1362
1363 if constexpr (DIM == 3) {
1364 t_diff(2, 2, 2, 2) = 1;
1365
1366 t_diff(2, 0, 2, 0) = 0.5;
1367 t_diff(2, 0, 0, 2) = 0.5;
1368 t_diff(0, 2, 0, 2) = 0.5;
1369 t_diff(0, 2, 2, 0) = 0.5;
1370
1371 t_diff(2, 1, 2, 1) = 0.5;
1372 t_diff(2, 1, 1, 2) = 0.5;
1373 t_diff(1, 2, 1, 2) = 0.5;
1374 t_diff(1, 2, 2, 1) = 0.5;
1375 }
1376
1377 return t_diff;
1378};
1379
1380template <int SPACE_DIM>
1391
1392 auto t_diff_symm = diff_symmetrize(FTensor::Number<SPACE_DIM>());
1393
1394 // get element volume
1395 const double vol = OP::getMeasure();
1396 // get integration weights
1397 auto t_w = OP::getFTensor0IntegrationWeight();
1398 // get Jacobian values
1399 auto t_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(jac));
1400 // get diff Jacobian values
1401 auto t_diff_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(diffJac));
1402 // get cofactor values
1403 auto t_cof = getFTensor0FromVec(*(cofVals));
1404 // get base function gradient on rows
1405 auto t_row_grad = row_data.getFTensor1DiffN<SPACE_DIM>();
1406 // get fradient of the field
1407 auto t_grad_u =
1408 getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(commPtr->matGradPtr));
1409 // get field gradient values
1410 auto t_cauchy_stress =
1411 getFTensor2SymmetricFromMat<SPACE_DIM>(*(commPtr->getMatCauchyStress()));
1412 // material stiffness tensor
1413 auto t_D =
1414 getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*(commPtr->matDPtr));
1415 // loop over integration points
1416 for (int gg = 0; gg != OP::nbIntegrationPts; gg++) {
1417 // take into account Jacobian
1418 const double alpha = t_w * vol;
1419
1420 auto t_det = determinantTensor(t_jac);
1422 CHKERR invertTensor(t_jac, t_det, t_inv_jac);
1423
1424 // Calculate the variation of the gradient due to geometry change
1426 t_diff_inv_jac(i, j) =
1427 -(t_inv_jac(i, I) * t_diff_jac(I, J)) * t_inv_jac(J, j);
1429 t_diff_grad(i, j) = t_grad_u(i, k) * t_diff_inv_jac(k, j);
1430
1431 // Calculate the variation of the strain tensor
1433 t_diff_strain(i, j) = t_diff_symm(i, j, k, l) * t_diff_grad(k, l);
1434
1435 // Calculate the variation of the stress tensor
1437 t_diff_stress(i, j) = t_D(i, j, k, l) * t_diff_strain(k, l);
1438
1439 // get rhs vector
1440 auto t_nf = OP::template getNf<SPACE_DIM>();
1441 // loop over rows base functions
1442 int rr = 0;
1443 for (; rr != OP::nbRows / SPACE_DIM; rr++) {
1444
1446 t_diff_row_grad(k) = t_row_grad(j) * t_diff_inv_jac(j, k);
1447
1448 // Variation due to change in geometry (diff base)
1449 t_nf(j) += alpha * t_diff_row_grad(i) * t_cauchy_stress(i, j);
1450
1451 // Variation due to change in domain (cofactor)
1452 t_nf(j) += (alpha * t_cof) * t_row_grad(i) * t_cauchy_stress(i, j);
1453
1454 // Variation due to change in stress (diff stress)
1455 t_nf(j) += alpha * t_row_grad(i) * t_diff_stress(i, j);
1456
1457 ++t_row_grad;
1458 ++t_nf;
1459 }
1460 for (; rr < OP::nbRowBaseFunctions; ++rr) {
1461 ++t_row_grad;
1462 }
1463
1464 ++t_grad_u;
1465 ++t_cauchy_stress;
1466 ++t_jac;
1467 ++t_diff_jac;
1468 ++t_cof;
1469 ++t_w; // move to another integration weight
1470 }
1472}
1473
1474//J = 1/2 * sigma : epsilon
1475// calculate dJ/du = (dJ/dsigma * dsigma/depsilon + dJ/depsilon) depsilon / dgrad_u * dgrad_u / du
1477{
1480 boost::shared_ptr<ObjectiveFunctionData> python_ptr,
1481 boost::shared_ptr<HookeOps::CommonData> comm_ptr,
1482 boost::shared_ptr<MatrixDouble>u_ptr)
1483 : OP(field_name, field_name, DomainEleOp::OPROW), pythonPtr(python_ptr), commPtr(comm_ptr),
1484 uPtr(u_ptr){}
1485
1487 {
1493
1494 constexpr auto symm_size = (SPACE_DIM * (SPACE_DIM + 1)) / 2;
1495 auto nb_gauss_pts = getGaussPts().size2();
1496
1497 auto objective_dstress = boost::make_shared<MatrixDouble>(symm_size, nb_gauss_pts);
1498 auto objective_dstrain = boost::make_shared<MatrixDouble>(symm_size, nb_gauss_pts);
1499 auto objective_du = boost::make_shared<MatrixDouble>(SPACE_DIM, nb_gauss_pts);
1500
1501
1502
1503 auto evaluate_python = [&]()
1504 {
1506 auto &coords = OP::getCoordsAtGaussPts();
1507 CHKERR pythonPtr->evalInteriorObjectiveGradientStress(
1508 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1509 objective_dstress);
1510 CHKERR pythonPtr->evalInteriorObjectiveGradientStrain(
1511 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1512 objective_dstrain);
1513 CHKERR pythonPtr->evalInteriorObjectiveGradientU(
1514 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1515 objective_du);
1516
1517 auto vol = OP::getMeasure();
1518 auto t_w = OP::getFTensor0IntegrationWeight();
1519
1520 auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*(commPtr->matDPtr));
1521 auto t_row_grad = data.getFTensor1DiffN<SPACE_DIM>();
1522 auto t_row_base = data.getFTensor0N();
1523
1524 auto t_obj_dstress = getFTensor2SymmetricFromMat<SPACE_DIM>(*objective_dstress);
1525 auto t_obj_dstrain = getFTensor2SymmetricFromMat<SPACE_DIM>(*objective_dstrain);
1526 auto t_obj_du = getFTensor1FromMat<SPACE_DIM>(*objective_du);
1527
1528
1529 for(int gg = 0; gg != nb_gauss_pts; ++gg)
1530 {
1531 const double alpha = t_w * vol;
1533 t_adjoint_stress(i,j) = t_D(i, j, k, l) * t_obj_dstress(k, l) + t_obj_dstrain(i, j);
1534
1535 auto t_nf = OP::template getNf<SPACE_DIM>();
1536 int rr = 0;
1537 for (; rr != OP::nbRows / SPACE_DIM; rr++)
1538 {
1539 t_nf(j) += alpha * t_row_grad(i) * t_adjoint_stress(i, j);
1540 t_nf(j) += alpha * t_row_base * t_obj_du(j);
1541
1542 ++t_row_grad;
1543 ++t_row_base;
1544 ++t_nf;
1545 }
1546
1547 for (; rr < OP::nbRowBaseFunctions; ++rr)
1548 {
1549 ++t_row_grad;
1550 ++t_row_base;
1551 }
1552 ++t_obj_dstrain;
1553 ++t_obj_dstress;
1554 ++t_obj_du;
1555 ++t_w;
1556
1557 }
1559
1560
1561 };
1562 CHKERR evaluate_python();
1564 }
1565
1566
1567 private:
1568 boost::shared_ptr<ObjectiveFunctionData> pythonPtr;
1569 boost::shared_ptr<HookeOps::CommonData> commPtr;
1570 boost::shared_ptr<MatrixDouble> uPtr;
1571
1572};
1575
1576 OpAdJointObjective(boost::shared_ptr<ObjectiveFunctionData> python_ptr,
1577 boost::shared_ptr<HookeOps::CommonData> comm_ptr,
1578 boost::shared_ptr<MatrixDouble> jac_ptr,
1579 boost::shared_ptr<MatrixDouble> diff_jac,
1580 boost::shared_ptr<VectorDouble> cof_vals,
1581 boost::shared_ptr<MatrixDouble> d_grad_ptr,
1582 boost::shared_ptr<MatrixDouble> d_u_ptr,
1583 boost::shared_ptr<MatrixDouble> u_ptr,
1584 boost::shared_ptr<double> glob_objective_ptr,
1585 boost::shared_ptr<double> glob_objective_grad_ptr)
1586 : OP(NOSPACE, OP::OPSPACE), pythonPtr(python_ptr), commPtr(comm_ptr),
1587 jacPtr(jac_ptr), diffJacPtr(diff_jac), cofVals(cof_vals),
1588 dGradPtr(d_grad_ptr), dUPtr(d_u_ptr), uPtr(u_ptr),
1589 globObjectivePtr(glob_objective_ptr),
1590 globObjectiveGradPtr(glob_objective_grad_ptr) {}
1591
1592 /**
1593 * @brief Compute objective function contributions at element level
1594 *
1595 * Evaluates Python objective function with current displacement and stress
1596 * state, and accumulates global objective value and gradients.
1597 */
1598 MoFEMErrorCode doWork(int side, EntityType type,
1601
1602 // Define tensor indices for calculations
1607
1610
1611 constexpr auto symm_size = (SPACE_DIM * (SPACE_DIM + 1)) / 2; // size of symmetric tensor in Voigt notation
1612
1613 auto t_diff_symm = diff_symmetrize(FTensor::Number<SPACE_DIM>()); // fourth order tensor for symmetrization to convert from gradient to strain
1614
1615 auto nb_gauss_pts = getGaussPts().size2(); // number of gauss points in the element
1616 auto objective_ptr = boost::make_shared<VectorDouble>(nb_gauss_pts); // objective function values at gauss points
1617 auto objective_dstress =
1618 boost::make_shared<MatrixDouble>(symm_size, nb_gauss_pts); // objective function gradient w.r.t. stress at gauss points
1619 auto objective_dstrain =
1620 boost::make_shared<MatrixDouble>(symm_size, nb_gauss_pts); // objective function gradient w.r.t. strain at gauss points
1621 auto objective_du = boost::make_shared<MatrixDouble>(SPACE_DIM, nb_gauss_pts);
1622// We have dJ/dsigma and dJ/depsilon at gauss points, we need to convert it to dJ/du
1623
1624 auto evaluate_python = [&]() {
1626 auto &coords = OP::getCoordsAtGaussPts();
1627 CHKERR pythonPtr->evalInteriorObjectiveFunction(
1628 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1629 objective_ptr);
1630 CHKERR pythonPtr->evalInteriorObjectiveGradientStress(
1631 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1632 objective_dstress);
1633 CHKERR pythonPtr->evalInteriorObjectiveGradientStrain(
1634 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1635 objective_dstrain);
1636 CHKERR pythonPtr->evalInteriorObjectiveGradientU(
1637 coords, uPtr, commPtr->getMatCauchyStress(), commPtr->getMatStrain(),
1638 objective_du);
1639
1640 auto t_grad_u =
1641 getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(commPtr->matGradPtr));
1642 auto t_D =
1643 getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*(commPtr->matDPtr));
1644 auto t_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(jacPtr));
1645 auto t_diff_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(diffJacPtr));
1646 auto t_cof = getFTensor0FromVec(*(cofVals));
1647 auto t_d_grad = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(dGradPtr));
1648
1649 auto t_obj = getFTensor0FromVec(*objective_ptr);
1650 auto t_obj_dstress =
1651 getFTensor2SymmetricFromMat<SPACE_DIM>(*objective_dstress);
1652 auto t_obj_dstrain =
1653 getFTensor2SymmetricFromMat<SPACE_DIM>(*objective_dstrain);
1654 auto t_obj_du = getFTensor1FromMat<SPACE_DIM>(*objective_du);
1655 auto t_d_u = getFTensor1FromMat<SPACE_DIM>(*dUPtr);
1656
1657 auto vol = OP::getMeasure();
1658 auto t_w = getFTensor0IntegrationWeight();
1659 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
1660
1661 auto t_det = determinantTensor(t_jac);
1663 CHKERR invertTensor(t_jac, t_det, t_inv_jac);
1664
1666 t_diff_inv_jac(i, j) =
1667 -(t_inv_jac(i, I) * t_diff_jac(I, J)) * t_inv_jac(J, j);
1669 t_diff_grad(i, j) = t_grad_u(i, k) * t_diff_inv_jac(k, j);
1670
1672 t_d_strain(i, j) = t_diff_symm(i, j, k, l) * (
1673
1674 t_d_grad(k, l)
1675
1676 +
1677
1678 t_diff_grad(k, l)
1679
1680 );
1681
1682 auto alpha = t_w * vol;
1683
1684 (*globObjectivePtr) += alpha * t_obj;
1685 (*globObjectiveGradPtr) +=
1686 alpha *
1687 (
1688
1689 t_obj_dstress(i, j) * (t_D(i, j, k, l) * t_d_strain(k, l))
1690
1691 +
1692
1693 t_obj_dstrain(i, j) * t_d_strain(i, j)
1694
1695 +
1696
1697 t_obj_du(i) * t_d_u(i)
1698
1699 +
1700
1701 t_obj * t_cof
1702
1703 );
1704
1705 ++t_w;
1706 ++t_jac;
1707 ++t_diff_jac;
1708 ++t_cof;
1709
1710 ++t_obj;
1711 ++t_obj_dstress;
1712 ++t_obj_dstrain;
1713 ++t_obj_du;
1714
1715 ++t_grad_u;
1716 ++t_d_grad;
1717 ++t_d_u;
1718 }
1720 };
1721
1722 CHKERR evaluate_python();
1723
1725 }
1726
1727private:
1728 boost::shared_ptr<ObjectiveFunctionData> pythonPtr;
1729 boost::shared_ptr<HookeOps::CommonData> commPtr;
1730 boost::shared_ptr<MatrixDouble> jacPtr;
1731 boost::shared_ptr<MatrixDouble> diffJacPtr;
1732 boost::shared_ptr<VectorDouble> cofVals;
1733 boost::shared_ptr<MatrixDouble> dGradPtr;
1734 boost::shared_ptr<MatrixDouble> dUPtr;
1735 boost::shared_ptr<MatrixDouble> uPtr;
1736
1737 boost::shared_ptr<double> globObjectivePtr;
1738 boost::shared_ptr<double> globObjectiveGradPtr;
1739};
1740
1741MoFEMErrorCode Example::calculateGradient(PetscReal *objective_function_value,
1742 Vec objective_function_gradient,
1743 Vec adjoint_vector) {
1744 MOFEM_LOG_CHANNEL("WORLD");
1746 auto simple = mField.getInterface<Simple>();
1747
1748 auto ents = get_range_from_block(mField, "OPTIMISE", SPACE_DIM - 1);
1749
1750 auto get_essential_fe = [this]() {
1751 auto post_proc_rhs = boost::make_shared<FEMethod>();
1752 auto get_post_proc_hook_rhs = [this, post_proc_rhs]() {
1754
1756 post_proc_rhs, 0)();
1758 };
1759 post_proc_rhs->postProcessHook = get_post_proc_hook_rhs;
1760 return post_proc_rhs;
1761 };
1762
1763 auto get_fd_direvative_fe = [&]() {
1764 auto fe = boost::make_shared<DomainEle>(mField);
1765 fe->getRuleHook = [](int, int, int p_data) {
1766 return 2 * p_data + p_data - 1;
1767 };
1768 auto &pip = fe->getOpPtrVector();
1770 // Add RHS operators for internal forces
1771 HookeOps::opFactoryDomainRhs<SPACE_DIM, A, I, DomainEleOp>(
1772 mField, pip, "U", "MAT_ELASTIC", Sev::noisy);
1773
1774 return fe;
1775 };
1776
1777 auto calulate_fd_residual = [&](auto eps, auto diff_vec, auto fd_vec) {
1779
1780 constexpr bool debug = false;
1781
1782 auto geom_norm = [](MoFEM::Interface &mField) {
1784 auto field_blas = mField.getInterface<FieldBlas>();
1785 double nrm2 = 0.0;
1786 auto norm2_field = [&](const double val) {
1787 nrm2 += val * val;
1788 return val;
1789 };
1790 CHKERR field_blas->fieldLambdaOnValues(norm2_field, "GEOMETRY");
1791 MPI_Allreduce(MPI_IN_PLACE, &nrm2, 1, MPI_DOUBLE, MPI_SUM,
1792 mField.get_comm());
1793 MOFEM_LOG("WORLD", Sev::inform) << "Geometry norm: " << sqrt(nrm2);
1795 };
1796
1797 if constexpr (debug)
1798 CHKERR geom_norm(mField);
1799
1800 auto initial_current_geometry = createDMVector(adjointDM);
1801 CHKERR mField.getInterface<VecManager>()->setOtherLocalGhostVector(
1802 "ADJOINT", "ADJOINT_FIELD", "GEOMETRY", RowColData::ROW,
1803 initial_current_geometry, INSERT_VALUES, SCATTER_FORWARD);
1804 CHKERR VecAssemblyBegin(initial_current_geometry);
1805 CHKERR VecAssemblyEnd(initial_current_geometry);
1806
1807 if constexpr (debug)
1808 CHKERR geom_norm(mField);
1809
1810 auto perturb_geometry = [&](auto eps, auto diff_vec) {
1812 auto current_geometry = vectorDuplicate(initial_current_geometry);
1813 CHKERR VecCopy(initial_current_geometry, current_geometry);
1814 CHKERR VecAXPY(current_geometry, eps, diff_vec);
1815 CHKERR mField.getInterface<VecManager>()->setOtherLocalGhostVector(
1816 "ADJOINT", "ADJOINT_FIELD", "GEOMETRY", RowColData::ROW,
1817 current_geometry, INSERT_VALUES, SCATTER_REVERSE);
1819 };
1820
1821 auto fe = get_fd_direvative_fe();
1822 auto fp = vectorDuplicate(diff_vec);
1823 auto fm = vectorDuplicate(diff_vec);
1824 auto calc_impl = [&](auto f, auto eps) { // is this finite difference!
1826 CHKERR VecZeroEntries(f);
1827 fe->f = f;
1828 CHKERR perturb_geometry(eps, diff_vec);
1830 simple->getDomainFEName(), fe);
1831 CHKERR VecAssemblyBegin(f);
1832 CHKERR VecAssemblyEnd(f);
1833 CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
1834 CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
1835 auto post_proc_rhs = get_essential_fe();
1836 post_proc_rhs->f = f;
1838 post_proc_rhs.get());
1840 };
1841 CHKERR calc_impl(fp, eps);
1842 CHKERR calc_impl(fm, -eps);
1843 CHKERR VecWAXPY(fd_vec, -1.0, fm, fp);
1844 CHKERR VecScale(fd_vec, 1.0 / (2.0 * eps));
1845
1846 CHKERR mField.getInterface<VecManager>()->setOtherLocalGhostVector(
1847 "ADJOINT", "ADJOINT_FIELD", "GEOMETRY", RowColData::ROW,
1848 initial_current_geometry, INSERT_VALUES, SCATTER_REVERSE);
1849
1850 if constexpr (debug)
1851 CHKERR geom_norm(mField);
1852
1854 };
1855 // here starts the preperation for the derivative
1856 auto get_direvative_fe = [&](auto diff_vec) {
1857 auto fe_adjoint = boost::make_shared<DomainEle>(mField);
1858 fe_adjoint->getRuleHook = [](int, int, int p_data) {
1859 return 2 * p_data + p_data - 1;
1860 };
1861 auto &pip = fe_adjoint->getOpPtrVector();
1862
1863 auto jac_ptr = boost::make_shared<MatrixDouble>();
1864 auto det_ptr = boost::make_shared<VectorDouble>();
1865 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
1866 auto diff_jac_ptr = boost::make_shared<MatrixDouble>();
1867 auto cof_ptr = boost::make_shared<VectorDouble>();
1868
1869 using OpCoFactor =
1870 AdJoint<DomainEleOp>::Integration<GAUSS>::OpGetCoFactor<SPACE_DIM>;
1871
1873
1875 "GEOMETRY", jac_ptr));
1877 "U", diff_jac_ptr, diff_vec));
1878
1879 pip.push_back(new OpCoFactor(jac_ptr, diff_jac_ptr, cof_ptr));
1880
1881 // Add RHS operators for internal forces
1882 auto common_ptr = HookeOps::commonDataFactory<SPACE_DIM, I, DomainEleOp>(
1883 mField, pip, "U", "MAT_ELASTIC", Sev::noisy);
1885 "U", common_ptr, jac_ptr, diff_jac_ptr, cof_ptr));
1886
1887 return fe_adjoint;
1888 };
1889
1890 auto get_objective_fe = [&](auto diff_vec, auto grad_vec,
1891 auto glob_objective_ptr,
1892 auto glob_objective_grad_ptr) {
1893 auto fe_adjoint = boost::make_shared<DomainEle>(mField);
1894 fe_adjoint->getRuleHook = [](int, int, int p_data) {
1895 return 2 * p_data + p_data - 1;
1896 };
1897 auto &pip = fe_adjoint->getOpPtrVector();
1899
1900 auto jac_ptr = boost::make_shared<MatrixDouble>();
1901 auto det_ptr = boost::make_shared<VectorDouble>();
1902 auto inv_jac_ptr = boost::make_shared<MatrixDouble>();
1903 auto diff_jac_ptr = boost::make_shared<MatrixDouble>();
1904 auto cof_ptr = boost::make_shared<VectorDouble>();
1905 auto d_grad_ptr = boost::make_shared<MatrixDouble>();
1906 auto d_u_ptr = boost::make_shared<MatrixDouble>();
1907 auto u_ptr = boost::make_shared<MatrixDouble>();
1908
1909 using OpCoFactor =
1910 AdJoint<DomainEleOp>::Integration<GAUSS>::OpGetCoFactor<SPACE_DIM>;
1912 "GEOMETRY", jac_ptr));
1914 "U", diff_jac_ptr,
1915 diff_vec)); // Note: that vector is stored on displacemnt vector, that
1916 // why is used here
1917 pip.push_back(new OpCoFactor(jac_ptr, diff_jac_ptr, cof_ptr));
1919 "U", d_grad_ptr, grad_vec));
1920 pip.push_back(
1921 new OpCalculateVectorFieldValues<SPACE_DIM>("U", d_u_ptr, grad_vec));
1922 pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
1923
1924 auto common_ptr = HookeOps::commonDataFactory<SPACE_DIM, I, DomainEleOp>(
1925 mField, pip, "U", "MAT_ELASTIC", Sev::noisy);
1926 pip.push_back(new OpAdJointObjective(
1927 pythonPtr, common_ptr, jac_ptr, diff_jac_ptr, cof_ptr, d_grad_ptr,
1928 d_u_ptr, u_ptr, glob_objective_ptr, glob_objective_grad_ptr));
1929
1930 return fe_adjoint;
1931 };
1932
1933 auto dm = simple->getDM();
1934 auto f = createDMVector(dm);
1935 auto d = vectorDuplicate(f);
1936 auto dm_diff_vec = vectorDuplicate(d);
1937 auto zero_diff_vec = vectorDuplicate(dm_diff_vec);
1938 auto zero_state_sensitivity = vectorDuplicate(d);
1939 CHKERR VecZeroEntries(zero_diff_vec);
1940 CHKERR VecZeroEntries(zero_state_sensitivity);
1941 CHKERR VecGhostUpdateBegin(zero_diff_vec, INSERT_VALUES, SCATTER_FORWARD);
1942 CHKERR VecGhostUpdateEnd(zero_diff_vec, INSERT_VALUES, SCATTER_FORWARD);
1943 CHKERR VecGhostUpdateBegin(zero_state_sensitivity, INSERT_VALUES,
1944 SCATTER_FORWARD);
1945 CHKERR VecGhostUpdateEnd(zero_state_sensitivity, INSERT_VALUES,
1946 SCATTER_FORWARD);
1947
1948 auto adjoint_fe = get_direvative_fe(dm_diff_vec);
1949 auto objective_ptr_direct = boost::make_shared<double>(0.0);
1950 auto objective_grad_ptr_direct = boost::make_shared<double>(0.0);
1951 auto objective_fe_direct = get_objective_fe(
1952 dm_diff_vec, d, objective_ptr_direct, objective_grad_ptr_direct);
1953 auto objective_ptr_explicit = boost::make_shared<double>(0.0);
1954 auto objective_grad_ptr_explicit = boost::make_shared<double>(0.0);
1955 auto objective_fe_explicit =
1956 get_objective_fe(dm_diff_vec, zero_state_sensitivity,
1957 objective_ptr_explicit, objective_grad_ptr_explicit);
1958 auto objective_ptr_value = boost::make_shared<double>(0.0);
1959 auto objective_grad_ptr_value = boost::make_shared<double>(0.0);
1960 auto objective_fe_value =
1961 get_objective_fe(zero_diff_vec, zero_state_sensitivity,
1962 objective_ptr_value, objective_grad_ptr_value);
1963
1964 auto set_variance_of_geometry =
1965 [&](auto mode, auto mod_vec) { // think of mod_vec as X(tau) = X + tau*v_h
1966 // take the mod_vec and write it into the adjoint DM then copy that
1967 // field to dm_diff_vec using the same layout as U
1969 CHKERR DMoFEMMeshToLocalVector(adjointDM, mod_vec, INSERT_VALUES,
1970 SCATTER_REVERSE);
1971 CHKERR mField.getInterface<VecManager>()->setOtherLocalGhostVector(
1972 simple->getProblemName(), "U", "ADJOINT_FIELD", RowColData::ROW,
1973 dm_diff_vec, INSERT_VALUES, SCATTER_FORWARD);
1974 CHKERR VecGhostUpdateBegin(dm_diff_vec, INSERT_VALUES, SCATTER_FORWARD);
1975 CHKERR VecGhostUpdateEnd(dm_diff_vec, INSERT_VALUES, SCATTER_FORWARD);
1977 };
1978
1979 auto calculate_variance_internal_forces = [&](auto mode, auto mod_vec) {
1981 CHKERR VecZeroEntries(f);
1982 CHKERR VecGhostUpdateBegin(f, INSERT_VALUES, SCATTER_FORWARD);
1983 CHKERR VecGhostUpdateEnd(f, INSERT_VALUES, SCATTER_FORWARD);
1984 adjoint_fe->f = f;
1985 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(), adjoint_fe);
1986 CHKERR VecAssemblyBegin(f);
1987 CHKERR VecAssemblyEnd(f);
1988 CHKERR VecGhostUpdateBegin(f, ADD_VALUES, SCATTER_REVERSE);
1989 CHKERR VecGhostUpdateEnd(f, ADD_VALUES, SCATTER_REVERSE);
1990 auto post_proc_rhs = get_essential_fe();
1991 post_proc_rhs->f = f;
1993 post_proc_rhs.get());
1994 CHKERR VecScale(f, -1.0);
1995
1996#ifndef NDEBUG
1997 constexpr bool debug = false;
1998 if constexpr (debug) {
1999 double norm0;
2000 CHKERR VecNorm(f, NORM_2, &norm0);
2001 auto fd_check = vectorDuplicate(f);
2002 double eps = 1e-5;
2003 CHKERR calulate_fd_residual(eps, dm_diff_vec, fd_check);
2004 double nrm;
2005 CHKERR VecAXPY(fd_check, -1.0, f);
2006 CHKERR VecNorm(fd_check, NORM_2, &nrm);
2007 MOFEM_LOG("WORLD", Sev::inform)
2008 << " FD check for internal forces [ " << mode << " ]: " << nrm
2009 << " / " << norm0 << " ( " << (nrm / norm0) << " )";
2010 }
2011#endif
2012
2014 };
2015
2016 auto calculate_variance_of_displacement = [&](auto mode, auto mod_vec) {
2018 CHKERR KSPSolve(kspElastic, f,
2019 d);
2020
2021 CHKERR VecGhostUpdateBegin(d, INSERT_VALUES, SCATTER_FORWARD);
2022 CHKERR VecGhostUpdateEnd(d, INSERT_VALUES, SCATTER_FORWARD);
2024 };
2025
2026 auto evaluate_objective_terms =
2027 [&](auto objective_fe, auto objective_ptr, auto objective_grad_ptr,
2028 double &objective_value, double &objective_gradient) {
2030 *objective_ptr = 0.0;
2031 *objective_grad_ptr = 0.0;
2032 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(), objective_fe);
2033 std::array<double, 2> array = {*objective_ptr, *objective_grad_ptr};
2034 MPI_Allreduce(MPI_IN_PLACE, array.data(), 2, MPI_DOUBLE, MPI_SUM,
2035 mField.get_comm());
2036 objective_value = array[0];
2037 objective_gradient = array[1];
2039 };
2040
2041 auto calculate_objective_value = [&]() {
2043 double objective_value = 0;
2044 double objective_gradient = 0;
2045 CHKERR evaluate_objective_terms(objective_fe_value, objective_ptr_value,
2046 objective_grad_ptr_value, objective_value,
2047 objective_gradient);
2048 *objective_function_value = objective_value;
2050 };
2051
2052 auto calculate_variance_of_objective_function_dJ_du = [&](Vec dJ_du) {
2054
2055 auto fe = boost::make_shared<DomainEle>(mField);
2056 fe->getRuleHook = [](int, int, int p_data) {
2057 return 2 * p_data + p_data - 1;
2058 };
2059 auto &pip = fe->getOpPtrVector();
2061
2062 auto u_ptr = boost::make_shared<MatrixDouble>();
2063 pip.push_back(new OpCalculateVectorFieldValues<SPACE_DIM>("U", u_ptr));
2064
2065 auto common_ptr = HookeOps::commonDataFactory<SPACE_DIM, I, DomainEleOp>(
2066 mField, pip, "U", "MAT_ELASTIC", Sev::noisy);
2067 pip.push_back(new OpStateSensitivity("U", pythonPtr, common_ptr, u_ptr));
2068
2069 CHKERR VecZeroEntries(dJ_du);
2070 fe->f = dJ_du;
2071 CHKERR DMoFEMLoopFiniteElements(dm, simple->getDomainFEName(), fe);
2072 CHKERR VecAssemblyBegin(dJ_du);
2073 CHKERR VecAssemblyEnd(dJ_du);
2074
2075 auto post_proc_rhs = get_essential_fe();
2076 post_proc_rhs->f = dJ_du;
2078 post_proc_rhs.get());
2079
2081
2082 };
2083
2084 auto calculate_adjoint_lambda = [&](auto lambda, auto dJ_du) {
2086
2087 MOFEM_LOG("WORLD", Sev::inform) << "Solving for adjoint variable lambda";
2088 CHKERR VecZeroEntries(lambda);
2089 CHKERR KSPSolveTranspose(kspElastic, dJ_du, lambda);
2090 CHKERR VecGhostUpdateBegin(lambda, INSERT_VALUES, SCATTER_FORWARD);
2091 CHKERR VecGhostUpdateEnd(lambda, INSERT_VALUES, SCATTER_FORWARD);
2093 };
2094
2095 CHKERR VecZeroEntries(objective_function_gradient);
2096 CHKERR VecZeroEntries(adjoint_vector);
2097
2098 CHKERR calculate_objective_value();
2099 MOFEM_LOG("WORLD", Sev::verbose)
2100 << "Objective function: " << *objective_function_value;
2101
2102 auto direct = [&]() {
2104 int mode = 0;
2105 for (auto mod_vec : modeVecs) {
2106 CHKERR set_variance_of_geometry(mode, mod_vec);
2107 CHKERR calculate_variance_internal_forces(mode, mod_vec);
2108 CHKERR calculate_variance_of_displacement(mode, mod_vec);
2109 double objective_value = 0;
2110 double objective_gradient = 0;
2111 CHKERR evaluate_objective_terms(objective_fe_direct, objective_ptr_direct,
2112 objective_grad_ptr_direct, objective_value,
2113 objective_gradient);
2114 CHKERR VecSetValue(objective_function_gradient, mode, objective_gradient,
2115 INSERT_VALUES);
2116 CHKERR VecAXPY(adjoint_vector, objective_gradient, dm_diff_vec);
2117 ++mode;
2118 }
2120 };
2121
2122 auto lambda = vectorDuplicate(f);
2123 auto dJ_du = vectorDuplicate(f);
2124
2125 auto adjoint = [&]() {
2127
2128 CHKERR calculate_variance_of_objective_function_dJ_du(dJ_du);
2129 CHKERR calculate_adjoint_lambda(lambda, dJ_du);
2130
2131 int mode = 0;
2132
2133 for (auto mod_vec : modeVecs) {
2134
2135 CHKERR set_variance_of_geometry(mode, mod_vec);
2136 CHKERR calculate_variance_internal_forces(mode, mod_vec);
2137 double objective_value = 0;
2138 double objective_gradient_explicit = 0;
2139 CHKERR evaluate_objective_terms(
2140 objective_fe_explicit, objective_ptr_explicit,
2141 objective_grad_ptr_explicit, objective_value,
2142 objective_gradient_explicit);
2143 double lambda_dot_residual_variation = 0;
2144 CHKERR VecDot(lambda, f, &lambda_dot_residual_variation);
2145 const double dJ_dp =
2146 objective_gradient_explicit + lambda_dot_residual_variation;
2147 CHKERR VecSetValue(objective_function_gradient, mode, dJ_dp,
2148 INSERT_VALUES);
2149 CHKERR VecAXPY(adjoint_vector, dJ_dp, dm_diff_vec);
2150 ++mode;
2151 }
2152 CHKERR VecAssemblyBegin(objective_function_gradient);
2153 CHKERR VecAssemblyEnd(objective_function_gradient);
2155 };
2156
2157 switch (derivative_type) {
2158 case DIRECT:
2159 MOFEM_LOG("WORLD", Sev::inform) << "Running Direct Sensitivity...";
2160 CHKERR direct();
2161 break;
2162
2163 case ADJOINT:
2164 MOFEM_LOG("WORLD", Sev::inform) << "Running Adjoint Sensitivity...";
2165 CHKERR adjoint();
2166 break;
2167
2168 default:
2169 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
2170 "Wrong sensitivity type selected");
2171 }
2172
2173 CHKERR VecAssemblyBegin(objective_function_gradient);
2174 CHKERR VecAssemblyEnd(objective_function_gradient);
2175
2176 CHKERR VecAssemblyBegin(adjoint_vector);
2177 CHKERR VecAssemblyEnd(adjoint_vector);
2178 CHKERR VecGhostUpdateBegin(adjoint_vector, INSERT_VALUES, SCATTER_FORWARD);
2179 CHKERR VecGhostUpdateEnd(adjoint_vector, INSERT_VALUES, SCATTER_FORWARD);
2180
2182}
2183//! [calculateGradient]
2184
2185static char help[] = "...\n\n";
2186
2187/**
2188 * @brief Main function for topology optimization tutorial using adjoint method
2189 *
2190 * This tutorial demonstrates structural topology optimization using:
2191 * - MoFEM finite element library for elasticity analysis
2192 * - Adjoint method for efficient gradient computation
2193 * - TAO optimization library for design optimization
2194 * - Python interface for flexible objective function definition
2195 *
2196 * Workflow:
2197 * 1. Initialize MoFEM, PETSc and Python environments
2198 * 2. Read mesh and setup finite element problem
2199 * 3. Define objective function via Python interface
2200 * 4. Compute topology optimization modes as design variables
2201 * 5. Run gradient-based optimization using adjoint sensitivities
2202 * 6. Post-process optimized design
2203 *
2204 * The adjoint method enables efficient gradient computation with cost
2205 * independent of the number of design variables, making it suitable
2206 * for large-scale topology optimization problems.
2207 *
2208 * Required input files:
2209 * - Mesh file (.h5m format from CUBIT)
2210 * - Parameter file (param_file.petsc)
2211 * - Objective function (objective_function.py)
2212 *
2213 * @param argc Command line argument count
2214 * @param argv Command line argument values
2215 * @return int Exit code (0 for success)
2216 */
2217int main(int argc, char *argv[]) {
2218
2219 // Initialize Python environment for objective function interface
2220 Py_Initialize();
2221 np::initialize();
2222
2223 // Initialize MoFEM/PETSc and MOAB data structures
2224 const char param_file[] = "param_file.petsc";
2225 MoFEM::Core::Initialize(&argc, &argv, param_file, help);
2226
2227 auto core_log = logging::core::get();
2228 core_log->add_sink(
2230
2231 core_log->add_sink(
2232 LogManager::createSink(LogManager::getStrmSync(), "FieldEvaluator"));
2233 LogManager::setLog("FieldEvaluator");
2234 MOFEM_LOG_TAG("FieldEvaluator", "field_eval");
2235
2236 try {
2237
2238 //! [Register MoFEM discrete manager in PETSc]
2239 DMType dm_name = "DMMOFEM";
2240 CHKERR DMRegister_MoFEM(dm_name);
2241 DMType dm_name_mg = "DMMOFEM_MG";
2243 //! [Register MoFEM discrete manager in PETSc
2244
2245 //! [Create MoAB]
2246 moab::Core mb_instance; ///< mesh database
2247 moab::Interface &moab = mb_instance; ///< mesh database interface
2248 //! [Create MoAB]
2249
2250 //! [Create MoFEM]
2251 MoFEM::Core core(moab); ///< finite element database
2252 MoFEM::Interface &m_field = core; ///< finite element database interface
2253 //! [Create MoFEM]
2254
2255 //! [Example]
2256
2257 Example ex(m_field);
2258 CHKERR ex.runProblem();
2259 //! [Example]
2260 }
2262
2264
2265 if (Py_FinalizeEx() < 0) {
2266 exit(120);
2267 }
2268}
2269
2271 const std::string block_name, int dim) {
2272 Range r;
2273
2274 auto mesh_mng = m_field.getInterface<MeshsetsManager>();
2275 auto bcs = mesh_mng->getCubitMeshsetPtr(
2276
2277 std::regex((boost::format("%s(.*)") % block_name).str())
2278
2279 );
2280
2281 for (auto bc : bcs) {
2282 Range faces;
2283 CHK_MOAB_THROW(bc->getMeshsetIdEntitiesByDimension(m_field.get_moab(), dim,
2284 faces, true),
2285 "get meshset ents");
2286 r.merge(faces);
2287 }
2288
2289 for (auto dd = dim - 1; dd >= 0; --dd) {
2290 if (dd >= 0) {
2291 Range ents;
2292 CHK_MOAB_THROW(m_field.get_moab().get_adjacencies(r, dd, false, ents,
2293 moab::Interface::UNION),
2294 "get adjs");
2295 r.merge(ents);
2296 } else {
2297 Range verts;
2298 CHK_MOAB_THROW(m_field.get_moab().get_connectivity(r, verts),
2299 "get verts");
2300 r.merge(verts);
2301 }
2303 m_field.getInterface<CommInterface>()->synchroniseEntities(r), "comm");
2304 }
2305
2306 return r;
2307};
2308
2309MoFEMErrorCode save_range(moab::Interface &moab, const std::string name,
2310 const Range r) {
2312 auto out_meshset = get_temp_meshset_ptr(moab);
2313 CHKERR moab.add_entities(*out_meshset, r);
2314 CHKERR moab.write_file(name.c_str(), "VTK", "", out_meshset->get_ptr(), 1);
2316};
#define MOFEM_LOG_SYNCHRONISE(comm)
Synchronise "SYNC" channel.
Interface for Python-based objective function evaluation in topology optimization.
#define FTENSOR_INDEX(DIM, I)
void simple(double P1[], double P2[], double P3[], double c[], const int N)
Definition acoustic.cpp:69
auto diff_symmetrize(FTensor::Number< DIM >)
Definition adjoint.cpp:1344
static char help[]
[calculateGradient]
Definition adjoint.cpp:2185
PetscBool is_plane_strain
Definition adjoint.cpp:42
constexpr int SPACE_DIM
[Define dimension]
Definition adjoint.cpp:28
SensitivityMethod derivative_type
Definition adjoint.cpp:95
constexpr double poisson_ratio
Poisson's ratio ν
Definition adjoint.cpp:38
constexpr int BASE_DIM
[Constants and material properties]
Definition adjoint.cpp:25
constexpr double shear_modulus_G
Shear modulus G = E/(2(1+ν))
Definition adjoint.cpp:40
constexpr IntegrationType I
Use Gauss quadrature for integration.
Definition adjoint.cpp:33
constexpr double bulk_modulus_K
Bulk modulus K = E/(3(1-2ν))
Definition adjoint.cpp:39
SensitivityMethod
Definition adjoint.cpp:90
@ DIRECT
Definition adjoint.cpp:91
@ ADJOINT
Definition adjoint.cpp:92
FormsIntegrators< DomainEleOp >::Assembly< A >::OpBase DomainBaseOp
[Postprocess results]
Definition adjoint.cpp:1318
constexpr AssemblyType A
[Define dimension]
Definition adjoint.cpp:32
Range get_range_from_block(MoFEM::Interface &m_field, const std::string block_name, int dim)
Definition adjoint.cpp:2270
constexpr double young_modulus
[Material properties for linear elasticity]
Definition adjoint.cpp:37
int main()
constexpr double a
constexpr int SPACE_DIM
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
@ QUIET
@ ROW
#define CATCH_ERRORS
Catch errors.
@ MF_EXIST
FieldApproximationBase
approximation base
Definition definitions.h:58
@ LASTBASE
Definition definitions.h:69
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition definitions.h:60
@ DEMKOWICZ_JACOBI_BASE
Definition definitions.h:66
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ H1
continuous field
Definition definitions.h:85
@ NOSPACE
Definition definitions.h:83
#define MYPCOMM_INDEX
default communicator number PCOMM
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
@ MOFEM_NOT_FOUND
Definition definitions.h:33
@ 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.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
PostProcEleByDim< SPACE_DIM >::PostProcEleDomain PostProcEleDomain
PostProcEleByDim< SPACE_DIM >::PostProcEleBdy PostProcEleBdy
@ F
auto integration_rule
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
Definition DMMoFEM.cpp:1113
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
Definition DMMoFEM.cpp:215
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
Definition DMMoFEM.cpp:488
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
Definition DMMoFEM.cpp:450
PetscErrorCode DMMoFEMCreateMoFEM(DM dm, MoFEM::Interface *m_field_ptr, const char problem_name[], const MoFEM::BitRefLevel bit_level, const MoFEM::BitRefLevel bit_mask=MoFEM::BitRefLevel().set())
Must be called by user to set MoFEM data structures.
Definition DMMoFEM.cpp:114
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition DMMoFEM.cpp:546
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
Definition DMMoFEM.cpp:238
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
Definition DMMoFEM.cpp:43
MoFEMErrorCode DMRegister_MGViaApproxOrders(const char sname[])
Register DM for Multi-Grid via approximation orders.
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:576
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1234
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
Definition DMMoFEM.cpp:280
auto createDMMatrix(DM dm)
Get smart matrix from DM.
Definition DMMoFEM.hpp:1191
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
Definition DMMoFEM.cpp:536
virtual MoFEMErrorCode add_ents_to_finite_element_by_dim(const EntityHandle entities, const int dim, const std::string name, const bool recursive=true)=0
add entities to finite element
virtual MoFEMErrorCode add_finite_element(const std::string &fe_name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
add finite element
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
virtual MoFEMErrorCode modify_finite_element_add_field_col(const std::string &fe_name, const std::string name_row)=0
set field col which finite element use
virtual MoFEMErrorCode modify_finite_element_add_field_row(const std::string &fe_name, const std::string name_row)=0
set field row which finite element use
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_field)=0
set finite element field data
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
virtual MoFEMErrorCode add_ents_to_field_by_dim(const Range &ents, const int dim, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
virtual MoFEMErrorCode set_field_order(const EntityHandle meshset, const EntityType type, const std::string &name, const ApproximationOrder order, int verb=DEFAULT_VERBOSITY)=0
Set order approximation of the entities in the field.
IntegrationType
Form integrator integration types.
AssemblyType
[Storage and set boundary conditions]
static LoggerType & setLog(const std::string channel)
Set ans resset chanel logger.
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_TAG(channel, tag)
Tag channel.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
virtual MoFEMErrorCode loop_dofs(const Problem *problem_ptr, const std::string &field_name, RowColData rc, DofMethod &method, int lower_rank, int upper_rank, int verb=DEFAULT_VERBOSITY)=0
Make a loop over dofs.
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
FTensor::Index< 'i', SPACE_DIM > i
static double lambda
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
const FTensor::Tensor2< T, Dim, Dim > Vec
MoFEMErrorCode postProcessElasticResults(MoFEM::Interface &mField, SmartPetscObj< DM > dm, const std::string &domain_fe_name, const std::string &out_file_name, SmartPetscObj< Vec > extra_vector=nullptr, const std::string &extra_vector_name="", const Sev hooke_ops_sev=Sev::verbose)
constexpr double eps
Definition HenckyOps.hpp:13
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
auto createKSP(MPI_Comm comm)
PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem)
Definition DMMoFEM.cpp:434
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data)
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
static const bool debug
PetscErrorCode PetscOptionsGetBool(PetscOptions *, const char pre[], const char name[], PetscBool *bval, PetscBool *set)
MoFEM::OpBrokenTopoBase OP
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
PetscErrorCode PetscOptionsGetRealArray(PetscOptions *, const char pre[], const char name[], PetscReal dval[], PetscInt *nmax, PetscBool *set)
auto createVectorMPI(MPI_Comm comm, PetscInt n, PetscInt N)
Create MPI Vector.
auto getDMKspCtx(DM dm)
Get KSP context data structure used by DM.
Definition DMMoFEM.hpp:1248
MoFEM::OpBrokenTopoBase GAUSS
Gaussian quadrature integration.
static MoFEMErrorCode invertTensor(FTensor::Tensor2< T1, DIM, DIM > &t, T2 &det, FTensor::Tensor2< T3, DIM, DIM > &inv_t)
static auto determinantTensor(FTensor::Tensor2< T, DIM, DIM > &t)
Calculate the determinant of a tensor of rank DIM.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
PetscErrorCode PetscOptionsGetEList(PetscOptions *, const char pre[], const char name[], const char *const *list, PetscInt next, PetscInt *value, PetscBool *set)
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
auto get_temp_meshset_ptr(moab::Interface &moab)
Create smart pointer to temporary meshset.
PetscErrorCode TaoSetObjectiveAndGradient(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx)
Sets the objective function value and gradient for a TAO optimization solver.
Definition TaoCtx.cpp:178
auto createDM(MPI_Comm comm, const std::string dm_type_name)
Creates smart DM object.
MoFEMErrorCode VecSetValues(Vec V, const EntitiesFieldData::EntData &data, const double *ptr, InsertMode iora)
Assemble PETSc vector.
auto createTao(MPI_Comm comm)
boost::shared_ptr< ObjectiveFunctionData > create_python_objective_function(std::string py_file)
Factory function to create Python-integrated objective function interface.
constexpr IntegrationType I
constexpr AssemblyType A
OpPostProcMapInMoab< SPACE_DIM, SPACE_DIM > OpPPMap
constexpr auto field_name
static constexpr int approx_order
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpMass< 1, SPACE_DIM > OpMass
[Only used with Hooke equation (linear material model)]
Definition seepage.cpp:56
PetscBool is_plane_strain
Definition seepage.cpp:177
constexpr double g
FTensor::Index< 'm', 3 > m
Boundary conditions marker.
Definition elastic.cpp:39
[Define entities]
Definition elastic.cpp:38
[Example]
Definition plastic.cpp:216
MoFEMErrorCode boundaryCondition()
Apply essential boundary conditions.
MoFEMErrorCode assembleSystem()
Setup operators in finite element pipeline.
MoFEMErrorCode readMesh()
Read mesh from file and setup meshsets.
SmartPetscObj< KSP > kspElastic
Linear solver for elastic problem.
Definition adjoint.cpp:147
std::vector< SmartPetscObj< Vec > > modeVecs
Topology mode vectors (design variables)
Definition adjoint.cpp:152
FieldApproximationBase base
Choice of finite element basis functions
Definition plot_base.cpp:68
std::vector< std::array< double, 3 > > modeCentroids
Centroids of optimization blocks
Definition adjoint.cpp:153
MoFEMErrorCode topologyModes()
Compute topology optimization modes.
Definition adjoint.cpp:702
SmartPetscObj< Vec > initialGeometry
Initial geometry field.
Definition adjoint.cpp:155
int fieldOrder
Polynomial order for approximation.
Definition adjoint.cpp:144
Example(MoFEM::Interface &m_field)
Definition adjoint.cpp:115
SmartPetscObj< Mat > M
MoFEMErrorCode runProblem()
Main driver function for the optimization process.
MoFEMErrorCode calculateGradient(PetscReal *objective_function_value, Vec objective_function_gradient, Vec adjoint_vector)
Calculate objective function gradient using adjoint method.
Definition adjoint.cpp:1741
MoFEMErrorCode setupAdJoint()
Setup adjoint fields and finite elements
Definition adjoint.cpp:579
MoFEM::Interface & mField
Reference to MoFEM interface.
Definition plastic.cpp:226
std::vector< std::array< double, 6 > > modeBBoxes
Bounding boxes of optimization blocks.
Definition adjoint.cpp:154
boost::shared_ptr< ObjectiveFunctionData > pythonPtr
Interface to Python objective function.
Definition adjoint.cpp:149
SmartPetscObj< DM > adjointDM
Data manager for adjoint problem.
Definition adjoint.cpp:148
MoFEMErrorCode setupProblem()
Setup fields, approximation spaces and DOFs.
MoFEMErrorCode postprocessElastic(int iter, SmartPetscObj< Vec > adjoint_vector=nullptr)
Post-process and output results.
Definition adjoint.cpp:1303
MoFEMErrorCode solveElastic()
Solve forward elastic problem.
Definition adjoint.cpp:1175
boost::shared_ptr< MatrixDouble > vectorFieldPtr
Field values at evaluation points.
Definition adjoint.cpp:123
Add operators pushing bases from local to physical configuration.
Boundary condition manager for finite element problem setup.
Managing BitRefLevels.
MoFEMErrorCode synchroniseEntities(Range &ent, std::map< int, Range > *received_ents, int verb=DEFAULT_VERBOSITY)
synchronize entity range on processors (collective)
virtual moab::Interface & get_moab()=0
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
virtual MoFEMErrorCode add_field(const std::string name, const FieldSpace space, const FieldApproximationBase base, const FieldCoefficientsNumber nb_of_coefficients, const TagType tag_type=MB_TAG_SPARSE, const enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add field.
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
Core (interface) class.
Definition Core.hpp:82
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
Definition Core.cpp:72
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition Core.cpp:118
Deprecated interface functions.
Definition of the displacement bc data structure.
Definition BCData.hpp:72
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1DiffN(const FieldApproximationBase base)
Get derivatives of base functions.
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
Class (Function) to enforce essential constrains on the left hand side diagonal.
Definition Essential.hpp:33
Class (Function) to enforce essential constrains on the right hand side diagonal.
Definition Essential.hpp:41
Class (Function) to enforce essential constrains.
Definition Essential.hpp:25
Basic algebra on fields.
Definition FieldBlas.hpp:21
MoFEMErrorCode fieldLambdaOnValues(OneFieldFunctionOnValues lambda, const std::string field_name, Range *ents_ptr=nullptr)
field lambda
Definition FieldBlas.cpp:21
Field evaluator interface.
static boost::shared_ptr< SinkType > createSink(boost::shared_ptr< std::ostream > stream_ptr, std::string comm_filter)
Create a sink object.
static boost::shared_ptr< std::ostream > getStrmWorld()
Get the strm world object.
static boost::shared_ptr< std::ostream > getStrmSync()
Get the strm sync object.
Interface for managing meshsets containing materials and boundary conditions.
Assembly methods.
Definition Natural.hpp:65
Get field gradients at integration pts for scalar field rank 0, i.e. vector field.
Specialization for double precision vector field values calculation.
Post post-proc data at points from hash maps.
Template struct for dimension-specific finite element types.
PipelineManager interface.
MoFEMErrorCode setDomainRhsIntegrationRule(RuleHookFun rule)
Set integration rule for domain right-hand side finite element.
Problem manager is used to build and partition problems.
Projection of edge entities with one mid-node on hierarchical basis.
Simple interface for fast problem set-up.
Definition Simple.hpp:27
MoFEMErrorCode getOptions()
get options
Definition Simple.cpp:180
MoFEMErrorCode getDM(DM *dm)
Get DM.
Definition Simple.cpp:799
intrusive_ptr for managing petsc objects
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Vector manager is used to create vectors \mofem_vectors.
OpAdJointGradTimesSymTensor(const std::string field_name, boost::shared_ptr< HookeOps::CommonData > comm_ptr, boost::shared_ptr< MatrixDouble > jac, boost::shared_ptr< MatrixDouble > diff_jac, boost::shared_ptr< VectorDouble > cof_vals)
Definition adjoint.cpp:1328
Forward declaration of operator for gradient times symmetric tensor operations.
Definition adjoint.cpp:88
boost::shared_ptr< MatrixDouble > dGradPtr
Definition adjoint.cpp:1733
boost::shared_ptr< double > globObjectiveGradPtr
Definition adjoint.cpp:1738
boost::shared_ptr< ObjectiveFunctionData > pythonPtr
Definition adjoint.cpp:1728
boost::shared_ptr< double > globObjectivePtr
Definition adjoint.cpp:1737
ForcesAndSourcesCore::UserDataOperator OP
Definition adjoint.cpp:1574
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Compute objective function contributions at element level.
Definition adjoint.cpp:1598
boost::shared_ptr< HookeOps::CommonData > commPtr
Definition adjoint.cpp:1729
boost::shared_ptr< VectorDouble > cofVals
Definition adjoint.cpp:1732
boost::shared_ptr< MatrixDouble > uPtr
Definition adjoint.cpp:1735
OpAdJointObjective(boost::shared_ptr< ObjectiveFunctionData > python_ptr, boost::shared_ptr< HookeOps::CommonData > comm_ptr, boost::shared_ptr< MatrixDouble > jac_ptr, boost::shared_ptr< MatrixDouble > diff_jac, boost::shared_ptr< VectorDouble > cof_vals, boost::shared_ptr< MatrixDouble > d_grad_ptr, boost::shared_ptr< MatrixDouble > d_u_ptr, boost::shared_ptr< MatrixDouble > u_ptr, boost::shared_ptr< double > glob_objective_ptr, boost::shared_ptr< double > glob_objective_grad_ptr)
Definition adjoint.cpp:1576
boost::shared_ptr< MatrixDouble > diffJacPtr
Definition adjoint.cpp:1731
boost::shared_ptr< MatrixDouble > jacPtr
Definition adjoint.cpp:1730
boost::shared_ptr< MatrixDouble > dUPtr
Definition adjoint.cpp:1734
OpStateSensitivity(const std::string field_name, boost::shared_ptr< ObjectiveFunctionData > python_ptr, boost::shared_ptr< HookeOps::CommonData > comm_ptr, boost::shared_ptr< MatrixDouble >u_ptr)
Definition adjoint.cpp:1479
boost::shared_ptr< ObjectiveFunctionData > pythonPtr
Definition adjoint.cpp:1568
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &data)
Definition adjoint.cpp:1486
boost::shared_ptr< MatrixDouble > uPtr
Definition adjoint.cpp:1570
boost::shared_ptr< HookeOps::CommonData > commPtr
Definition adjoint.cpp:1569
PipelineManager::ElementsAndOpsByDim< 2 >::FaceSideEle SideEle
PipelineManager::ElementsAndOpsByDim< 3 >::FaceSideEle SideEle
#define EXECUTABLE_DIMENSION
Definition plastic.cpp:13
PetscBool do_eval_field
Evaluate field.
Definition plastic.cpp:119
ElementsAndOps< SPACE_DIM >::SideEle SideEle
Definition plastic.cpp:61
auto save_range
constexpr int SPACE_DIM
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::BoundaryEle BoundaryEle
PostProcEleByDim< SPACE_DIM >::SideEle SideEle
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::DomainEle DomainEle