v0.14.0
Tools.hpp
Go to the documentation of this file.
1 /** \file Tools.hpp
2  * \brief Tools interface
3  *
4  * Implementation of some useful and very often used methods * in MoFEM.
5  *
6  * \ingroup mofem_tools
7  */
8 
9 #ifndef __TOOLS_HPP__
10 #define __TOOLS_HPP__
11 
12 namespace MoFEM {
13 
14 /**
15  * \brief Auxiliary tools
16  * \nosubgrouping
17  * \ingroup mofem_tools
18  */
19 struct Tools : public UnknownInterface {
20 
21  MoFEMErrorCode query_interface(boost::typeindex::type_index type_index,
22  UnknownInterface **iface) const;
23 
25  Tools(const MoFEM::Core &core) : cOre(const_cast<MoFEM::Core &>(core)) {}
26 
27  /** \name Computational */
28 
29  /**
30  * \brief Calculate tetrahedron volume length quality
31  * @param coords tet coordinates
32  * @return Volume-length quality
33  */
34  static double volumeLengthQuality(const double *coords);
35 
36  /**
37  * @brief Calculate volume of tetrahedron
38  *
39  * @param coords
40  * @return double volume
41  */
42  static double tetVolume(const double *coords);
43 
44  /**
45  * \brief calculate minimal quality of tetrahedra in range
46  * @param tets range
47  * @param min_quality minimal quality
48  * @return error code
49  */
51  const Range &tets, double &min_quality, Tag th = nullptr,
52  boost::function<double(double, double)> f =
53  [](double a, double b) -> double { return std::min(a, b); });
54 
55  static constexpr double shapeFunMBEDGE0At00 = N_MBEDGE0(0);
56  static constexpr double shapeFunMBEDGE1At00 = N_MBEDGE1(0);
57 
58  /**
59  * @brief Array of shape function at zero local point on reference element
60  *
61  */
62  static constexpr std::array<double, 2> shapeFunMBEDGEAt00 = {
64 
65  static constexpr double diffN_MBEDGE0x = diffN_MBEDGE0;
66  static constexpr double diffN_MBEDGE1x = diffN_MBEDGE1;
67 
68  static constexpr std::array<double, 2> diffShapeFunMBEDGE = {diffN_MBEDGE0x,
70 
71  static inline double shapeFunMBEDGE0(const double x);
72 
73  static inline double shapeFunMBEDGE1(const double x);
74 
75  /**
76  * @brief Calculate shape functions on edge
77  *
78  * \note Template parameter is leading dimension of point coordinate arrays,
79  * such that \f$ksi_{n+1} = ksi[n + LDB]\f$
80  *
81  * @tparam 1
82  * @param shape shape functions
83  * @param ksi pointer to first local coordinates
84  * @param nb number of points
85  * @return MoFEMErrorCode
86  */
87  template <int LDB = 1>
88  static MoFEMErrorCode shapeFunMBEDGE(double *shape, const double *ksi,
89  const int nb);
90 
91  static constexpr double diffShapeFunMBTRI0x =
92  diffN_MBTRI0x; ///< derivative of triangle shape function
93  static constexpr double diffShapeFunMBTRI0y =
94  diffN_MBTRI0y; ///< derivative of triangle shape function
95  static constexpr double diffShapeFunMBTRI1x =
96  diffN_MBTRI1x; ///< derivative of triangle shape function
97  static constexpr double diffShapeFunMBTRI1y =
98  diffN_MBTRI1y; ///< derivative of triangle shape function
99  static constexpr double diffShapeFunMBTRI2x =
100  diffN_MBTRI2x; ///< derivative of triangle shape function
101  static constexpr double diffShapeFunMBTRI2y =
102  diffN_MBTRI2y; ///< derivative of triangle shape function
103 
104  static constexpr std::array<double, 6> diffShapeFunMBTRI = {
105 
107 
109 
111 
112  static constexpr double shapeFunMBTRI0At00 = N_MBTRI0(0, 0);
113  static constexpr double shapeFunMBTRI1At00 = N_MBTRI1(0, 0);
114  static constexpr double shapeFunMBTRI2At00 = N_MBTRI2(0, 0);
115 
116  /**
117  * @brief Array of shape function at zero local point on reference element
118  *
119  */
120  static constexpr std::array<double, 3> shapeFunMBTRIAt00 = {
122 
123  static inline double shapeFunMBTRI0(const double x, const double y);
124 
125  static inline double shapeFunMBTRI1(const double x, const double y);
126 
127  static inline double shapeFunMBTRI2(const double x, const double y);
128 
129  /**
130  * @brief Calculate shape functions on triangle
131  *
132  * \note Template parameter is leading dimension of point coordinate arrays,
133  * such that \f$ksi_{n+1} = ksi[n + LDB]\f$
134  *
135  * @tparam 1
136  * @param shape shape functions
137  * @param ksi pointer to first local coordinates
138  * @param eta pointer to second local coordinates
139  * @param nb number of points
140  * @return MoFEMErrorCode
141  */
142  template <int LDB = 1>
143  static MoFEMErrorCode shapeFunMBTRI(double *shape, const double *ksi,
144  const double *eta, const int nb);
145 
146  static constexpr double diffShapeFunMBQUADAtCenter0x =
147  diffN_MBQUAD0x(0.5); ///< derivative of quad shape function
148  static constexpr double diffShapeFunMBQUADAtCenter0y =
149  diffN_MBQUAD0y(0.5); ///< derivative of quad shape function
150  static constexpr double diffShapeFunMBQUADAtCenter1x =
151  diffN_MBQUAD1x(0.5); ///< derivative of quad shape function
152  static constexpr double diffShapeFunMBQUADAtCenter1y =
153  diffN_MBQUAD1y(0.5); ///< derivative of quad shape function
154  static constexpr double diffShapeFunMBQUADAtCenter2x =
155  diffN_MBQUAD2x(0.5); ///< derivative of quad shape function
156  static constexpr double diffShapeFunMBQUADAtCenter2y =
157  diffN_MBQUAD2y(0.5); ///< derivative of quad shape function
158  static constexpr double diffShapeFunMBQUADAtCenter3x =
159  diffN_MBQUAD3x(0.5); ///< derivative of quad shape function
160  static constexpr double diffShapeFunMBQUADAtCenter3y =
161  diffN_MBQUAD3y(0.5); ///< derivative of quad shape function
162 
163  static constexpr double diffShapeFunMBHEXAtCenter0x =
164  diffN_MBHEX0x(0.5, 0.5); ///< derivative of HEX shape function
165  static constexpr double diffShapeFunMBHEXAtCenter0y =
166  diffN_MBHEX0y(0.5, 0.5); ///< derivative of HEX shape function
167  static constexpr double diffShapeFunMBHEXAtCenter0z =
168  diffN_MBHEX0z(0.5, 0.5); ///< derivative of HEX shape function
169  static constexpr double diffShapeFunMBHEXAtCenter1x =
170  diffN_MBHEX1x(0.5, 0.5); ///< derivative of HEX shape function
171  static constexpr double diffShapeFunMBHEXAtCenter1y =
172  diffN_MBHEX1y(0.5, 0.5); ///< derivative of HEX shape function
173  static constexpr double diffShapeFunMBHEXAtCenter1z =
174  diffN_MBHEX1z(0.5, 0.5); ///< derivative of HEX shape function
175  static constexpr double diffShapeFunMBHEXAtCenter2x =
176  diffN_MBHEX2x(0.5, 0.5); ///< derivative of HEX shape function
177  static constexpr double diffShapeFunMBHEXAtCenter2y =
178  diffN_MBHEX2y(0.5, 0.5); ///< derivative of HEX shape function
179  static constexpr double diffShapeFunMBHEXAtCenter2z =
180  diffN_MBHEX2z(0.5, 0.5); ///< derivative of HEX shape function
181  static constexpr double diffShapeFunMBHEXAtCenter3x =
182  diffN_MBHEX3x(0.5, 0.5); ///< derivative of HEX shape function
183  static constexpr double diffShapeFunMBHEXAtCenter3y =
184  diffN_MBHEX3y(0.5, 0.5); ///< derivative of quad shape function
185  static constexpr double diffShapeFunMBHEXAtCenter3z =
186  diffN_MBHEX3z(0.5, 0.5); ///< derivative of quad shape function
187  static constexpr double diffShapeFunMBHEXAtCenter4x =
188  diffN_MBHEX4x(0.5, 0.5); ///< derivative of HEX shape function
189  static constexpr double diffShapeFunMBHEXAtCenter4y =
190  diffN_MBHEX4y(0.5, 0.5); ///< derivative of quad shape function
191  static constexpr double diffShapeFunMBHEXAtCenter4z =
192  diffN_MBHEX4z(0.5, 0.5); ///< derivative of quad shape function
193  static constexpr double diffShapeFunMBHEXAtCenter5x =
194  diffN_MBHEX5x(0.5, 0.5); ///< derivative of HEX shape function
195  static constexpr double diffShapeFunMBHEXAtCenter5y =
196  diffN_MBHEX5y(0.5, 0.5); ///< derivative of quad shape function
197  static constexpr double diffShapeFunMBHEXAtCenter5z =
198  diffN_MBHEX5z(0.5, 0.5); ///< derivative of quad shape function
199  static constexpr double diffShapeFunMBHEXAtCenter6x =
200  diffN_MBHEX6x(0.5, 0.5); ///< derivative of HEX shape function
201  static constexpr double diffShapeFunMBHEXAtCenter6y =
202  diffN_MBHEX6y(0.5, 0.5); ///< derivative of quad shape function
203  static constexpr double diffShapeFunMBHEXAtCenter6z =
204  diffN_MBHEX6z(0.5, 0.5); ///< derivative of quad shape function
205  static constexpr double diffShapeFunMBHEXAtCenter7x =
206  diffN_MBHEX7x(0.5, 0.5); ///< derivative of HEX shape function
207  static constexpr double diffShapeFunMBHEXAtCenter7y =
208  diffN_MBHEX7y(0.5, 0.5); ///< derivative of quad shape function
209  static constexpr double diffShapeFunMBHEXAtCenter7z =
210  diffN_MBHEX7z(0.5, 0.5); ///< derivative of quad shape function
211 
212  static constexpr std::array<double, 8> diffShapeFunMBQUADAtCenter = {
217 
218  static constexpr std::array<double, 24> diffShapeFunMBHEXAtCenter = {
219 
222 
225 
228 
231 
234 
237 
240 
243 
244  };
245 
246  static constexpr double diffShapeFunMBTET0x =
247  diffN_MBTET0x; ///< derivative of tetrahedral shape function
248  static constexpr double diffShapeFunMBTET0y =
249  diffN_MBTET0y; ///< derivative of tetrahedral shape function
250  static constexpr double diffShapeFunMBTET0z =
251  diffN_MBTET0z; ///< derivative of tetrahedral shape function
252  static constexpr double diffShapeFunMBTET1x =
253  diffN_MBTET1x; ///< derivative of tetrahedral shape function
254  static constexpr double diffShapeFunMBTET1y =
255  diffN_MBTET1y; ///< derivative of tetrahedral shape function
256  static constexpr double diffShapeFunMBTET1z =
257  diffN_MBTET1z; ///< derivative of tetrahedral shape function
258  static constexpr double diffShapeFunMBTET2x =
259  diffN_MBTET2x; ///< derivative of tetrahedral shape function
260  static constexpr double diffShapeFunMBTET2y =
261  diffN_MBTET2y; ///< derivative of tetrahedral shape function
262  static constexpr double diffShapeFunMBTET2z =
263  diffN_MBTET2z; ///< derivative of tetrahedral shape function
264  static constexpr double diffShapeFunMBTET3x =
265  diffN_MBTET3x; ///< derivative of tetrahedral shape function
266  static constexpr double diffShapeFunMBTET3y =
267  diffN_MBTET3y; ///< derivative of tetrahedral shape function
268  static constexpr double diffShapeFunMBTET3z =
269  diffN_MBTET3z; ///< derivative of tetrahedral shape function
270 
271  static constexpr std::array<double, 12> diffShapeFunMBTET = {
272 
274 
276 
278 
280 
281  static inline double shapeFunMBTET0(const double x, const double y,
282  const double z);
283 
284  static inline double shapeFunMBTET1(const double x, const double y,
285  const double z);
286 
287  static inline double shapeFunMBTET2(const double x, const double y,
288  const double z);
289 
290  static inline double shapeFunMBTET3(const double x, const double y,
291  const double z);
292 
293  static constexpr double shapeFunMBTET0At000 = N_MBTET0(0, 0, 0);
294  static constexpr double shapeFunMBTET1At000 = N_MBTET1(0, 0, 0);
295  static constexpr double shapeFunMBTET2At000 = N_MBTET2(0, 0, 0);
296  static constexpr double shapeFunMBTET3At000 = N_MBTET3(0, 0, 0);
297 
298  static constexpr double shapeFunMBTET0AtOneThird =
299  N_MBTET0(1. / 3., 1. / 3., 1. / 3.);
300  static constexpr double shapeFunMBTET1AtOneThird =
301  N_MBTET1(1. / 3., 1. / 3., 1. / 3.);
302  static constexpr double shapeFunMBTET2AtOneThird =
303  N_MBTET2(1. / 3., 1. / 3., 1. / 3.);
304  static constexpr double shapeFunMBTET3AtOneThird =
305  N_MBTET3(1. / 3., 1. / 3., 1. / 3.);
306 
307  /**
308  * @brief Calculate shape functions on tetrahedron
309  *
310  * \note Template parameter is leading dimension of point coordinate arrays,
311  * such that \f$ksi_{n+1} = ksi[n + LDB]\f$
312  *
313  * @tparam 1
314  * @param shape shape functions
315  * @param ksi pointer to first local coordinates
316  * @param eta pointer to second local coordinates
317  * @param zeta pointer to first third coordinates
318  * @param nb number of points
319  * @return MoFEMErrorCode
320  */
321  template <int LDB = 1>
322  static MoFEMErrorCode shapeFunMBTET(double *shape, const double *ksi,
323  const double *eta, const double *zeta,
324  const double nb);
325 
326  /**
327  * @brief Array of shape function at zero local point on reference element
328  *
329  */
330  static constexpr std::array<double, 4> shapeFunMBTETAt000 = {
333 
334  /**
335  * @brief Array of shape function at center on reference element
336  *
337  */
338  static constexpr std::array<double, 4> shapeFunMBTETAtOneThird = {
341 
342  /**
343  * @brief Get the local coordinates on reference four node tet object
344  *
345  * \code
346  * MatrixDouble elem_coords(4, 3);
347  * // Set nodal coordinates
348  * MatrixDouble global_coords(5, 3);
349  * // Set global coordinates
350  * MatrixDouble local_coords(global_coords.size1(), 3);
351  * CHKERR Tools::getLocalCoordinatesOnReferenceFourNodeTet(
352  * &elem_coords(0, 0), &global_coords(0, 0), global_coords.size1(),
353  * &local_coords(0, 0))
354  * \endcode
355  *
356  * @param elem_coords Global element node coordinates
357  * @param glob_coords Global coordinates
358  * @param nb_nodes Number of points
359  * @param local_coords Result
360  * @return MoFEMErrorCode
361  */
363  const double *elem_coords, const double *glob_coords, const int nb_nodes,
364  double *local_coords);
365 
366  /**
367  * @brief Get the local coordinates on reference three node tri object
368  *
369  * \code
370  * MatrixDouble elem_coords(4, 3);
371  * // Set nodal coordinates
372  * MatrixDouble global_coords(5, 3);
373  * // Set global coordinates
374  * MatrixDouble local_coords(global_coords.size1(), 3);
375  * CHKERR Tools::getLocalCoordinatesOnReferenceFourNodeTet(
376  * &elem_coords(0, 0), &global_coords(0, 0), global_coords.size1(),
377  * &local_coords(0, 0))
378  * \endcode
379  *
380  * @param elem_coords Global element node coordinates
381  * @param glob_coords Global coordinates
382  * @param nb_nodes Number of points
383  * @param local_coords Result
384  * @param d_elem_coords Derivative of local coordinates
385  * @param d_global_coords
386  * @return MoFEMErrorCode
387  */
389  const double *elem_coords, const double *glob_coords, const int nb_nodes,
390  double *local_coords);
391 
392  /** @copydoc getLocalCoordinatesOnReferenceThreeNodeTri */
394  const double *elem_coords, const std::complex<double> *glob_coords,
395  const int nb_nodes, std::complex<double> *local_coords);
396 
397 /** @deprecated use getLocalCoordinatesOnReferenceThreeNodeTri */
399  const double *elem_coords, const double *glob_coords, const int nb_nodes,
400  double *local_coords) {
401  return getLocalCoordinatesOnReferenceThreeNodeTri(elem_coords, glob_coords,
402  nb_nodes, local_coords);
403  }
404 
405  /**
406  * @brief Get the local coordinates on reference four node tet object
407  *
408  * \code
409  * MatrixDouble elem_coords(4, 3);
410  * // Set nodal coordinates
411  * MatrixDouble global_coords(5, 3);
412  * // Set global coordinates
413  * MatrixDouble local_coords(global_coords.size1(), 3);
414  * CHKERR Tools::getLocalCoordinatesOnReferenceEdgeNodeEdge(
415  * &elem_coords(0, 0), &global_coords(0, 0), global_coords.size1(),
416  * &local_coords(0, 0))
417  * \endcode
418  *
419  * @param elem_coords Global element node coordinates
420  * @param glob_coords Global coordinates
421  * @param nb_nodes Number of points
422  * @param local_coords Result
423  * @return MoFEMErrorCode
424  */
426  const double *elem_coords, const double *glob_coords, const int nb_nodes,
427  double *local_coords);
428 
429  /**
430  * @brief Get the Tets With Quality
431  *
432  * @param out_tets
433  * @param tets
434  * @param th
435  * @param f
436  * @return MoFEMErrorCode
437  */
439  Range &out_tets, const Range &tets, Tag th = nullptr,
440  boost::function<bool(double)> f = [](double q) -> bool {
441  if (q <= 0)
442  return true;
443  else
444  return false;
445  });
446  /**
447  * @brief Write file with tetrahedral of given quality
448  *
449  * @param file_name
450  * @param file_type
451  * @param options
452  * @param tets
453  * @param th
454  * @param f
455  * @return MoFEMErrorCode
456  */
458  const char *file_name, const char *file_type, const char *options,
459  const Range &tets, Tag th = nullptr,
460  boost::function<bool(double)> f = [](double q) -> bool {
461  if (q <= 0)
462  return true;
463  else
464  return false;
465  });
466 
467  /**
468  * @brief Check of point is in tetrahedral
469  *
470  * @param tet_coords
471  * @param global_coord
472  * @param tol
473  * @param result
474  * @return MoFEMErrorCode
475  */
476  static MoFEMErrorCode checkIfPointIsInTet(const double tet_coords[],
477  const double global_coord[],
478  const double tol, bool &result);
479 
480  /**
481  * @brief Get the Tri Normal objectGet triangle normal
482  *
483  * @param coords
484  * @param normal
485  * @param d_normal derbiative, if pointer is null, derbiative is not
486  * calculated
487  * @return MoFEMErrorCode
488  */
489  static MoFEMErrorCode getTriNormal(const double *coords, double *normal,
490  double *d_normal = nullptr);
491 
492  /**
493  * @brief Get triangle normal
494  *
495  * @param tri
496  * @param normal
497  * @return MoFEMErrorCode
498  */
499  MoFEMErrorCode getTriNormal(const EntityHandle tri, double *normal) const;
500 
501  /**
502  * @brief Get triangle area
503  *
504  * @param tri
505  * @return double
506  */
507  double getTriArea(const EntityHandle tri) const;
508 
509  /**
510  * @brief Get edge length
511  *
512  * @param edge_coords
513  * @return double
514  */
515  static double getEdgeLength(const double *edge_coords);
516 
517  /**
518  * @brief Get edge length
519  *
520  * @param edge
521  * @return double
522  */
523  double getEdgeLength(const EntityHandle edge);
524 
531  };
532 
533  /**
534  * @brief Find closet point on the segment from the point
535  *
536  * @param w_ptr segment first vertex coordinate
537  * @param v_ptr segment second vertex coordinate
538  * @param p_ptr coordinate of point
539  * @param t_ptr distance on the segment
540  *
541  * \note If t is outside bounds [ 0,-1 ] point is on the line point
542  * beyond segment.
543  *
544  * \code
545  *
546  * double w[] = {-1, 0, 0};
547  * double v[] = {1, 0, 0};
548  * double p[] = {0, 1, 0};
549  * double t;
550  * CHKERR Toolas::minDistancePointFromOnSegment(w, v, p, &t);
551  * double point_on_segment[3];
552  * for (int i = 0; i != 3; ++i)
553  * point_on_segment[i] = w[i] + t * (v[i] - w[i]);
554  *
555  * \endcode
556  *
557  * @return SEGMENT_MIN_DISTANCE
558  */
559  static SEGMENT_MIN_DISTANCE
560  minDistancePointFromOnSegment(const double *w_ptr, const double *v_ptr,
561  const double *p_ptr,
562  double *const t_ptr = nullptr);
563  /**
564  * @brief Find points on two segments in closest distance
565  *
566  * @param w_ptr
567  * @param v_ptr
568  * @param k_ptr
569  * @param l_ptr
570  * @param tvw_ptr
571  * @param tlk_ptr
572  * @return SEGMENT_MIN_DISTANCE
573  *
574  * \note If tvwk or tlk are outside bound [0,-1], it means that points
575  * are on the lines beyond segments, respectively for segment vw and
576  * lk.
577  *
578  */
579  static SEGMENT_MIN_DISTANCE
580  minDistanceFromSegments(const double *w_ptr, const double *v_ptr,
581  const double *k_ptr, const double *l_ptr,
582  double *const tvw_ptr = nullptr,
583  double *const tlk_ptr = nullptr);
584 
585  /**
586  * @brief Find minimal distance to edges
587  *
588  * \note Finding only edges with have smaller distance than distance
589  * set on the input by min_dist_ptr
590  *
591  * @param v_ptr point coordinates
592  * @param nb nb points
593  * @param edges range of edges
594  * @param min_dist_ptr on return minimal distance, on input starting
595  * distance
596  * @param o_ptr coordinates of the point on edge
597  * @param o_segments closest segments
598  * @return MoFEMErrorCode
599  */
601  findMinDistanceFromTheEdges(const double *v_ptr, const int nb, Range edges,
602  double *min_dist_ptr, double *o_ptr = nullptr,
603  EntityHandle *o_segments = nullptr) const;
604 
605  /** \name Debugging */
606 
607  /**@{*/
608 
609  /** \brief Print all DOFs for which element of vector is not a number
610  *
611  */
613  const RowColData row_or_col, Vec v);
614 
615  /**@}*/
616 
617  static MoFEMErrorCode
619  const int edge1);
620 
621  static MoFEMErrorCode
623  const int edge1, const int edge2);
624 
625  /** \name Mesh refinement */
626 
627  /**@{*/
628 
629  static constexpr std::array<int, 12> uniformTriangleRefineTriangles = {
630 
631  0, 3, 5, // 0
632  3, 1, 4, // 1
633  5, 4, 2, // 2
634  5, 3, 4 // 3
635 
636  };
637 
638  using RefineTrianglesReturn =
639  std::tuple<std::vector<double>, std::vector<int>, std::vector<int>>;
640 
641  /**
642  * @brief create uniform triangle mesh of refined elements
643  *
644  * @param nb_levels
645  * @return RefineTrianglesReturn
646  */
648 
649  /**
650  * @brief generate integration points for refined triangle mesh for last level
651  *
652  * @param pts
653  * @param refined
654  * @return MatrixDouble
655  */
656  static MatrixDouble
658 
659  /**
660  * @brief generate integration points for refined triangle mesh for last level
661  *
662  * @param rule Gauss integration rule
663  * @param refined
664  * @return MatrixDouble
665  */
666  static MatrixDouble
668 
669  /**@}*/
670 
671 };
672 
673 double Tools::shapeFunMBEDGE0(const double x) {
674  return N_MBEDGE0(x);
675 }
676 
677 double Tools::shapeFunMBEDGE1(const double x) {
678  return N_MBEDGE1(x);
679 }
680 
681 template <int LDB>
682 MoFEMErrorCode Tools::shapeFunMBEDGE(double *shape, const double *ksi,
683  const int nb) {
685  for (int n = 0; n != nb; ++n) {
686  shape[0] = shapeFunMBEDGE0(*ksi);
687  shape[1] = shapeFunMBEDGE1(*ksi);
688  shape += 2;
689  ksi += LDB;
690  }
692 }
693 
694 double Tools::shapeFunMBTRI0(const double x, const double y) {
695  return N_MBTRI0(x, y);
696 }
697 
698 double Tools::shapeFunMBTRI1(const double x, const double y) {
699  return N_MBTRI1(x, y);
700 }
701 
702 double Tools::shapeFunMBTRI2(const double x, const double y) {
703  return N_MBTRI2(x, y);
704 }
705 
706 template <int LDB>
707 MoFEMErrorCode Tools::shapeFunMBTRI(double *shape, const double *ksi,
708  const double *eta, const int nb) {
710  for (int n = 0; n != nb; ++n) {
711  shape[0] = shapeFunMBTRI0(*ksi, *eta);
712  shape[1] = shapeFunMBTRI1(*ksi, *eta);
713  shape[2] = shapeFunMBTRI2(*ksi, *eta);
714  shape += 3;
715  ksi += LDB;
716  eta += LDB;
717  }
719 }
720 
721 double Tools::shapeFunMBTET0(const double x, const double y, const double z) {
722  return N_MBTET0(x, y, z);
723 }
724 
725 double Tools::shapeFunMBTET1(const double x, const double y, const double z) {
726  return N_MBTET1(x, y, z);
727 }
728 
729 double Tools::shapeFunMBTET2(const double x, const double y, const double z) {
730  return N_MBTET2(x, y, z);
731 }
732 
733 double Tools::shapeFunMBTET3(const double x, const double y, const double z) {
734  return N_MBTET3(x, y, z);
735 };
736 
737 template <int LDB>
738 MoFEMErrorCode Tools::shapeFunMBTET(double *shape, const double *ksi,
739  const double *eta, const double *zeta,
740  const double nb) {
742  for (int n = 0; n != nb; ++n) {
743  shape[0] = shapeFunMBTET0(*ksi, *eta, *zeta);
744  shape[1] = shapeFunMBTET1(*ksi, *eta, *zeta);
745  shape[2] = shapeFunMBTET2(*ksi, *eta, *zeta);
746  shape[3] = shapeFunMBTET3(*ksi, *eta, *zeta);
747  shape += 4;
748  ksi += LDB;
749  eta += LDB;
750  zeta += LDB;
751  }
753 }
754 
755 
756 } // namespace MoFEM
757 
758 #endif // __TOOLS_HPP__a
759 
760 /**
761  * \defgroup mofem_tools Tools interface
762  * \brief Interface for tools
763  *
764  * \ingroup mofem
765  */
UBlasMatrix< double >
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
diffN_MBEDGE0
#define diffN_MBEDGE0
derivative of edge shape function
Definition: fem_tools.h:107
MoFEM::Tools::diffShapeFunMBTET0z
static constexpr double diffShapeFunMBTET0z
derivative of tetrahedral shape function
Definition: Tools.hpp:250
MoFEM::Tools::shapeFunMBTET1At000
static constexpr double shapeFunMBTET1At000
Definition: Tools.hpp:294
MoFEM::Tools::diffShapeFunMBHEXAtCenter0y
static constexpr double diffShapeFunMBHEXAtCenter0y
derivative of HEX shape function
Definition: Tools.hpp:165
MoFEM::Tools::query_interface
MoFEMErrorCode query_interface(boost::typeindex::type_index type_index, UnknownInterface **iface) const
Definition: Tools.cpp:9
MoFEM::Tools::cOre
MoFEM::Core & cOre
Definition: Tools.hpp:24
diffN_MBTET3z
#define diffN_MBTET3z
derivative of tetrahedral shape function
Definition: fem_tools.h:43
MoFEM::Tools::diffShapeFunMBEDGE
static constexpr std::array< double, 2 > diffShapeFunMBEDGE
Definition: Tools.hpp:68
diffN_MBTRI0y
#define diffN_MBTRI0y
derivative of triangle shape function
Definition: fem_tools.h:50
diffN_MBHEX0y
#define diffN_MBHEX0y(x, z)
Definition: fem_tools.h:87
MoFEM::Tools::diffShapeFunMBQUADAtCenter1y
static constexpr double diffShapeFunMBQUADAtCenter1y
derivative of quad shape function
Definition: Tools.hpp:152
MoFEM::Tools::minDistancePointFromOnSegment
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.
Definition: Tools.cpp:446
MoFEM::CoreTmp< 0 >
Core (interface) class.
Definition: Core.hpp:82
DEPRECATED
#define DEPRECATED
Definition: definitions.h:17
MoFEM::Tools::diffShapeFunMBTRI0x
static constexpr double diffShapeFunMBTRI0x
derivative of triangle shape function
Definition: Tools.hpp:91
MoFEM::Tools::diffShapeFunMBHEXAtCenter0z
static constexpr double diffShapeFunMBHEXAtCenter0z
derivative of HEX shape function
Definition: Tools.hpp:167
MoFEM::Tools::SEGMENT_MIN_DISTANCE
SEGMENT_MIN_DISTANCE
Definition: Tools.hpp:525
diffN_MBQUAD1x
#define diffN_MBQUAD1x(y)
Definition: fem_tools.h:63
MoFEM::Tools::diffN_MBEDGE0x
static constexpr double diffN_MBEDGE0x
Definition: Tools.hpp:65
diffN_MBTRI2y
#define diffN_MBTRI2y
derivative of triangle shape function
Definition: fem_tools.h:54
MoFEM::Tools::shapeFunMBTET3AtOneThird
static constexpr double shapeFunMBTET3AtOneThird
Definition: Tools.hpp:304
MoFEM::Tools::diffShapeFunMBHEXAtCenter7y
static constexpr double diffShapeFunMBHEXAtCenter7y
derivative of quad shape function
Definition: Tools.hpp:207
EntityHandle
MoFEM::Tools::NO_SOLUTION
@ NO_SOLUTION
Definition: Tools.hpp:530
MoFEM::Tools::shapeFunMBEDGE0At00
static constexpr double shapeFunMBEDGE0At00
Definition: Tools.hpp:55
diffN_MBHEX1z
#define diffN_MBHEX1z(x, y)
Definition: fem_tools.h:96
MoFEM::Tools::diffN_MBEDGE1x
static constexpr double diffN_MBEDGE1x
Definition: Tools.hpp:66
MoFEM::Tools::diffShapeFunMBQUADAtCenter3y
static constexpr double diffShapeFunMBQUADAtCenter3y
derivative of quad shape function
Definition: Tools.hpp:160
MoFEM::Tools::shapeFunMBTET2
static double shapeFunMBTET2(const double x, const double y, const double z)
Definition: Tools.hpp:729
MoFEM::Tools
Auxiliary tools.
Definition: Tools.hpp:19
MoFEM::Tools::diffShapeFunMBHEXAtCenter3x
static constexpr double diffShapeFunMBHEXAtCenter3x
derivative of HEX shape function
Definition: Tools.hpp:181
MoFEM::Tools::diffShapeFunMBQUADAtCenter1x
static constexpr double diffShapeFunMBQUADAtCenter1x
derivative of quad shape function
Definition: Tools.hpp:150
MoFEM::Tools::diffShapeFunMBTET2x
static constexpr double diffShapeFunMBTET2x
derivative of tetrahedral shape function
Definition: Tools.hpp:258
diffN_MBHEX3x
#define diffN_MBHEX3x(y, z)
Definition: fem_tools.h:82
MoFEM::Tools::diffShapeFunMBTET2z
static constexpr double diffShapeFunMBTET2z
derivative of tetrahedral shape function
Definition: Tools.hpp:262
diffN_MBQUAD1y
#define diffN_MBQUAD1y(x)
Definition: fem_tools.h:64
MoFEM::Tools::checkIfPointIsInTet
static MoFEMErrorCode checkIfPointIsInTet(const double tet_coords[], const double global_coord[], const double tol, bool &result)
Check of point is in tetrahedral.
Definition: Tools.cpp:287
MoFEM::Tools::diffShapeFunMBTET3x
static constexpr double diffShapeFunMBTET3x
derivative of tetrahedral shape function
Definition: Tools.hpp:264
MoFEM::Tools::getTetsWithQuality
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.
Definition: Tools.cpp:244
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
MoFEM::Tools::diffShapeFunMBHEXAtCenter5z
static constexpr double diffShapeFunMBHEXAtCenter5z
derivative of quad shape function
Definition: Tools.hpp:197
MoFEM::Tools::shapeFunMBTET3
static double shapeFunMBTET3(const double x, const double y, const double z)
Definition: Tools.hpp:733
diffN_MBTET2z
#define diffN_MBTET2z
derivative of tetrahedral shape function
Definition: fem_tools.h:40
nb_levels
constexpr int nb_levels
Definition: level_set.cpp:58
MoFEM::Tools::shapeFunMBTRI1At00
static constexpr double shapeFunMBTRI1At00
Definition: Tools.hpp:113
MoFEM::Tools::diffShapeFunMBTET0y
static constexpr double diffShapeFunMBTET0y
derivative of tetrahedral shape function
Definition: Tools.hpp:248
MoFEM::Tools::diffShapeFunMBHEXAtCenter
static constexpr std::array< double, 24 > diffShapeFunMBHEXAtCenter
Definition: Tools.hpp:218
MoFEM::th
Tag th
Definition: Projection10NodeCoordsOnField.cpp:122
diffN_MBQUAD2x
#define diffN_MBQUAD2x(y)
Definition: fem_tools.h:65
diffN_MBTET2y
#define diffN_MBTET2y
derivative of tetrahedral shape function
Definition: fem_tools.h:39
MoFEM::Tools::diffShapeFunMBTET
static constexpr std::array< double, 12 > diffShapeFunMBTET
Definition: Tools.hpp:271
MoFEM::Tools::shapeFunMBEDGEAt00
static constexpr std::array< double, 2 > shapeFunMBEDGEAt00
Array of shape function at zero local point on reference element.
Definition: Tools.hpp:62
N_MBEDGE0
#define N_MBEDGE0(x)
edge shape function
Definition: fem_tools.h:105
MoFEM::Tools::diffShapeFunMBQUADAtCenter0y
static constexpr double diffShapeFunMBQUADAtCenter0y
derivative of quad shape function
Definition: Tools.hpp:148
MoFEM::Tools::diffShapeFunMBHEXAtCenter0x
static constexpr double diffShapeFunMBHEXAtCenter0x
derivative of HEX shape function
Definition: Tools.hpp:163
MoFEM::Tools::diffShapeFunMBTET2y
static constexpr double diffShapeFunMBTET2y
derivative of tetrahedral shape function
Definition: Tools.hpp:260
MoFEM::Tools::diffShapeFunMBTET1y
static constexpr double diffShapeFunMBTET1y
derivative of tetrahedral shape function
Definition: Tools.hpp:254
MoFEM::Tools::getTriArea
double getTriArea(const EntityHandle tri) const
Get triangle area.
Definition: Tools.cpp:408
MoFEM::Tools::outerProductOfEdgeIntegrationPtsForHex
static MoFEMErrorCode outerProductOfEdgeIntegrationPtsForHex(MatrixDouble &pts, const int edge0, const int edge1, const int edge2)
Definition: Tools.cpp:659
MoFEM::Tools::shapeFunMBTRI0
static double shapeFunMBTRI0(const double x, const double y)
Definition: Tools.hpp:694
MoFEM::Tools::shapeFunMBEDGE1
static double shapeFunMBEDGE1(const double x)
Definition: Tools.hpp:677
MoFEM::Tools::RefineTrianglesReturn
std::tuple< std::vector< double >, std::vector< int >, std::vector< int > > RefineTrianglesReturn
Definition: Tools.hpp:639
diffN_MBHEX5x
#define diffN_MBHEX5x(y, z)
Definition: fem_tools.h:84
MoFEM::Tools::diffShapeFunMBHEXAtCenter2z
static constexpr double diffShapeFunMBHEXAtCenter2z
derivative of HEX shape function
Definition: Tools.hpp:179
MoFEM::Tools::shapeFunMBEDGE0
static double shapeFunMBEDGE0(const double x)
Definition: Tools.hpp:673
MoFEM::Tools::diffShapeFunMBHEXAtCenter4y
static constexpr double diffShapeFunMBHEXAtCenter4y
derivative of quad shape function
Definition: Tools.hpp:189
diffN_MBHEX2y
#define diffN_MBHEX2y(x, z)
Definition: fem_tools.h:89
N_MBTET0
#define N_MBTET0(x, y, z)
tetrahedral shape function
Definition: fem_tools.h:28
MoFEM::Tools::diffShapeFunMBHEXAtCenter7z
static constexpr double diffShapeFunMBHEXAtCenter7z
derivative of quad shape function
Definition: Tools.hpp:209
zeta
double zeta
Viscous hardening.
Definition: plastic.cpp:126
MoFEM::Tools::getLocalCoordinatesOnReferenceEdgeNodeEdge
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.
Definition: Tools.cpp:203
MoFEM::Tools::minDistanceFromSegments
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.
Definition: Tools.cpp:469
MoFEM::Tools::refineTriangleIntegrationPts
static MatrixDouble refineTriangleIntegrationPts(MatrixDouble pts, RefineTrianglesReturn refined)
generate integration points for refined triangle mesh for last level
MoFEM::Tools::diffShapeFunMBHEXAtCenter2y
static constexpr double diffShapeFunMBHEXAtCenter2y
derivative of HEX shape function
Definition: Tools.hpp:177
diffN_MBHEX5y
#define diffN_MBHEX5y(x, z)
Definition: fem_tools.h:92
MoFEM::Tools::diffShapeFunMBQUADAtCenter2x
static constexpr double diffShapeFunMBQUADAtCenter2x
derivative of quad shape function
Definition: Tools.hpp:154
MoFEM::Tools::getLocalCoordinatesOnReferenceThreeNodeTri
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.
Definition: Tools.cpp:188
MoFEM::Tools::diffShapeFunMBQUADAtCenter3x
static constexpr double diffShapeFunMBQUADAtCenter3x
derivative of quad shape function
Definition: Tools.hpp:158
N_MBTET3
#define N_MBTET3(x, y, z)
tetrahedral shape function
Definition: fem_tools.h:31
MoFEM::Tools::shapeFunMBTET1
static double shapeFunMBTET1(const double x, const double y, const double z)
Definition: Tools.hpp:725
N_MBTET2
#define N_MBTET2(x, y, z)
tetrahedral shape function
Definition: fem_tools.h:30
diffN_MBTET2x
#define diffN_MBTET2x
derivative of tetrahedral shape function
Definition: fem_tools.h:38
MoFEM::Tools::shapeFunMBTRI
static MoFEMErrorCode shapeFunMBTRI(double *shape, const double *ksi, const double *eta, const int nb)
Calculate shape functions on triangle.
Definition: Tools.hpp:707
MoFEM::Tools::uniformTriangleRefineTriangles
static constexpr std::array< int, 12 > uniformTriangleRefineTriangles
Definition: Tools.hpp:629
diffN_MBHEX1x
#define diffN_MBHEX1x(y, z)
Definition: fem_tools.h:80
N_MBEDGE1
#define N_MBEDGE1(x)
edge shape function
Definition: fem_tools.h:106
diffN_MBQUAD3x
#define diffN_MBQUAD3x(y)
Definition: fem_tools.h:67
MoFEM
implementation of Data Operators for Forces and Sources
Definition: Common.hpp:10
a
constexpr double a
Definition: approx_sphere.cpp:30
diffN_MBHEX4z
#define diffN_MBHEX4z(x, y)
Definition: fem_tools.h:99
MoFEM::Tools::diffShapeFunMBQUADAtCenter2y
static constexpr double diffShapeFunMBQUADAtCenter2y
derivative of quad shape function
Definition: Tools.hpp:156
diffN_MBHEX2x
#define diffN_MBHEX2x(y, z)
Definition: fem_tools.h:81
RowColData
RowColData
RowColData.
Definition: definitions.h:136
diffN_MBTRI0x
#define diffN_MBTRI0x
derivative of triangle shape function
Definition: fem_tools.h:49
MoFEM::Tools::SEGMENT_ONE_IS_POINT
@ SEGMENT_ONE_IS_POINT
Definition: Tools.hpp:527
MoFEM::Tools::diffShapeFunMBHEXAtCenter5y
static constexpr double diffShapeFunMBHEXAtCenter5y
derivative of quad shape function
Definition: Tools.hpp:195
MoFEM::Tools::findMinDistanceFromTheEdges
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.
Definition: Tools.cpp:535
MoFEM::Tools::diffShapeFunMBTRI0y
static constexpr double diffShapeFunMBTRI0y
derivative of triangle shape function
Definition: Tools.hpp:93
diffN_MBHEX3z
#define diffN_MBHEX3z(x, y)
Definition: fem_tools.h:98
MoFEM::Tools::diffShapeFunMBQUADAtCenter0x
static constexpr double diffShapeFunMBQUADAtCenter0x
derivative of quad shape function
Definition: Tools.hpp:146
eta
double eta
Definition: free_surface.cpp:170
MoFEM::Tools::diffShapeFunMBTET0x
static constexpr double diffShapeFunMBTET0x
derivative of tetrahedral shape function
Definition: Tools.hpp:246
diffN_MBHEX2z
#define diffN_MBHEX2z(x, y)
Definition: fem_tools.h:97
MoFEM::Tools::shapeFunMBTETAtOneThird
static constexpr std::array< double, 4 > shapeFunMBTETAtOneThird
Array of shape function at center on reference element.
Definition: Tools.hpp:338
MoFEM::Tools::diffShapeFunMBHEXAtCenter4z
static constexpr double diffShapeFunMBHEXAtCenter4z
derivative of quad shape function
Definition: Tools.hpp:191
diffN_MBHEX7x
#define diffN_MBHEX7x(y, z)
Definition: fem_tools.h:86
MoFEM::Tools::diffShapeFunMBTET3y
static constexpr double diffShapeFunMBTET3y
derivative of tetrahedral shape function
Definition: Tools.hpp:266
diffN_MBTRI2x
#define diffN_MBTRI2x
derivative of triangle shape function
Definition: fem_tools.h:53
MoFEM::Tools::diffShapeFunMBTRI1x
static constexpr double diffShapeFunMBTRI1x
derivative of triangle shape function
Definition: Tools.hpp:95
diffN_MBQUAD2y
#define diffN_MBQUAD2y(x)
Definition: fem_tools.h:66
diffN_MBTET3y
#define diffN_MBTET3y
derivative of tetrahedral shape function
Definition: fem_tools.h:42
diffN_MBHEX7y
#define diffN_MBHEX7y(x, z)
Definition: fem_tools.h:94
MoFEM::Tools::diffShapeFunMBQUADAtCenter
static constexpr std::array< double, 8 > diffShapeFunMBQUADAtCenter
Definition: Tools.hpp:212
MoFEM::Tools::diffShapeFunMBTRI
static constexpr std::array< double, 6 > diffShapeFunMBTRI
Definition: Tools.hpp:104
MoFEM::Tools::diffShapeFunMBTRI2y
static constexpr double diffShapeFunMBTRI2y
derivative of triangle shape function
Definition: Tools.hpp:101
diffN_MBHEX5z
#define diffN_MBHEX5z(x, y)
Definition: fem_tools.h:100
MoFEM::Tools::shapeFunMBTRI0At00
static constexpr double shapeFunMBTRI0At00
Definition: Tools.hpp:112
MoFEM::Tools::diffShapeFunMBHEXAtCenter1x
static constexpr double diffShapeFunMBHEXAtCenter1x
derivative of HEX shape function
Definition: Tools.hpp:169
MoFEM::Tools::getLocalCoordinatesOnReferenceFourNodeTet
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.
Definition: Tools.cpp:88
diffN_MBHEX3y
#define diffN_MBHEX3y(x, z)
Definition: fem_tools.h:90
MoFEM::Tools::shapeFunMBTET2At000
static constexpr double shapeFunMBTET2At000
Definition: Tools.hpp:295
MoFEM::Tools::shapeFunMBEDGE
static MoFEMErrorCode shapeFunMBEDGE(double *shape, const double *ksi, const int nb)
Calculate shape functions on edge.
Definition: Tools.hpp:682
MoFEM::Tools::SEGMENT_TWO_IS_POINT
@ SEGMENT_TWO_IS_POINT
Definition: Tools.hpp:528
diffN_MBQUAD0x
#define diffN_MBQUAD0x(y)
Definition: fem_tools.h:61
diffN_MBHEX1y
#define diffN_MBHEX1y(x, z)
Definition: fem_tools.h:88
diffN_MBTET1z
#define diffN_MBTET1z
derivative of tetrahedral shape function
Definition: fem_tools.h:37
diffN_MBQUAD0y
#define diffN_MBQUAD0y(x)
Definition: fem_tools.h:62
diffN_MBTET0x
#define diffN_MBTET0x
derivative of tetrahedral shape function
Definition: fem_tools.h:32
diffN_MBQUAD3y
#define diffN_MBQUAD3y(x)
Definition: fem_tools.h:68
diffN_MBEDGE1
#define diffN_MBEDGE1
derivative of edge shape function
Definition: fem_tools.h:108
MoFEM::Tools::shapeFunMBTET2AtOneThird
static constexpr double shapeFunMBTET2AtOneThird
Definition: Tools.hpp:302
MoFEM::Tools::shapeFunMBTET0AtOneThird
static constexpr double shapeFunMBTET0AtOneThird
Definition: Tools.hpp:298
convert.n
n
Definition: convert.py:82
MoFEM::Tools::diffShapeFunMBHEXAtCenter7x
static constexpr double diffShapeFunMBHEXAtCenter7x
derivative of HEX shape function
Definition: Tools.hpp:205
MoFEM::Tools::diffShapeFunMBHEXAtCenter5x
static constexpr double diffShapeFunMBHEXAtCenter5x
derivative of HEX shape function
Definition: Tools.hpp:193
N_MBTRI0
#define N_MBTRI0(x, y)
triangle shape function
Definition: fem_tools.h:46
diffN_MBHEX4y
#define diffN_MBHEX4y(x, z)
Definition: fem_tools.h:91
MoFEM::Tools::diffShapeFunMBTET1x
static constexpr double diffShapeFunMBTET1x
derivative of tetrahedral shape function
Definition: Tools.hpp:252
MoFEM::UnknownInterface
base class for all interface classes
Definition: UnknownInterface.hpp:34
diffN_MBTET0y
#define diffN_MBTET0y
derivative of tetrahedral shape function
Definition: fem_tools.h:33
diffN_MBHEX6x
#define diffN_MBHEX6x(y, z)
Definition: fem_tools.h:85
MoFEM::Tools::SEGMENT_TWO_AND_TWO_ARE_POINT
@ SEGMENT_TWO_AND_TWO_ARE_POINT
Definition: Tools.hpp:529
edge_coords
static const double edge_coords[6][6]
Definition: forces_and_sources_h1_continuity_check.cpp:18
MoFEM::Tools::diffShapeFunMBHEXAtCenter6y
static constexpr double diffShapeFunMBHEXAtCenter6y
derivative of quad shape function
Definition: Tools.hpp:201
MoFEM::Tools::shapeFunMBTRI2
static double shapeFunMBTRI2(const double x, const double y)
Definition: Tools.hpp:702
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
MoFEM::Tools::shapeFunMBTET
static MoFEMErrorCode shapeFunMBTET(double *shape, const double *ksi, const double *eta, const double *zeta, const double nb)
Calculate shape functions on tetrahedron.
Definition: Tools.hpp:738
MoFEM::Tools::minTetsQuality
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
Definition: Tools.cpp:49
diffN_MBHEX4x
#define diffN_MBHEX4x(y, z)
Definition: fem_tools.h:83
MoFEM::Tools::tetVolume
static double tetVolume(const double *coords)
Calculate volume of tetrahedron.
Definition: Tools.cpp:30
Range
diffN_MBTET1x
#define diffN_MBTET1x
derivative of tetrahedral shape function
Definition: fem_tools.h:35
MoFEM::Tools::shapeFunMBTRI1
static double shapeFunMBTRI1(const double x, const double y)
Definition: Tools.hpp:698
MoFEM::Tools::diffShapeFunMBTRI2x
static constexpr double diffShapeFunMBTRI2x
derivative of triangle shape function
Definition: Tools.hpp:99
MoFEM::Tools::shapeFunMBTET0
static double shapeFunMBTET0(const double x, const double y, const double z)
Definition: Tools.hpp:721
MoFEM::Tools::Tools
Tools(const MoFEM::Core &core)
Definition: Tools.hpp:25
diffN_MBTRI1y
#define diffN_MBTRI1y
derivative of triangle shape function
Definition: fem_tools.h:52
MoFEM::Tools::diffShapeFunMBHEXAtCenter1z
static constexpr double diffShapeFunMBHEXAtCenter1z
derivative of HEX shape function
Definition: Tools.hpp:173
N_MBTET1
#define N_MBTET1(x, y, z)
tetrahedral shape function
Definition: fem_tools.h:29
diffN_MBHEX6y
#define diffN_MBHEX6y(x, z)
Definition: fem_tools.h:93
diffN_MBTET1y
#define diffN_MBTET1y
derivative of tetrahedral shape function
Definition: fem_tools.h:36
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
MoFEM::Tools::refineTriangle
static RefineTrianglesReturn refineTriangle(int nb_levels)
create uniform triangle mesh of refined elements
Definition: Tools.cpp:721
N_MBTRI1
#define N_MBTRI1(x, y)
triangle shape function
Definition: fem_tools.h:47
MoFEM::Tools::checkVectorForNotANumber
MoFEMErrorCode checkVectorForNotANumber(const Problem *prb_ptr, const RowColData row_or_col, Vec v)
Print all DOFs for which element of vector is not a number.
Definition: Tools.cpp:307
MoFEM::Tools::writeTetsWithQuality
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.
Definition: Tools.cpp:269
MoFEM::Tools::diffShapeFunMBTET1z
static constexpr double diffShapeFunMBTET1z
derivative of tetrahedral shape function
Definition: Tools.hpp:256
EigenMatrix::Vec
const FTensor::Tensor2< T, Dim, Dim > Vec
Definition: MatrixFunction.hpp:66
diffN_MBHEX0z
#define diffN_MBHEX0z(x, y)
Definition: fem_tools.h:95
MoFEM::Tools::diffShapeFunMBHEXAtCenter3y
static constexpr double diffShapeFunMBHEXAtCenter3y
derivative of quad shape function
Definition: Tools.hpp:183
MoFEM::Tools::shapeFunMBEDGE1At00
static constexpr double shapeFunMBEDGE1At00
Definition: Tools.hpp:56
diffN_MBTRI1x
#define diffN_MBTRI1x
derivative of triangle shape function
Definition: fem_tools.h:51
MoFEM::Tools::volumeLengthQuality
static double volumeLengthQuality(const double *coords)
Calculate tetrahedron volume length quality.
Definition: Tools.cpp:15
MoFEM::Tools::getEdgeLength
static double getEdgeLength(const double *edge_coords)
Get edge length.
Definition: Tools.cpp:415
diffN_MBHEX6z
#define diffN_MBHEX6z(x, y)
Definition: fem_tools.h:101
MoFEM::Tools::diffShapeFunMBTRI1y
static constexpr double diffShapeFunMBTRI1y
derivative of triangle shape function
Definition: Tools.hpp:97
MoFEM::Tools::diffShapeFunMBTET3z
static constexpr double diffShapeFunMBTET3z
derivative of tetrahedral shape function
Definition: Tools.hpp:268
MoFEM::Tools::diffShapeFunMBHEXAtCenter6x
static constexpr double diffShapeFunMBHEXAtCenter6x
derivative of HEX shape function
Definition: Tools.hpp:199
MoFEM::Tools::diffShapeFunMBHEXAtCenter1y
static constexpr double diffShapeFunMBHEXAtCenter1y
derivative of HEX shape function
Definition: Tools.hpp:171
MoFEMFunctionBeginHot
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:453
MoFEM::Tools::shapeFunMBTRI2At00
static constexpr double shapeFunMBTRI2At00
Definition: Tools.hpp:114
MoFEM::Tools::outerProductOfEdgeIntegrationPtsForQuad
static MoFEMErrorCode outerProductOfEdgeIntegrationPtsForQuad(MatrixDouble &pts, const int edge0, const int edge1)
Definition: Tools.cpp:610
MoFEM::Problem
keeps basic data about problem
Definition: ProblemsMultiIndices.hpp:54
diffN_MBTET0z
#define diffN_MBTET0z
derivative of tetrahedral shape function
Definition: fem_tools.h:34
MoFEM::Tools::shapeFunMBTETAt000
static constexpr std::array< double, 4 > shapeFunMBTETAt000
Array of shape function at zero local point on reference element.
Definition: Tools.hpp:330
MoFEM::Tools::diffShapeFunMBHEXAtCenter6z
static constexpr double diffShapeFunMBHEXAtCenter6z
derivative of quad shape function
Definition: Tools.hpp:203
MoFEM::Tools::shapeFunMBTRIAt00
static constexpr std::array< double, 3 > shapeFunMBTRIAt00
Array of shape function at zero local point on reference element.
Definition: Tools.hpp:120
MoFEM::Tools::diffShapeFunMBHEXAtCenter2x
static constexpr double diffShapeFunMBHEXAtCenter2x
derivative of HEX shape function
Definition: Tools.hpp:175
MoFEM::Tools::getLocalCoordinatesOnReferenceTriNodeTri
static DEPRECATED MoFEMErrorCode getLocalCoordinatesOnReferenceTriNodeTri(const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords)
Definition: Tools.hpp:398
MoFEM::Tools::shapeFunMBTET3At000
static constexpr double shapeFunMBTET3At000
Definition: Tools.hpp:296
MoFEM::Tools::getTriNormal
static MoFEMErrorCode getTriNormal(const double *coords, double *normal, double *d_normal=nullptr)
Get the Tri Normal objectGet triangle normal.
Definition: Tools.cpp:353
MoFEM::Tools::SOLUTION_EXIST
@ SOLUTION_EXIST
Definition: Tools.hpp:526
diffN_MBHEX7z
#define diffN_MBHEX7z(x, y)
Definition: fem_tools.h:102
MoFEM::Tools::diffShapeFunMBHEXAtCenter3z
static constexpr double diffShapeFunMBHEXAtCenter3z
derivative of quad shape function
Definition: Tools.hpp:185
MoFEM::Tools::diffShapeFunMBHEXAtCenter4x
static constexpr double diffShapeFunMBHEXAtCenter4x
derivative of HEX shape function
Definition: Tools.hpp:187
diffN_MBHEX0x
#define diffN_MBHEX0x(y, z)
Definition: fem_tools.h:79
N_MBTRI2
#define N_MBTRI2(x, y)
triangle shape function
Definition: fem_tools.h:48
diffN_MBTET3x
#define diffN_MBTET3x
derivative of tetrahedral shape function
Definition: fem_tools.h:41
tol
double tol
Definition: mesh_smoothing.cpp:26
MoFEM::Tools::shapeFunMBTET0At000
static constexpr double shapeFunMBTET0At000
Definition: Tools.hpp:293
MoFEM::Tools::shapeFunMBTET1AtOneThird
static constexpr double shapeFunMBTET1AtOneThird
Definition: Tools.hpp:300