DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Elements.hpp
Go to the documentation of this file.
1#pragma once
2// Elements.hpp -- Finite element definitions for the Geom module.
3//
4// Architecture
5// ------------
6// Element metadata is split across three layers:
7//
8// ElemEnum.hpp ElemType / ParamSpace enumerations
9// ElementTraits.hpp Per-element compile-time data (nodes, faces, coords,
10// elevation spans, bisection topology, VTK mapping, ...)
11// Elements.hpp (this file) Runtime API consumed by the rest of the
12// codebase: thin wrappers that dispatch to ElementTraits
13// via DispatchElementType(), the Element struct, shape
14// function dispatch, and geometry helpers.
15//
16// Shape functions are auto-generated by tools/gen_shape_functions/ and live
17// in per-element headers under Geom/Elements/<ElemName>.hpp. Each header
18// provides a full specialization of ShapeFuncImpl<ElemType> with static
19// Diff0..Diff3 methods. ShapeFunc_DiNj() resolves the runtime ElemType to
20// the correct specialization via DispatchElementType().
21//
22// Remaining runtime-only logic
23// ----------------------------
24// GetO2ElemBisectVariant() -- chooses the shortest-diagonal bisection variant
25// for Tet10 (3 variants) and Pyramid14 (2 variants) based on physical
26// coordinates. This cannot be made constexpr because it depends on runtime
27// geometry.
28
29#include <type_traits>
30#include <cstdint>
31#include <array>
32#include <tuple>
33#include "DNDS/HardEigen.hpp"
34
35#include "DNDS/Defines.hpp"
36#include "Geometric.hpp"
37#include "BaseFunction.hpp"
38#include "ElemEnum.hpp"
39#include "ElementTraits.hpp"
40
41// Auto-generated per-element shape function headers
56
57namespace DNDS::Geom::Elem
58{
59
60 // ================================================================
61 // Thin wrappers that delegate to ElementTraits via DispatchElementType.
62 // These preserve the old free-function signatures for backward compat.
63 // ================================================================
64
66 {
67 return DispatchElementType(t, [](auto tr) -> ParamSpace
68 { return decltype(tr)::paramSpace; });
69 }
70
72 {
73 return DispatchElementType(t_v, [iFace](auto tr) -> ElemType
74 { return decltype(tr)::GetFaceType(iFace); });
75 }
76
78 {
79 return Elem::ParamSpaceVolume(ps);
80 }
81
83 {
84 return DispatchElementType(t, [](auto tr) -> int
85 { return decltype(tr)::numElevNodes; });
86 }
87
89 {
90 return DispatchElementType(t, [](auto tr) -> ElemType
91 { return decltype(tr)::elevatedType; });
92 }
93
95 {
96 return DispatchElementType(t, [ine](auto tr) -> ElemType
97 {
98 using T = decltype(tr);
99 if constexpr (T::numElevNodes > 0)
100 return T::elevNodeSpanTypes[ine];
101 else
102 return UnknownElem;
103 });
104 }
105
107 {
108 return ParamSpaceO1Elem(ps);
109 }
110
112 {
113 return DispatchElementType(t, [](auto tr) -> t_index
114 {
115 using T = decltype(tr);
116 if constexpr (T::order == 2)
117 return T::numBisect;
118 else
119 return 0;
120 });
121 }
122
124 {
125 return DispatchElementType(t, [i](auto tr) -> ElemType
126 {
127 using T = decltype(tr);
128 if constexpr (T::order == 2)
129 return T::GetBisectElemType(i);
130 else
131 return UnknownElem;
132 });
133 }
134
135 DNDS_DEVICE_CALLABLE static constexpr real __iipow(real x, int y)
136 {
137 if (y == 0)
138 return 1.;
139 if (y == 1)
140 return x;
141 if (y > 1)
142 return x * __iipow(x, y - 1);
143 if (y == -1)
144 return 1. / x;
145 if (y < -1)
146 return 1. / x * __iipow(x, y + 1);
147 return UnInitReal;
148 }
149
150 /// Primary template for per-element shape function implementations.
151 /// Full specializations are auto-generated in Geom/Elements/<Elem>.hpp.
152 template <ElemType>
154
155 /**
156 * @brief Dispatch to per-element shape function code via ShapeFuncImpl<T>.
157 *
158 * Each ShapeFuncImpl specialization provides static Diff0..Diff3 methods.
159 * DispatchElementType (from ElementTraits.hpp) resolves the runtime ElemType
160 * to the correct compile-time specialization.
161 */
162 template <int diffOrder, class TPoint, class TArray>
163 DNDS_DEVICE_CALLABLE void ShapeFunc_DiNj(ElemType t, const TPoint &p, TArray &&v)
164 {
165 static_assert(diffOrder == 0 || diffOrder == 1 || diffOrder == 2 || diffOrder == 3);
166
167 DispatchElementType(t, [&](auto traits)
168 {
169 constexpr ElemType eType = decltype(traits)::elemType;
170 if constexpr (diffOrder == 0)
171 ShapeFuncImpl<eType>::Diff0(p, std::forward<TArray>(v));
172 else if constexpr (diffOrder == 1)
173 ShapeFuncImpl<eType>::Diff1(p, std::forward<TArray>(v));
174 else if constexpr (diffOrder == 2)
175 ShapeFuncImpl<eType>::Diff2(p, std::forward<TArray>(v));
176 else if constexpr (diffOrder == 3)
177 ShapeFuncImpl<eType>::Diff3(p, std::forward<TArray>(v));
178 });
179 }
180
181 using tNj = Eigen::RowVector<t_real, Eigen::Dynamic>;
182 using tD1Nj = Eigen::Matrix<t_real, 3, Eigen::Dynamic>;
183 using tD01Nj = Eigen::Matrix<t_real, 4, Eigen::Dynamic>;
184 using tDiNj = Eigen::Matrix<t_real, Eigen::Dynamic, Eigen::Dynamic>;
185
186 struct Element
187 {
189
191
192 DNDS_DEVICE_CALLABLE [[nodiscard]] constexpr ParamSpace GetParamSpace() const
193 {
195 }
196
197 DNDS_DEVICE_CALLABLE [[nodiscard]] constexpr t_index GetDim() const
198 {
199 return DispatchElementType(type, [](auto tr) -> t_index
200 { return decltype(tr)::dim; });
201 }
202
203 DNDS_DEVICE_CALLABLE [[nodiscard]] constexpr t_index GetOrder() const
204 {
205 return DispatchElementType(type, [](auto tr) -> t_index
206 { return decltype(tr)::order; });
207 }
208
209 DNDS_DEVICE_CALLABLE [[nodiscard]] constexpr t_index GetNumVertices() const
210 {
211 return DispatchElementType(type, [](auto tr) -> t_index
212 { return decltype(tr)::numVertices; });
213 }
214
215 DNDS_DEVICE_CALLABLE [[nodiscard]] constexpr t_index GetNumNodes() const
216 {
217 return DispatchElementType(type, [](auto tr) -> t_index
218 { return decltype(tr)::numNodes; });
219 }
220
221 DNDS_DEVICE_CALLABLE [[nodiscard]] constexpr t_index GetNumFaces() const
222 {
223 return DispatchElementType(type, [](auto tr) -> t_index
224 { return decltype(tr)::numFaces; });
225 }
226
227 DNDS_DEVICE_CALLABLE [[nodiscard]] constexpr t_index GetNumElev_O1O2() const
228 {
230 }
231
233 {
234 return GetO2ElemBisectNum(type);
235 }
236
237 DNDS_DEVICE_CALLABLE [[nodiscard]] constexpr Element ObtainFace(t_index iFace) const
238 {
239 DNDS_assert(iFace < this->GetNumFaces());
240 return Element{GetFaceType(type, iFace)};
241 }
242
243 /**
244 * @warning assuming Out has correct size
245 */
246 template <class TIn, class TOut>
247 void ExtractFaceNodes(t_index iFace, const TIn &nodes, TOut &faceNodesOut)
248 {
249 DNDS_assert(iFace < this->GetNumFaces());
250 DispatchElementType(type, [&](auto tr)
251 {
252 using T = decltype(tr);
253 if constexpr (T::numFaces > 0)
254 {
255 for (t_index i = 0; i < ObtainFace(iFace).GetNumNodes(); i++)
256 faceNodesOut[i] = nodes[T::faceNodes[iFace][i]];
257 }
258 });
259 }
260
261 [[nodiscard]] constexpr Element ObtainElevNodeSpan(t_index iNodeElev) const
262 {
263 DNDS_assert(iNodeElev < this->GetNumElev_O1O2());
265 }
266
267 [[nodiscard]] constexpr Element ObtainElevatedElem() const
268 {
270 }
271
272 [[nodiscard]] constexpr Element ObtainO1Elem() const
273 {
275 }
276
277 [[nodiscard]] Element ObtainO2BisectElem(t_index iSubElem) const
278 {
279 DNDS_assert(iSubElem < this->GetO2NumBisect());
280 return Element{GetO2ElemBisectElem(type, iSubElem)};
281 }
282
283 template <class TIn, class TOut>
284 void ExtractElevNodeSpanNodes(t_index iNodeElev, const TIn &nodes, TOut &spanNodes)
285 {
286 DNDS_assert(iNodeElev < this->GetNumElev_O1O2());
287 DispatchElementType(type, [&](auto tr)
288 {
289 using T = decltype(tr);
290 if constexpr (T::numElevNodes > 0)
291 {
292 auto spanElem = Element{T::elevNodeSpanTypes[iNodeElev]};
293 for (t_index i = 0; i < spanElem.GetNumNodes(); i++)
294 spanNodes[i] = nodes[T::elevSpans[iNodeElev][i]];
295 }
296 });
297 }
298
299 template <class TIn, class TOut>
300 void ExtractO2BisectElemNodes(t_index iSubElem, t_index iVariant, const TIn &nodes, TOut &subNodes)
301 {
302 DNDS_assert(iSubElem < this->GetO2NumBisect());
303 DispatchElementType(type, [&](auto tr)
304 {
305 using T = decltype(tr);
306 if constexpr (T::order == 2)
307 {
308 DNDS_assert(iVariant < T::numBisectVariants);
309 t_index offset = iSubElem + T::numBisect * iVariant;
310 auto subElem = Element{T::GetBisectElemType(iSubElem)};
311 for (t_index i = 0; i < subElem.GetNumNodes(); i++)
312 subNodes[i] = nodes[T::bisectElements[offset][i]];
313 }
314 });
315 }
316
317 /**
318 * @warning Nj resized within
319 */
320 DNDS_DEVICE_CALLABLE void GetNj(const tPoint &pParam, tNj &Nj) const
321 {
322 Nj.setZero(1, this->GetNumNodes());
323 ShapeFunc_DiNj<0>(type, pParam, Nj);
324 }
325
326 /**
327 * @warning D1Nj resized within
328 */
329 DNDS_DEVICE_CALLABLE void GetD1Nj(const tPoint &pParam, tD1Nj &D1Nj) const
330 {
331 D1Nj.setZero(3, this->GetNumNodes());
332 ShapeFunc_DiNj<1>(type, pParam, D1Nj);
333 }
334
335 /**
336 * @warning DiNj resized within
337 */
338 void GetD01Nj(const tPoint &pParam, tD01Nj &D01Nj) const
339 {
340 tNj Nj;
341 tD1Nj D1Nj;
342 this->GetNj(pParam, Nj);
343 this->GetD1Nj(pParam, D1Nj);
344
345 D01Nj.resize(4, this->GetNumNodes());
346 D01Nj(0, EigenAll) = Nj;
347 D01Nj({1, 2, 3}, EigenAll) = D1Nj;
348 }
349 /**
350 * @warning DiNj resized within
351 */
352 void GetDiNj(const tPoint &pParam, tDiNj &DiNj, int maxOrder) const
353 {
354 if (this->GetDim() == 2)
355 {
356 int nDiffCur = Base::ndiffSizC2D.at(maxOrder);
357 DiNj.setZero(nDiffCur, this->GetNumNodes());
358 for (; maxOrder >= 0; maxOrder--)
359 {
360 if (maxOrder == 3)
361 ShapeFunc_DiNj<3>(
362 type, pParam,
363 DiNj(
364 Eigen::seq(Base::ndiffSizC2D.at(2), Base::ndiffSizC2D.at(3) - 1),
365 EigenAll));
366 if (maxOrder == 2)
367 ShapeFunc_DiNj<2>(
368 type, pParam,
369 DiNj(
370 Eigen::seq(Base::ndiffSizC2D.at(1), Base::ndiffSizC2D.at(2) - 1),
371 EigenAll));
372 if (maxOrder == 1)
373 ShapeFunc_DiNj<1>(
374 type, pParam,
375 DiNj(
376 Eigen::seq(Base::ndiffSizC2D.at(0), Base::ndiffSizC2D.at(1) - 1),
377 EigenAll));
378 if (maxOrder == 0)
379 ShapeFunc_DiNj<0>(
380 type, pParam,
381 DiNj(
382 0,
383 EigenAll));
384 }
385 }
386 else if (this->GetDim() == 3)
387 {
388 int nDiffCur = Base::ndiffSizC.at(maxOrder);
389 DiNj.setZero(nDiffCur, this->GetNumNodes());
390 for (; maxOrder >= 0; maxOrder--)
391 {
392 if (maxOrder == 3)
393 ShapeFunc_DiNj<3>(
394 type, pParam,
395 DiNj(
396 Eigen::seq(Base::ndiffSizC.at(2), Base::ndiffSizC.at(3) - 1),
397 EigenAll));
398 if (maxOrder == 2)
399 ShapeFunc_DiNj<2>(
400 type, pParam,
401 DiNj(
402 Eigen::seq(Base::ndiffSizC.at(1), Base::ndiffSizC.at(2) - 1),
403 EigenAll));
404 if (maxOrder == 1)
405 ShapeFunc_DiNj<1>(
406 type, pParam,
407 DiNj(
408 Eigen::seq(Base::ndiffSizC.at(0), Base::ndiffSizC.at(1) - 1),
409 EigenAll));
410 if (maxOrder == 0)
411 ShapeFunc_DiNj<0>(
412 type, pParam,
413 DiNj(
414 0,
415 EigenAll));
416 }
417 }
418 else if (this->GetDim() == 1)
419 {
420 int nDiffCur = maxOrder + 1;
421 DiNj.setZero(nDiffCur, this->GetNumNodes());
422 for (; maxOrder >= 0; maxOrder--)
423 {
424 if (maxOrder == 3)
425 ShapeFunc_DiNj<3>(
426 type, pParam,
427 DiNj(
428 3,
429 EigenAll));
430 if (maxOrder == 2)
431 ShapeFunc_DiNj<2>(
432 type, pParam,
433 DiNj(
434 2,
435 EigenAll));
436 if (maxOrder == 1)
437 ShapeFunc_DiNj<1>(
438 type, pParam,
439 DiNj(
440 1,
441 EigenAll));
442 if (maxOrder == 0)
443 ShapeFunc_DiNj<0>(
444 type, pParam,
445 DiNj(
446 0,
447 EigenAll));
448 }
449 }
450 }
451 };
452
453 Eigen::Matrix<t_real, 3, Eigen::Dynamic> GetStandardCoord(ElemType t);
454}
455
456namespace DNDS::Geom::Elem
457{
459 const std::vector<DNDS::index> &nodes_A,
460 const std::vector<DNDS::index> &nodes_B,
461 Element eA,
462 Element eB)
463 {
464 DNDS_assert(nodes_A.size() >= eA.GetNumNodes());
465 DNDS_assert(nodes_B.size() >= eB.GetNumVertices());
466 std::vector<DNDS::index> nodes_B_Vert{nodes_B.begin(), nodes_B.begin() + eB.GetNumVertices()};
467 std::sort(nodes_B_Vert.begin(), nodes_B_Vert.end());
468 for (int iF = 0; iF < eA.GetNumFaces(); iF++)
469 {
470 auto eF = eA.ObtainFace(iF);
471 std::vector<DNDS::index> fNodes(eF.GetNumNodes());
472 eA.ExtractFaceNodes(iF, nodes_A, fNodes);
473 std::sort(fNodes.begin(), fNodes.begin() + eF.GetNumVertices());
474 if (std::includes(
475 nodes_B_Vert.begin(), nodes_B_Vert.end(),
476 fNodes.begin(), fNodes.begin() + eF.GetNumVertices()))
477 return true;
478 }
479 return false;
480 }
481
483 const std::vector<DNDS::index> &verts_A,
484 const std::vector<DNDS::index> &nodes_B)
485 {
486
487 return false;
488 }
489
490 template <class tCoordsIn>
491 tJacobi ShapeJacobianCoordD01Nj(const tCoordsIn &cs, Eigen::Ref<const tD01Nj> DiNj)
492 {
493 return cs * DiNj({1, 2, 3}, EigenAll).transpose();
494 }
495
496 template <class tCoordsIn>
497 tPoint PPhysicsCoordD01Nj(const tCoordsIn &cs, Eigen::Ref<const tD01Nj> DiNj)
498 {
499 return cs * DiNj(0, EigenAll).transpose();
500 }
501
502 template <int dim>
503 real JacobiDetFace(Eigen::Ref<const tJacobi> J)
504 {
505 static_assert(dim == 2 || dim == 3);
506 if constexpr (dim == 2)
507 return J(EigenAll, 0).stableNorm();
508 else
509 return J(EigenAll, 0).cross(J(EigenAll, 1)).stableNorm();
510 }
511
512 // TODO: add a integration-based counterpart
514 {
515 tPoint c = coords.rowwise().mean();
516 tSmallCoords coordsC = coords.colwise() - c;
517 tPoint ve0 = coordsC(EigenAll, 1) - coordsC(EigenAll, 0);
518 tJacobi inertia = coordsC * coordsC.transpose();
519 real cond01 = HardEigen::Eigen3x3RealSymEigenDecompositionGetCond01(inertia); // ratio of 2 largest eigenvalues
520 if (cond01 < 1 + 1e-6)
521 inertia += ve0 * ve0.transpose() * 1e-4; // first edge gets priority
523 coordsC = decRet.transpose() * coordsC;
524 tPoint ret = coordsC.rowwise().maxCoeff() - coordsC.rowwise().minCoeff();
525 std::sort(ret.begin(), ret.end(), std::greater_equal<real>());
526 return ret;
527 }
528
530 {
531 DNDS_assert(e.GetOrder() == 2);
532 //! if serendipity ?
533 switch (e.type)
534 {
535 case Tet10:
536 {
537 Eigen::Vector<real, 3> lengths{
538 (coords(EigenAll, 5 - 1) - coords(EigenAll, 10 - 1)).norm(),
539 (coords(EigenAll, 6 - 1) - coords(EigenAll, 8 - 1)).norm(),
540 (coords(EigenAll, 7 - 1) - coords(EigenAll, 9 - 1)).norm(),
541 };
542 Eigen::Index iMin{-1};
543 real minLen = lengths.minCoeff(&iMin);
544 (void)minLen;
545 return t_index(iMin);
546 }
547 break;
548 case Pyramid14:
549 {
550 Eigen::Vector<real, 2> lengths{
551 (coords(EigenAll, 10 - 1) - coords(EigenAll, 12 - 1)).norm(),
552 (coords(EigenAll, 11 - 1) - coords(EigenAll, 13 - 1)).norm(),
553 };
554 Eigen::Index iMin{-1};
555 real minLen = lengths.minCoeff(&iMin);
556 (void)minLen;
557 return t_index(iMin);
558 }
559 break;
560 default:
561 return 0;
562 }
563 }
564
565 template <class TIn>
566 std::pair<int, std::vector<index>> ToVTKVertsAndData(Element e, const TIn &vin)
567 {
568 return DispatchElementType(e.type, [&](auto tr) -> std::pair<int, std::vector<index>>
569 {
570 using T = decltype(tr);
571 constexpr int nOut = static_cast<int>(
572 std::tuple_size<decltype(T::vtkNodeOrder)>::value);
573 std::vector<index> v(nOut);
574 for (int i = 0; i < nOut; i++)
575 v[i] = vin[T::vtkNodeOrder[i]];
576 return {T::vtkCellType, std::move(v)};
577 });
578 }
579}
Core type aliases, constants, and metaprogramming utilities for the DNDS framework.
#define DNDS_DEVICE_TRIVIAL_COPY_DEFINE(T, T_Self)
Definition Defines.hpp:83
#define DNDS_DEVICE_CALLABLE
Definition Defines.hpp:76
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:108
Robust linear-algebra primitives for small / moderately-sized Eigen matrices that are more numericall...
Eigen::Matrix< t_real, Eigen::Dynamic, Eigen::Dynamic > tDiNj
Definition Elements.hpp:184
Eigen::Matrix< t_real, 3, Eigen::Dynamic > tD1Nj
Definition Elements.hpp:182
tJacobi ShapeJacobianCoordD01Nj(const tCoordsIn &cs, Eigen::Ref< const tD01Nj > DiNj)
Definition Elements.hpp:491
DNDS_DEVICE_CALLABLE constexpr ParamSpace ElemType_to_ParamSpace(const ElemType t)
Definition Elements.hpp:65
DNDS_DEVICE_CALLABLE constexpr ElemType GetElemO1(ParamSpace ps)
Definition Elements.hpp:106
DNDS_DEVICE_CALLABLE constexpr decltype(auto) DispatchElementType(ElemType t, Func &&func)
Static dispatch over element types at compile time.
DNDS_DEVICE_CALLABLE constexpr int GetElemElevation_O1O2_NumNode(ElemType t)
Definition Elements.hpp:82
DNDS_DEVICE_CALLABLE constexpr t_index GetO2ElemBisectNum(ElemType t)
Definition Elements.hpp:111
DNDS_DEVICE_CALLABLE constexpr ElemType GetO2ElemBisectElem(ElemType t, t_index i)
Definition Elements.hpp:123
DNDS_DEVICE_CALLABLE constexpr ElemType GetElemElevation_O1O2_ElevatedType(ElemType t)
Definition Elements.hpp:88
DNDS_DEVICE_CALLABLE constexpr t_real ParamSpaceVolume(ParamSpace ps)
Get the volume of a parametric space.
Eigen::Matrix< t_real, 4, Eigen::Dynamic > tD01Nj
Definition Elements.hpp:183
bool cellsAreFaceConnected(const std::vector< DNDS::index > &nodes_A, const std::vector< DNDS::index > &nodes_B, Element eA, Element eB)
Definition Elements.hpp:458
std::pair< int, std::vector< index > > ToVTKVertsAndData(Element e, const TIn &vin)
Definition Elements.hpp:566
DNDS_DEVICE_CALLABLE void ShapeFunc_DiNj(ElemType t, const TPoint &p, TArray &&v)
Dispatch to per-element shape function code via ShapeFuncImpl<T>.
Definition Elements.hpp:163
tPoint GetElemNodeMajorSpan(const tSmallCoords &coords)
Definition Elements.hpp:513
real JacobiDetFace(Eigen::Ref< const tJacobi > J)
Definition Elements.hpp:503
DNDS_DEVICE_CALLABLE constexpr t_real ParamSpaceVol(ParamSpace ps)
Definition Elements.hpp:77
DNDS_DEVICE_CALLABLE constexpr ElemType GetElemElevation_O1O2_NodeSpanType(ElemType t, t_index ine)
Definition Elements.hpp:94
DNDS_DEVICE_CALLABLE constexpr ElemType GetFaceType(ElemType t_v, t_index iFace)
Definition Elements.hpp:71
DNDS_DEVICE_CALLABLE constexpr ElemType ParamSpaceO1Elem(ParamSpace ps)
Get the order-1 (linear) element type for a parametric space.
Eigen::Matrix< t_real, 3, Eigen::Dynamic > GetStandardCoord(ElemType t)
Definition Elements.cpp:6
Eigen::RowVector< t_real, Eigen::Dynamic > tNj
Definition Elements.hpp:181
tPoint PPhysicsCoordD01Nj(const tCoordsIn &cs, Eigen::Ref< const tD01Nj > DiNj)
Definition Elements.hpp:497
t_index GetO2ElemBisectVariant(Element e, const tSmallCoords &coords)
Definition Elements.hpp:529
int32_t t_index
Definition Geometric.hpp:6
Eigen::Matrix3d tJacobi
Definition Geometric.hpp:10
Eigen::Vector3d tPoint
Definition Geometric.hpp:9
Eigen::Matrix< real, 3, Eigen::Dynamic > tSmallCoords
Definition Geometric.hpp:50
double t_real
Definition Geometric.hpp:8
real Eigen3x3RealSymEigenDecompositionGetCond01(const Eigen::Matrix3d &A)
Like Eigen3x3RealSymEigenDecompositionGetCond but returns lambda0 / lambda1 only (ignores the smalles...
Eigen::Matrix3d Eigen3x3RealSymEigenDecompositionNormalized(const Eigen::Matrix3d &A)
Eigen-decomposition with eigenvector columns normalised to unit length.
DNDS_CONSTANT const real UnInitReal
Sentinel "not initialised" real value (NaN). Cheap to detect with std::isnan or IsUnInitReal; survive...
Definition Defines.hpp:174
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:105
constexpr Element ObtainO1Elem() const
Definition Elements.hpp:272
void ExtractO2BisectElemNodes(t_index iSubElem, t_index iVariant, const TIn &nodes, TOut &subNodes)
Definition Elements.hpp:300
DNDS_DEVICE_CALLABLE constexpr ParamSpace GetParamSpace() const
Definition Elements.hpp:192
constexpr Element ObtainElevNodeSpan(t_index iNodeElev) const
Definition Elements.hpp:261
DNDS_DEVICE_CALLABLE void GetNj(const tPoint &pParam, tNj &Nj) const
Definition Elements.hpp:320
DNDS_DEVICE_CALLABLE constexpr Element ObtainFace(t_index iFace) const
Definition Elements.hpp:237
DNDS_DEVICE_CALLABLE t_index GetO2NumBisect() const
Definition Elements.hpp:232
constexpr Element ObtainElevatedElem() const
Definition Elements.hpp:267
DNDS_DEVICE_CALLABLE constexpr t_index GetNumElev_O1O2() const
Definition Elements.hpp:227
DNDS_DEVICE_CALLABLE constexpr t_index GetNumVertices() const
Definition Elements.hpp:209
DNDS_DEVICE_CALLABLE constexpr t_index GetNumFaces() const
Definition Elements.hpp:221
DNDS_DEVICE_CALLABLE void GetD1Nj(const tPoint &pParam, tD1Nj &D1Nj) const
Definition Elements.hpp:329
void GetD01Nj(const tPoint &pParam, tD01Nj &D01Nj) const
Definition Elements.hpp:338
void ExtractFaceNodes(t_index iFace, const TIn &nodes, TOut &faceNodesOut)
Definition Elements.hpp:247
void GetDiNj(const tPoint &pParam, tDiNj &DiNj, int maxOrder) const
Definition Elements.hpp:352
Element ObtainO2BisectElem(t_index iSubElem) const
Definition Elements.hpp:277
DNDS_DEVICE_CALLABLE constexpr t_index GetDim() const
Definition Elements.hpp:197
void ExtractElevNodeSpanNodes(t_index iNodeElev, const TIn &nodes, TOut &spanNodes)
Definition Elements.hpp:284
DNDS_DEVICE_CALLABLE constexpr t_index GetOrder() const
Definition Elements.hpp:203
DNDS_DEVICE_CALLABLE constexpr t_index GetNumNodes() const
Definition Elements.hpp:215
Eigen::Matrix< real, 5, 1 > v
tVec x(NCells)
double order
Definition test_ODE.cpp:257