14static char help[] =
"...\n\n";
20static constexpr std::array<double, 3>
d3 = {0, 0, 0};
21static constexpr std::array<double, 3>
d4 = {0, 0,
delta};
38typedef multi_index_container<
42 hashed_unique<composite_key<
44 member<CoordsAndHandle, int, &CoordsAndHandle::y>,
45 member<CoordsAndHandle, int, &CoordsAndHandle::z>>>
78template <
typename OP>
struct Op :
public OP {
96 int getRule(
int order_row,
int order_col,
int order_data);
110 int getRule(
int order_row,
int order_col,
int order_data);
123 int getRule(
int order_row,
int order_col,
int order_data);
132int main(
int argc,
char *argv[]) {
138 moab::Core mb_instance;
139 moab::Interface &moab = mb_instance;
141 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
144 PetscBool flg = PETSC_TRUE;
148 if (flg != PETSC_TRUE)
150 "error -my_file (MESH FILE NEEDED)");
158 CHKERR moab.get_entities_by_type(0, MBTRI, tris,
false);
160 CHKERR moab.get_connectivity(tris, tris_verts);
162 CHKERR moab.get_coords(tris_verts, &tri_coords(0, 0));
167 std::array<Range, 3> edge_block;
168 for (
auto b : {1, 2, 3})
179 Range add_prims_layer;
180 Range extrude_prisms = prisms;
185 extrude_prisms,
false, add_prims_layer);
186 prisms_from_surface_interface->
setThickness(add_prims_layer,
d3.data(),
188 extrude_prisms = add_prims_layer;
189 prisms.merge(add_prims_layer);
190 add_prims_layer.clear();
194 CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER, meshset);
195 CHKERR moab.add_entities(meshset, prisms);
199 CHKERR moab.get_connectivity(prisms, verts);
201 CHKERR moab.get_coords(verts, &coords(0, 0));
203 for (
size_t v = 0;
v != verts.size(); ++
v)
207 CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER,
210 std::array<double, 18> one_prism_coords = {0, 0, 0, 1, 0, 0, 0, 1, 0,
211 0, 0, 1, 1, 0, 1, 0, 1, 1};
212 std::array<EntityHandle, 6> one_prism_nodes;
213 for (
int n = 0;
n != 6; ++
n)
214 CHKERR moab.create_vertex(&one_prism_coords[3 *
n], one_prism_nodes[
n]);
216 CHKERR m_field.
get_moab().create_element(MBPRISM, one_prism_nodes.data(), 6,
218 Range one_prism_range;
219 one_prism_range.insert(one_prism);
220 CHKERR moab.add_entities(one_prism_meshset, one_prism_range);
221 Range one_prism_verts;
222 CHKERR moab.get_connectivity(one_prism_range, one_prism_verts);
223 CHKERR moab.add_entities(one_prism_meshset, one_prism_verts);
224 Range one_prism_adj_ents;
225 for (
int d = 1; d != 3; ++d)
226 CHKERR moab.get_adjacencies(one_prism_range, d,
true, one_prism_adj_ents,
227 moab::Interface::UNION);
228 CHKERR moab.add_entities(one_prism_meshset, one_prism_adj_ents);
233 one_prism_range, bit_level0);
310 PrismFE fe_prism(m_field, tri_coords);
314 EdgeFE fe_edge(m_field, edge_block, one_prism);
317 moab, map_coords, one_prism));
320 TriFE fe_tri(m_field, tri_coords, one_prism);
323 moab, map_coords, one_prism));
326 QuadFE fe_quad(m_field, edge_block, one_prism);
329 moab, map_coords, one_prism));
332 CHKERR moab.write_file(
"prism_mesh.vtk",
"VTK",
"", &meshset, 1);
333 CHKERR moab.write_file(
"one_prism_mesh.vtk",
"VTK",
"", &one_prism_meshset,
345 postProc(post_proc), mapCoords(map_coords) {}
349 constexpr double def_val[] = {0, 0, 0};
361 if (type == MBTRI && (side != 3 && side != 4))
363 if (type == MBQUAD && (side == 3 || side == 4))
367 for (
int dd = 0; dd != nb_dofs; ++dd)
370 if (type == MBVERTEX) {
375 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
390 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
391 constexpr double eps = 1e-12;
392 for (
auto d : {0, 1, 2})
395 "Inconsistency between node coords and integration "
398 for (
auto d : {0, 1, 2})
401 "Inconsistency between node coords and integration "
407 MB_TAG_CREAT | MB_TAG_DENSE, def_val);
412 auto to_str = [](
auto i) {
return boost::lexical_cast<std::string>(
i); };
413 std::string tag_name_base =
414 "PrismType" + to_str(type) +
"Side" + to_str(side);
415 std::cout << tag_name_base << endl;
421 for (
size_t rr = 0; rr != trans_base.size1(); ++rr) {
422 auto tag_name = tag_name_base +
"Base" + to_str(rr);
425 MB_TAG_CREAT | MB_TAG_DENSE, def_val);
443 for (
int gg = 0; gg !=
triCoords.size1(); ++gg)
444 for (
int dd = 0; dd != 2; ++dd)
471 int side, sense, offset;
473 std::array<int, 2> swap = {0, 1};
475 swap = std::array<int, 2>{1, 0};
479 for (
int gg = 0; gg !=
triCoords.size1(); ++gg)
480 for (
int dd = 0; dd != 2; ++dd)
498 int side, sense, offset;
506 constexpr double normal[3][2] = {{-1, 0}, {-0.5, 0.5}, {0, 1}};
507 constexpr double origin[3][2] = {{1, 0}, {1, 0}, {0, 0}};
508 constexpr int swap[3][2] = {{0, 1}, {0, 1}, {1, 0}};
511 for (
size_t rr = 0; rr != edge_verts.size(); ++rr) {
513 const double y =
edge_coords(rr, 1) - origin[side][1];
514 const double edge_dist =
x * normal[side][0] + y * normal[side][1];
516 gaussPts(swap[side][0], gg) = edge_dist;
526 :
EdgeEle(m_field), edgeBlocks(edge_blocks),
536 int side, sense, offset;
539 if (side >= 3 && side <= 5) {
549 constexpr double normal[3][2] = {{1, 0}, {-0.5, 0.5}, {0, 1}};
550 constexpr double origin[3][2] = {{0, 1}, {1, 0}, {0, 0}};
551 constexpr int side_map[9] = {0, 1, 2, -1, -1, -1, 0, 1, 2};
558 &*coords.data().begin());
559 side = side_map[side];
565 gaussPts.resize(2, edge_verts.size(),
false);
567 for (
size_t gg = 0; gg != edge_verts.size(); ++gg) {
569 const double y =
edge_coords(gg, 1) - origin[side][1];
570 const double edge_dist =
x * normal[side][0] + y * normal[side][1];
578template <
typename OP>
582 postProc(post_proc), mapCoords(map_coords), prism(prism) {}
584template <
typename OP>
587 constexpr double def_val[] = {0, 0, 0};
600 for (
int dd = 0; dd != nb_dofs; ++dd)
602 std::cerr <<
"Indices: " << data.
getIndices() << std::endl;
604 std::cerr <<
"Data: " << data.
getFieldData() << std::endl;
606 "Indicices/data inconsistency %3.1f != %d",
612 int side_prism, sense, offset;
613 if (type == MBVERTEX) {
614 CHKERR postProc.side_number(prism, fe_ent, side_prism, sense, offset);
616 CHKERR postProc.side_number(prism, ent, side_prism, sense, offset);
618 if (type == MBVERTEX) {
619 auto &coords_at_pts = OP::getCoordsAtGaussPts();
620 const size_t nb_gauss_pts = coords_at_pts.size1();
622 nodeHandles.reserve(nb_gauss_pts);
624 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
629 auto it = mapCoords.find(
t);
630 if (it != mapCoords.end())
631 nodeHandles.emplace_back(it->node);
637 auto to_str = [](
auto i) {
return boost::lexical_cast<std::string>(
i); };
638 std::string tag_name_base =
639 "FEType" + to_str(OP::getNumeredEntFiniteElementPtr()->getEntType()) +
640 "Type" + to_str(type) +
"Side" + to_str(side_prism);
642 std::string tag_prism_name_base =
643 "PrismType" + to_str(type) +
"Side" + to_str(side_prism);
646 MatrixDouble prism_base(trans_base.size1(), trans_base.size2());
647 if (trans_base.size2() != nodeHandles.size())
649 trans_base.size2(), nodeHandles.size());
651 for (
size_t rr = 0; rr != trans_base.size1(); ++rr) {
653 std::string tag_name = tag_name_base +
"Base";
654 if (type == MBVERTEX) {
657 CHKERR postProc.side_number(prism, node, prism_mode_side, sense, offset);
658 tag_name += to_str(prism_mode_side);
660 tag_name += to_str(rr);
663 std::cout << tag_name << endl;
666 CHKERR postProc.tag_get_handle(tag_name.c_str(), 1, MB_TYPE_DOUBLE,
th,
667 MB_TAG_CREAT | MB_TAG_DENSE, def_val);
668 CHKERR postProc.tag_set_data(
th, &nodeHandles[0], nodeHandles.size(),
671 if (type != MBVERTEX) {
672 auto tag_prism_name = tag_prism_name_base +
"Base" + to_str(rr);
674 CHKERR postProc.tag_get_handle(tag_prism_name.c_str(), th_prism);
675 CHKERR postProc.tag_get_data(th_prism, &nodeHandles[0],
676 nodeHandles.size(), &prism_base(rr, 0));
682 for (
unsigned int ii = 0; ii <
m.size1(); ii++) {
683 for (
unsigned int jj = 0; jj <
m.size2(); jj++) {
684 s += std::abs(
m(ii, jj));
690 if (type != MBVERTEX) {
691 prism_base -= trans_base;
693 constexpr double eps = 1e-6;
695 if (std::abs(sum) >
eps)
696 cout <<
"Inconsistent base " << tag_prism_name_base <<
" "
697 << tag_name_base <<
" sum " << sum << endl;
699 if (std::abs(sum) >
eps)
701 "Inconsistent base %s sum %6.4e", tag_prism_name_base.c_str(),
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)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
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)
CoordsAndHandle(const double *coords, EntityHandle v)
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.
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
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