28 const std::string &root_path,
29 const bool non_conforming)
34 if (j_mesh[
"extract"].get<std::string>() !=
"volume")
42 if (j_mesh[
"advanced"][
"normalize_mesh"])
48 mesh->bounding_box(bbox[0], bbox[1]);
50 const std::string unit = j_mesh[
"unit"];
51 double unit_scale = 1;
60 j_mesh[
"transformation"],
61 (bbox[1] - bbox[0]).cwiseAbs().transpose(),
63 mesh->apply_affine_transformation(A, b);
66 mesh->bounding_box(bbox[0], bbox[1]);
70 const int n_refs = j_mesh[
"n_refs"];
71 const double refinement_location = j_mesh[
"advanced"][
"refinement_location"];
83 assert(mesh->n_elements() > 0);
84 const int uniform_value = mesh->get_body_id(0);
85 for (
int i = 1; i < mesh->n_elements(); ++i)
86 if (mesh->get_body_id(i) != uniform_value)
87 log_and_throw_error(
"Unable to apply stored nonuniform volume_selection because n_refs={} > 0!", n_refs);
89 logger().info(
"Performing global h-refinement with {} refinements", n_refs);
90 mesh->refine(n_refs, refinement_location);
91 mesh->set_body_ids(std::vector<int>(mesh->n_elements(), uniform_value));
96 if (j_mesh[
"advanced"][
"min_component"].get<int>() != -1)
128 if (j_mesh[
"advanced"][
"force_linear_geometry"].get<bool>())
142 const std::vector<std::shared_ptr<Selection>> node_selections =
145 if (!node_selections.empty())
147 mesh->compute_node_ids([&](
const size_t n_id,
const RowVectorNd &p,
bool is_boundary) {
151 const std::vector<int> tmp = {int(n_id)};
152 for (
const auto &selection : node_selections)
154 if (selection->inside(n_id, tmp, p))
155 return selection->id(n_id, tmp, p);
157 return std::numeric_limits<int>::max();
161 if (!j_mesh[
"curve_selection"].is_null())
166 std::vector<std::shared_ptr<Selection>> surface_selections =
169 if (!surface_selections.empty())
171 mesh->compute_boundary_ids([&](
const size_t p_id,
const std::vector<int> &vs,
const RowVectorNd &p,
bool is_boundary) {
175 for (
const auto &selection : surface_selections)
177 if (selection->inside(p_id, vs, p))
178 return selection->id(p_id, vs, p);
180 return std::numeric_limits<int>::max();
187 const json volume_selection = j_mesh[
"volume_selection"];
188 if (volume_selection.is_object()
189 && volume_selection.size() == 1
190 && volume_selection.contains(
"id_offset"))
192 const int id_offset = volume_selection[
"id_offset"].get<
int>();
195 const int n_body_ids = mesh->n_elements();
196 std::vector<int> body_ids(n_body_ids);
197 for (
int i = 0; i < n_body_ids; ++i)
198 body_ids[i] = mesh->get_body_id(i) + id_offset;
199 mesh->set_body_ids(body_ids);
205 std::vector<std::shared_ptr<Selection>> volume_selections =
209 if (mesh->has_body_ids())
210 volume_selections.push_back(std::make_shared<SpecifiedSelection>(mesh->get_body_ids()));
212 mesh->compute_body_ids([&](
const size_t cell_id,
const std::vector<int> &vs,
const RowVectorNd &p) ->
int {
213 for (
const auto &selection : volume_selections)
216 if (selection->inside(cell_id, vs, p))
217 return selection->id(cell_id, vs, p);
232 const json &geometry,
233 const std::string &root_path,
234 const std::vector<std::string> &_names,
235 const std::vector<Eigen::MatrixXd> &_vertices,
236 const std::vector<Eigen::MatrixXi> &_cells,
237 const bool non_conforming)
257 assert(_names.empty());
258 assert(_vertices.empty());
259 assert(_cells.empty());
263 if (geometry.empty())
270 std::unique_ptr<Mesh> mesh =
nullptr;
272 for (
const json &geometry : geometries)
274 if (!geometry[
"enabled"].get<bool>() || geometry[
"is_obstacle"].get<
bool>())
277 if (geometry[
"type"] !=
"mesh" && geometry[
"type"] !=
"mesh_array")
280 const std::unique_ptr<Mesh> tmp_mesh =
read_fem_mesh(units, geometry, root_path, non_conforming);
283 mesh = tmp_mesh->copy();
285 mesh->append(tmp_mesh);
287 if (geometry[
"type"] ==
"mesh_array")
290 tmp_mesh->bounding_box(bbox[0], bbox[1]);
292 const long dim = tmp_mesh->dimension();
293 const bool is_offset_relative = geometry[
"array"][
"relative"];
294 const double offset = geometry[
"array"][
"offset"];
295 const VectorNd dimensions = (bbox[1] - bbox[0]);
296 const VectorNi size = geometry[
"array"][
"size"];
298 for (
int i = 0; i < size[0]; ++i)
300 for (
int j = 0; j < size[1]; ++j)
302 for (
int k = 0; k < (size.size() > 2 ? size[2] : 1); ++k)
304 if (i == 0 && j == 0 && k == 0)
307 RowVectorNd translation = offset * Eigen::RowVector3d(i, j, k).head(dim);
308 if (is_offset_relative)
309 translation.array() *= dimensions.array();
311 const std::unique_ptr<Mesh> copy_mesh = tmp_mesh->copy();
312 copy_mesh->apply_affine_transformation(MatrixNd::Identity(dim, dim), translation);
313 mesh->append(copy_mesh);
330 const std::string &root_path,
332 Eigen::MatrixXd &vertices,
333 Eigen::VectorXi &codim_vertices,
334 Eigen::MatrixXi &codim_edges,
335 Eigen::MatrixXi &
faces)
340 const std::string mesh_path =
resolve_path(j_mesh[
"mesh"], root_path);
343 mesh_path, vertices, codim_vertices, codim_edges,
faces);
347 throw std::runtime_error(fmt::format(
"Unable to read mesh: {}", mesh_path));
349 const int prev_dim = vertices.cols();
350 vertices.conservativeResize(vertices.rows(), dim);
352 vertices.rightCols(dim - prev_dim).setZero();
357 const std::string unit = j_mesh[
"unit"];
358 double unit_scale = 1;
362 const VectorNd mesh_dimensions = (vertices.colwise().maxCoeff() - vertices.colwise().minCoeff()).cwiseAbs();
366 vertices = vertices * A.transpose();
367 vertices.rowwise() += b.transpose();
370 std::string extract = j_mesh[
"extract"];
372 if (extract ==
"volume")
375 if (extract ==
"points")
378 codim_edges.resize(0, 0);
380 codim_vertices.resize(vertices.rows());
381 for (
int i = 0; i < codim_vertices.size(); ++i)
382 codim_vertices[i] = i;
384 else if (extract ==
"edges" &&
faces.size() != 0)
387 Eigen::MatrixXi edges;
388 igl::edges(
faces, edges);
390 codim_edges.conservativeResize(codim_edges.rows() + edges.rows(), 2);
391 codim_edges.bottomRows(edges.rows()) = edges;
393 else if (extract ==
"surface" && dim == 2 &&
faces.size() != 0)
396 Eigen::MatrixXi boundary_edges;
397 igl::boundary_facets(
faces, boundary_edges);
398 codim_edges.conservativeResize(codim_edges.rows() + boundary_edges.rows(), 2);
399 codim_edges.bottomRows(boundary_edges.rows()) = boundary_edges;
404 else if (extract ==
"volume")
410 if (j_mesh[
"n_refs"].get<int>() != 0)
412 if (
faces.size() != 0)
415 const int n_refs = j_mesh[
"n_refs"];
416 const double refinement_location = j_mesh[
"advanced"][
"refinement_location"];
417 for (
int i = 0; i < n_refs; i++)
419 const size_t n_vertices = vertices.rows();
420 const size_t n_edges = codim_edges.rows();
421 vertices.conservativeResize(n_vertices + n_edges, vertices.cols());
422 codim_edges.conservativeResize(2 * n_edges, codim_edges.cols());
423 for (
size_t ei = 0; ei < n_edges; ei++)
425 const int v0i = codim_edges(ei, 0);
426 const int v1i = codim_edges(ei, 1);
427 const int v2i = n_vertices + ei;
428 vertices.row(v2i) = (vertices.row(v1i) - vertices.row(v0i)) * refinement_location + vertices.row(v0i);
429 codim_edges.row(ei) << v0i, v2i;
430 codim_edges.row(n_edges + ei) << v2i, v1i;
440 const json &geometry,
441 const std::vector<json> &displacements,
442 const std::vector<json> &dirichlets,
443 const std::string &root_path,
445 const std::vector<std::string> &_names,
446 const std::vector<Eigen::MatrixXd> &_vertices,
447 const std::vector<Eigen::MatrixXi> &_cells,
448 const bool non_conforming)
468 assert(_names.empty());
469 assert(_vertices.empty());
470 assert(_cells.empty());
474 if (geometry.empty())
479 for (
const json &geometry : geometries)
482 if (!geometry[
"is_obstacle"].get<bool>())
485 if (!geometry[
"enabled"].get<bool>())
488 if (geometry[
"type"] ==
"mesh" || geometry[
"type"] ==
"mesh_array")
490 Eigen::MatrixXd vertices;
491 Eigen::VectorXi codim_vertices;
492 Eigen::MatrixXi codim_edges;
493 Eigen::MatrixXi
faces;
495 geometry, root_path, dim, vertices, codim_vertices,
498 if (geometry[
"type"] ==
"mesh_array")
500 const Selection::BBox bbox{{vertices.colwise().minCoeff(), vertices.colwise().maxCoeff()}};
502 const bool is_offset_relative = geometry[
"array"][
"relative"];
503 const double offset = geometry[
"array"][
"offset"];
504 const VectorNd dimensions = (bbox[1] - bbox[0]);
505 const VectorNi size = geometry[
"array"][
"size"];
507 const int N = size.head(dim).prod();
508 const int nV = vertices.rows(), nCV = codim_vertices.rows(), nCE = codim_edges.rows(), nF =
faces.rows();
510 vertices.conservativeResize(N * nV, Eigen::NoChange);
511 codim_vertices.conservativeResize(N * nCV, Eigen::NoChange);
512 codim_edges.conservativeResize(N * nCE, Eigen::NoChange);
513 faces.conservativeResize(N * nF, Eigen::NoChange);
515 for (
int i = 0; i < size[0]; ++i)
517 for (
int j = 0; j < size[1]; ++j)
519 for (
int k = 0; k < (size.size() > 2 ? size[2] : 1); ++k)
521 RowVectorNd translation = offset * Eigen::RowVector3d(i, j, k).head(vertices.cols());
522 if (is_offset_relative)
523 translation.array() *= dimensions.array();
525 int n = i * size[1] + j;
531 vertices.middleRows(n * nV, nV) = vertices.topRows(nV).rowwise() + translation;
533 codim_vertices.segment(n * nV, nV) = codim_vertices.head(nV).array() + n * nV;
535 codim_edges.middleRows(n * nCE, nCE) = codim_edges.topRows(nCE).array() + n * nV;
537 faces.middleRows(n * nF, nF) =
faces.topRows(nF).array() + n * nV;
543 json displacement =
"{\"value\":[0, 0, 0]}"_json;
546 if (!geometry[
"surface_selection"].is_number())
549 const int id = geometry[
"surface_selection"];
550 for (
const json &disp : dirichlets)
552 if ((disp[
"id"].is_string() && disp[
"id"].get<std::string>() ==
"all")
553 || (disp[
"id"].is_number_integer() && disp[
"id"].get<
int>() ==
id))
558 else if (disp[
"id"].is_array())
560 for (
const json &disp_id : disp[
"id"])
562 assert(disp_id.is_number_integer());
563 if (disp_id.get<
int>() ==
id)
571 for (
const json &disp : displacements)
573 if ((disp[
"id"].is_string() && disp[
"id"].get<std::string>() ==
"all")
574 || (disp[
"id"].is_number_integer() && disp[
"id"].get<
int>() ==
id))
579 else if (disp[
"id"].is_array())
581 for (
const json &disp_id : disp[
"id"])
583 assert(disp_id.is_number_integer());
584 if (disp_id.get<
int>() ==
id)
595 vertices, codim_vertices, codim_edges,
faces, displacement);
597 else if (geometry[
"type"] ==
"plane")
599 obstacle.
append_plane(geometry[
"point"], geometry[
"normal"]);
601 else if (geometry[
"type"] ==
"ground")
603 VectorNd gravity = VectorNd::Zero(dim);
605 const double height = geometry[
"height"];
606 assert(gravity.norm() != 0);
607 const VectorNd normal = -gravity.normalized();
608 const VectorNd point = height * normal;
611 else if (geometry[
"type"] ==
"mesh_sequence")
613 namespace fs = std::filesystem;
614 std::vector<fs::path> mesh_files;
615 if (geometry[
"mesh_sequence"].is_array())
617 mesh_files = geometry[
"mesh_sequence"].get<std::vector<fs::path>>();
621 assert(geometry[
"mesh_sequence"].is_string());
622 const fs::path meshes(
resolve_path(geometry[
"mesh_sequence"], root_path));
624 if (fs::is_directory(meshes))
626 for (
const auto &entry : std::filesystem::directory_iterator(meshes))
628 if (entry.is_regular_file())
629 mesh_files.push_back(entry.path());
634 mesh_files = glob::rglob(meshes.string());
637 std::sort(mesh_files.begin(), mesh_files.end(), [](
const fs::path &p1,
const fs::path &p2) {
638 return strnatcmp(p1.string().c_str(), p2.string().c_str()) < 0;
642 std::vector<Eigen::MatrixXd> vertices(mesh_files.size());
643 Eigen::VectorXi codim_vertices;
644 Eigen::MatrixXi codim_edges;
645 Eigen::MatrixXi
faces;
647 for (
int i = 0; i < mesh_files.size(); ++i)
649 json jmesh = geometry;
650 jmesh[
"mesh"] = mesh_files[i];
653 Eigen::VectorXi tmp_codim_vertices;
654 Eigen::MatrixXi tmp_codim_edges;
655 Eigen::MatrixXi tmp_faces;
657 jmesh, root_path, dim, vertices[i],
658 tmp_codim_vertices, tmp_codim_edges, tmp_faces);
661 codim_vertices = tmp_codim_vertices;
662 codim_edges = tmp_codim_edges;
667 assert((codim_vertices.array() == tmp_codim_vertices.array()).all());
668 assert((codim_edges.array() == tmp_codim_edges.array()).all());
669 assert((
faces.array() == tmp_faces.array()).all());
674 vertices, codim_vertices, codim_edges,
faces, geometry[
"fps"]);
689 const double unit_scale,
690 const json &transform,
695 const int dim = mesh_dimensions.size();
702 if (transform[
"dimensions"].is_array())
705 (mesh_dimensions.array() == 0).select(1, mesh_dimensions);
707 scale = transform[
"dimensions"];
708 const int scale_size = scale.size();
709 scale.conservativeResize(dim);
710 if (scale_size < dim)
711 scale.tail(dim - scale_size).setZero();
713 scale.array() /= modified_dimensions.array();
715 else if (transform[
"scale"].is_number())
717 scale.setConstant(dim, transform[
"scale"].get<double>());
721 assert(transform[
"scale"].is_array());
722 scale = transform[
"scale"];
723 const int scale_size = scale.size();
724 scale.conservativeResize(dim);
725 if (scale_size < dim)
726 scale.tail(dim - scale_size).setZero();
732 A = (unit_scale * scale).asDiagonal();
740 MatrixNd R = MatrixNd::Identity(dim, dim);
741 if (!transform[
"rotation"].is_null())
745 if (transform[
"rotation"].is_number())
746 R = Eigen::Rotation2Dd(
deg2rad(transform[
"rotation"].get<double>()))
748 else if (!transform[
"rotation"].is_array() || !transform[
"rotation"].empty())
763 b = transform[
"translation"];
764 const int translation_size = b.size();
765 b.conservativeResize(dim);
766 if (translation_size < dim)
767 b.tail(dim - translation_size).setZero();
Obstacle read_obstacle_geometry(const Units &units, const json &geometry, const std::vector< json > &displacements, const std::vector< json > &dirichlets, const std::string &root_path, const int dim, const std::vector< std::string > &_names, const std::vector< Eigen::MatrixXd > &_vertices, const std::vector< Eigen::MatrixXi > &_cells, const bool non_conforming)
read a FEM mesh from a geometry JSON