v0.14.0
prism_elements_from_surface.cpp

Adding prims on the surface and checking conformity between quads, triangles and prism

/** \file prism_elements_from_surface.cpp
\example prism_elements_from_surface.cpp
\brief Adding prims on the surface and checking conformity between quads,
triangles and prism
*/
#include <MoFEM.hpp>
using namespace MoFEM;
static char help[] = "...\n\n";
static constexpr int precision_exponent = 5;
static constexpr int number_of_prisms_layers = 18;
static constexpr double delta =
1. / static_cast<double>(number_of_prisms_layers);
static constexpr std::array<double, 3> d3 = {0, 0, 0};
static constexpr std::array<double, 3> d4 = {0, 0, delta};
inline static double getArg(double x) {
return std::round(x * pow(10., precision_exponent - 1));
};
int x, y, z;
CoordsAndHandle(const double *coords, EntityHandle v)
: x(getArg(coords[0])), y(getArg(coords[1])), z(getArg(coords[2])),
node(v) {}
};
typedef multi_index_container<
indexed_by<
hashed_unique<composite_key<
CoordsAndHandle, member<CoordsAndHandle, int, &CoordsAndHandle::x>,
member<CoordsAndHandle, int, &CoordsAndHandle::y>,
member<CoordsAndHandle, int, &CoordsAndHandle::z>>>
>>
PrismOp(moab::Interface &post_proc, MapCoords &map_coords);
MoFEMErrorCode doWork(int side, EntityType type,
public:
moab::Interface &postProc;
MapCoords &mapCoords;
std::vector<EntityHandle> nodeHandles;
};
PrismFE(MoFEM::Interface &m_field, MatrixDouble &tri_coords);
int getRuleTrianglesOnly(int order);
int getRuleThroughThickness(int order);
MoFEMErrorCode setGaussPtsTrianglesOnly(int order_triangles_only);
MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness);
private:
MatrixDouble &triCoords;
EntityHandle prims;
};
template <typename OP> struct Op : public OP {
Op(moab::Interface &post_proc, MapCoords &map_coords, EntityHandle prism);
MoFEMErrorCode doWork(int side, EntityType type,
public:
moab::Interface &postProc;
MapCoords &mapCoords;
EntityHandle prism;
std::vector<EntityHandle> nodeHandles;
};
TriFE(MoFEM::Interface &m_field, MatrixDouble &tri_coords,
EntityHandle prims);
int getRule(int order_row, int order_col, int order_data);
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data);
private:
MatrixDouble &triCoords;
EntityHandle prism;
};
QuadFE(MoFEM::Interface &m_field, std::array<Range, 3> &edges_blocks,
EntityHandle prims);
int getRule(int order_row, int order_col, int order_data);
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data);
private:
std::array<Range, 3> &edgeBlocks;
EntityHandle prism;
};
struct EdgeFE : public EdgeEle {
EdgeFE(MoFEM::Interface &m_field, std::array<Range, 3> &edges_blocks,
EntityHandle prims);
int getRule(int order_row, int order_col, int order_data);
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data);
private:
std::array<Range, 3> &edgeBlocks;
EntityHandle prism;
};
int main(int argc, char *argv[]) {
MoFEM::Core::Initialize(&argc, &argv, (char *)0, help);
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
int rank;
MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
// Read parameters from line command
PetscBool flg = PETSC_TRUE;
char mesh_file_name[255];
CHKERR PetscOptionsGetString(PETSC_NULL, "", "-my_file", mesh_file_name,
255, &flg);
if (flg != PETSC_TRUE)
SETERRQ(PETSC_COMM_SELF, MOFEM_INVALID_DATA,
"error -my_file (MESH FILE NEEDED)");
// Read mesh to MOAB
const char *option;
option = ""; //"PARALLEL=BCAST;";//;DEBUG_IO";
CHKERR moab.load_file(mesh_file_name, 0, option);
Range tris;
CHKERR moab.get_entities_by_type(0, MBTRI, tris, false);
Range tris_verts;
CHKERR moab.get_connectivity(tris, tris_verts);
MatrixDouble tri_coords(tris_verts.size(), 3);
CHKERR moab.get_coords(tris_verts, &tri_coords(0, 0));
MoFEM::Core core(moab);
MoFEM::Interface &m_field = core;
std::array<Range, 3> edge_block;
for (auto b : {1, 2, 3})
b, BLOCKSET, 1, edge_block[b - 1]);
PrismsFromSurfaceInterface *prisms_from_surface_interface;
CHKERR m_field.getInterface(prisms_from_surface_interface);
Range prisms;
CHKERR prisms_from_surface_interface->createPrisms(
prisms_from_surface_interface->setThickness(prisms, d3.data(), d4.data());
Range add_prims_layer;
Range extrude_prisms = prisms;
for (int ll = 1; ll != number_of_prisms_layers; ++ll) {
prisms_from_surface_interface->createdVertices.clear();
CHKERR prisms_from_surface_interface->createPrismsFromPrisms(
extrude_prisms, false, add_prims_layer);
prisms_from_surface_interface->setThickness(add_prims_layer, d3.data(),
d4.data());
extrude_prisms = add_prims_layer;
prisms.merge(add_prims_layer);
add_prims_layer.clear();
}
EntityHandle meshset;
CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER, meshset);
CHKERR moab.add_entities(meshset, prisms);
MapCoords map_coords;
Range verts;
CHKERR moab.get_connectivity(prisms, verts);
MatrixDouble coords(verts.size(), 3);
CHKERR moab.get_coords(verts, &coords(0, 0));
for (size_t v = 0; v != verts.size(); ++v)
map_coords.insert(CoordsAndHandle(&coords(v, 0), verts[v]));
EntityHandle one_prism_meshset;
CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER,
one_prism_meshset);
std::array<double, 18> one_prism_coords = {0, 0, 0, 1, 0, 0, 0, 1, 0,
0, 0, 1, 1, 0, 1, 0, 1, 1};
std::array<EntityHandle, 6> one_prism_nodes;
for (int n = 0; n != 6; ++n)
CHKERR moab.create_vertex(&one_prism_coords[3 * n], one_prism_nodes[n]);
EntityHandle one_prism;
CHKERR m_field.get_moab().create_element(MBPRISM, one_prism_nodes.data(), 6,
one_prism);
Range one_prism_range;
one_prism_range.insert(one_prism);
CHKERR moab.add_entities(one_prism_meshset, one_prism_range);
Range one_prism_verts;
CHKERR moab.get_connectivity(one_prism_range, one_prism_verts);
CHKERR moab.add_entities(one_prism_meshset, one_prism_verts);
Range one_prism_adj_ents;
for (int d = 1; d != 3; ++d)
CHKERR moab.get_adjacencies(one_prism_range, d, true, one_prism_adj_ents,
moab::Interface::UNION);
CHKERR moab.add_entities(one_prism_meshset, one_prism_adj_ents);
BitRefLevel bit_level0;
bit_level0.set(0);
CHKERR m_field.getInterface<BitRefManager>()->setEntitiesBitRefLevel(
one_prism_range, bit_level0);
CHKERR prisms_from_surface_interface->seedPrismsEntities(one_prism_range,
bit_level0);
// Fields
CHKERR m_field.add_field("FIELD1", H1, AINSWORTH_LEGENDRE_BASE, 1);
CHKERR m_field.add_ents_to_field_by_type(one_prism_meshset, MBPRISM,
"FIELD1", VERY_NOISY);
CHKERR m_field.set_field_order(one_prism_meshset, MBVERTEX, "FIELD1", 1);
CHKERR m_field.set_field_order(one_prism_meshset, MBEDGE, "FIELD1", 5,
CHKERR m_field.set_field_order(one_prism_meshset, MBTRI, "FIELD1", 5);
CHKERR m_field.set_field_order(one_prism_meshset, MBQUAD, "FIELD1", 5,
CHKERR m_field.set_field_order(one_prism_meshset, MBPRISM, "FIELD1", 7,
// FE
CHKERR m_field.add_finite_element("PRISM");
CHKERR m_field.add_finite_element("EDGE");
CHKERR m_field.add_finite_element("TRI");
CHKERR m_field.add_finite_element("QUAD");
// Define rows/cols and element data
CHKERR m_field.modify_finite_element_add_field_row("PRISM", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_col("PRISM", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_data("PRISM", "FIELD1");
// Define rows/cols and element data
CHKERR m_field.modify_finite_element_add_field_row("EDGE", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_col("EDGE", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_data("EDGE", "FIELD1");
// Define rows/cols and element data
CHKERR m_field.modify_finite_element_add_field_row("TRI", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_col("TRI", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_data("TRI", "FIELD1");
// Define rows/cols and element data
CHKERR m_field.modify_finite_element_add_field_row("QUAD", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_col("QUAD", "FIELD1");
CHKERR m_field.modify_finite_element_add_field_data("QUAD", "FIELD1");
CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBEDGE,
"EDGE");
CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBTRI,
"TRI");
CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset, MBQUAD,
"QUAD");
CHKERR m_field.add_ents_to_finite_element_by_type(one_prism_meshset,
MBPRISM, "PRISM");
// build finite elemnts
// //build adjacencies
CHKERR m_field.build_adjacencies(bit_level0);
// Problem
CHKERR m_field.add_problem("TEST_PROBLEM");
// set finite elements for problem
CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "PRISM");
CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "EDGE");
CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "TRI");
CHKERR m_field.modify_problem_add_finite_element("TEST_PROBLEM", "QUAD");
// set refinement level for problem
CHKERR m_field.modify_problem_ref_level_add_bit("TEST_PROBLEM", bit_level0);
// build problem
ProblemsManager *prb_mng_ptr;
CHKERR m_field.getInterface(prb_mng_ptr);
CHKERR prb_mng_ptr->buildProblem("TEST_PROBLEM", true);
// partition
CHKERR prb_mng_ptr->partitionSimpleProblem("TEST_PROBLEM");
CHKERR prb_mng_ptr->partitionFiniteElements("TEST_PROBLEM");
// what are ghost nodes, see Petsc Manual
CHKERR prb_mng_ptr->partitionGhostDofs("TEST_PROBLEM");
PrismFE fe_prism(m_field, tri_coords);
fe_prism.getOpPtrVector().push_back(new PrismOp(moab, map_coords));
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "PRISM", fe_prism);
EdgeFE fe_edge(m_field, edge_block, one_prism);
fe_edge.getOpPtrVector().push_back(
moab, map_coords, one_prism));
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "EDGE", fe_edge);
TriFE fe_tri(m_field, tri_coords, one_prism);
fe_tri.getOpPtrVector().push_back(
moab, map_coords, one_prism));
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "TRI", fe_tri);
QuadFE fe_quad(m_field, edge_block, one_prism);
fe_quad.getOpPtrVector().push_back(
moab, map_coords, one_prism));
CHKERR m_field.loop_finite_elements("TEST_PROBLEM", "QUAD", fe_quad);
CHKERR moab.write_file("prism_mesh.vtk", "VTK", "", &meshset, 1);
CHKERR moab.write_file("one_prism_mesh.vtk", "VTK", "", &one_prism_meshset,
1);
}
return 0;
}
PrismOp::PrismOp(moab::Interface &post_proc, MapCoords &map_coords)
"FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROW),
postProc(post_proc), mapCoords(map_coords) {}
MoFEMErrorCode PrismOp::doWork(int side, EntityType type,
constexpr double def_val[] = {0, 0, 0};
switch (type) {
case MBVERTEX:
case MBEDGE:
case MBTRI:
case MBQUAD:
case MBPRISM:
break;
default:
}
if (type == MBTRI && (side != 3 && side != 4))
if (type == MBQUAD && (side == 3 || side == 4))
const int nb_dofs = data.getIndices().size();
for (int dd = 0; dd != nb_dofs; ++dd)
data.getFieldDofs()[dd]->getFieldData() = data.getIndices()[dd];
if (type == MBVERTEX) {
const size_t nb_gauss_pts = getGaussPts().size2();
auto &coords_at_pts = getGaussPts();
nodeHandles.reserve(nb_gauss_pts);
nodeHandles.clear();
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
auto t = boost::make_tuple(CoordsAndHandle::getArg(coords_at_pts(0, gg)),
CoordsAndHandle::getArg(coords_at_pts(1, gg)),
CoordsAndHandle::getArg(coords_at_pts(2, gg)));
auto it = mapCoords.find(t);
if (it != mapCoords.end())
nodeHandles.emplace_back(it->node);
else
SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Vertex not found");
}
MatrixDouble node_coords(nb_gauss_pts, 3);
CHKERR postProc.get_coords(&nodeHandles[0], nodeHandles.size(),
&node_coords(0, 0));
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
constexpr double eps = 1e-12;
for (auto d : {0, 1, 2})
if (std::abs(node_coords(gg, d) - getGaussPts()(d, gg)) > eps) {
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistency between node coords and integration "
"points");
}
for (auto d : {0, 1, 2})
if (std::abs(node_coords(gg, d) - getCoordsAtGaussPts()(gg, d)) > eps)
SETERRQ(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistency between node coords and integration "
"points");
}
Tag th;
CHKERR postProc.tag_get_handle("Coords", 3, MB_TYPE_DOUBLE, th,
MB_TAG_CREAT | MB_TAG_DENSE, def_val);
CHKERR postProc.tag_set_data(th, &nodeHandles[0], nodeHandles.size(),
}
auto to_str = [](auto i) { return boost::lexical_cast<std::string>(i); };
std::string tag_name_base =
"PrismType" + to_str(type) + "Side" + to_str(side);
std::cout << tag_name_base << endl;
MatrixDouble trans_base = trans(data.getN());
if (trans_base.size2() != nodeHandles.size())
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "wrong size %d != %d",
trans_base.size2(), nodeHandles.size());
for (size_t rr = 0; rr != trans_base.size1(); ++rr) {
auto tag_name = tag_name_base + "Base" + to_str(rr);
Tag th;
CHKERR postProc.tag_get_handle(tag_name.c_str(), 1, MB_TYPE_DOUBLE, th,
MB_TAG_CREAT | MB_TAG_DENSE, def_val);
CHKERR postProc.tag_set_data(th, &nodeHandles[0], nodeHandles.size(),
&trans_base(rr, 0));
}
}
: FatPrismElementForcesAndSourcesCore(m_field), triCoords(tri_coords) {}
int PrismFE::getRuleTrianglesOnly(int order) { return -1; };
int PrismFE::getRuleThroughThickness(int order) { return -1; };
gaussPtsTrianglesOnly.resize(3, triCoords.size1(), false);
for (int gg = 0; gg != triCoords.size1(); ++gg)
for (int dd = 0; dd != 2; ++dd)
}
for (int gg = 0; gg != number_of_prisms_layers + 1; ++gg)
}
EntityHandle prism)
: FaceElementForcesAndSourcesCore(m_field), triCoords(tri_coords),
prism(prism) {}
int TriFE::getRule(int order_row, int order_col, int order_data) { return -1; }
MoFEMErrorCode TriFE::setGaussPts(int order_row, int order_col,
int order_data) {
int side, sense, offset;
CHKERR mField.get_moab().side_number(prism, ent, side, sense, offset);
std::array<int, 2> swap = {0, 1};
if (side == 3)
swap = std::array<int, 2>{1, 0};
gaussPts.resize(3, triCoords.size1(), false);
gaussPts.clear();
for (int gg = 0; gg != triCoords.size1(); ++gg)
for (int dd = 0; dd != 2; ++dd)
gaussPts(dd, gg) = triCoords(gg, swap[dd]);
}
QuadFE::QuadFE(MoFEM::Interface &m_field, std::array<Range, 3> &edge_blocks,
EntityHandle prism)
: FaceElementForcesAndSourcesCore(m_field), edgeBlocks(edge_blocks),
prism(prism) {}
int QuadFE::getRule(int order_row, int order_col, int order_data) { return -1; }
MoFEMErrorCode QuadFE::setGaussPts(int order_row, int order_col,
int order_data) {
int side, sense, offset;
CHKERR mField.get_moab().side_number(prism, ent, side, sense, offset);
Range edge_verts;
CHKERR mField.get_moab().get_connectivity(edgeBlocks[side], edge_verts);
MatrixDouble edge_coords(edge_verts.size(), 3);
CHKERR mField.get_moab().get_coords(edge_verts, &edge_coords(0, 0));
constexpr double normal[3][2] = {{-1, 0}, {-0.5, 0.5}, {0, 1}};
constexpr double origin[3][2] = {{1, 0}, {1, 0}, {0, 0}};
constexpr int swap[3][2] = {{0, 1}, {0, 1}, {1, 0}};
gaussPts.resize(3, edge_verts.size() * (number_of_prisms_layers + 1), false);
int gg = 0;
for (size_t rr = 0; rr != edge_verts.size(); ++rr) {
const double x = edge_coords(rr, 0) - origin[side][0];
const double y = edge_coords(rr, 1) - origin[side][1];
const double edge_dist = x * normal[side][0] + y * normal[side][1];
for (size_t cc = 0; cc != number_of_prisms_layers + 1; ++cc, ++gg) {
gaussPts(swap[side][0], gg) = edge_dist;
gaussPts(swap[side][1], gg) = delta * cc;
}
}
}
EdgeFE::EdgeFE(MoFEM::Interface &m_field, std::array<Range, 3> &edge_blocks,
EntityHandle prism)
: EdgeEle(m_field), edgeBlocks(edge_blocks),
prism(prism) {}
int EdgeFE::getRule(int order_row, int order_col, int order_data) { return -1; }
MoFEMErrorCode EdgeFE::setGaussPts(int order_row, int order_col,
int order_data) {
int side, sense, offset;
CHKERR mField.get_moab().side_number(prism, ent, side, sense, offset);
if (side >= 3 && side <= 5) {
gaussPts.resize(2, number_of_prisms_layers + 1, false);
for (size_t gg = 0; gg != number_of_prisms_layers + 1; ++gg)
gaussPts(0, gg) = delta * gg;
// cerr << gaussPts << endl;
} else {
constexpr double normal[3][2] = {{1, 0}, {-0.5, 0.5}, {0, 1}};
constexpr double origin[3][2] = {{0, 1}, {1, 0}, {0, 0}};
constexpr int side_map[9] = {0, 1, 2, -1, -1, -1, 0, 1, 2};
int num_nodes;
const EntityHandle *conn;
CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
MatrixDouble coords(num_nodes, 3);
CHKERR mField.get_moab().get_coords(conn, num_nodes,
&*coords.data().begin());
side = side_map[side];
Range edge_verts;
CHKERR mField.get_moab().get_connectivity(edgeBlocks[side], edge_verts);
MatrixDouble edge_coords(edge_verts.size(), 3);
CHKERR mField.get_moab().get_coords(edge_verts, &edge_coords(0, 0));
gaussPts.resize(2, edge_verts.size(), false);
for (size_t gg = 0; gg != edge_verts.size(); ++gg) {
const double x = edge_coords(gg, 0) - origin[side][0];
const double y = edge_coords(gg, 1) - origin[side][1];
const double edge_dist = x * normal[side][0] + y * normal[side][1];
gaussPts(0, gg) = edge_dist;
}
}
}
template <typename OP>
Op<OP>::Op(moab::Interface &post_proc, MapCoords &map_coords,
EntityHandle prism)
: OP("FIELD1", "FIELD1", ForcesAndSourcesCore::UserDataOperator::OPROW),
postProc(post_proc), mapCoords(map_coords), prism(prism) {}
template <typename OP>
MoFEMErrorCode Op<OP>::doWork(int side, EntityType type,
constexpr double def_val[] = {0, 0, 0};
switch (type) {
case MBVERTEX:
case MBEDGE:
case MBTRI:
case MBQUAD:
break;
default:
}
const int nb_dofs = data.getIndices().size();
for (int dd = 0; dd != nb_dofs; ++dd)
if (data.getFieldData()[dd] != data.getIndices()[dd]) {
std::cerr << "Indices: " << data.getIndices() << std::endl;
std::cerr << "Local indices: " << data.getLocalIndices() << std::endl;
std::cerr << "Data: " << data.getFieldData() << std::endl;
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Indicices/data inconsistency %3.1f != %d",
data.getFieldData()[dd], data.getIndices()[dd]);
}
const EntityHandle fe_ent = OP::getFEEntityHandle();
const EntityHandle ent = OP::getSideEntity(side, type);
int side_prism, sense, offset;
if (type == MBVERTEX) {
CHKERR postProc.side_number(prism, fe_ent, side_prism, sense, offset);
} else
CHKERR postProc.side_number(prism, ent, side_prism, sense, offset);
if (type == MBVERTEX) {
auto &coords_at_pts = OP::getCoordsAtGaussPts();
const size_t nb_gauss_pts = coords_at_pts.size1();
nodeHandles.reserve(nb_gauss_pts);
nodeHandles.clear();
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
auto t = boost::make_tuple(CoordsAndHandle::getArg(coords_at_pts(gg, 0)),
CoordsAndHandle::getArg(coords_at_pts(gg, 1)),
CoordsAndHandle::getArg(coords_at_pts(gg, 2)));
auto it = mapCoords.find(t);
if (it != mapCoords.end())
nodeHandles.emplace_back(it->node);
else
SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Vertex not found");
}
}
auto to_str = [](auto i) { return boost::lexical_cast<std::string>(i); };
std::string tag_name_base =
"FEType" + to_str(OP::getNumeredEntFiniteElementPtr()->getEntType()) +
"Type" + to_str(type) + "Side" + to_str(side_prism);
std::string tag_prism_name_base =
"PrismType" + to_str(type) + "Side" + to_str(side_prism);
MatrixDouble trans_base = trans(data.getN());
MatrixDouble prism_base(trans_base.size1(), trans_base.size2());
if (trans_base.size2() != nodeHandles.size())
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID, "wrong size %d != %d",
trans_base.size2(), nodeHandles.size());
for (size_t rr = 0; rr != trans_base.size1(); ++rr) {
std::string tag_name = tag_name_base + "Base";
if (type == MBVERTEX) {
EntityHandle node = data.getFieldDofs()[rr]->getEnt();
int prism_mode_side;
CHKERR postProc.side_number(prism, node, prism_mode_side, sense, offset);
tag_name += to_str(prism_mode_side);
} else {
tag_name += to_str(rr);
}
std::cout << tag_name << endl;
Tag th;
CHKERR postProc.tag_get_handle(tag_name.c_str(), 1, MB_TYPE_DOUBLE, th,
MB_TAG_CREAT | MB_TAG_DENSE, def_val);
CHKERR postProc.tag_set_data(th, &nodeHandles[0], nodeHandles.size(),
&trans_base(rr, 0));
if (type != MBVERTEX) {
auto tag_prism_name = tag_prism_name_base + "Base" + to_str(rr);
Tag th_prism;
CHKERR postProc.tag_get_handle(tag_prism_name.c_str(), th_prism);
CHKERR postProc.tag_get_data(th_prism, &nodeHandles[0],
nodeHandles.size(), &prism_base(rr, 0));
}
}
auto sum_matrix = [](MatrixDouble &m) {
double s = 0;
for (unsigned int ii = 0; ii < m.size1(); ii++) {
for (unsigned int jj = 0; jj < m.size2(); jj++) {
s += std::abs(m(ii, jj));
}
}
return s;
};
if (type != MBVERTEX) {
prism_base -= trans_base;
double sum = sum_matrix(prism_base);
constexpr double eps = 1e-6;
if (std::abs(sum) > eps)
cout << "Inconsistent base " << tag_prism_name_base << " "
<< tag_name_base << " sum " << sum << endl;
if (std::abs(sum) > eps)
SETERRQ2(PETSC_COMM_SELF, MOFEM_ATOM_TEST_INVALID,
"Inconsistent base %s sum %6.4e", tag_prism_name_base.c_str(),
sum);
}
}
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:447
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Definition: UnknownInterface.hpp:93
main
int main(int argc, char *argv[])
Definition: prism_elements_from_surface.cpp:132
MoFEM::EdgeElementForcesAndSourcesCore
Edge finite element.
Definition: EdgeElementForcesAndSourcesCore.hpp:30
EdgeFE::prism
EntityHandle prism
Definition: prism_elements_from_surface.cpp:129
QuadFE::edgeBlocks
std::array< Range, 3 > & edgeBlocks
Definition: prism_elements_from_surface.cpp:115
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:127
MoFEM::CoreInterface::loop_finite_elements
virtual MoFEMErrorCode loop_finite_elements(const std::string problem_name, const std::string &fe_name, FEMethod &method, boost::shared_ptr< NumeredEntFiniteElement_multiIndex > fe_ptr=nullptr, MoFEMTypes bh=MF_EXIST, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr(), int verb=DEFAULT_VERBOSITY)=0
Make a loop over finite elements.
MoFEM::PrismsFromSurfaceInterface
merge node from two bit levels
Definition: PrismsFromSurfaceInterface.hpp:15
MoFEM::ProblemsManager::buildProblem
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
Definition: ProblemsManager.cpp:87
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
TriFE::setGaussPts
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
Definition: prism_elements_from_surface.cpp:466
H1
@ H1
continuous field
Definition: definitions.h:85
PrismOp::PrismOp
PrismOp(moab::Interface &post_proc, MapCoords &map_coords)
Definition: prism_elements_from_surface.cpp:342
MoFEM::PetscData::x
Vec x
Definition: LoopMethods.hpp:51
MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator
default operator for Flat Prism element
Definition: FatPrismElementForcesAndSourcesCore.hpp:54
sdf_hertz.d
float d
Definition: sdf_hertz.py:5
EntityHandle
MoFEM::ProblemsManager
Problem manager is used to build and partition problems.
Definition: ProblemsManager.hpp:21
PrismFE::setGaussPtsThroughThickness
MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness)
Definition: prism_elements_from_surface.cpp:449
PrismFE::setGaussPtsTrianglesOnly
MoFEMErrorCode setGaussPtsTrianglesOnly(int order_triangles_only)
Definition: prism_elements_from_surface.cpp:439
d3
static constexpr std::array< double, 3 > d3
Definition: prism_elements_from_surface.cpp:20
MoFEM::CoreInterface::modify_finite_element_add_field_row
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
MoFEM::EntitiesFieldData::EntData::getLocalIndices
const VectorInt & getLocalIndices() const
get local indices of dofs on entity
Definition: EntitiesFieldData.hpp:1216
PrismOp::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition: prism_elements_from_surface.cpp:347
PrismFE::getRuleTrianglesOnly
int getRuleTrianglesOnly(int order)
Definition: prism_elements_from_surface.cpp:436
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
PrismOp
Definition: prism_elements_from_surface.cpp:50
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
VERY_NOISY
@ VERY_NOISY
Definition: definitions.h:212
MoFEM::th
Tag th
Definition: Projection10NodeCoordsOnField.cpp:122
MoFEM::EntitiesFieldData::EntData::getFieldData
const VectorDouble & getFieldData() const
get dofs values
Definition: EntitiesFieldData.hpp:1241
MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsTrianglesOnly
MatrixDouble gaussPtsTrianglesOnly
Definition: FatPrismElementForcesAndSourcesCore.hpp:189
EdgeFE::EdgeFE
EdgeFE(MoFEM::Interface &m_field, std::array< Range, 3 > &edges_blocks, EntityHandle prims)
Definition: prism_elements_from_surface.cpp:524
EdgeFE::getRule
int getRule(int order_row, int order_col, int order_data)
another variant of getRule
Definition: prism_elements_from_surface.cpp:529
MoFEM.hpp
EdgeFE
Definition: prism_elements_from_surface.cpp:118
TriFE
Definition: prism_elements_from_surface.cpp:91
MoFEM::CoreTmp< 0 >::Finalize
static MoFEMErrorCode Finalize()
Checks for options to be called at the conclusion of the program.
Definition: Core.cpp:112
PrismOp::postProc
moab::Interface & postProc
Definition: prism_elements_from_surface.cpp:57
TriFE::TriFE
TriFE(MoFEM::Interface &m_field, MatrixDouble &tri_coords, EntityHandle prims)
Definition: prism_elements_from_surface.cpp:459
QuadFE::prism
EntityHandle prism
Definition: prism_elements_from_surface.cpp:116
MoFEM::CoreInterface::add_ents_to_field_by_type
virtual MoFEMErrorCode add_ents_to_field_by_type(const Range &ents, const EntityType type, const std::string &name, int verb=DEFAULT_VERBOSITY)=0
Add entities to field meshset.
order
constexpr int order
Definition: dg_projection.cpp:18
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
CoordsAndHandle
Definition: prism_elements_from_surface.cpp:25
QuadFE::getRule
int getRule(int order_row, int order_col, int order_data)
another variant of getRule
Definition: prism_elements_from_surface.cpp:491
MoFEM::Interface
DeprecatedCoreInterface Interface
Definition: Interface.hpp:1975
Op::doWork
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Definition: prism_elements_from_surface.cpp:585
TriFE::getRule
int getRule(int order_row, int order_col, int order_data)
another variant of getRule
Definition: prism_elements_from_surface.cpp:464
Op::Op
Op(moab::Interface &post_proc, MapCoords &map_coords, EntityHandle prism)
Definition: prism_elements_from_surface.cpp:579
MoFEM::CoreInterface::add_ents_to_finite_element_by_type
virtual MoFEMErrorCode add_ents_to_finite_element_by_type(const EntityHandle entities, const EntityType type, const std::string &name, const bool recursive=true)=0
add entities to finite element
TriFE::prism
EntityHandle prism
Definition: prism_elements_from_surface.cpp:102
CoordsAndHandle::getArg
static double getArg(double x)
Definition: prism_elements_from_surface.cpp:27
number_of_prisms_layers
static constexpr int number_of_prisms_layers
Definition: prism_elements_from_surface.cpp:17
PrismFE::triCoords
MatrixDouble & triCoords
Definition: prism_elements_from_surface.cpp:74
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:535
MoFEM::CoreInterface::add_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
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
QuadFE::QuadFE
QuadFE(MoFEM::Interface &m_field, std::array< Range, 3 > &edges_blocks, EntityHandle prims)
Definition: prism_elements_from_surface.cpp:486
PrismFE::getRuleThroughThickness
int getRuleThroughThickness(int order)
Definition: prism_elements_from_surface.cpp:437
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
MoFEM::CoreInterface::modify_finite_element_add_field_col
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
delta
static constexpr double delta
Definition: prism_elements_from_surface.cpp:18
MoFEM::CoreInterface::build_finite_elements
virtual MoFEMErrorCode build_finite_elements(int verb=DEFAULT_VERBOSITY)=0
Build finite elements.
MoFEM::ForcesAndSourcesCore::UserDataOperator
Definition: ForcesAndSourcesCore.hpp:549
MoFEM::MeshsetsManager::getEntitiesByDimension
MoFEMErrorCode getEntitiesByDimension(const int ms_id, const unsigned int cubit_bc_type, const int dimension, Range &entities, const bool recursive=true) const
get entities from CUBIT/meshset of a particular entity dimension
Definition: MeshsetsManager.cpp:666
MoFEM::PrismsFromSurfaceInterface::seedPrismsEntities
MoFEMErrorCode seedPrismsEntities(Range &prisms, const BitRefLevel &bit, int verb=-1)
Seed prism entities by bit level.
Definition: PrismsFromSurfaceInterface.cpp:110
MoFEM::CoreInterface::add_field
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.
convert.type
type
Definition: convert.py:64
MoFEM::PrismsFromSurfaceInterface::createPrismsFromPrisms
MoFEMErrorCode createPrismsFromPrisms(const Range &prisms, bool from_down, Range &out_prisms, int verb=-1)
Make prisms by extruding top or bottom prisms.
Definition: PrismsFromSurfaceInterface.cpp:143
EdgeFE::setGaussPts
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
Definition: prism_elements_from_surface.cpp:531
MoFEM::EntitiesFieldData::EntData::getIndices
const VectorInt & getIndices() const
Get global indices of dofs on entity.
Definition: EntitiesFieldData.hpp:1201
help
static char help[]
Definition: prism_elements_from_surface.cpp:14
precision_exponent
static constexpr int precision_exponent
Definition: prism_elements_from_surface.cpp:16
MoFEM::FaceElementForcesAndSourcesCore
Face finite element.
Definition: FaceElementForcesAndSourcesCore.hpp:23
PrismFE
Definition: prism_elements_from_surface.cpp:62
MoFEM::EntitiesFieldData::EntData::getFieldDofs
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
Definition: EntitiesFieldData.hpp:1256
MoFEM::FatPrismElementForcesAndSourcesCore
FatPrism finite element.
Definition: FatPrismElementForcesAndSourcesCore.hpp:31
t
constexpr double t
plate stiffness
Definition: plate.cpp:59
TriFE::triCoords
MatrixDouble & triCoords
Definition: prism_elements_from_surface.cpp:101
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
MoFEM::ProblemsManager::partitionFiniteElements
MoFEMErrorCode partitionFiniteElements(const std::string name, bool part_from_moab=false, int low_proc=-1, int hi_proc=-1, int verb=VERBOSE)
partition finite elements
Definition: ProblemsManager.cpp:2167
MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getCoordsAtGaussPts
MatrixDouble & getCoordsAtGaussPts()
get coordinates at Gauss pts.
Definition: FatPrismElementForcesAndSourcesCore.hpp:260
MoFEM::CoreInterface::modify_problem_ref_level_add_bit
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
mesh_file_name
char mesh_file_name[255]
Definition: mesh_smoothing.cpp:23
convert.n
n
Definition: convert.py:82
MoFEM::FatPrismElementForcesAndSourcesCore::gaussPtsThroughThickness
MatrixDouble gaussPtsThroughThickness
Definition: FatPrismElementForcesAndSourcesCore.hpp:191
MoFEM::ForcesAndSourcesCore
structure to get information form mofem into EntitiesFieldData
Definition: ForcesAndSourcesCore.hpp:22
edge_coords
static const double edge_coords[6][6]
Definition: forces_and_sources_h1_continuity_check.cpp:18
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
Range
MoFEM::CoreTmp< 0 >::Initialize
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
FTensor::dd
const Tensor2_symmetric_Expr< const ddTensor0< T, Dim, i, j >, typename promote< T, double >::V, Dim, i, j > dd(const Tensor0< T * > &a, const Index< i, Dim > index1, const Index< j, Dim > index2, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
Definition: ddTensor0.hpp:33
CATCH_ERRORS
#define CATCH_ERRORS
Catch errors.
Definition: definitions.h:372
PrismOp::nodeHandles
std::vector< EntityHandle > nodeHandles
Definition: prism_elements_from_surface.cpp:59
MoFEM::EntitiesFieldData::EntData::getN
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
Definition: EntitiesFieldData.hpp:1305
PrismFE::PrismFE
PrismFE(MoFEM::Interface &m_field, MatrixDouble &tri_coords)
Definition: prism_elements_from_surface.cpp:433
MoFEM::Core
CoreTmp< 0 > Core
Definition: Core.hpp:1094
Op
Definition: prism_elements_from_surface.cpp:78
MoFEM::ProblemsManager::partitionSimpleProblem
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
Definition: ProblemsManager.cpp:1537
BLOCKSET
@ BLOCKSET
Definition: definitions.h:148
eps
static const double eps
Definition: check_base_functions_derivatives_on_tet.cpp:11
AINSWORTH_LEGENDRE_BASE
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
Definition: definitions.h:60
MoFEM::PrismsFromSurfaceInterface::NO_SWAP
@ NO_SWAP
Definition: PrismsFromSurfaceInterface.hpp:38
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
MoFEM::ForcesAndSourcesCore::gaussPts
MatrixDouble gaussPts
Matrix of integration points.
Definition: ForcesAndSourcesCore.hpp:109
MoFEM::PetscOptionsGetString
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
Definition: DeprecatedPetsc.hpp:172
MoFEM::MeshsetsManager
Interface for managing meshsets containing materials and boundary conditions.
Definition: MeshsetsManager.hpp:104
QuadFE
Definition: prism_elements_from_surface.cpp:105
MoFEM::CoreInterface::build_fields
virtual MoFEMErrorCode build_fields(int verb=DEFAULT_VERBOSITY)=0
MoFEM::CoreInterface::modify_finite_element_add_field_data
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
MoFEM::PrismsFromSurfaceInterface::createdVertices
std::map< EntityHandle, EntityHandle > createdVertices
Definition: PrismsFromSurfaceInterface.hpp:24
d4
static constexpr std::array< double, 3 > d4
Definition: prism_elements_from_surface.cpp:21
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
MoFEM::BitRefManager
Managing BitRefLevels.
Definition: BitRefManager.hpp:21
MoFEM::CoreInterface::modify_problem_add_finite_element
virtual MoFEMErrorCode modify_problem_add_finite_element(const std::string name_problem, const std::string &fe_name)=0
add finite element to problem, this add entities assigned to finite element to a particular problem
MoFEM::Types::BitRefLevel
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
Definition: Types.hpp:40
MOFEM_ATOM_TEST_INVALID
@ MOFEM_ATOM_TEST_INVALID
Definition: definitions.h:40
MoFEM::CoreInterface::build_adjacencies
virtual MoFEMErrorCode build_adjacencies(const Range &ents, int verb=DEFAULT_VERBOSITY)=0
build adjacencies
QuadFE::setGaussPts
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
Definition: prism_elements_from_surface.cpp:493
MoFEM::FEMethod::numeredEntFiniteElementPtr
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
Definition: LoopMethods.hpp:383
MoFEM::PrismsFromSurfaceInterface::createPrisms
MoFEMErrorCode createPrisms(const Range &ents, const SwapType swap_type, Range &prisms, int verb=-1)
Make prisms from triangles.
Definition: PrismsFromSurfaceInterface.cpp:22
MoFEM::CoreInterface::set_field_order
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.
sum_matrix
static double sum_matrix(MatrixDouble &m)
Definition: edge_and_bubble_shape_functions_on_quad.cpp:12
PrismOp::mapCoords
MapCoords & mapCoords
Definition: prism_elements_from_surface.cpp:58
MoFEM::ProblemsManager::partitionGhostDofs
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
Definition: ProblemsManager.cpp:2339
MoFEM::CoreInterface::add_problem
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
OP
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
EdgeFE::edgeBlocks
std::array< Range, 3 > & edgeBlocks
Definition: prism_elements_from_surface.cpp:128
MoFEM::FatPrismElementForcesAndSourcesCore::UserDataOperator::getGaussPts
MatrixDouble & getGaussPts()
get Gauss pts. in the prism
Definition: FatPrismElementForcesAndSourcesCore.hpp:243
MoFEM::PrismsFromSurfaceInterface::setThickness
MoFEMErrorCode setThickness(const Range &prisms, const double director3[], const double director4[])
Definition: PrismsFromSurfaceInterface.cpp:161
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
MOFEM_INVALID_DATA
@ MOFEM_INVALID_DATA
Definition: definitions.h:36
MoFEM::ForcesAndSourcesCore::mField
Interface & mField
Definition: ForcesAndSourcesCore.hpp:24
MapCoords
multi_index_container< CoordsAndHandle, indexed_by< hashed_unique< composite_key< CoordsAndHandle, member< CoordsAndHandle, int, &CoordsAndHandle::x >, member< CoordsAndHandle, int, &CoordsAndHandle::y >, member< CoordsAndHandle, int, &CoordsAndHandle::z > > > > > MapCoords
Definition: prism_elements_from_surface.cpp:48