DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Quadrature.hpp
Go to the documentation of this file.
1#pragma once
2
3// =============================================================================
4// Quadrature.hpp - High-level quadrature interface
5// =============================================================================
6// This file provides the Quadrature class and related utilities for numerical
7// integration over finite elements. The actual quadrature data and low-level
8// access functions are in the Quadratures/ subdirectory.
9//
10// Structure:
11// - Quadratures/QuadratureBase.hpp: Constants and type aliases
12// - Quadratures/*.hpp: Individual parametric space quadrature data
13// - QuadratureHub.hpp: Unified interface for scheme selection and point access
14// - Quadrature.hpp: This file - Quadrature class and integration utilities
15// =============================================================================
16
17#include "QuadratureHub.hpp"
18
19namespace DNDS::Geom::Elem
20{
21 // =========================================================================
22 // SummationNoOp - Helper class for no-op accumulation
23 // =========================================================================
25 {
26 public:
28 {
29 }
30
32 {
33 return {};
34 }
35
37 {
38 return {};
39 }
40
42 {
43 return {};
44 }
45 };
46
47 // =========================================================================
48 // Precomputed Shape Function Buffers at Quadrature Points
49 // =========================================================================
50 // This static structure caches shape function values (D01Nj) at all
51 // quadrature points for all element types and integration orders.
52 // This avoids recomputing shape functions during integration loops.
53 // =========================================================================
54 static struct __TNBufferAtQuadrature
55 {
56 std::array<std::array<std::vector<tD01Nj>, INT_ORDER_MAX + 1>, ElemType_NUM> buf;
57
58 __TNBufferAtQuadrature()
59 {
60 for (t_index i = 1; i < ElemType_NUM; i++)
61 {
62 Element c_elem{ElemType(i)};
63 for (int order = 0; order <= INT_ORDER_MAX; order++)
64 {
65 auto int_scheme = GetQuadratureScheme(c_elem.GetParamSpace(), order);
66 buf.at(i).at(order).resize(int_scheme);
67 for (auto &m : buf.at(i).at(order))
68 m.resize(4, c_elem.GetNumNodes());
69 for (int iG = 0; iG < int_scheme; iG++)
70 {
71 tPoint pParam{0, 0, 0};
72 t_real w;
73 GetQuadraturePoint(c_elem.GetParamSpace(), int_scheme, iG, pParam, w);
74 c_elem.GetD01Nj(pParam, buf.at(i).at(order).at(iG));
75 }
76 }
77 }
78 }
79
80 } __NBufferAtQuadrature{};
81
82 // =========================================================================
83 // Quadrature Class
84 // =========================================================================
85 // Main interface for numerical integration over elements.
86 //
87 // Usage:
88 // Element elem{Quad4};
89 // Quadrature quad(elem, 3); // Order 3 integration
90 //
91 // real sum = 0;
92 // quad.Integration(sum, [](real& acc, int iG, tPoint p, tD01Nj D01Nj) {
93 // acc = function_value(p);
94 // });
95 // =========================================================================
97 {
102
104 : elem(n_elem), int_order(n_int_order), ps(elem.GetParamSpace())
105 {
107 }
108
109 /**
110 * @brief General integration method with shape function access
111 * @param f Callback: f(TAcc& acc, int iG, tPoint pParam, tD01Nj D01Nj)
112 */
113 template <class TAcc, class TFunc>
114 DNDS_DEVICE_CALLABLE void Integration(TAcc &buf, TFunc &&f)
115 {
116 for (t_index iG = 0; iG < int_scheme; iG++)
117 {
118 tPoint pParam{0, 0, 0};
119 t_real w;
120 GetQuadraturePoint(ps, int_scheme, iG, pParam, w);
121 TAcc acc;
122 f(acc, iG, pParam, __NBufferAtQuadrature.buf.at(elem.type).at(int_order).at(iG));
123 buf += acc * w;
124 }
125 }
126
127 /**
128 * @brief Simple integration without shape function access
129 * @param f Callback: f(TAcc& acc, int iG)
130 */
131 template <class TAcc, class TFunc>
132 DNDS_DEVICE_CALLABLE std::enable_if_t<std::is_invocable_v<TFunc, TAcc &, int>> IntegrationSimple(TAcc &buf, TFunc &&f)
133 {
134 for (t_index iG = 0; iG < int_scheme; iG++)
135 {
136 tPoint pParam{0, 0, 0};
137 t_real w;
138 GetQuadraturePoint(ps, int_scheme, iG, pParam, w);
139 TAcc acc;
140 f(acc, iG);
141 buf += acc * w;
142 }
143 }
144
145 /**
146 * @brief Simple integration with weight access
147 * @param f Callback: f(TAcc& acc, int iG, real w)
148 */
149 template <class TAcc, class TFunc>
150 DNDS_DEVICE_CALLABLE std::enable_if_t<std::is_invocable_v<TFunc, TAcc &, int, real>> IntegrationSimple(TAcc &buf, TFunc &&f)
151 {
152 for (t_index iG = 0; iG < int_scheme; iG++)
153 {
154 tPoint pParam{0, 0, 0};
155 t_real w;
156 GetQuadraturePoint(ps, int_scheme, iG, pParam, w);
157 TAcc acc;
158 f(acc, iG, w);
159 buf += acc * w;
160 }
161 }
162
163 DNDS_DEVICE_CALLABLE [[nodiscard]] auto GetQuadraturePointInfo(int iG) const
164 {
165 tPoint pParam{0, 0, 0};
166 t_real w;
167 GetQuadraturePoint(ps, int_scheme, iG, pParam, w);
168 return std::make_tuple(pParam, w);
169 }
170
171 DNDS_DEVICE_CALLABLE [[nodiscard]] real GetWeight(int iG) const
172 {
173 return std::get<1>(GetQuadraturePointInfo(iG));
174 }
175
176 DNDS_DEVICE_CALLABLE [[nodiscard]] t_index GetNumPoints() const { return int_scheme; }
177 };
178
179 // =========================================================================
180 // Quadrature Patch Generation (for visualization/post-processing)
181 // =========================================================================
182 // Generates triangle patches for visualizing quadrature point distributions.
183 // Used for debugging and visualization purposes.
184 // =========================================================================
186 {
187 std::vector<std::array<t_index, 3>> ret;
188 using tArr = std::array<t_index, 3>;
189 ret.reserve(40);
191 {
192 if (q.elem.GetOrder() == 1 || q.elem.GetOrder() == 2) // omitting the middle point
193 {
194 ret.emplace_back(tArr{1, -1, 0});
195 for (int i = 1; i < q.GetNumPoints(); i++)
196 ret.emplace_back(tArr{-i, -1 - i, 0});
197 ret.emplace_back(tArr{-q.GetNumPoints(), 2, 0});
198 }
199 }
201 {
202 if (q.elem.GetOrder() == 1 || q.elem.GetOrder() == 2) // omitting the middle point
203 {
204 t_index nGL = static_cast<t_index>(std::round(std::sqrt(q.GetNumPoints())));
205 auto ijG = [=](int i, int j)
206 { return nGL * j + i; };
207 t_index iMax = nGL - 1;
208
209 if (q.elem.GetOrder() == 1)
210 {
211 ret.emplace_back(tArr{1, 2, -1 - ijG(0, 0)});
212 for (int i = 0; i < iMax; i++)
213 ret.emplace_back(tArr{-1 - ijG(i + 1, 0), -1 - ijG(i, 0), 2});
214 ret.emplace_back(tArr{3, 4, -1 - ijG(iMax, iMax)});
215 for (int i = 0; i < iMax; i++)
216 ret.emplace_back(tArr{-1 - ijG(i, iMax), -1 - ijG(i + 1, iMax), 4});
217 ret.emplace_back(tArr{2, 3, -1 - ijG(iMax, 0)});
218 for (int i = 0; i < iMax; i++)
219 ret.emplace_back(tArr{-1 - ijG(iMax, i + 1), -1 - ijG(iMax, i), 3});
220 ret.emplace_back(tArr{4, 1, -1 - ijG(0, iMax)});
221 for (int i = 0; i < iMax; i++)
222 ret.emplace_back(tArr{-1 - ijG(0, i), -1 - ijG(0, i + 1), 1});
223 }
224 else if (q.elem.GetOrder() == 2)
225 {
226 ret.emplace_back(tArr{1, 5, -1 - ijG(0, 0)});
227 ret.emplace_back(tArr{5, 2, -1 - ijG(iMax, 0)});
228 for (int i = 0; i < iMax; i++)
229 ret.emplace_back(tArr{-1 - ijG(i + 1, 0), -1 - ijG(i, 0), 5});
230 ret.emplace_back(tArr{3, 7, -1 - ijG(iMax, iMax)});
231 ret.emplace_back(tArr{7, 4, -1 - ijG(0, iMax)});
232 for (int i = 0; i < iMax; i++)
233 ret.emplace_back(tArr{-1 - ijG(i, iMax), -1 - ijG(i + 1, iMax), 7});
234 ret.emplace_back(tArr{2, 6, -1 - ijG(iMax, 0)});
235 ret.emplace_back(tArr{6, 3, -1 - ijG(iMax, iMax)});
236 for (int i = 0; i < iMax; i++)
237 ret.emplace_back(tArr{-1 - ijG(iMax, i + 1), -1 - ijG(iMax, i), 6});
238 ret.emplace_back(tArr{4, 8, -1 - ijG(0, iMax)});
239 ret.emplace_back(tArr{8, 1, -1 - ijG(0, 0)});
240 for (int i = 0; i < iMax; i++)
241 ret.emplace_back(tArr{-1 - ijG(0, i), -1 - ijG(0, i + 1), 8});
242 }
243
244 for (int i = 0; i < iMax; i++)
245 for (int j = 0; j < iMax; j++)
246 {
247 ret.emplace_back(tArr{-1 - ijG(i, j), -1 - ijG(i + 1, j), -1 - ijG(i, j + 1)});
248 ret.emplace_back(tArr{-1 - ijG(i + 1, j + 1), -1 - ijG(i, j + 1), -1 - ijG(i + 1, j)});
249 }
250 }
251 }
252 else if (q.elem.GetParamSpace() == ParamSpace::TriSpace)
253 {
254 if (q.GetNumPoints() == 1)
255 {
256 if (q.elem.GetOrder() == 1)
257 {
258 ret.emplace_back(tArr{1, 2, -1});
259 ret.emplace_back(tArr{2, 3, -1});
260 ret.emplace_back(tArr{3, 1, -1});
261 }
262 else if (q.elem.GetOrder() == 2)
263 {
264 ret.emplace_back(tArr{1, 4, -1});
265 ret.emplace_back(tArr{4, 2, -1});
266 ret.emplace_back(tArr{2, 5, -1});
267 ret.emplace_back(tArr{5, 3, -1});
268 ret.emplace_back(tArr{3, 6, -1});
269 ret.emplace_back(tArr{6, 1, -1});
270 }
271 }
272 else if (q.GetNumPoints() == 3)
273 {
274 if (q.elem.GetOrder() == 1)
275 {
276 ret.emplace_back(tArr{1, 2, -3});
277 ret.emplace_back(tArr{2, -1, -3});
278 ret.emplace_back(tArr{2, 3, -1});
279 ret.emplace_back(tArr{3, -2, -1});
280 ret.emplace_back(tArr{3, 1, -2});
281 ret.emplace_back(tArr{1, -3, -2});
282
283 ret.emplace_back(tArr{-1, -2, -3});
284 }
285 else if (q.elem.GetOrder() == 2)
286 {
287 ret.emplace_back(tArr{1, 4, -3});
288 ret.emplace_back(tArr{4, 2, -1});
289 ret.emplace_back(tArr{2, 5, -1});
290 ret.emplace_back(tArr{5, 3, -2});
291 ret.emplace_back(tArr{3, 6, -2});
292 ret.emplace_back(tArr{6, 1, -3});
293
294 ret.emplace_back(tArr{-1, -3, 4});
295 ret.emplace_back(tArr{-2, -1, 5});
296 ret.emplace_back(tArr{-3, -2, 6});
297
298 ret.emplace_back(tArr{-1, -2, -3});
299 }
300 }
301 else if (q.GetNumPoints() == 6)
302 {
303 if (q.elem.GetOrder() == 1)
304 {
305 ret.emplace_back(tArr{1, 2, -6});
306 ret.emplace_back(tArr{2, 3, -4});
307 ret.emplace_back(tArr{3, 1, -5});
308
309 ret.emplace_back(tArr{-6, -1, 1});
310 ret.emplace_back(tArr{-2, -6, 2});
311 ret.emplace_back(tArr{-4, -2, 2});
312 ret.emplace_back(tArr{-3, -4, 3});
313 ret.emplace_back(tArr{-5, -3, 3});
314 ret.emplace_back(tArr{-1, -5, 1});
315
316 ret.emplace_back(tArr{-4, -5, -6});
317 ret.emplace_back(tArr{-1, -6, -5});
318 ret.emplace_back(tArr{-6, -2, -4});
319 ret.emplace_back(tArr{-5, -4, -3});
320 }
321 else if (q.elem.GetOrder() == 2)
322 {
323 ret.emplace_back(tArr{1, 4, -6});
324 ret.emplace_back(tArr{4, 2, -6});
325 ret.emplace_back(tArr{2, 5, -4});
326 ret.emplace_back(tArr{5, 3, -4});
327 ret.emplace_back(tArr{3, 6, -5});
328 ret.emplace_back(tArr{6, 1, -5});
329
330 ret.emplace_back(tArr{-6, -1, 1});
331 ret.emplace_back(tArr{-2, -6, 2});
332 ret.emplace_back(tArr{-4, -2, 2});
333 ret.emplace_back(tArr{-3, -4, 3});
334 ret.emplace_back(tArr{-5, -3, 3});
335 ret.emplace_back(tArr{-1, -5, 1});
336
337 ret.emplace_back(tArr{-4, -5, -6});
338 ret.emplace_back(tArr{-1, -6, -5});
339 ret.emplace_back(tArr{-6, -2, -4});
340 ret.emplace_back(tArr{-5, -4, -3});
341 }
342 }
343 else if (q.GetNumPoints() == 7)
344 {
345 static const std::array<int, 7> idxTri6ToTri7{
346 -1, 5, 6, 7, 2, 3, 4};
347 if (q.elem.GetOrder() == 1)
348 {
349 ret.emplace_back(tArr{1, 2, -idxTri6ToTri7[6]});
350 ret.emplace_back(tArr{2, 3, -idxTri6ToTri7[4]});
351 ret.emplace_back(tArr{3, 1, -idxTri6ToTri7[5]});
352
353 ret.emplace_back(tArr{-idxTri6ToTri7[6], -idxTri6ToTri7[1], 1});
354 ret.emplace_back(tArr{-idxTri6ToTri7[2], -idxTri6ToTri7[6], 2});
355 ret.emplace_back(tArr{-idxTri6ToTri7[4], -idxTri6ToTri7[2], 2});
356 ret.emplace_back(tArr{-idxTri6ToTri7[3], -idxTri6ToTri7[4], 3});
357 ret.emplace_back(tArr{-idxTri6ToTri7[5], -idxTri6ToTri7[3], 3});
358 ret.emplace_back(tArr{-idxTri6ToTri7[1], -idxTri6ToTri7[5], 1});
359
360 ret.emplace_back(tArr{-idxTri6ToTri7[4], -idxTri6ToTri7[5], -1});
361 ret.emplace_back(tArr{-idxTri6ToTri7[5], -idxTri6ToTri7[6], -1});
362 ret.emplace_back(tArr{-idxTri6ToTri7[6], -idxTri6ToTri7[4], -1});
363
364 ret.emplace_back(tArr{-idxTri6ToTri7[1], -idxTri6ToTri7[6], -idxTri6ToTri7[5]});
365 ret.emplace_back(tArr{-idxTri6ToTri7[6], -idxTri6ToTri7[2], -idxTri6ToTri7[4]});
366 ret.emplace_back(tArr{-idxTri6ToTri7[5], -idxTri6ToTri7[4], -idxTri6ToTri7[3]});
367 }
368 else if (q.elem.GetOrder() == 2)
369 {
370 ret.emplace_back(tArr{1, 4, -idxTri6ToTri7[6]});
371 ret.emplace_back(tArr{4, 2, -idxTri6ToTri7[6]});
372 ret.emplace_back(tArr{2, 5, -idxTri6ToTri7[4]});
373 ret.emplace_back(tArr{5, 3, -idxTri6ToTri7[4]});
374 ret.emplace_back(tArr{3, 6, -idxTri6ToTri7[5]});
375 ret.emplace_back(tArr{6, 1, -idxTri6ToTri7[5]});
376
377 ret.emplace_back(tArr{-idxTri6ToTri7[6], -idxTri6ToTri7[1], 1});
378 ret.emplace_back(tArr{-idxTri6ToTri7[2], -idxTri6ToTri7[6], 2});
379 ret.emplace_back(tArr{-idxTri6ToTri7[4], -idxTri6ToTri7[2], 2});
380 ret.emplace_back(tArr{-idxTri6ToTri7[3], -idxTri6ToTri7[4], 3});
381 ret.emplace_back(tArr{-idxTri6ToTri7[5], -idxTri6ToTri7[3], 3});
382 ret.emplace_back(tArr{-idxTri6ToTri7[1], -idxTri6ToTri7[5], 1});
383
384 ret.emplace_back(tArr{-idxTri6ToTri7[4], -idxTri6ToTri7[5], -1});
385 ret.emplace_back(tArr{-idxTri6ToTri7[5], -idxTri6ToTri7[6], -1});
386 ret.emplace_back(tArr{-idxTri6ToTri7[6], -idxTri6ToTri7[4], -1});
387
388 ret.emplace_back(tArr{-idxTri6ToTri7[1], -idxTri6ToTri7[6], -idxTri6ToTri7[5]});
389 ret.emplace_back(tArr{-idxTri6ToTri7[6], -idxTri6ToTri7[2], -idxTri6ToTri7[4]});
390 ret.emplace_back(tArr{-idxTri6ToTri7[5], -idxTri6ToTri7[4], -idxTri6ToTri7[3]});
391 }
392 }
393 else if (q.GetNumPoints() == 12)
394 {
395 if (q.elem.GetOrder() == 1)
396 {
397 ret.emplace_back(tArr{1, 2, -4});
398 ret.emplace_back(tArr{2, 3, -5});
399 ret.emplace_back(tArr{3, 1, -6});
400
401 ret.emplace_back(tArr{-7, -4, 2});
402 ret.emplace_back(tArr{-9, -7, 2});
403 ret.emplace_back(tArr{-5, -9, 2});
404
405 ret.emplace_back(tArr{-11, -5, 3});
406 ret.emplace_back(tArr{-12, -11, 3});
407 ret.emplace_back(tArr{-6, -12, 3});
408
409 ret.emplace_back(tArr{-10, -6, 1});
410 ret.emplace_back(tArr{-8, -10, 1});
411 ret.emplace_back(tArr{-4, -8, 1});
412
413 ret.emplace_back(tArr{-4, -7, -8});
414 ret.emplace_back(tArr{-9, -5, -11});
415 ret.emplace_back(tArr{-10, -12, -6});
416 ret.emplace_back(tArr{-8, -7, -1});
417 ret.emplace_back(tArr{-2, -9, -11});
418 ret.emplace_back(tArr{-3, -12, -10});
419
420 ret.emplace_back(tArr{-7, -9, -2});
421 ret.emplace_back(tArr{-7, -2, -1});
422 ret.emplace_back(tArr{-11, -12, -3});
423 ret.emplace_back(tArr{-11, -3, -2});
424 ret.emplace_back(tArr{-10, -8, -1});
425 ret.emplace_back(tArr{-10, -1, -3});
426
427 ret.emplace_back(tArr{-1, -2, -3});
428 }
429 else if (q.elem.GetOrder() == 2)
430 {
431 ret.emplace_back(tArr{1, 4, -4});
432 ret.emplace_back(tArr{4, 2, -5});
433 ret.emplace_back(tArr{2, 5, -5});
434 ret.emplace_back(tArr{5, 3, -6});
435 ret.emplace_back(tArr{3, 6, -6});
436 ret.emplace_back(tArr{6, 1, -4});
437
438 ret.emplace_back(tArr{-7, -4, 4});
439 ret.emplace_back(tArr{-9, -7, 4});
440 ret.emplace_back(tArr{-5, -9, 4});
441
442 ret.emplace_back(tArr{-11, -5, 5});
443 ret.emplace_back(tArr{-12, -11, 5});
444 ret.emplace_back(tArr{-6, -12, 5});
445
446 ret.emplace_back(tArr{-10, -6, 6});
447 ret.emplace_back(tArr{-8, -10, 6});
448 ret.emplace_back(tArr{-4, -8, 6});
449
450 ret.emplace_back(tArr{-4, -7, -8});
451 ret.emplace_back(tArr{-9, -5, -11});
452 ret.emplace_back(tArr{-10, -12, -6});
453 ret.emplace_back(tArr{-8, -7, -1});
454 ret.emplace_back(tArr{-2, -9, -11});
455 ret.emplace_back(tArr{-3, -12, -10});
456
457 ret.emplace_back(tArr{-7, -9, -2});
458 ret.emplace_back(tArr{-7, -2, -1});
459 ret.emplace_back(tArr{-11, -12, -3});
460 ret.emplace_back(tArr{-11, -3, -2});
461 ret.emplace_back(tArr{-10, -8, -1});
462 ret.emplace_back(tArr{-10, -1, -3});
463
464 ret.emplace_back(tArr{-1, -2, -3});
465 }
466 }
467 else
468 DNDS_assert(false);
469 }
470 else
471 DNDS_assert(false);
472 return ret;
473 }
474
475 // =========================================================================
476 // Jacobian Determinant Helpers
477 // =========================================================================
478 template <int dim>
480 {
481 using namespace Geom;
482 real JDet{0};
483 tJacobi J = Elem::ShapeJacobianCoordD01Nj(coordsCell, DiNj);
484 if constexpr (dim == 2)
485 JDet = J(EigenAll, 0).cross(J(EigenAll, 1)).stableNorm();
486 else
487 JDet = J.fullPivLu().determinant();
488 return JDet;
489 }
490
491 inline real CellJacobianDet(int dim, const Geom::tSmallCoords &coordsCell, const Geom::Elem::tD01Nj &DiNj)
492 {
493 if (dim == 2)
494 return CellJacobianDet<2>(coordsCell, DiNj);
495 else
496 return CellJacobianDet<3>(coordsCell, DiNj);
497 }
498
499 template <int dim>
501 {
502 using namespace Geom;
503 real JDet{0};
504 tJacobi J = Elem::ShapeJacobianCoordD01Nj(coords, DiNj);
505 if constexpr (dim == 2)
506 JDet = J(EigenAll, 0).stableNorm();
507 else
508 JDet = J(EigenAll, 0).cross(J(EigenAll, 1)).stableNorm();
509 return JDet;
510 }
511
512 inline real FaceJacobianDet(int dim, const Geom::tSmallCoords &coords, const Geom::Elem::tD01Nj &DiNj)
513 {
514 if (dim == 2)
515 return FaceJacobianDet<2>(coords, DiNj);
516 else
517 return FaceJacobianDet<3>(coords, DiNj);
518 }
519
520} // namespace DNDS::Geom::Elem
#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
SummationNoOp operator*(const SummationNoOp &R) const
friend SummationNoOp operator*(const SummationNoOp &L, t_real R)
void operator+=(const SummationNoOp &R)
friend SummationNoOp operator*(t_real L, const SummationNoOp &R)
real FaceJacobianDet(const Geom::tSmallCoords &coords, const Geom::Elem::tD01Nj &DiNj)
constexpr t_index GetQuadratureScheme(ParamSpace ps, int int_order)
void GetQuadraturePoint(ParamSpace ps, t_index scheme, int iG, TPoint &pParam, t_real &w)
Eigen::Matrix< t_real, 4, Eigen::Dynamic > tD01Nj
Definition Elements.hpp:183
real CellJacobianDet(const Geom::tSmallCoords &coordsCell, const Geom::Elem::tD01Nj &DiNj)
auto GetQuadPatches(Quadrature &q)
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
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:105
DNDS_DEVICE_CALLABLE constexpr ParamSpace GetParamSpace() const
Definition Elements.hpp:192
DNDS_DEVICE_CALLABLE constexpr t_index GetOrder() const
Definition Elements.hpp:203
DNDS_DEVICE_CALLABLE std::enable_if_t< std::is_invocable_v< TFunc, TAcc &, int > > IntegrationSimple(TAcc &buf, TFunc &&f)
Simple integration without shape function access.
DNDS_DEVICE_CALLABLE t_index GetNumPoints() const
DNDS_DEVICE_CALLABLE Quadrature(Element n_elem=Element{UnknownElem}, int n_int_order=0)
DNDS_DEVICE_CALLABLE auto GetQuadraturePointInfo(int iG) const
DNDS_DEVICE_CALLABLE std::enable_if_t< std::is_invocable_v< TFunc, TAcc &, int, real > > IntegrationSimple(TAcc &buf, TFunc &&f)
Simple integration with weight access.
DNDS_DEVICE_CALLABLE real GetWeight(int iG) const
DNDS_DEVICE_CALLABLE void Integration(TAcc &buf, TFunc &&f)
General integration method with shape function access.
double order
Definition test_ODE.cpp:257