Adding prims on the surface and checking conformity between quads, triangles and prism.
Adding prims on the surface and checking conformity between quads, triangles and prism
static char help[] =
"...\n\n";
static constexpr double delta =
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) {
};
};
typedef multi_index_container<
indexed_by<
hashed_unique<composite_key<
member<CoordsAndHandle, int, &CoordsAndHandle::y>,
member<CoordsAndHandle, int, &CoordsAndHandle::z>>>
>>
public:
};
private:
};
template <
typename OP>
struct Op :
public OP {
public:
};
int getRule(
int order_row,
int order_col,
int order_data);
private:
};
int getRule(
int order_row,
int order_col,
int order_data);
private:
};
int getRule(
int order_row,
int order_col,
int order_data);
private:
};
int main(
int argc,
char *argv[]) {
try {
moab::Core mb_instance;
moab::Interface &moab = mb_instance;
int rank;
MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
PetscBool flg = PETSC_TRUE;
255, &flg);
if (flg != PETSC_TRUE)
"error -my_file (MESH FILE NEEDED)");
const char *option;
option = "";
CHKERR moab.get_entities_by_type(0, MBTRI, tris,
false);
CHKERR moab.get_connectivity(tris, tris_verts);
CHKERR moab.get_coords(tris_verts, &tri_coords(0, 0));
std::array<Range, 3> edge_block;
for (auto b : {1, 2, 3})
tris, PrismsFromSurfaceInterface::NO_SWAP, prisms);
Range extrude_prisms = prisms;
extrude_prisms, false, add_prims_layer);
prisms_from_surface_interface->
setThickness(add_prims_layer,
d3.data(),
extrude_prisms = add_prims_layer;
prisms.merge(add_prims_layer);
add_prims_layer.clear();
}
CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER, meshset);
CHKERR moab.add_entities(meshset, prisms);
CHKERR moab.get_connectivity(prisms, verts);
CHKERR moab.get_coords(verts, &coords(0, 0));
for (
size_t v = 0;
v != verts.size(); ++
v)
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]);
CHKERR m_field.
get_moab().create_element(MBPRISM, one_prism_nodes.data(), 6,
one_prism);
one_prism_range.insert(one_prism);
CHKERR moab.add_entities(one_prism_meshset, one_prism_range);
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);
bit_level0.set(0);
one_prism_range, bit_level0);
bit_level0);
"EDGE");
"TRI");
"QUAD");
MBPRISM, "PRISM");
PrismFE fe_prism(m_field, tri_coords);
fe_prism.getOpPtrVector().push_back(
new PrismOp(moab, map_coords));
EdgeFE fe_edge(m_field, edge_block, one_prism);
fe_edge.getOpPtrVector().push_back(
moab, map_coords, one_prism));
TriFE fe_tri(m_field, tri_coords, one_prism);
fe_tri.getOpPtrVector().push_back(
moab, map_coords, one_prism));
QuadFE fe_quad(m_field, edge_block, one_prism);
fe_quad.getOpPtrVector().push_back(
moab, map_coords, one_prism));
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;
}
postProc(post_proc), mapCoords(map_coords) {}
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))
for (
int dd = 0;
dd != nb_dofs; ++
dd)
if (type == MBVERTEX) {
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
else
}
&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})
"Inconsistency between node coords and integration "
"points");
}
for (auto d : {0, 1, 2})
"Inconsistency between node coords and integration "
"points");
}
MB_TAG_CREAT | MB_TAG_DENSE, def_val);
}
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;
for (size_t rr = 0; rr != trans_base.size1(); ++rr) {
auto tag_name = tag_name_base + "Base" + to_str(rr);
CHKERR postProc.tag_get_handle(tag_name.c_str(), 1, MB_TYPE_DOUBLE, th,
MB_TAG_CREAT | MB_TAG_DENSE, def_val);
&trans_base(rr, 0));
}
}
for (
int gg = 0; gg !=
triCoords.size1(); ++gg)
for (
int dd = 0;
dd != 2; ++
dd)
}
}
prism(prism) {}
int TriFE::getRule(
int order_row,
int order_col,
int order_data) {
return -1; }
int order_data) {
int side, sense, offset;
std::array<int, 2> swap = {0, 1};
if (side == 3)
swap = std::array<int, 2>{1, 0};
for (
int gg = 0; gg !=
triCoords.size1(); ++gg)
for (
int dd = 0;
dd != 2; ++
dd)
}
prism(prism) {}
int QuadFE::getRule(
int order_row,
int order_col,
int order_data) {
return -1; }
int order_data) {
int side, sense, offset;
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}};
int gg = 0;
for (size_t rr = 0; rr != edge_verts.size(); ++rr) {
const double edge_dist =
x * normal[side][0] + y * normal[side][1];
gaussPts(swap[side][0], gg) = edge_dist;
}
}
}
:
EdgeEle(m_field), edgeBlocks(edge_blocks),
prism(prism) {}
int EdgeFE::getRule(
int order_row,
int order_col,
int order_data) {
return -1; }
int order_data) {
int side, sense, offset;
if (side >= 3 && side <= 5) {
} 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;
&*coords.data().begin());
side = side_map[side];
gaussPts.resize(2, edge_verts.size(),
false);
for (size_t gg = 0; gg != edge_verts.size(); ++gg) {
const double edge_dist =
x * normal[side][0] + y * normal[side][1];
}
}
}
template <typename OP>
postProc(post_proc), mapCoords(map_coords), prism(prism) {}
template <typename OP>
constexpr double def_val[] = {0, 0, 0};
switch (type) {
case MBVERTEX:
case MBEDGE:
case MBTRI:
case MBQUAD:
break;
default:
}
for (
int dd = 0;
dd != nb_dofs; ++
dd)
std::cerr <<
"Indices: " << data.
getIndices() << std::endl;
"Indicices/data inconsistency %3.1f != %d",
}
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 it = mapCoords.find(
t);
if (it != mapCoords.end())
nodeHandles.emplace_back(it->node);
else
}
}
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 prism_base(trans_base.size1(), trans_base.size2());
if (trans_base.size2() != nodeHandles.size())
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) {
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;
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));
}
}
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;
constexpr double eps = 1e-6;
cout << "Inconsistent base " << tag_prism_name_base << " "
<< tag_name_base << " sum " << sum << endl;
"Inconsistent base %s sum %6.4e", tag_prism_name_base.c_str(),
sum);
}
}
ForcesAndSourcesCore::UserDataOperator UserDataOperator
#define CATCH_ERRORS
Catch errors.
@ AINSWORTH_LEGENDRE_BASE
Ainsworth Cole (Legendre) approx. base .
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_ATOM_TEST_INVALID
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
static double sum_matrix(MatrixDouble &m)
FTensor::Index< 'n', SPACE_DIM > n
FTensor::Index< 'm', SPACE_DIM > m
static const double edge_coords[6][6]
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 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
virtual MoFEMErrorCode modify_finite_element_add_field_data(const std::string &fe_name, const std::string name_filed)=0
set finite element field data
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 build_fields(int verb=DEFAULT_VERBOSITY)=0
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.
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.
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.
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
MoFEMErrorCode partitionGhostDofs(const std::string name, int verb=VERBOSE)
determine ghost nodes
MoFEMErrorCode partitionSimpleProblem(const std::string name, int verb=VERBOSE)
partition problem dofs
MoFEMErrorCode buildProblem(const std::string name, const bool square_matrix, int verb=VERBOSE)
build problem data structures
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
virtual MoFEMErrorCode add_problem(const std::string &name, enum MoFEMTypes bh=MF_EXCL, int verb=DEFAULT_VERBOSITY)=0
Add problem.
virtual MoFEMErrorCode modify_problem_ref_level_add_bit(const std::string &name_problem, const BitRefLevel &bit)=0
add ref level to problem
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
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
const Tensor1_Expr< const dTensor0< T, Dim, i >, typename promote< T, double >::V, Dim, i > d(const Tensor0< T * > &a, const Index< i, Dim > index, const Tensor1< int, Dim > &d_ijk, const Tensor1< double, Dim > &d_xyz)
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)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
UBlasMatrix< double > MatrixDouble
std::bitset< BITREFLEVEL_SIZE > BitRefLevel
Bit structure attached to each entity identifying to what mesh entity is attached.
implementation of Data Operators for Forces and Sources
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
constexpr double t
plate stiffness
static constexpr int number_of_prisms_layers
static constexpr int precision_exponent
static constexpr std::array< double, 3 > d4
static constexpr double delta
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
static constexpr std::array< double, 3 > d3
static double getArg(double x)
std::array< Range, 3 > & edgeBlocks
EdgeFE(MoFEM::Interface &m_field, std::array< Range, 3 > &edges_blocks, EntityHandle prims)
int getRule(int order_row, int order_col, int order_data)
another variant of getRule
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
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.
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.
Data on single entity (This is passed as argument to DataOperator::doWork)
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
const VectorDouble & getFieldData() const
get dofs values
const VectorInt & getLocalIndices() const
get local indices of dofs on entity
const VectorDofs & getFieldDofs() const
get dofs data stature FEDofEntity
const VectorInt & getIndices() const
Get global indices of dofs on entity.
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
default operator for Flat Prism element
MatrixDouble & getGaussPts()
get Gauss pts. in the prism
MatrixDouble & getCoordsAtGaussPts()
get coordinates at Gauss pts.
MatrixDouble gaussPtsTrianglesOnly
MatrixDouble gaussPtsThroughThickness
structure to get information form mofem into EntitiesFieldData
MatrixDouble gaussPts
Matrix of integration points.
Interface for managing meshsets containing materials and boundary conditions.
merge node from two bit levels
std::map< EntityHandle, EntityHandle > createdVertices
MoFEMErrorCode createPrisms(const Range &ents, const SwapType swap_type, Range &prisms, int verb=-1)
Make prisms from triangles.
MoFEMErrorCode setThickness(const Range &prisms, const double director3[], const double director4[])
MoFEMErrorCode seedPrismsEntities(Range &prisms, const BitRefLevel &bit, int verb=-1)
Seed prism entities by bit level.
MoFEMErrorCode createPrismsFromPrisms(const Range &prisms, bool from_down, Range &out_prisms, int verb=-1)
Make prisms by extruding top or bottom prisms.
Problem manager is used to build and partition problems.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface refernce to pointer of interface.
Op(moab::Interface &post_proc, MapCoords &map_coords, EntityHandle prism)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
moab::Interface & postProc
std::vector< EntityHandle > nodeHandles
int getRuleThroughThickness(int order)
PrismFE(MoFEM::Interface &m_field, MatrixDouble &tri_coords)
MoFEMErrorCode setGaussPtsTrianglesOnly(int order_triangles_only)
MoFEMErrorCode setGaussPtsThroughThickness(int order_thickness)
int getRuleTrianglesOnly(int order)
std::vector< EntityHandle > nodeHandles
PrismOp(moab::Interface &post_proc, MapCoords &map_coords)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
moab::Interface & postProc
std::array< Range, 3 > & edgeBlocks
int getRule(int order_row, int order_col, int order_data)
another variant of getRule
QuadFE(MoFEM::Interface &m_field, std::array< Range, 3 > &edges_blocks, EntityHandle prims)
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
set user specific integration rule
TriFE(MoFEM::Interface &m_field, MatrixDouble &tri_coords, EntityHandle prims)
int getRule(int order_row, int order_col, int order_data)
another variant of getRule