|
| v0.14.0
|
Auxiliary tools.
More...
#include <src/interfaces/Tools.hpp>
|
enum | SEGMENT_MIN_DISTANCE {
SOLUTION_EXIST,
SEGMENT_ONE_IS_POINT,
SEGMENT_TWO_IS_POINT,
SEGMENT_TWO_AND_TWO_ARE_POINT,
NO_SOLUTION
} |
|
static constexpr double | shapeFunMBEDGE0At00 = N_MBEDGE0(0) |
|
static constexpr double | shapeFunMBEDGE1At00 = N_MBEDGE1(0) |
|
static constexpr std::array< double, 2 > | shapeFunMBEDGEAt00 |
| Array of shape function at zero local point on reference element. More...
|
|
static constexpr double | diffN_MBEDGE0x = diffN_MBEDGE0 |
|
static constexpr double | diffN_MBEDGE1x = diffN_MBEDGE1 |
|
static constexpr std::array< double, 2 > | diffShapeFunMBEDGE |
|
static constexpr double | diffShapeFunMBTRI0x |
| derivative of triangle shape function More...
|
|
static constexpr double | diffShapeFunMBTRI0y |
| derivative of triangle shape function More...
|
|
static constexpr double | diffShapeFunMBTRI1x |
| derivative of triangle shape function More...
|
|
static constexpr double | diffShapeFunMBTRI1y |
| derivative of triangle shape function More...
|
|
static constexpr double | diffShapeFunMBTRI2x |
| derivative of triangle shape function More...
|
|
static constexpr double | diffShapeFunMBTRI2y |
| derivative of triangle shape function More...
|
|
static constexpr std::array< double, 6 > | diffShapeFunMBTRI |
|
static constexpr double | shapeFunMBTRI0At00 = N_MBTRI0(0, 0) |
|
static constexpr double | shapeFunMBTRI1At00 = N_MBTRI1(0, 0) |
|
static constexpr double | shapeFunMBTRI2At00 = N_MBTRI2(0, 0) |
|
static constexpr std::array< double, 3 > | shapeFunMBTRIAt00 |
| Array of shape function at zero local point on reference element. More...
|
|
static constexpr double | diffShapeFunMBQUADAtCenter0x |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBQUADAtCenter0y |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBQUADAtCenter1x |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBQUADAtCenter1y |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBQUADAtCenter2x |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBQUADAtCenter2y |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBQUADAtCenter3x |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBQUADAtCenter3y |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter0x |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter0y |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter0z |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter1x |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter1y |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter1z |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter2x |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter2y |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter2z |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter3x |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter3y |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter3z |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter4x |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter4y |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter4z |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter5x |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter5y |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter5z |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter6x |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter6y |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter6z |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter7x |
| derivative of HEX shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter7y |
| derivative of quad shape function More...
|
|
static constexpr double | diffShapeFunMBHEXAtCenter7z |
| derivative of quad shape function More...
|
|
static constexpr std::array< double, 8 > | diffShapeFunMBQUADAtCenter |
|
static constexpr std::array< double, 24 > | diffShapeFunMBHEXAtCenter |
|
static constexpr double | diffShapeFunMBTET0x |
| derivative of tetrahedral shape function More...
|
|
static constexpr double | diffShapeFunMBTET0y |
| derivative of tetrahedral shape function More...
|
|
static constexpr double | diffShapeFunMBTET0z |
| derivative of tetrahedral shape function More...
|
|
static constexpr double | diffShapeFunMBTET1x |
| derivative of tetrahedral shape function More...
|
|
static constexpr double | diffShapeFunMBTET1y |
| derivative of tetrahedral shape function More...
|
|
static constexpr double | diffShapeFunMBTET1z |
| derivative of tetrahedral shape function More...
|
|
static constexpr double | diffShapeFunMBTET2x |
| derivative of tetrahedral shape function More...
|
|
static constexpr double | diffShapeFunMBTET2y |
| derivative of tetrahedral shape function More...
|
|
static constexpr double | diffShapeFunMBTET2z |
| derivative of tetrahedral shape function More...
|
|
static constexpr double | diffShapeFunMBTET3x |
| derivative of tetrahedral shape function More...
|
|
static constexpr double | diffShapeFunMBTET3y |
| derivative of tetrahedral shape function More...
|
|
static constexpr double | diffShapeFunMBTET3z |
| derivative of tetrahedral shape function More...
|
|
static constexpr std::array< double, 12 > | diffShapeFunMBTET |
|
static constexpr double | shapeFunMBTET0At000 = N_MBTET0(0, 0, 0) |
|
static constexpr double | shapeFunMBTET1At000 = N_MBTET1(0, 0, 0) |
|
static constexpr double | shapeFunMBTET2At000 = N_MBTET2(0, 0, 0) |
|
static constexpr double | shapeFunMBTET3At000 = N_MBTET3(0, 0, 0) |
|
static constexpr double | shapeFunMBTET0AtOneThird |
|
static constexpr double | shapeFunMBTET1AtOneThird |
|
static constexpr double | shapeFunMBTET2AtOneThird |
|
static constexpr double | shapeFunMBTET3AtOneThird |
|
static constexpr std::array< double, 4 > | shapeFunMBTETAt000 |
| Array of shape function at zero local point on reference element. More...
|
|
static constexpr std::array< double, 4 > | shapeFunMBTETAtOneThird |
| Array of shape function at center on reference element. More...
|
|
MoFEMErrorCode | minTetsQuality (const Range &tets, double &min_quality, Tag th=nullptr, boost::function< double(double, double)> f=[](double a, double b) -> double { return std::min(a, b);}) |
| calculate minimal quality of tetrahedra in range More...
|
|
MoFEMErrorCode | getTetsWithQuality (Range &out_tets, const Range &tets, Tag th=nullptr, boost::function< bool(double)> f=[](double q) -> bool { if(q<=0) return true;else return false;}) |
| Get the Tets With Quality. More...
|
|
MoFEMErrorCode | writeTetsWithQuality (const char *file_name, const char *file_type, const char *options, const Range &tets, Tag th=nullptr, boost::function< bool(double)> f=[](double q) -> bool { if(q<=0) return true;else return false;}) |
| Write file with tetrahedral of given quality. More...
|
|
MoFEMErrorCode | getTriNormal (const EntityHandle tri, double *normal) const |
| Get triangle normal. More...
|
|
double | getTriArea (const EntityHandle tri) const |
| Get triangle area. More...
|
|
double | getEdgeLength (const EntityHandle edge) |
| Get edge length. More...
|
|
MoFEMErrorCode | findMinDistanceFromTheEdges (const double *v_ptr, const int nb, Range edges, double *min_dist_ptr, double *o_ptr=nullptr, EntityHandle *o_segments=nullptr) const |
| Find minimal distance to edges. More...
|
|
static double | volumeLengthQuality (const double *coords) |
| Calculate tetrahedron volume length quality. More...
|
|
static double | tetVolume (const double *coords) |
| Calculate volume of tetrahedron. More...
|
|
static double | shapeFunMBEDGE0 (const double x) |
|
static double | shapeFunMBEDGE1 (const double x) |
|
template<int LDB = 1> |
static MoFEMErrorCode | shapeFunMBEDGE (double *shape, const double *ksi, const int nb) |
| Calculate shape functions on edge. More...
|
|
static double | shapeFunMBTRI0 (const double x, const double y) |
|
static double | shapeFunMBTRI1 (const double x, const double y) |
|
static double | shapeFunMBTRI2 (const double x, const double y) |
|
template<int LDB = 1> |
static MoFEMErrorCode | shapeFunMBTRI (double *shape, const double *ksi, const double *eta, const int nb) |
| Calculate shape functions on triangle. More...
|
|
static double | shapeFunMBTET0 (const double x, const double y, const double z) |
|
static double | shapeFunMBTET1 (const double x, const double y, const double z) |
|
static double | shapeFunMBTET2 (const double x, const double y, const double z) |
|
static double | shapeFunMBTET3 (const double x, const double y, const double z) |
|
template<int LDB = 1> |
static MoFEMErrorCode | shapeFunMBTET (double *shape, const double *ksi, const double *eta, const double *zeta, const double nb) |
| Calculate shape functions on tetrahedron. More...
|
|
static MoFEMErrorCode | getLocalCoordinatesOnReferenceFourNodeTet (const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords) |
| Get the local coordinates on reference four node tet object. More...
|
|
static MoFEMErrorCode | getLocalCoordinatesOnReferenceThreeNodeTri (const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords) |
| Get the local coordinates on reference three node tri object. More...
|
|
static MoFEMErrorCode | getLocalCoordinatesOnReferenceThreeNodeTri (const double *elem_coords, const std::complex< double > *glob_coords, const int nb_nodes, std::complex< double > *local_coords) |
| Get the local coordinates on reference three node tri object. More...
|
|
static DEPRECATED MoFEMErrorCode | getLocalCoordinatesOnReferenceTriNodeTri (const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords) |
|
static MoFEMErrorCode | getLocalCoordinatesOnReferenceEdgeNodeEdge (const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords) |
| Get the local coordinates on reference four node tet object. More...
|
|
static MoFEMErrorCode | checkIfPointIsInTet (const double tet_coords[], const double global_coord[], const double tol, bool &result) |
| Check of point is in tetrahedral. More...
|
|
static MoFEMErrorCode | getTriNormal (const double *coords, double *normal, double *d_normal=nullptr) |
| Get the Tri Normal objectGet triangle normal. More...
|
|
static double | getEdgeLength (const double *edge_coords) |
| Get edge length. More...
|
|
static SEGMENT_MIN_DISTANCE | minDistancePointFromOnSegment (const double *w_ptr, const double *v_ptr, const double *p_ptr, double *const t_ptr=nullptr) |
| Find closet point on the segment from the point. More...
|
|
static SEGMENT_MIN_DISTANCE | minDistanceFromSegments (const double *w_ptr, const double *v_ptr, const double *k_ptr, const double *l_ptr, double *const tvw_ptr=nullptr, double *const tlk_ptr=nullptr) |
| Find points on two segments in closest distance. More...
|
|
Auxiliary tools.
Definition at line 19 of file Tools.hpp.
◆ RefineTrianglesReturn
◆ SEGMENT_MIN_DISTANCE
Enumerator |
---|
SOLUTION_EXIST | |
SEGMENT_ONE_IS_POINT | |
SEGMENT_TWO_IS_POINT | |
SEGMENT_TWO_AND_TWO_ARE_POINT | |
NO_SOLUTION | |
Definition at line 525 of file Tools.hpp.
◆ Tools()
◆ checkIfPointIsInTet()
Check of point is in tetrahedral.
- Parameters
-
tet_coords | |
global_coord | |
tol | |
result | |
- Returns
- MoFEMErrorCode
Definition at line 287 of file Tools.cpp.
290 double loc_coord[] = {0, 0, 0};
291 double N[4], diffN[12];
298 for (
int n = 0;
n != 4; ++
n) {
◆ checkVectorForNotANumber()
Print all DOFs for which element of vector is not a number.
Definition at line 307 of file Tools.cpp.
312 CHKERR VecGetLocalSize(
v, &loc_size);
313 int prb_loc_size = 0;
314 boost::shared_ptr<NumeredDofEntity_multiIndex> prb_dofs;
315 switch (row_or_col) {
317 prb_loc_size = prb_ptr->getNbLocalDofsRow();
318 prb_dofs = prb_ptr->getNumeredRowDofsPtr();
321 prb_loc_size = prb_ptr->getNbLocalDofsCol();
322 prb_dofs = prb_ptr->getNumeredColDofsPtr();
327 "Wrong argument, row_or_col should be row or column");
329 if (loc_size != prb_loc_size) {
331 "Inconsistent size of vector and problem %d != %d", loc_size,
336 MPI_Comm comm = PetscObjectComm((PetscObject)
v);
337 for (
int ii = 0; ii != loc_size; ++ii) {
338 if (!boost::math::isfinite(
a[ii])) {
339 NumeredDofEntityByLocalIdx::iterator dit =
340 prb_dofs->get<PetscLocalIdx_mi_tag>().find(ii);
341 std::ostringstream ss;
342 ss <<
"Not a number " <<
a[ii] <<
" on dof: " << endl
345 PetscSynchronizedPrintf(comm,
"%s", ss.str().c_str());
349 PetscSynchronizedFlush(comm, PETSC_STDOUT);
◆ findMinDistanceFromTheEdges()
Find minimal distance to edges.
- Note
- Finding only edges with have smaller distance than distance set on the input by min_dist_ptr
- Parameters
-
v_ptr | point coordinates |
nb | nb points |
edges | range of edges |
min_dist_ptr | on return minimal distance, on input starting distance |
o_ptr | coordinates of the point on edge |
o_segments | closest segments |
- Returns
- MoFEMErrorCode
Definition at line 535 of file Tools.cpp.
544 auto get_point = [
i](
auto &t_w,
auto &t_delta,
auto t) {
546 t = std::max(0., std::min(1.,
t));
547 t_p(
i) = t_w(
i) +
t * t_delta(
i);
551 auto get_distance = [
i](
auto &t_p,
auto &t_n) {
553 t_dist_vector(
i) = t_p(
i) - t_n(
i);
554 return sqrt(t_dist_vector(
i) * t_dist_vector(
i));
557 for (
auto e : edges) {
561 CHKERR moab.get_connectivity(e, conn_fixed, num_nodes,
true);
563 CHKERR moab.get_coords(conn_fixed, num_nodes, &coords_fixed[0]);
565 &coords_fixed[0], &coords_fixed[1], &coords_fixed[2]);
567 &coords_fixed[3], &coords_fixed[4], &coords_fixed[5]);
570 t_edge_delta(
i) = t_f1(
i) - t_f0(
i);
573 v_ptr, v_ptr + 1, v_ptr + 2);
575 o_ptr, o_ptr + 1, o_ptr + 2);
580 colsest_segment_it = o_segments;
582 for (
int n = 0;
n != nb; ++
n) {
587 auto t_p = get_point(t_f0, t_edge_delta,
t);
588 auto dist_n = get_distance(t_p, t_n);
589 if (dist_n < t_min_dist || t_min_dist < 0) {
592 t_min_coords(
i) = t_p(
i);
594 *colsest_segment_it = e;
603 ++colsest_segment_it;
◆ getEdgeLength() [1/2]
double MoFEM::Tools::getEdgeLength |
( |
const double * |
edge_coords | ) |
|
|
static |
Get edge length.
- Parameters
-
- Returns
- double
Definition at line 415 of file Tools.cpp.
421 t_coords_n0(
i) -= t_coords_n1(
i);
422 return sqrt(t_coords_n0(
i) * t_coords_n0(
i));
◆ getEdgeLength() [2/2]
Get edge length.
- Parameters
-
- Returns
- double
Definition at line 425 of file Tools.cpp.
428 auto get_edge_coords = [edge, &moab](
double *
const coords) {
435 CHKERR moab.get_connectivity(edge, conn, num_nodes,
true);
436 CHKERR moab.get_coords(conn, 2, coords);
440 ierr = get_edge_coords(coords);
441 CHKERRABORT(PETSC_COMM_SELF,
ierr);
◆ getLocalCoordinatesOnReferenceEdgeNodeEdge()
MoFEMErrorCode MoFEM::Tools::getLocalCoordinatesOnReferenceEdgeNodeEdge |
( |
const double * |
elem_coords, |
|
|
const double * |
glob_coords, |
|
|
const int |
nb_nodes, |
|
|
double * |
local_coords |
|
) |
| |
|
static |
Get the local coordinates on reference four node tet object.
&elem_coords(0, 0), &global_coords(0, 0), global_coords.size1(),
&local_coords(0, 0))
- Parameters
-
elem_coords | Global element node coordinates |
glob_coords | Global coordinates |
nb_nodes | Number of points |
local_coords | Result |
- Returns
- MoFEMErrorCode
Definition at line 203 of file Tools.cpp.
211 &elem_coords[0], &elem_coords[3], &elem_coords[6]};
217 for (
auto ii : {0, 1, 2}) {
226 &global_coords[0], &global_coords[1], &global_coords[2]};
230 const double b = t_a(i3) * t_a(i3);
231 const double inv_b = 1 / b;
234 for (
int ii = 0; ii != nb_nodes; ++ii) {
236 inv_b * (t_a(i3) * (t_global_coords(i3) - t_coords_at_0(i3)));
◆ getLocalCoordinatesOnReferenceFourNodeTet()
MoFEMErrorCode MoFEM::Tools::getLocalCoordinatesOnReferenceFourNodeTet |
( |
const double * |
elem_coords, |
|
|
const double * |
glob_coords, |
|
|
const int |
nb_nodes, |
|
|
double * |
local_coords |
|
) |
| |
|
static |
Get the local coordinates on reference four node tet object.
&elem_coords(0, 0), &global_coords(0, 0), global_coords.size1(),
&local_coords(0, 0))
- Parameters
-
elem_coords | Global element node coordinates |
glob_coords | Global coordinates |
nb_nodes | Number of points |
local_coords | Result |
- Returns
- MoFEMErrorCode
Definition at line 88 of file Tools.cpp.
95 &elem_coords[0], &elem_coords[3], &elem_coords[6], &elem_coords[9]};
106 for (
auto ii : {0, 1, 2}) {
110 for (
auto jj : {0, 1, 2}) {
111 a(jj, ii) = t_diff(
i) * t_elem_coords(
i);
114 t_coords_at_0(ii) = t_n(
i) * t_elem_coords(
i);
119 &global_coords[0], &global_coords[1], &global_coords[2]};
121 &local_coords[0], &local_coords[1], &local_coords[2]};
125 for (
int ii = 0; ii != nb_nodes; ++ii) {
126 t_local_coords(
j) = t_global_coords(
j) - t_coords_at_0(
j);
133 int info =
lapack_dgesv(3, nb_nodes, &
a(0, 0), 3, IPIV, local_coords, 3);
135 SETERRQ1(PETSC_COMM_SELF, 1,
"info == %d", info);
◆ getLocalCoordinatesOnReferenceThreeNodeTri() [1/2]
MoFEMErrorCode MoFEM::Tools::getLocalCoordinatesOnReferenceThreeNodeTri |
( |
const double * |
elem_coords, |
|
|
const double * |
glob_coords, |
|
|
const int |
nb_nodes, |
|
|
double * |
local_coords |
|
) |
| |
|
static |
Get the local coordinates on reference three node tri object.
&elem_coords(0, 0), &global_coords(0, 0), global_coords.size1(),
&local_coords(0, 0))
- Parameters
-
elem_coords | Global element node coordinates |
glob_coords | Global coordinates |
nb_nodes | Number of points |
local_coords | Result |
d_elem_coords | Derivative of local coordinates |
d_global_coords | |
- Returns
- MoFEMErrorCode
Definition at line 188 of file Tools.cpp.
191 return getLocalCoordinatesOnReferenceThreeNodeTriImpl<double, double>(
192 elem_coords, glob_coords, nb_nodes, local_coords);
◆ getLocalCoordinatesOnReferenceThreeNodeTri() [2/2]
MoFEMErrorCode MoFEM::Tools::getLocalCoordinatesOnReferenceThreeNodeTri |
( |
const double * |
elem_coords, |
|
|
const std::complex< double > * |
glob_coords, |
|
|
const int |
nb_nodes, |
|
|
std::complex< double > * |
local_coords |
|
) |
| |
|
static |
Get the local coordinates on reference three node tri object.
&elem_coords(0, 0), &global_coords(0, 0), global_coords.size1(),
&local_coords(0, 0))
- Parameters
-
elem_coords | Global element node coordinates |
glob_coords | Global coordinates |
nb_nodes | Number of points |
local_coords | Result |
d_elem_coords | Derivative of local coordinates |
d_global_coords | |
- Returns
- MoFEMErrorCode
Definition at line 195 of file Tools.cpp.
199 std::complex<double>>(
200 elem_coords, glob_coords, nb_nodes, local_coords);
◆ getLocalCoordinatesOnReferenceTriNodeTri()
- Deprecated:
- use getLocalCoordinatesOnReferenceThreeNodeTri
Definition at line 398 of file Tools.hpp.
402 nb_nodes, local_coords);
◆ getTetsWithQuality()
Get the Tets With Quality.
- Parameters
-
- Returns
- MoFEMErrorCode
Definition at line 244 of file Tools.cpp.
254 for (
auto tet : tets) {
255 CHKERR m_field.
get_moab().get_connectivity(tet, conn, num_nodes,
true);
257 CHKERR moab.tag_get_data(
th, conn, num_nodes, coords);
259 CHKERR moab.get_coords(conn, num_nodes, coords);
263 out_tets.insert(tet);
◆ getTriArea()
Get triangle area.
- Parameters
-
- Returns
- double
Definition at line 408 of file Tools.cpp.
412 return sqrt(t_normal(
i) * t_normal(
i)) * 0.5;
◆ getTriNormal() [1/2]
Get the Tri Normal objectGet triangle normal.
- Parameters
-
coords | |
normal | |
d_normal | derbiative, if pointer is null, derbiative is not calculated |
- Returns
- MoFEMErrorCode
Definition at line 353 of file Tools.cpp.
366 auto t_diff_tensor = getFTensor2FromPtr<3, 2>(
const_cast<double *
>(diff_ptr));
367 auto t_coords = getFTensor2FromPtr<3, 3>(
const_cast<double *
>(coords));
369 t_tangent(
i,
J) = t_coords(
n,
i) * t_diff_tensor(
n,
J);
370 auto t_normal = getFTensor1FromPtr<3>(normal);
378 t_d_tangent(
i,
k,
l,
J) = t_d_coords(
n,
i,
k,
l) * t_diff_tensor(
n,
J);
381 t_d_tangent(
k,
m,
n, N0)
386 t_d_tangent(
i,
m,
n, N1);
◆ getTriNormal() [2/2]
Get triangle normal.
- Parameters
-
- Returns
- MoFEMErrorCode
Definition at line 391 of file Tools.cpp.
402 CHKERR moab.get_connectivity(tri, conn, num_nodes,
true);
403 CHKERR moab.get_coords(conn, num_nodes, coords);
◆ minDistanceFromSegments()
Find points on two segments in closest distance.
- Parameters
-
w_ptr | |
v_ptr | |
k_ptr | |
l_ptr | |
tvw_ptr | |
tlk_ptr | |
- Returns
- SEGMENT_MIN_DISTANCE
- Note
- If tvwk or tlk are outside bound [0,-1], it means that points are on the lines beyond segments, respectively for segment vw and lk.
Definition at line 469 of file Tools.cpp.
481 t_vw(
i) = t_v(
i) - t_w(
i);
482 double dot_vw = t_vw(
i) * t_vw(
i);
483 if (std::abs(dot_vw) < std::numeric_limits<double>::epsilon()) {
495 t_lk(
i) = t_l(
i) - t_k(
i);
496 double dot_lk = t_lk(
i) * t_lk(
i);
497 if (std::abs(dot_lk) < std::numeric_limits<double>::epsilon()) {
507 const double a = t_vw(
i) * t_vw(
i);
508 const double b = -t_vw(
i) * t_lk(
i);
509 const double c = t_lk(
i) * t_lk(
i);
511 const double det =
a *
c - b * b;
512 if (std::abs(det) < std::numeric_limits<double>::epsilon()) {
519 t_wk(
i) = t_w(
i) - t_k(
i);
521 const double ft0 = t_vw(
i) * t_wk(
i);
522 const double ft1 = -t_lk(
i) * t_wk(
i);
523 const double t0 = (ft1 * b - ft0 *
c) / det;
524 const double t1 = (ft0 * b - ft1 *
a) / det;
◆ minDistancePointFromOnSegment()
Find closet point on the segment from the point.
- Parameters
-
w_ptr | segment first vertex coordinate |
v_ptr | segment second vertex coordinate |
p_ptr | coordinate of point |
t_ptr | distance on the segment |
- Note
- If t is outside bounds [ 0,-1 ] point is on the line point beyond segment.
double p[] = {0, 1, 0};
CHKERR Toolas::minDistancePointFromOnSegment(
w,
v, p, &
t);
double point_on_segment[3];
for (
int i = 0;
i != 3; ++
i)
point_on_segment[
i] =
w[
i] +
t * (
v[
i] -
w[
i]);
- Returns
- SEGMENT_MIN_DISTANCE
Definition at line 446 of file Tools.cpp.
453 t_vw(
i) = t_v(
i) - t_w(
i);
454 const double dot_vw = t_vw(
i) * t_vw(
i);
455 if (std::abs(dot_vw) < std::numeric_limits<double>::epsilon()) {
461 t_pw(
i) = t_p(
i) - t_w(
i);
462 const double t = t_pw(
i) * t_vw(
i) / dot_vw;
◆ minTetsQuality()
calculate minimal quality of tetrahedra in range
- Parameters
-
tets | range |
min_quality | minimal quality |
- Returns
- error code
Definition at line 49 of file Tools.cpp.
57 for (
auto tet : tets) {
58 CHKERR m_field.
get_moab().get_connectivity(tet, conn, num_nodes,
true);
60 CHKERR moab.tag_get_data(
th, conn, num_nodes, coords);
62 CHKERR moab.get_coords(conn, num_nodes, coords);
65 if (!std::isnormal(q))
67 min_quality =
f(q, min_quality);
◆ outerProductOfEdgeIntegrationPtsForHex()
MoFEMErrorCode MoFEM::Tools::outerProductOfEdgeIntegrationPtsForHex |
( |
MatrixDouble & |
pts, |
|
|
const int |
edge0, |
|
|
const int |
edge1, |
|
|
const int |
edge2 |
|
) |
| |
|
static |
- Examples
- gauss_points_on_outer_product.cpp.
Definition at line 659 of file Tools.cpp.
664 auto check_rule_edge = [](
int rule) {
668 "Wrong integration rule: %d", rule);
684 CHKERR check_rule_edge(rule_ksi);
685 CHKERR check_rule_edge(rule_eta);
686 CHKERR check_rule_edge(rule_zeta);
691 gauss_pts.resize(4, nb_gauss_pts_ksi * nb_gauss_pts_eta * nb_gauss_pts_zeta,
695 for (
size_t i = 0;
i != nb_gauss_pts_ksi; ++
i) {
698 for (
size_t j = 0;
j != nb_gauss_pts_eta; ++
j) {
701 for (
size_t k = 0;
k != nb_gauss_pts_zeta; ++
k, ++gg) {
704 gauss_pts(0, gg) = ksi;
705 gauss_pts(1, gg) =
eta;
706 gauss_pts(2, gg) =
zeta;
707 gauss_pts(3, gg) = wk;
712 if (gg != gauss_pts.size2())
◆ outerProductOfEdgeIntegrationPtsForQuad()
- Examples
- gauss_points_on_outer_product.cpp.
Definition at line 610 of file Tools.cpp.
614 auto check_rule_edge = [](
int rule) {
618 "Wrong integration rule: %d", rule);
634 CHKERR check_rule_edge(rule_ksi);
635 CHKERR check_rule_edge(rule_eta);
639 gauss_pts.resize(3, nb_gauss_pts_ksi * nb_gauss_pts_eta,
false);
642 for (
size_t i = 0;
i != nb_gauss_pts_ksi; ++
i) {
645 for (
size_t j = 0;
j != nb_gauss_pts_eta; ++
j, ++gg) {
648 gauss_pts(0, gg) = ksi;
649 gauss_pts(1, gg) =
eta;
650 gauss_pts(2, gg) = wk;
653 if (gg != gauss_pts.size2())
◆ query_interface()
◆ refineTriangle()
std::tuple< std::vector< double >, std::vector< int >, std::vector< int > > MoFEM::Tools::refineTriangle |
( |
int |
nb_levels | ) |
|
|
static |
create uniform triangle mesh of refined elements
- Parameters
-
- Returns
- RefineTrianglesReturn
- Examples
- gauss_points_on_outer_product.cpp.
Definition at line 721 of file Tools.cpp.
723 std::vector<int> triangles{0, 1, 2, 3, 4, 5};
724 std::vector<double> nodes{
734 std::map<std::pair<int, int>,
int> edges{
735 {{0, 1}, 3}, {{1, 2}, 4}, {{0, 2}, 5}};
737 auto add_edge = [&](
auto a,
auto b) {
741 auto it = edges.find(std::make_pair(
a, b));
742 if (it == edges.end()) {
743 int e = 3 + edges.size();
744 edges[std::make_pair(
a, b)] = e;
745 for (
auto n : {0, 1}) {
746 nodes.push_back((nodes[2 *
a +
n] + nodes[2 * b +
n]) / 2);
749 if (e != nodes.size() / 2 - 1)
757 auto add_triangle = [&](
auto t) {
758 for (
auto tt : {0, 1, 2, 3}) {
759 auto last = triangles.size() / 6;
760 for (
auto n : {0, 1, 2}) {
766 auto cycle = std::array<int, 4>{0, 1, 2, 0};
767 for (
auto e : {0, 1, 2}) {
768 triangles.push_back(add_edge(triangles[6 * last + cycle[e]],
769 triangles[6 * last + cycle[e + 1]]));
774 std::vector<int> level_index{0, 1};
777 auto first_tet = level_index[
l];
778 auto nb_last_level_test = level_index.back() - level_index[
l];
779 for (
auto t = first_tet;
t != (first_tet + nb_last_level_test); ++
t) {
782 level_index.push_back(triangles.size() / 6);
785 return std::make_tuple(nodes, triangles, level_index);
◆ refineTriangleIntegrationPts() [1/2]
generate integration points for refined triangle mesh for last level
- Parameters
-
rule | Gauss integration rule |
refined | |
- Returns
- MatrixDouble
Definition at line 853 of file Tools.cpp.
858 gauss_pts.resize(3, nb_gauss_pts,
false);
859 cblas_dcopy(nb_gauss_pts, &
QUAD_2D_TABLE[rule]->points[1], 3,
860 &gauss_pts(0, 0), 1);
861 cblas_dcopy(nb_gauss_pts, &
QUAD_2D_TABLE[rule]->points[2], 3,
862 &gauss_pts(1, 0), 1);
863 cblas_dcopy(nb_gauss_pts,
QUAD_2D_TABLE[rule]->weights, 1, &gauss_pts(2, 0),
◆ refineTriangleIntegrationPts() [2/2]
◆ shapeFunMBEDGE()
Calculate shape functions on edge.
- Note
- Template parameter is leading dimension of point coordinate arrays, such that \(ksi_{n+1} = ksi[n + LDB]\)
- Template Parameters
-
- Parameters
-
shape | shape functions |
ksi | pointer to first local coordinates |
nb | number of points |
- Returns
- MoFEMErrorCode
Definition at line 682 of file Tools.hpp.
685 for (
int n = 0;
n != nb; ++
n) {
◆ shapeFunMBEDGE0()
◆ shapeFunMBEDGE1()
◆ shapeFunMBTET()
Calculate shape functions on tetrahedron.
- Note
- Template parameter is leading dimension of point coordinate arrays, such that \(ksi_{n+1} = ksi[n + LDB]\)
- Template Parameters
-
- Parameters
-
shape | shape functions |
ksi | pointer to first local coordinates |
eta | pointer to second local coordinates |
zeta | pointer to first third coordinates |
nb | number of points |
- Returns
- MoFEMErrorCode
Definition at line 738 of file Tools.hpp.
742 for (
int n = 0;
n != nb; ++
n) {
◆ shapeFunMBTET0()
◆ shapeFunMBTET1()
◆ shapeFunMBTET2()
◆ shapeFunMBTET3()
◆ shapeFunMBTRI()
Calculate shape functions on triangle.
- Note
- Template parameter is leading dimension of point coordinate arrays, such that \(ksi_{n+1} = ksi[n + LDB]\)
- Template Parameters
-
- Parameters
-
shape | shape functions |
ksi | pointer to first local coordinates |
eta | pointer to second local coordinates |
nb | number of points |
- Returns
- MoFEMErrorCode
Definition at line 707 of file Tools.hpp.
710 for (
int n = 0;
n != nb; ++
n) {
◆ shapeFunMBTRI0()
◆ shapeFunMBTRI1()
◆ shapeFunMBTRI2()
◆ tetVolume()
Calculate volume of tetrahedron.
- Parameters
-
- Returns
- double volume
- Examples
- EshelbianPlasticity.cpp.
Definition at line 30 of file Tools.cpp.
40 for (
int nn = 0; nn != 4; nn++) {
41 jac(
i,
j) += t_coords(
i) * t_diff_n(
j);
◆ volumeLengthQuality()
double MoFEM::Tools::volumeLengthQuality |
( |
const double * |
coords | ) |
|
|
static |
Calculate tetrahedron volume length quality.
- Parameters
-
- Returns
- Volume-length quality
- Examples
- mesh_smoothing.cpp.
Definition at line 15 of file Tools.cpp.
17 for (
int dd = 0;
dd != 3;
dd++) {
18 lrms += pow(coords[0 * 3 +
dd] - coords[1 * 3 +
dd], 2) +
19 pow(coords[0 * 3 +
dd] - coords[2 * 3 +
dd], 2) +
20 pow(coords[0 * 3 +
dd] - coords[3 * 3 +
dd], 2) +
21 pow(coords[1 * 3 +
dd] - coords[2 * 3 +
dd], 2) +
22 pow(coords[1 * 3 +
dd] - coords[3 * 3 +
dd], 2) +
23 pow(coords[2 * 3 +
dd] - coords[3 * 3 +
dd], 2);
25 lrms = sqrt((1. / 6.) * lrms);
27 return 6. * sqrt(2.) * volume / pow(lrms, 3);
◆ writeTetsWithQuality()
MoFEMErrorCode MoFEM::Tools::writeTetsWithQuality |
( |
const char * |
file_name, |
|
|
const char * |
file_type, |
|
|
const char * |
options, |
|
|
const Range & |
tets, |
|
|
Tag |
th = nullptr , |
|
|
boost::function< bool(double)> |
f = [](double q) -> bool { if (q <= 0) return true; else return false; } |
|
) |
| |
Write file with tetrahedral of given quality.
- Parameters
-
file_name | |
file_type | |
options | |
tets | |
th | |
f | |
- Returns
- MoFEMErrorCode
Definition at line 269 of file Tools.cpp.
280 CHKERR moab.create_meshset(MESHSET_SET, meshset);
281 CHKERR moab.add_entities(meshset, out_tets);
282 CHKERR moab.write_file(file_name, file_type, options, &meshset, 1);
283 CHKERR moab.delete_entities(&meshset, 1);
◆ cOre
◆ diffN_MBEDGE0x
◆ diffN_MBEDGE1x
◆ diffShapeFunMBEDGE
constexpr std::array< double, 2 > MoFEM::Tools::diffShapeFunMBEDGE |
|
staticconstexpr |
◆ diffShapeFunMBHEXAtCenter
constexpr std::array< double, 24 > MoFEM::Tools::diffShapeFunMBHEXAtCenter |
|
staticconstexpr |
◆ diffShapeFunMBHEXAtCenter0x
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter0x |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 163 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter0y
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter0y |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 165 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter0z
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter0z |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 167 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter1x
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter1x |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 169 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter1y
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter1y |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 171 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter1z
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter1z |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 173 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter2x
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter2x |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 175 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter2y
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter2y |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 177 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter2z
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter2z |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 179 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter3x
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter3x |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 181 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter3y
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter3y |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 183 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter3z
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter3z |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 185 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter4x
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter4x |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 187 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter4y
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter4y |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 189 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter4z
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter4z |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 191 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter5x
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter5x |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 193 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter5y
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter5y |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 195 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter5z
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter5z |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 197 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter6x
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter6x |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 199 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter6y
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter6y |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 201 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter6z
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter6z |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 203 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter7x
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter7x |
|
staticconstexpr |
Initial value:
derivative of HEX shape function
Definition at line 205 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter7y
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter7y |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 207 of file Tools.hpp.
◆ diffShapeFunMBHEXAtCenter7z
constexpr double MoFEM::Tools::diffShapeFunMBHEXAtCenter7z |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 209 of file Tools.hpp.
◆ diffShapeFunMBQUADAtCenter
constexpr std::array< double, 8 > MoFEM::Tools::diffShapeFunMBQUADAtCenter |
|
staticconstexpr |
◆ diffShapeFunMBQUADAtCenter0x
constexpr double MoFEM::Tools::diffShapeFunMBQUADAtCenter0x |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 146 of file Tools.hpp.
◆ diffShapeFunMBQUADAtCenter0y
constexpr double MoFEM::Tools::diffShapeFunMBQUADAtCenter0y |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 148 of file Tools.hpp.
◆ diffShapeFunMBQUADAtCenter1x
constexpr double MoFEM::Tools::diffShapeFunMBQUADAtCenter1x |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 150 of file Tools.hpp.
◆ diffShapeFunMBQUADAtCenter1y
constexpr double MoFEM::Tools::diffShapeFunMBQUADAtCenter1y |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 152 of file Tools.hpp.
◆ diffShapeFunMBQUADAtCenter2x
constexpr double MoFEM::Tools::diffShapeFunMBQUADAtCenter2x |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 154 of file Tools.hpp.
◆ diffShapeFunMBQUADAtCenter2y
constexpr double MoFEM::Tools::diffShapeFunMBQUADAtCenter2y |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 156 of file Tools.hpp.
◆ diffShapeFunMBQUADAtCenter3x
constexpr double MoFEM::Tools::diffShapeFunMBQUADAtCenter3x |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 158 of file Tools.hpp.
◆ diffShapeFunMBQUADAtCenter3y
constexpr double MoFEM::Tools::diffShapeFunMBQUADAtCenter3y |
|
staticconstexpr |
Initial value:
derivative of quad shape function
Definition at line 160 of file Tools.hpp.
◆ diffShapeFunMBTET
constexpr std::array< double, 12 > MoFEM::Tools::diffShapeFunMBTET |
|
staticconstexpr |
◆ diffShapeFunMBTET0x
constexpr double MoFEM::Tools::diffShapeFunMBTET0x |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 246 of file Tools.hpp.
◆ diffShapeFunMBTET0y
constexpr double MoFEM::Tools::diffShapeFunMBTET0y |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 248 of file Tools.hpp.
◆ diffShapeFunMBTET0z
constexpr double MoFEM::Tools::diffShapeFunMBTET0z |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 250 of file Tools.hpp.
◆ diffShapeFunMBTET1x
constexpr double MoFEM::Tools::diffShapeFunMBTET1x |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 252 of file Tools.hpp.
◆ diffShapeFunMBTET1y
constexpr double MoFEM::Tools::diffShapeFunMBTET1y |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 254 of file Tools.hpp.
◆ diffShapeFunMBTET1z
constexpr double MoFEM::Tools::diffShapeFunMBTET1z |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 256 of file Tools.hpp.
◆ diffShapeFunMBTET2x
constexpr double MoFEM::Tools::diffShapeFunMBTET2x |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 258 of file Tools.hpp.
◆ diffShapeFunMBTET2y
constexpr double MoFEM::Tools::diffShapeFunMBTET2y |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 260 of file Tools.hpp.
◆ diffShapeFunMBTET2z
constexpr double MoFEM::Tools::diffShapeFunMBTET2z |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 262 of file Tools.hpp.
◆ diffShapeFunMBTET3x
constexpr double MoFEM::Tools::diffShapeFunMBTET3x |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 264 of file Tools.hpp.
◆ diffShapeFunMBTET3y
constexpr double MoFEM::Tools::diffShapeFunMBTET3y |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 266 of file Tools.hpp.
◆ diffShapeFunMBTET3z
constexpr double MoFEM::Tools::diffShapeFunMBTET3z |
|
staticconstexpr |
Initial value:
derivative of tetrahedral shape function
Definition at line 268 of file Tools.hpp.
◆ diffShapeFunMBTRI
constexpr std::array< double, 6 > MoFEM::Tools::diffShapeFunMBTRI |
|
staticconstexpr |
◆ diffShapeFunMBTRI0x
constexpr double MoFEM::Tools::diffShapeFunMBTRI0x |
|
staticconstexpr |
Initial value:
derivative of triangle shape function
Definition at line 91 of file Tools.hpp.
◆ diffShapeFunMBTRI0y
constexpr double MoFEM::Tools::diffShapeFunMBTRI0y |
|
staticconstexpr |
Initial value:
derivative of triangle shape function
Definition at line 93 of file Tools.hpp.
◆ diffShapeFunMBTRI1x
constexpr double MoFEM::Tools::diffShapeFunMBTRI1x |
|
staticconstexpr |
Initial value:
derivative of triangle shape function
Definition at line 95 of file Tools.hpp.
◆ diffShapeFunMBTRI1y
constexpr double MoFEM::Tools::diffShapeFunMBTRI1y |
|
staticconstexpr |
Initial value:
derivative of triangle shape function
Definition at line 97 of file Tools.hpp.
◆ diffShapeFunMBTRI2x
constexpr double MoFEM::Tools::diffShapeFunMBTRI2x |
|
staticconstexpr |
Initial value:
derivative of triangle shape function
Definition at line 99 of file Tools.hpp.
◆ diffShapeFunMBTRI2y
constexpr double MoFEM::Tools::diffShapeFunMBTRI2y |
|
staticconstexpr |
Initial value:
derivative of triangle shape function
Definition at line 101 of file Tools.hpp.
◆ shapeFunMBEDGE0At00
◆ shapeFunMBEDGE1At00
◆ shapeFunMBEDGEAt00
constexpr std::array< double, 2 > MoFEM::Tools::shapeFunMBEDGEAt00 |
|
staticconstexpr |
Initial value:
Array of shape function at zero local point on reference element.
Definition at line 62 of file Tools.hpp.
◆ shapeFunMBTET0At000
constexpr double MoFEM::Tools::shapeFunMBTET0At000 = N_MBTET0(0, 0, 0) |
|
staticconstexpr |
◆ shapeFunMBTET0AtOneThird
constexpr double MoFEM::Tools::shapeFunMBTET0AtOneThird |
|
staticconstexpr |
◆ shapeFunMBTET1At000
constexpr double MoFEM::Tools::shapeFunMBTET1At000 = N_MBTET1(0, 0, 0) |
|
staticconstexpr |
◆ shapeFunMBTET1AtOneThird
constexpr double MoFEM::Tools::shapeFunMBTET1AtOneThird |
|
staticconstexpr |
◆ shapeFunMBTET2At000
constexpr double MoFEM::Tools::shapeFunMBTET2At000 = N_MBTET2(0, 0, 0) |
|
staticconstexpr |
◆ shapeFunMBTET2AtOneThird
constexpr double MoFEM::Tools::shapeFunMBTET2AtOneThird |
|
staticconstexpr |
◆ shapeFunMBTET3At000
constexpr double MoFEM::Tools::shapeFunMBTET3At000 = N_MBTET3(0, 0, 0) |
|
staticconstexpr |
◆ shapeFunMBTET3AtOneThird
constexpr double MoFEM::Tools::shapeFunMBTET3AtOneThird |
|
staticconstexpr |
◆ shapeFunMBTETAt000
constexpr std::array< double, 4 > MoFEM::Tools::shapeFunMBTETAt000 |
|
staticconstexpr |
Initial value:
Array of shape function at zero local point on reference element.
Definition at line 330 of file Tools.hpp.
◆ shapeFunMBTETAtOneThird
constexpr std::array<double, 4> MoFEM::Tools::shapeFunMBTETAtOneThird |
|
staticconstexpr |
Initial value:
Array of shape function at center on reference element.
Definition at line 338 of file Tools.hpp.
◆ shapeFunMBTRI0At00
◆ shapeFunMBTRI1At00
◆ shapeFunMBTRI2At00
◆ shapeFunMBTRIAt00
constexpr std::array< double, 3 > MoFEM::Tools::shapeFunMBTRIAt00 |
|
staticconstexpr |
Initial value:
Array of shape function at zero local point on reference element.
Definition at line 120 of file Tools.hpp.
◆ uniformTriangleRefineTriangles
constexpr std::array< int, 12 > MoFEM::Tools::uniformTriangleRefineTriangles |
|
staticconstexpr |
Initial value:= {
0, 3, 5,
3, 1, 4,
5, 4, 2,
5, 3, 4
}
Definition at line 629 of file Tools.hpp.
The documentation for this struct was generated from the following files:
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
VectorBoundedArray< double, 6 > VectorDouble6
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
static QUAD *const QUAD_2D_TABLE[]
static __CLPK_integer lapack_dgesv(__CLPK_integer n, __CLPK_integer nrhs, __CLPK_doublereal *a, __CLPK_integer lda, __CLPK_integer *ipiv, __CLPK_doublereal *b, __CLPK_integer ldb)
UBlasMatrix< double > MatrixDouble
MoFEMErrorCode getLocalCoordinatesOnReferenceThreeNodeTriImpl(const T1 *elem_coords, const T2 *global_coords, const int nb_nodes, typename FTensor::promote< T1, T2 >::V *local_coords)
FTensor::Index< 'J', DIM1 > J
constexpr std::enable_if<(Dim0<=2 &&Dim1<=2), Tensor2_Expr< Levi_Civita< T >, T, Dim0, Dim1, i, j > >::type levi_civita(const Index< i, Dim0 > &, const Index< j, Dim1 > &)
levi_civita functions to make for easy adhoc use
double zeta
Viscous hardening.
Deprecated interface functions.
DeprecatedCoreInterface Interface
const double c
speed of light (cm/ns)
static QUAD *const QUAD_1D_TABLE[]
#define CHKERR
Inline error check.
virtual moab::Interface & get_moab()=0
auto type_from_handle(const EntityHandle h)
get type from entity handle
constexpr double t
plate stiffness
FTensor::Index< 'i', SPACE_DIM > i
static const double edge_coords[6][6]
static auto determinantTensor3by3(T &t)
Calculate the determinant of a 3x3 matrix or a tensor of rank 2.
const double v
phase velocity of light in medium (cm/ns)
#define QUAD_1D_TABLE_SIZE
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)
FTensor::Index< 'j', 3 > j
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
@ MOFEM_DATA_INCONSISTENCY
MatrixBoundedArray< double, 9 > MatrixDouble3by3
FTensor::Index< 'm', 3 > m
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
FTensor::Index< 'k', 3 > k
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FTensor::Tensor3< FTensor::PackPtr< double *, 27 >, 3, 3, 3 > getFTensor3FromPtr< 3, 3, 3 >(double *ptr)
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
FTensor::Index< 'l', 3 > l