v0.15.0
Loading...
Searching...
No Matches
VolumeElementForcesAndSourcesCoreOnSide.cpp
Go to the documentation of this file.
1/** \file VolumeElementForcesAndSourcesCoreOnSide.cpp
2
3\brief Implementation of volume element on side
4
5*/
6
7namespace MoFEM {
8
10
13
14 if (!sidePtrFE)
15 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Side element not set");
16
17 if (!operatorImplCalled) {
19 }
20
21 const auto type = numeredEntFiniteElementPtr->getEntType();
22 auto face_ptr_fe =
24
25 const int nb_gauss_pts = face_ptr_fe->gaussPts.size2();
26 gaussPts.resize(4, nb_gauss_pts, false);
27 gaussPts.clear();
28
29 constexpr double tet_coords[12] = {0, 0, 0, 1, 0, 0, 0, 1, 0, 0, 0, 1};
30 constexpr double hex_coords[24] = {0, 0, 0, 1, 0, 0, 1, 1, 0, 0, 1, 0,
31 0, 0, 1, 1, 0, 1, 1, 1, 1, 0, 1, 1};
32
33 auto set_gauss = [&](auto &coords) {
35 FTensor::Index<'i', 3> i;
36 auto &data = *face_ptr_fe->dataOnElement[H1];
37 auto t_base = data.dataOnEntities[MBVERTEX][0].getFTensor0N(NOBASE);
39 &gaussPts(0, 0), &gaussPts(1, 0), &gaussPts(2, 0)};
40 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
41 for (int bb = 0; bb != nbNodesOnFace; ++bb) {
42 const int shift = 3 * faceConnMap[bb];
44 coords[shift + 0], coords[shift + 1], coords[shift + 2]};
45 t_gauss_coords(i) += t_base * t_coords(i);
46 ++t_base;
47 }
48 ++t_gauss_coords;
49 gaussPts(3, gg) = face_ptr_fe->gaussPts(2, gg);
50 }
52 };
53
54 switch (type) {
55 case MBTET:
56 CHKERR set_gauss(tet_coords);
57 break;
58 case MBHEX:
59 CHKERR set_gauss(hex_coords);
60 break;
61 default:
62 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
63 "Element type not implemented: %d", type);
64 }
65
67}
68
70VolumeElementForcesAndSourcesCoreOnSide::UserDataOperator::setPtrFE(
73 if (!(ptrFE = dynamic_cast<VolumeElementForcesAndSourcesCoreOnSide *>(ptr)))
74 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
75 "User operator and finite element do not work together");
77}
78
88
91
92 if (!sidePtrFE)
93 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Side element not set");
94
95 const auto type = numeredEntFiniteElementPtr->getEntType();
96 const auto nb_nodes_on_ele = CN::VerticesPerEntity(type);
97
98 auto face_ptr_fe = static_cast<FaceElementForcesAndSourcesCore *>(sidePtrFE);
99
100 const auto face_entity = sidePtrFE->numeredEntFiniteElementPtr->getEnt();
101 auto &side_table = const_cast<SideNumber_multiIndex &>(
102 numeredEntFiniteElementPtr->getSideNumberTable());
103 auto sit = side_table.get<0>().find(face_entity);
104 if (sit == side_table.get<0>().end())
105 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
106 "Face can not be found on volume element");
107
108 nbNodesOnFace = CN::VerticesPerEntity((*sit)->getEntType());
109
110 const auto ent = numeredEntFiniteElementPtr->getEnt();
111 CHKERR mField.get_moab().get_connectivity(ent, conn, num_nodes, true);
112
113 faceSense = (*sit)->sense;
114 faceSideNumber = (*sit)->side_number;
115 fill(&tetConnMap[0], &tetConnMap[nb_nodes_on_ele], -1);
116 for (int nn = 0; nn != nbNodesOnFace; ++nn) {
117 faceConnMap[nn] = std::distance(
118 conn, find(conn, &conn[nb_nodes_on_ele], face_ptr_fe->conn[nn]));
119 tetConnMap[faceConnMap[nn]] = nn;
120 if (faceConnMap[nn] >= nb_nodes_on_ele)
121 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
122 "No common node on face and element can not be found");
123 }
124
125 oppositeNode = std::distance(
126 &tetConnMap[0], find(&tetConnMap[0], &tetConnMap[nb_nodes_on_ele], -1));
127
128
130}
131
132} // namespace MoFEM
@ NOBASE
Definition definitions.h:59
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ H1
continuous field
Definition definitions.h:85
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
constexpr int order
multi_index_container< boost::shared_ptr< SideNumber >, indexed_by< ordered_unique< member< SideNumber, EntityHandle, &SideNumber::ent > >, ordered_non_unique< composite_key< SideNumber, const_mem_fun< SideNumber, EntityType, &SideNumber::getEntType >, member< SideNumber, signed char, &SideNumber::side_number > > > > > SideNumber_multiIndex
SideNumber_multiIndex for SideNumber.
FTensor::Index< 'i', SPACE_DIM > i
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
virtual moab::Interface & get_moab()=0
boost::shared_ptr< const NumeredEntFiniteElement > numeredEntFiniteElementPtr
structure to get information form mofem into EntitiesFieldData
ForcesAndSourcesCore * sidePtrFE
Element to integrate on the sides.
MatrixDouble gaussPts
Matrix of integration points.
MoFEMErrorCode operator()()
Get side and sense and call operator from derived class.
MoFEMErrorCode operator()()
function is run for every finite element