24 using namespace MoFEM;
26 static char help[] =
"...\n\n";
32 static constexpr std::array<double, 3>
d3 = {0, 0, 0};
33 static constexpr std::array<double, 3>
d4 = {0, 0,
delta};
39 inline static double getArg(
double x) {
46 : x(getArg(coords[0])), y(getArg(coords[1])), z(getArg(coords[2])),
50 typedef multi_index_container<
54 hashed_unique<composite_key<
56 member<CoordsAndHandle, int, &CoordsAndHandle::y>,
57 member<CoordsAndHandle, int, &CoordsAndHandle::z>>>
78 int getRuleTrianglesOnly(
int order);
79 int getRuleThroughThickness(
int order);
90 template <
typename OP>
struct Op :
public OP {
108 int getRule(
int order_row,
int order_col,
int order_data);
110 MoFEMErrorCode setGaussPts(
int order_row,
int order_col,
int order_data);
122 int getRule(
int order_row,
int order_col,
int order_data);
124 MoFEMErrorCode setGaussPts(
int order_row,
int order_col,
int order_data);
135 int getRule(
int order_row,
int order_col,
int order_data);
137 MoFEMErrorCode setGaussPts(
int order_row,
int order_col,
int order_data);
144 int main(
int argc,
char *argv[]) {
153 MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
156 PetscBool flg = PETSC_TRUE;
160 if (flg != PETSC_TRUE)
162 "error -my_file (MESH FILE NEEDED)");
170 CHKERR moab.get_entities_by_type(0, MBTRI, tris,
false);
172 CHKERR moab.get_connectivity(tris, tris_verts);
174 CHKERR moab.get_coords(tris_verts, &tri_coords(0, 0));
179 std::array<Range, 3> edge_block;
180 for (
auto b : {1, 2, 3})
190 Range add_prims_layer;
191 Range extrude_prisms = prisms;
196 extrude_prisms,
false, add_prims_layer);
197 prisms_from_surface_interface->
setThickness(add_prims_layer,
d3.data(),
199 extrude_prisms = add_prims_layer;
200 prisms.merge(add_prims_layer);
201 add_prims_layer.clear();
205 CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER, meshset);
206 CHKERR moab.add_entities(meshset, prisms);
210 CHKERR moab.get_connectivity(prisms, verts);
212 CHKERR moab.get_coords(verts, &coords(0, 0));
214 for (
size_t v = 0;
v != verts.size(); ++
v)
218 CHKERR moab.create_meshset(MESHSET_SET | MESHSET_TRACK_OWNER,
221 std::array<double, 18> one_prism_coords = {0, 0, 0, 1, 0, 0, 0, 1, 0,
222 0, 0, 1, 1, 0, 1, 0, 1, 1};
223 std::array<EntityHandle, 6> one_prism_nodes;
224 for (
int n = 0;
n != 6; ++
n)
225 CHKERR moab.create_vertex(&one_prism_coords[3 *
n], one_prism_nodes[
n]);
227 CHKERR m_field.
get_moab().create_element(MBPRISM, one_prism_nodes.data(), 6,
229 Range one_prism_range;
230 one_prism_range.insert(one_prism);
231 CHKERR moab.add_entities(one_prism_meshset, one_prism_range);
232 Range one_prism_verts;
233 CHKERR moab.get_connectivity(one_prism_range, one_prism_verts);
234 CHKERR moab.add_entities(one_prism_meshset, one_prism_verts);
235 Range one_prism_adj_ents;
236 for (
int d = 1;
d != 3; ++
d)
237 CHKERR moab.get_adjacencies(one_prism_range,
d,
true, one_prism_adj_ents,
238 moab::Interface::UNION);
239 CHKERR moab.add_entities(one_prism_meshset, one_prism_adj_ents);
244 one_prism_range, bit_level0);
321 PrismFE fe_prism(m_field, tri_coords);
325 EdgeFE fe_edge(m_field, edge_block, one_prism);
328 moab, map_coords, one_prism));
331 TriFE fe_tri(m_field, tri_coords, one_prism);
332 fe_tri.getOpPtrVector().push_back(
334 moab, map_coords, one_prism));
337 QuadFE fe_quad(m_field, edge_block, one_prism);
338 fe_quad.getOpPtrVector().push_back(
340 moab, map_coords, one_prism));
343 CHKERR moab.write_file(
"prism_mesh.vtk",
"VTK",
"", &meshset, 1);
344 CHKERR moab.write_file(
"one_prism_mesh.vtk",
"VTK",
"", &one_prism_meshset,
356 postProc(post_proc), mapCoords(map_coords) {}
360 constexpr
double def_val[] = {0, 0, 0};
372 if (
type == MBTRI && (side != 3 && side != 4))
374 if (
type == MBQUAD && (side == 3 || side == 4))
378 for (
int dd = 0;
dd != nb_dofs; ++
dd)
381 if (
type == MBVERTEX) {
386 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
401 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
402 constexpr
double eps = 1e-12;
403 for (
auto d : {0, 1, 2})
406 "Inconsistency between node coords and integration "
409 for (
auto d : {0, 1, 2})
412 "Inconsistency between node coords and integration "
418 MB_TAG_CREAT | MB_TAG_DENSE, def_val);
423 auto to_str = [](
auto i) {
return boost::lexical_cast<std::string>(
i); };
424 std::string tag_name_base =
425 "PrismType" + to_str(
type) +
"Side" + to_str(side);
426 std::cout << tag_name_base << endl;
432 for (
size_t rr = 0; rr != trans_base.size1(); ++rr) {
433 auto tag_name = tag_name_base +
"Base" + to_str(rr);
436 MB_TAG_CREAT | MB_TAG_DENSE, def_val);
454 for (
int gg = 0; gg !=
triCoords.size1(); ++gg)
455 for (
int dd = 0;
dd != 2; ++
dd)
481 const EntityHandle ent = numeredEntFiniteElementPtr->getEnt();
482 int side, sense, offset;
483 CHKERR mField.get_moab().side_number(
prism, ent, side, sense, offset);
484 std::array<int, 2> swap = {0, 1};
486 swap = std::array<int, 2>{1, 0};
488 gaussPts.resize(3,
triCoords.size1(),
false);
490 for (
int gg = 0; gg !=
triCoords.size1(); ++gg)
491 for (
int dd = 0;
dd != 2; ++
dd)
508 const EntityHandle ent = numeredEntFiniteElementPtr->getEnt();
509 int side, sense, offset;
510 CHKERR mField.get_moab().side_number(
prism, ent, side, sense, offset);
517 constexpr
double normal[3][2] = {{-1, 0}, {-0.5, 0.5}, {0, 1}};
518 constexpr
double origin[3][2] = {{1, 0}, {1, 0}, {0, 0}};
519 constexpr
int swap[3][2] = {{0, 1}, {0, 1}, {1, 0}};
522 for (
size_t rr = 0; rr != edge_verts.size(); ++rr) {
523 const double x =
edge_coords(rr, 0) - origin[side][0];
524 const double y =
edge_coords(rr, 1) - origin[side][1];
525 const double edge_dist = x *
normal[side][0] + y *
normal[side][1];
527 gaussPts(swap[side][0], gg) = edge_dist;
528 gaussPts(swap[side][1], gg) =
delta * cc;
537 :
EdgeEle(m_field), edgeBlocks(edge_blocks),
547 int side, sense, offset;
550 if (side >= 3 && side <= 5) {
560 constexpr
double normal[3][2] = {{1, 0}, {-0.5, 0.5}, {0, 1}};
561 constexpr
double origin[3][2] = {{0, 1}, {1, 0}, {0, 0}};
562 constexpr
int side_map[9] = {0, 1, 2, -1, -1, -1, 0, 1, 2};
569 &*coords.data().begin());
570 side = side_map[side];
576 gaussPts.resize(2, edge_verts.size(),
false);
578 for (
size_t gg = 0; gg != edge_verts.size(); ++gg) {
580 const double y =
edge_coords(gg, 1) - origin[side][1];
581 const double edge_dist =
x *
normal[side][0] + y *
normal[side][1];
589 template <
typename OP>
593 postProc(post_proc), mapCoords(map_coords), prism(prism) {}
595 template <
typename OP>
598 constexpr
double def_val[] = {0, 0, 0};
611 for (
int dd = 0;
dd != nb_dofs; ++
dd)
613 std::cerr <<
"Indices: " << data.
getIndices() << std::endl;
615 std::cerr <<
"Data: " << data.
getFieldData() << std::endl;
617 "Indicices/data inconsistency %3.1f != %d",
623 int side_prism, sense, offset;
624 if (
type == MBVERTEX) {
625 CHKERR postProc.side_number(prism, fe_ent, side_prism, sense, offset);
627 CHKERR postProc.side_number(prism, ent, side_prism, sense, offset);
629 if (
type == MBVERTEX) {
630 auto &coords_at_pts = OP::getCoordsAtGaussPts();
631 const size_t nb_gauss_pts = coords_at_pts.size1();
633 nodeHandles.reserve(nb_gauss_pts);
635 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
640 auto it = mapCoords.find(
t);
641 if (it != mapCoords.end())
642 nodeHandles.emplace_back(it->node);
648 auto to_str = [](
auto i) {
return boost::lexical_cast<std::string>(
i); };
649 std::string tag_name_base =
650 "FEType" + to_str(OP::getNumeredEntFiniteElementPtr()->getEntType()) +
651 "Type" + to_str(
type) +
"Side" + to_str(side_prism);
653 std::string tag_prism_name_base =
654 "PrismType" + to_str(
type) +
"Side" + to_str(side_prism);
657 MatrixDouble prism_base(trans_base.size1(), trans_base.size2());
658 if (trans_base.size2() != nodeHandles.size())
660 trans_base.size2(), nodeHandles.size());
662 for (
size_t rr = 0; rr != trans_base.size1(); ++rr) {
664 std::string tag_name = tag_name_base +
"Base";
665 if (
type == MBVERTEX) {
668 CHKERR postProc.side_number(prism, node, prism_mode_side, sense, offset);
669 tag_name += to_str(prism_mode_side);
671 tag_name += to_str(rr);
674 std::cout << tag_name << endl;
677 CHKERR postProc.tag_get_handle(tag_name.c_str(), 1, MB_TYPE_DOUBLE,
th,
678 MB_TAG_CREAT | MB_TAG_DENSE, def_val);
679 CHKERR postProc.tag_set_data(
th, &nodeHandles[0], nodeHandles.size(),
682 if (
type != MBVERTEX) {
683 auto tag_prism_name = tag_prism_name_base +
"Base" + to_str(rr);
685 CHKERR postProc.tag_get_handle(tag_prism_name.c_str(), th_prism);
686 CHKERR postProc.tag_get_data(th_prism, &nodeHandles[0],
687 nodeHandles.size(), &prism_base(rr, 0));
693 for (
unsigned int ii = 0; ii <
m.size1(); ii++) {
694 for (
unsigned int jj = 0; jj <
m.size2(); jj++) {
695 s += std::abs(
m(ii, jj));
701 if (
type != MBVERTEX) {
702 prism_base -= trans_base;
704 constexpr
double eps = 1e-6;
706 if (std::abs(sum) >
eps)
707 cout <<
"Inconsistent base " << tag_prism_name_base <<
" "
708 << tag_name_base <<
" sum " << sum << endl;
710 if (std::abs(sum) >
eps)
712 "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)
static const double edge_coords[6][6]
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 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 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_col(const std::string &fe_name, const std::string &name_row)=0
set field col 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.
EdgeElementForcesAndSourcesCoreSwitch< 0 > EdgeElementForcesAndSourcesCore
Edge finite element default.
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 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
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
FTensor::Index< 'i', SPACE_DIM > i
double v
phase velocity of light in medium (cm/ns)
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)
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)
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)
DeprecatedCoreInterface Interface
constexpr double t
plate stiffness
int main(int argc, char *argv[])
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
FTensor::Index< 'm', 3 > m
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 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 moab::Interface & get_moab()=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.
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_vector< 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)
QuadFE(MoFEM::Interface &m_field, std::array< Range, 3 > &edges_blocks, EntityHandle prims)
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
MoFEMErrorCode setGaussPts(int order_row, int order_col, int order_data)
TriFE(MoFEM::Interface &m_field, MatrixDouble &tri_coords, EntityHandle prims)
int getRule(int order_row, int order_col, int order_data)