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;
34 IntegrationType::GAUSS;
83template <
int SPACE_DIM, IntegrationType I,
typename OpBase>
87template <
int SPACE_DIM, IntegrationType I,
typename OpBase>
97#include <ElasticSpring.hpp>
98#include <FluidLevel.hpp>
99#include <CalculateTraction.hpp>
100#include <NaturalDomainBC.hpp>
101#include <NaturalBoundaryBC.hpp>
102#include <HookeOps.hpp>
109 const std::string block_name,
int dim);
139 Vec objective_function_gradient,
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);
188 auto file_exists = [](std::string myfile) {
189 std::ifstream file(myfile.c_str());
195 if (!file_exists(objective_function_file_name)) {
196 MOFEM_LOG(
"WORLD", Sev::error) <<
"Objective function file NOT found: "
197 << objective_function_file_name;
204 char sensitivity_method_name[32] =
"adjoint";
206 "-sensitivity_method", sensitivity_method_name,
207 sizeof(sensitivity_method_name),
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") {
219 "Unknown -sensitivity_method. Use 'direct' or 'adjoint'.");
222 <<
"Sensitivity method: " << sensitivity_method;
240 auto create_vec_modes = [&](
auto block_name) {
245 std::regex((boost::format(
"%s(.*)") % block_name).str())
249 int nb_total_modes = 0;
250 for (
auto &bc : bcs) {
251 auto id = bc->getMeshsetId();
254 nb_total_modes += nb_modes;
258 <<
"Total number of modes to apply: " << nb_total_modes;
265 auto get_modes_bounding_boxes = [&](
auto block_name) {
270 std::regex((boost::format(
"%s(.*)") % block_name).str())
274 for (
auto &bc : bcs) {
275 auto meshset = bc->getMeshset();
280 std::vector<double> x(verts.size());
281 std::vector<double> y(verts.size());
282 std::vector<double> z(verts.size());
284 std::array<double, 3> centroid = {0, 0, 0};
285 for (
int i = 0;
i != verts.size(); ++
i) {
290 MPI_Allreduce(MPI_IN_PLACE, centroid.data(), 3, MPI_DOUBLE, MPI_SUM,
292 int nb_vertex = verts.size();
293 MPI_Allreduce(MPI_IN_PLACE, &nb_vertex, 1, MPI_INT, MPI_SUM,
296 centroid[0] /= nb_vertex;
297 centroid[1] /= nb_vertex;
298 centroid[2] /= nb_vertex;
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]);
310 MPI_Allreduce(MPI_IN_PLACE, &bbox[0], 3, MPI_DOUBLE, MPI_MIN,
312 MPI_Allreduce(MPI_IN_PLACE, &bbox[3], 3, MPI_DOUBLE, MPI_MAX,
316 <<
"Block: " << bc->getName() <<
" centroid: " << centroid[0] <<
" "
317 << centroid[1] <<
" " << centroid[2];
319 <<
"Block: " << bc->getName() <<
" bbox: " << bbox[0] <<
" "
320 << bbox[1] <<
" " << bbox[2] <<
" " << bbox[3] <<
" " << bbox[4]
330 auto eval_objective_and_gradient = [](Tao tao,
Vec x, PetscReal *
f,
Vec g,
331 void *ctx) -> PetscErrorCode {
336 CHKERR TaoGetIterationNumber(tao, &iter);
338 auto set_geometry = [&](
Vec x) {
343 CHKERR VecScatterCreateToAll(x, &ctx, &x_local);
345 CHKERR VecScatterBegin(ctx, x, x_local, INSERT_VALUES, SCATTER_FORWARD);
346 CHKERR VecScatterEnd(ctx, x, x_local, INSERT_VALUES, SCATTER_FORWARD);
348 CHKERR VecScatterDestroy(&ctx);
351 CHKERR VecCopy(ex_ptr->initialGeometry, current_geometry);
353 CHKERR VecGetArrayRead(x_local, &
a);
354 const double *coeff =
a;
355 for (
auto &mode_vec : ex_ptr->
modeVecs) {
357 <<
"Adding mode with coeff: " << *coeff;
358 CHKERR VecAXPY(current_geometry, (*coeff), mode_vec);
361 CHKERR VecRestoreArrayRead(x_local, &
a);
362 CHKERR VecGhostUpdateBegin(current_geometry, INSERT_VALUES,
364 CHKERR VecGhostUpdateEnd(current_geometry, INSERT_VALUES,
367 ->setOtherLocalGhostVector(
"ADJOINT",
"ADJOINT_FIELD",
"GEOMETRY",
369 INSERT_VALUES, SCATTER_REVERSE);
371 CHKERR VecDestroy(&x_local);
376 CHKERR KSPReset(ex_ptr->kspElastic);
377 CHKERR ex_ptr->solveElastic();
380 CHKERR ex_ptr->calculateGradient(
f,
g, adjoint_vector);
381 CHKERR ex_ptr->postprocessElastic(iter, adjoint_vector);
387 CHKERR create_vec_modes(
"OPTIMISE");
388 CHKERR get_modes_bounding_boxes(
"OPTIMISE");
401 CHKERR TaoSetType(tao, TAOLMVM);
413 INSERT_VALUES, SCATTER_FORWARD);
442 CHKERR TaoSetSolution(tao, x0);
443 CHKERR TaoSetFromOptions(tao);
447 MOFEM_LOG(
"WORLD", Sev::inform) <<
"Topology optimization completed";
500 enum bases { AINSWORTH, DEMKOWICZ, LASBASETOPT };
501 const char *list_bases[LASBASETOPT] = {
"ainsworth",
"demkowicz"};
502 PetscInt choice_base_value = AINSWORTH;
504 LASBASETOPT, &choice_base_value, PETSC_NULLPTR);
506 switch (choice_base_value) {
510 <<
"Set AINSWORTH_LEGENDRE_BASE for displacements";
515 <<
"Set DEMKOWICZ_JACOBI_BASE for displacements";
549 auto project_ho_geometry = [&]() {
553 CHKERR project_ho_geometry();
584 auto create_adjoint_dm = [&]() {
587 auto add_field = [&]() {
592 for (
auto tt = MBEDGE; tt <= moab::CN::TypeDimensionMap[
SPACE_DIM].second;
602 auto add_adjoint_fe_impl = [&]() {
627 auto block_name =
"OPTIMISE";
631 std::regex((boost::format(
"%s(.*)") % block_name).str())
635 for (
auto bc : bcs) {
637 bc->getMeshset(),
SPACE_DIM - 1,
"ADJOINT_BOUNDARY_FE");
643 simple->getBitRefLevelMask());
648 auto set_adjoint_dm_imp = [&]() {
652 simple->getBitRefLevelMask());
654 CHKERR DMSetFromOptions(adjoint_dm);
661 CHKERR DMSetUp(adjoint_dm);
686 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
"REMOVE_X",
688 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
"REMOVE_Y",
690 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
"REMOVE_Z",
692 CHKERR bc_mng->removeBlockDOFsOnEntities(
simple->getProblemName(),
693 "REMOVE_ALL",
"U", 0, 3);
695 simple->getProblemName(),
"U");
711 boost::make_shared<Range>(opt_ents));
713 boost::make_shared<Range>(opt_ents));
714 CHKERR DMSetUp(subset_dm_bdy);
722 CHKERR DMSetUp(subset_dm_domain);
725 auto remove_dofs = [&]() {
728 std::array<Range, 3> remove_dim_ents;
736 for (
int d = 0; d != 3; ++d) {
738 <<
"Removing topology modes on block OPT_REMOVE_" << (char)(
'X' + d)
739 <<
" with " << remove_dim_ents[d].size() <<
" entities";
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]);
751 ParallelComm *pcomm =
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) {
760 moab::Interface::UNION);
761 boundary_ents.merge(ents);
765 boundary_ents.merge(verts);
770 boundary_ents.merge(opt_ents);
772 "SUBSET_DOMAIN",
"ADJOINT_FIELD", boundary_ents);
773 for (
int d = 0; d != 3; ++d) {
775 "SUBSET_DOMAIN",
"ADJOINT_FIELD", remove_dim_ents[d], d, d);
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;
795 auto &pip = fe_lhs->getOpPtrVector();
800 pip.push_back(
new OpMass(
"ADJOINT_FIELD",
"ADJOINT_FIELD",
801 [](
double,
double,
double) {
return 1.; }));
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;
810 auto &pip = fe_rhs->getOpPtrVector();
817 auto block_name =
"OPTIMISE";
821 std::regex((boost::format(
"%s(.*)") % block_name).str())
830 struct OpMode :
public OP {
831 OpMode(
const std::string name,
832 boost::shared_ptr<ObjectiveFunctionData> python_ptr,
int id,
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) {}
846 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
854 auto nb_base_functions = data.
getN().size2();
857 CHKERR pythonPtr->blockModes(iD, OP::getCoordsAtGaussPts(),
858 modeCentroids[blockCounter],
859 modeBboxes[blockCounter], blockModes);
861 auto nb_integration_pts = getGaussPts().size2();
862 if (blockModes.size2() != 3 * nb_integration_pts) {
864 <<
"Number of modes does not match number of integration points: "
865 << blockModes.size2() <<
"!=" << 3 * nb_integration_pts;
871 int nb_modes = blockModes.size1();
872 for (
auto mode = 0; mode != nb_modes; ++mode) {
875 auto t_mode = getFTensor1FromPtr<3>(&blockModes(mode, 0));
877 const double vol = OP::getMeasure();
879 auto t_w = OP::getFTensor0IntegrationWeight();
883 for (
int gg = 0; gg != nb_integration_pts; gg++) {
886 const double alpha = t_w * vol;
888 auto t_nf = getFTensor1FromPtr<SPACE_DIM>(nf.data().data());
890 for (; rr != nb_rows /
SPACE_DIM; ++rr) {
891 t_nf(
i) += alpha * t_base * t_mode(
i);
895 for (; rr < nb_base_functions; ++rr)
900 Vec vec = modeVecs[modeCounter + mode];
902 auto *indices = data.
getIndices().data().data();
903 auto *nf_data = nf.data().data();
911 boost::shared_ptr<ObjectiveFunctionData> pythonPtr;
913 std::vector<std::array<double, 3>> modeCentroids;
914 std::vector<std::array<double, 6>> modeBboxes;
916 std::vector<SmartPetscObj<Vec>> modeVecs;
921 auto solve_bdy = [&]() {
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();
933 auto range = boost::make_shared<Range>(ents);
934 auto &pip_rhs = fe_rhs->getOpPtrVector();
937 mode_counter, range));
944 mode_counter += nb_modes;
946 <<
"Setting mode block block: " << bc->getName()
947 <<
" with ID: " << bc->getMeshsetId()
948 <<
" total modes: " << mode_counter;
954 CHKERR VecGhostUpdateBegin(
v, ADD_VALUES, SCATTER_REVERSE);
955 CHKERR VecGhostUpdateEnd(
v, ADD_VALUES, SCATTER_REVERSE);
962 CHKERR MatAssemblyBegin(
M, MAT_FINAL_ASSEMBLY);
963 CHKERR MatAssemblyEnd(
M, MAT_FINAL_ASSEMBLY);
966 CHKERR KSPSetOperators(solver,
M,
M);
967 CHKERR KSPSetFromOptions(solver);
976 CHKERR VecGhostUpdateBegin(
v, INSERT_VALUES, SCATTER_FORWARD);
977 CHKERR VecGhostUpdateEnd(
v, INSERT_VALUES, SCATTER_FORWARD);
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;
990 auto &pip = fe->getOpPtrVector();
993 CHKERR HookeOps::opFactoryDomainLhs<SPACE_DIM, A, I, DomainEleOp>(
994 mField, pip,
"ADJOINT_FIELD",
"MAT_ADJOINT", Sev::noisy);
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;
1003 auto &pip = fe->getOpPtrVector();
1006 CHKERR HookeOps::opFactoryDomainRhs<SPACE_DIM, A, I, DomainEleOp>(
1007 mField, pip,
"ADJOINT_FIELD",
"MAT_ADJOINT", Sev::noisy);
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,
1017 auto post_proc_end = boost::make_shared<PostProcBrokenMeshInMoabBaseEnd>(
1020 auto geom_vec = boost::make_shared<MatrixDouble>();
1023 boost::make_shared<PostProcEleDomain>(
mField, post_proc_mesh);
1025 post_proc_fe->getOpPtrVector(), {H1},
"GEOMETRY");
1026 post_proc_fe->getOpPtrVector().push_back(
1032 post_proc_fe->getOpPtrVector().push_back(
1036 post_proc_fe->getPostProcMesh(), post_proc_fe->getMapGaussPts(),
1040 {{
"MODE", geom_vec}},
1051 post_proc_begin->getFEMethod());
1055 post_proc_begin->getFEMethod());
1057 CHKERR post_proc_end->writeFile(
"mode_" + std::to_string(mode) +
".h5m");
1062 auto solve_domain = [&]() {
1064 auto fe_lhs = get_elastic_fe_lhs();
1065 auto fe_rhs = get_elastic_fe_rhs();
1074 CHKERR MatAssemblyBegin(M, MAT_FINAL_ASSEMBLY);
1075 CHKERR MatAssemblyEnd(M, MAT_FINAL_ASSEMBLY);
1077 auto solver =
createKSP(mField.get_comm());
1078 CHKERR KSPSetOperators(solver, M, M);
1079 CHKERR KSPSetFromOptions(solver);
1082 int mode_counter = 0;
1083 for (
auto &f : modeVecs) {
1084 CHKERR mField.getInterface<
FieldBlas>()->setField(0,
"ADJOINT_FIELD");
1092 CHKERR VecGhostUpdateBegin(
F, ADD_VALUES, SCATTER_REVERSE);
1093 CHKERR VecGhostUpdateEnd(
F, ADD_VALUES, SCATTER_REVERSE);
1095 CHKERR VecGhostUpdateBegin(
v, INSERT_VALUES, SCATTER_FORWARD);
1096 CHKERR VecGhostUpdateEnd(
v, INSERT_VALUES, SCATTER_FORWARD);
1110 for (
int i = 0;
i < modeVecs.size(); ++
i) {
1111 CHKERR adjoint_gradient_postprocess(
i);
1134 CHKERR pip->setBoundaryRhsIntegrationRule(integration_rule_bc);
1135 CHKERR pip->setBoundaryLhsIntegrationRule(integration_rule_bc);
1139 pip->getOpDomainLhsPipeline(), {H1},
"GEOMETRY");
1141 pip->getOpDomainRhsPipeline(), {H1},
"GEOMETRY");
1143 pip->getOpBoundaryRhsPipeline(), {NOSPACE},
"GEOMETRY");
1145 pip->getOpBoundaryLhsPipeline(), {NOSPACE},
"GEOMETRY");
1149 CHKERR HookeOps::opFactoryDomainLhs<SPACE_DIM, A, I, DomainEleOp>(
1150 mField, pip->getOpDomainLhsPipeline(),
"U",
"MAT_ELASTIC", Sev::noisy);
1154 CHKERR HookeOps::opFactoryDomainRhs<SPACE_DIM, A, I, DomainEleOp>(
1155 mField, pip->getOpDomainRhsPipeline(),
"U",
"MAT_ELASTIC", Sev::noisy);
1159 pip->getOpDomainRhsPipeline(),
mField,
"U", Sev::inform);
1165 pip->getOpBoundaryRhsPipeline(),
mField,
"U", -1, Sev::inform);
1168 pip->getOpBoundaryLhsPipeline(),
mField,
"U", Sev::noisy);
1181 CHKERR VecZeroEntries(d);
1184 auto set_essential_bc = [&]() {
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>();
1193 auto get_pre_proc_hook = [&]() {
1197 pre_proc_rhs->preProcessHook = get_pre_proc_hook();
1199 auto get_post_proc_hook_rhs = [
this, post_proc_rhs]() {
1203 post_proc_rhs, 1.)();
1207 auto get_post_proc_hook_lhs = [
this, post_proc_lhs]() {
1211 post_proc_lhs, 1.)();
1215 post_proc_rhs->postProcessHook = get_post_proc_hook_rhs;
1216 post_proc_lhs->postProcessHook = get_post_proc_hook_lhs;
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);
1224 auto setup_and_solve = [&](
auto solver) {
1226 BOOST_LOG_SCOPED_THREAD_ATTR(
"Timeline", attrs::timer());
1227 MOFEM_LOG(
"TIMER", Sev::noisy) <<
"KSPSetUp";
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";
1239 CHKERR set_essential_bc();
1242 CHKERR VecGhostUpdateBegin(d, INSERT_VALUES, SCATTER_FORWARD);
1243 CHKERR VecGhostUpdateEnd(d, INSERT_VALUES, SCATTER_FORWARD);
1246 auto evaluate_field_at_the_point = [&]() {
1250 std::array<double, 3> field_eval_coords{0.0, 0.0, 0.0};
1253 field_eval_coords.data(), &coords_dim,
1259 auto field_eval_data =
1263 ->buildTree<SPACE_DIM>(field_eval_data,
simple->getDomainFEName());
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;
1270 field_eval_fe_ptr->getOpPtrVector().push_back(
1274 ->evalFEAtThePoint<SPACE_DIM>(
1275 field_eval_coords.data(), 1e-12,
simple->getProblemName(),
1276 simple->getDomainFEName(), field_eval_data,
1283 MOFEM_LOG(
"FieldEvaluator", Sev::inform)
1284 <<
"U_X: " << t_disp(0) <<
" U_Y: " << t_disp(1);
1286 MOFEM_LOG(
"FieldEvaluator", Sev::inform)
1287 <<
"U_X: " << t_disp(0) <<
" U_Y: " << t_disp(1)
1288 <<
" U_Z: " << t_disp(2);
1296 CHKERR evaluate_field_at_the_point();
1310 "out_elastic_" + std::to_string(iter) +
".h5m", adjoint_vector,
1311 "ADJOINT", Sev::noisy);
1322template <
int SPACE_DIM>
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)
1334 jac(jac), diffJac(diff_jac), cofVals(cof_vals) {}
1338 boost::shared_ptr<MatrixDouble>
jac;
1353 t_diff(
i,
j,
k,
l) = 0;
1354 t_diff(0, 0, 0, 0) = 1;
1355 t_diff(1, 1, 1, 1) = 1;
1357 t_diff(1, 0, 1, 0) = 0.5;
1358 t_diff(1, 0, 0, 1) = 0.5;
1360 t_diff(0, 1, 0, 1) = 0.5;
1361 t_diff(0, 1, 1, 0) = 0.5;
1363 if constexpr (DIM == 3) {
1364 t_diff(2, 2, 2, 2) = 1;
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;
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;
1380template <
int SPACE_DIM>
1395 const double vol = OP::getMeasure();
1397 auto t_w = OP::getFTensor0IntegrationWeight();
1399 auto t_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(jac));
1401 auto t_diff_jac = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(diffJac));
1408 getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(commPtr->matGradPtr));
1410 auto t_cauchy_stress =
1411 getFTensor2SymmetricFromMat<SPACE_DIM>(*(commPtr->getMatCauchyStress()));
1414 getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*(commPtr->matDPtr));
1416 for (
int gg = 0; gg != OP::nbIntegrationPts; gg++) {
1418 const double alpha = t_w * vol;
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);
1433 t_diff_strain(
i,
j) = t_diff_symm(
i,
j,
k,
l) * t_diff_grad(
k,
l);
1437 t_diff_stress(
i,
j) = t_D(
i,
j,
k,
l) * t_diff_strain(
k,
l);
1440 auto t_nf = OP::template getNf<SPACE_DIM>();
1443 for (; rr != OP::nbRows /
SPACE_DIM; rr++) {
1446 t_diff_row_grad(
k) = t_row_grad(
j) * t_diff_inv_jac(
j,
k);
1449 t_nf(
j) += alpha * t_diff_row_grad(
i) * t_cauchy_stress(
i,
j);
1452 t_nf(
j) += (alpha * t_cof) * t_row_grad(
i) * t_cauchy_stress(
i,
j);
1455 t_nf(
j) += alpha * t_row_grad(
i) * t_diff_stress(
i,
j);
1460 for (; rr < OP::nbRowBaseFunctions; ++rr) {
1480 boost::shared_ptr<ObjectiveFunctionData> python_ptr,
1481 boost::shared_ptr<HookeOps::CommonData> comm_ptr,
1482 boost::shared_ptr<MatrixDouble>u_ptr)
1495 auto nb_gauss_pts = getGaussPts().size2();
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);
1503 auto evaluate_python = [&]()
1506 auto &coords = OP::getCoordsAtGaussPts();
1517 auto vol = OP::getMeasure();
1518 auto t_w = OP::getFTensor0IntegrationWeight();
1520 auto t_D = getFTensor4DdgFromMat<SPACE_DIM, SPACE_DIM, 0>(*(
commPtr->matDPtr));
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);
1529 for(
int gg = 0; gg != nb_gauss_pts; ++gg)
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);
1535 auto t_nf = OP::template getNf<SPACE_DIM>();
1537 for (; rr != OP::nbRows /
SPACE_DIM; rr++)
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);
1547 for (; rr < OP::nbRowBaseFunctions; ++rr)
1562 CHKERR evaluate_python();
1570 boost::shared_ptr<MatrixDouble>
uPtr;
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)
1615 auto nb_gauss_pts = getGaussPts().size2();
1616 auto objective_ptr = boost::make_shared<VectorDouble>(nb_gauss_pts);
1617 auto objective_dstress =
1618 boost::make_shared<MatrixDouble>(symm_size, nb_gauss_pts);
1619 auto objective_dstrain =
1620 boost::make_shared<MatrixDouble>(symm_size, nb_gauss_pts);
1621 auto objective_du = boost::make_shared<MatrixDouble>(
SPACE_DIM, nb_gauss_pts);
1624 auto evaluate_python = [&]() {
1626 auto &coords = OP::getCoordsAtGaussPts();
1641 getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(
commPtr->matGradPtr));
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));
1647 auto t_d_grad = getFTensor2FromMat<SPACE_DIM, SPACE_DIM>(*(
dGradPtr));
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);
1657 auto vol = OP::getMeasure();
1658 auto t_w = getFTensor0IntegrationWeight();
1659 for (
auto gg = 0; gg != nb_gauss_pts; ++gg) {
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);
1672 t_d_strain(
i,
j) = t_diff_symm(
i,
j,
k,
l) * (
1682 auto alpha = t_w * vol;
1684 (*globObjectivePtr) += alpha * t_obj;
1685 (*globObjectiveGradPtr) +=
1689 t_obj_dstress(
i,
j) * (t_D(
i,
j,
k,
l) * t_d_strain(
k,
l))
1693 t_obj_dstrain(
i,
j) * t_d_strain(
i,
j)
1697 t_obj_du(
i) * t_d_u(
i)
1722 CHKERR evaluate_python();
1735 boost::shared_ptr<MatrixDouble>
uPtr;
1742 Vec objective_function_gradient,
1743 Vec adjoint_vector) {
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]() {
1756 post_proc_rhs, 0)();
1759 post_proc_rhs->postProcessHook = get_post_proc_hook_rhs;
1760 return post_proc_rhs;
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;
1768 auto &pip = fe->getOpPtrVector();
1771 HookeOps::opFactoryDomainRhs<SPACE_DIM, A, I, DomainEleOp>(
1772 mField, pip,
"U",
"MAT_ELASTIC", Sev::noisy);
1777 auto calulate_fd_residual = [&](
auto eps,
auto diff_vec,
auto fd_vec) {
1780 constexpr bool debug =
false;
1786 auto norm2_field = [&](
const double val) {
1791 MPI_Allreduce(MPI_IN_PLACE, &nrm2, 1, MPI_DOUBLE, MPI_SUM,
1793 MOFEM_LOG(
"WORLD", Sev::inform) <<
"Geometry norm: " << sqrt(nrm2);
1797 if constexpr (
debug)
1803 initial_current_geometry, INSERT_VALUES, SCATTER_FORWARD);
1804 CHKERR VecAssemblyBegin(initial_current_geometry);
1805 CHKERR VecAssemblyEnd(initial_current_geometry);
1807 if constexpr (
debug)
1810 auto perturb_geometry = [&](
auto eps,
auto diff_vec) {
1813 CHKERR VecCopy(initial_current_geometry, current_geometry);
1814 CHKERR VecAXPY(current_geometry,
eps, diff_vec);
1817 current_geometry, INSERT_VALUES, SCATTER_REVERSE);
1821 auto fe = get_fd_direvative_fe();
1824 auto calc_impl = [&](
auto f,
auto eps) {
1830 simple->getDomainFEName(), fe);
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());
1843 CHKERR VecWAXPY(fd_vec, -1.0, fm, fp);
1844 CHKERR VecScale(fd_vec, 1.0 / (2.0 *
eps));
1848 initial_current_geometry, INSERT_VALUES, SCATTER_REVERSE);
1850 if constexpr (
debug)
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;
1861 auto &pip = fe_adjoint->getOpPtrVector();
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>();
1875 "GEOMETRY", jac_ptr));
1877 "U", diff_jac_ptr, diff_vec));
1879 pip.push_back(
new OpCoFactor(jac_ptr, diff_jac_ptr, cof_ptr));
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));
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;
1897 auto &pip = fe_adjoint->getOpPtrVector();
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>();
1912 "GEOMETRY", jac_ptr));
1917 pip.push_back(
new OpCoFactor(jac_ptr, diff_jac_ptr, cof_ptr));
1919 "U", d_grad_ptr, grad_vec));
1924 auto common_ptr = HookeOps::commonDataFactory<SPACE_DIM, I, DomainEleOp>(
1925 mField, pip,
"U",
"MAT_ELASTIC", Sev::noisy);
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));
1933 auto dm =
simple->getDM();
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,
1945 CHKERR VecGhostUpdateEnd(zero_state_sensitivity, INSERT_VALUES,
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);
1964 auto set_variance_of_geometry =
1965 [&](
auto mode,
auto mod_vec) {
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);
1979 auto calculate_variance_internal_forces = [&](
auto mode,
auto mod_vec) {
1982 CHKERR VecGhostUpdateBegin(
f, INSERT_VALUES, SCATTER_FORWARD);
1983 CHKERR VecGhostUpdateEnd(
f, INSERT_VALUES, SCATTER_FORWARD);
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());
1997 constexpr bool debug =
false;
1998 if constexpr (
debug) {
2000 CHKERR VecNorm(
f, NORM_2, &norm0);
2003 CHKERR calulate_fd_residual(
eps, dm_diff_vec, fd_check);
2005 CHKERR VecAXPY(fd_check, -1.0,
f);
2006 CHKERR VecNorm(fd_check, NORM_2, &nrm);
2008 <<
" FD check for internal forces [ " << mode <<
" ]: " << nrm
2009 <<
" / " << norm0 <<
" ( " << (nrm / norm0) <<
" )";
2016 auto calculate_variance_of_displacement = [&](
auto mode,
auto mod_vec) {
2021 CHKERR VecGhostUpdateBegin(d, INSERT_VALUES, SCATTER_FORWARD);
2022 CHKERR VecGhostUpdateEnd(d, INSERT_VALUES, SCATTER_FORWARD);
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;
2033 std::array<double, 2> array = {*objective_ptr, *objective_grad_ptr};
2034 MPI_Allreduce(MPI_IN_PLACE, array.data(), 2, MPI_DOUBLE, MPI_SUM,
2036 objective_value = array[0];
2037 objective_gradient = array[1];
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;
2052 auto calculate_variance_of_objective_function_dJ_du = [&](Vec dJ_du) {
2055 auto fe = boost::make_shared<DomainEle>(
mField);
2056 fe->getRuleHook = [](int, int,
int p_data) {
2057 return 2 * p_data + p_data - 1;
2059 auto &pip = fe->getOpPtrVector();
2062 auto u_ptr = boost::make_shared<MatrixDouble>();
2065 auto common_ptr = HookeOps::commonDataFactory<SPACE_DIM, I, DomainEleOp>(
2066 mField, pip,
"U",
"MAT_ELASTIC", Sev::noisy);
2069 CHKERR VecZeroEntries(dJ_du);
2072 CHKERR VecAssemblyBegin(dJ_du);
2073 CHKERR VecAssemblyEnd(dJ_du);
2075 auto post_proc_rhs = get_essential_fe();
2076 post_proc_rhs->f = dJ_du;
2078 post_proc_rhs.get());
2084 auto calculate_adjoint_lambda = [&](
auto lambda,
auto dJ_du) {
2087 MOFEM_LOG(
"WORLD", Sev::inform) <<
"Solving for adjoint variable lambda";
2090 CHKERR VecGhostUpdateBegin(
lambda, INSERT_VALUES, SCATTER_FORWARD);
2091 CHKERR VecGhostUpdateEnd(
lambda, INSERT_VALUES, SCATTER_FORWARD);
2095 CHKERR VecZeroEntries(objective_function_gradient);
2096 CHKERR VecZeroEntries(adjoint_vector);
2098 CHKERR calculate_objective_value();
2100 <<
"Objective function: " << *objective_function_value;
2102 auto direct = [&]() {
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,
2116 CHKERR VecAXPY(adjoint_vector, objective_gradient, dm_diff_vec);
2125 auto adjoint = [&]() {
2128 CHKERR calculate_variance_of_objective_function_dJ_du(dJ_du);
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;
2145 const double dJ_dp =
2146 objective_gradient_explicit + lambda_dot_residual_variation;
2147 CHKERR VecSetValue(objective_function_gradient, mode, dJ_dp,
2149 CHKERR VecAXPY(adjoint_vector, dJ_dp, dm_diff_vec);
2152 CHKERR VecAssemblyBegin(objective_function_gradient);
2153 CHKERR VecAssemblyEnd(objective_function_gradient);
2159 MOFEM_LOG(
"WORLD", Sev::inform) <<
"Running Direct Sensitivity...";
2164 MOFEM_LOG(
"WORLD", Sev::inform) <<
"Running Adjoint Sensitivity...";
2170 "Wrong sensitivity type selected");
2173 CHKERR VecAssemblyBegin(objective_function_gradient);
2174 CHKERR VecAssemblyEnd(objective_function_gradient);
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);
2224 const char param_file[] =
"param_file.petsc";
2227 auto core_log = logging::core::get();
2239 DMType dm_name =
"DMMOFEM";
2241 DMType dm_name_mg =
"DMMOFEM_MG";
2246 moab::Core mb_instance;
2247 moab::Interface &moab = mb_instance;
2265 if (Py_FinalizeEx() < 0) {
2271 const std::string block_name,
int dim) {
2277 std::regex((boost::format(
"%s(.*)") % block_name).str())
2281 for (
auto bc : bcs) {
2285 "get meshset ents");
2289 for (
auto dd = dim - 1; dd >= 0; --dd) {
2293 moab::Interface::UNION),
2313 CHKERR moab.add_entities(*out_meshset, r);
2314 CHKERR moab.write_file(name.c_str(),
"VTK",
"", out_meshset->get_ptr(), 1);
#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)
auto diff_symmetrize(FTensor::Number< DIM >)
static char help[]
[calculateGradient]
PetscBool is_plane_strain
constexpr int SPACE_DIM
[Define dimension]
SensitivityMethod derivative_type
constexpr double poisson_ratio
Poisson's ratio ν
constexpr int BASE_DIM
[Constants and material properties]
constexpr double shear_modulus_G
Shear modulus G = E/(2(1+ν))
constexpr IntegrationType I
Use Gauss quadrature for integration.
constexpr double bulk_modulus_K
Bulk modulus K = E/(3(1-2ν))
FormsIntegrators< DomainEleOp >::Assembly< A >::OpBase DomainBaseOp
[Postprocess results]
constexpr AssemblyType A
[Define dimension]
Range get_range_from_block(MoFEM::Interface &m_field, const std::string block_name, int dim)
constexpr double young_modulus
[Material properties for linear elasticity]
ElementsAndOps< SPACE_DIM >::DomainEle DomainEle
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
#define CATCH_ERRORS
Catch errors.
FieldApproximationBase
approximation base
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
#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()
#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_DATA_INCONSISTENCY
#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
PetscErrorCode DMMoFEMSetIsPartitioned(DM dm, PetscBool is_partitioned)
PetscErrorCode DMMoFEMCreateSubDM(DM subdm, DM dm, const char problem_name[])
Must be called by user to set Sub DM MoFEM data structures.
PetscErrorCode DMMoFEMAddElement(DM dm, std::string fe_name)
add element to dm
PetscErrorCode DMMoFEMSetSquareProblem(DM dm, PetscBool square_problem)
set squared problem
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.
PetscErrorCode DMoFEMPostProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
PetscErrorCode DMMoFEMAddSubFieldRow(DM dm, const char field_name[])
PetscErrorCode DMRegister_MoFEM(const char sname[])
Register MoFEM problem.
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.
auto createDMVector(DM dm)
Get smart vector from DM.
PetscErrorCode DMMoFEMAddSubFieldCol(DM dm, const char field_name[])
auto createDMMatrix(DM dm)
Get smart matrix from DM.
PetscErrorCode DMoFEMPreProcessFiniteElements(DM dm, MoFEM::FEMethod *method)
execute finite element method for each element in dm (problem)
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.
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
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'J', DIM1 > J
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)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
auto createKSP(MPI_Comm comm)
PetscErrorCode DMMoFEMSetDestroyProblem(DM dm, PetscBool destroy_problem)
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data)
PetscErrorCode PetscOptionsGetInt(PetscOptions *, const char pre[], const char name[], PetscInt *ivalue, PetscBool *set)
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.
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.
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
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)]
PetscBool is_plane_strain
FTensor::Index< 'm', 3 > m
Boundary conditions marker.
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.
std::vector< SmartPetscObj< Vec > > modeVecs
Topology mode vectors (design variables)
FieldApproximationBase base
Choice of finite element basis functions
std::vector< std::array< double, 3 > > modeCentroids
Centroids of optimization blocks
MoFEMErrorCode topologyModes()
Compute topology optimization modes.
SmartPetscObj< Vec > initialGeometry
Initial geometry field.
int fieldOrder
Polynomial order for approximation.
Example(MoFEM::Interface &m_field)
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.
MoFEMErrorCode setupAdJoint()
Setup adjoint fields and finite elements
MoFEM::Interface & mField
Reference to MoFEM interface.
std::vector< std::array< double, 6 > > modeBBoxes
Bounding boxes of optimization blocks.
boost::shared_ptr< ObjectiveFunctionData > pythonPtr
Interface to Python objective function.
SmartPetscObj< DM > adjointDM
Data manager for adjoint problem.
MoFEMErrorCode setupProblem()
Setup fields, approximation spaces and DOFs.
MoFEMErrorCode postprocessElastic(int iter, SmartPetscObj< Vec > adjoint_vector=nullptr)
Post-process and output results.
MoFEMErrorCode solveElastic()
Solve forward elastic problem.
boost::shared_ptr< MatrixDouble > vectorFieldPtr
Field values at evaluation points.
Add operators pushing bases from local to physical configuration.
Boundary condition manager for finite element problem setup.
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
static MoFEMErrorCode Initialize(int *argc, char ***args, const char file[], const char help[])
Initializes the MoFEM database PETSc, MOAB and MPI.
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Deprecated interface functions.
Definition of the displacement bc data structure.
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.
Class (Function) to enforce essential constrains on the right hand side diagonal.
Class (Function) to enforce essential constrains.
MoFEMErrorCode fieldLambdaOnValues(OneFieldFunctionOnValues lambda, const std::string field_name, Range *ents_ptr=nullptr)
field lambda
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.
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.
MoFEMErrorCode getOptions()
get options
MoFEMErrorCode getDM(DM *dm)
Get DM.
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.
boost::shared_ptr< MatrixDouble > diffJac
boost::shared_ptr< MatrixDouble > jac
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)
boost::shared_ptr< HookeOps::CommonData > commPtr
boost::shared_ptr< VectorDouble > cofVals
Forward declaration of operator for gradient times symmetric tensor operations.
boost::shared_ptr< MatrixDouble > dGradPtr
boost::shared_ptr< double > globObjectiveGradPtr
boost::shared_ptr< ObjectiveFunctionData > pythonPtr
boost::shared_ptr< double > globObjectivePtr
ForcesAndSourcesCore::UserDataOperator OP
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Compute objective function contributions at element level.
boost::shared_ptr< HookeOps::CommonData > commPtr
boost::shared_ptr< VectorDouble > cofVals
boost::shared_ptr< MatrixDouble > uPtr
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)
boost::shared_ptr< MatrixDouble > diffJacPtr
boost::shared_ptr< MatrixDouble > jacPtr
boost::shared_ptr< MatrixDouble > dUPtr
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)
boost::shared_ptr< ObjectiveFunctionData > pythonPtr
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &data)
boost::shared_ptr< MatrixDouble > uPtr
boost::shared_ptr< HookeOps::CommonData > commPtr
PipelineManager::ElementsAndOpsByDim< 2 >::FaceSideEle SideEle
PipelineManager::ElementsAndOpsByDim< 3 >::FaceSideEle SideEle
#define EXECUTABLE_DIMENSION
PetscBool do_eval_field
Evaluate field.
ElementsAndOps< SPACE_DIM >::SideEle SideEle
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::BoundaryEle BoundaryEle
PostProcEleByDim< SPACE_DIM >::SideEle SideEle
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::DomainEle DomainEle