DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
Tri.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "QuadratureBase.hpp"
4
5namespace DNDS::Geom::Elem
6{
7 // ===================================================================
8 // Hammer Quadrature on Reference Triangle
9 // ===================================================================
10 // Reference triangle: vertices at (0,0), (1,0), (0,1)
11 // Area = 1/2
12 //
13 // Data format: [3][n] array where:
14 // [0][i] = xi coordinate (barycentric coord for vertex 1)
15 // [1][i] = eta coordinate (barycentric coord for vertex 2)
16 // [2][i] = quadrature weight (wi)
17 //
18 // Note: The third barycentric coordinate is 1 - xi - eta
19 // ===================================================================
20
21 /// Hammer rule with 1 point: degree 1 exact
22 /// Point: centroid (1/3, 1/3), Weight: 1/2 (total area)
23 static const t_real __HammerTri_1[3][1]{
24 {1. / 3.},
25 {1. / 3.},
26 {1. / 2.}};
27
28 /// Hammer rule with 3 points: degree 2 exact
29 /// Points: (2/3, 1/6), (1/6, 2/3), (1/6, 1/6)
30 /// Weights: 1/6 each (sum = 1/2)
31 static const t_real __HammerTri_3[3][3]{
32 {2. / 3., 1. / 6., 1. / 6.},
33 {1. / 6., 2. / 3., 1. / 6.},
34 {1. / 6., 1. / 6., 1. / 6.}};
35
36 /// Hammer rule with 6 points: degree 4 exact
37 /// Points include permutations of (a, a), (b, b), (c, d)
38 static const t_real __HammerTri_6[3][6]{
39 {0.09157621350977102, 0.8168475729804585, 0.09157621350977052,
40 0.445948490915965, 0.10810301816807, 0.445948490915965},
41 {0.09157621350977102, 0.09157621350977052, 0.8168475729804585,
42 0.445948490915965, 0.445948490915965, 0.10810301816807},
43 {0.054975871827661, 0.054975871827661, 0.054975871827661,
44 0.1116907948390057, 0.1116907948390057, 0.1116907948390057}};
45
46 /// Hammer rule with 7 points: degree 5 exact
47 /// Includes centroid + 6 points in 3 pairs
48 static const t_real __HammerTri_7[3][7]{
49 {0.3333333333333335, 0.470142064105115, 0.05971587178977, 0.470142064105115,
50 0.1012865073234565, 0.7974269853530875, 0.1012865073234565},
51 {0.3333333333333335, 0.470142064105115, 0.470142064105115, 0.05971587178977,
52 0.1012865073234565, 0.1012865073234565, 0.7974269853530875},
53 {0.1125, 0.066197076394253, 0.066197076394253, 0.066197076394253,
54 0.0629695902724135, 0.0629695902724135, 0.0629695902724135}};
55
56 /// Hammer rule with 12 points: degree 6 exact
57 static const t_real __HammerTri_12[3][12]{
58 {0.2492867451709105, 0.501426509658179, 0.2492867451709105,
59 0.063089014491502, 0.8738219710169954, 0.063089014491502,
60 0.3103524510337845, 0.05314504984481699, 0.6365024991213986,
61 0.05314504984481699, 0.6365024991213986, 0.3103524510337845},
62 {0.2492867451709105, 0.2492867451709105, 0.501426509658179,
63 0.063089014491502, 0.063089014491502, 0.8738219710169954,
64 0.05314504984481699, 0.3103524510337845, 0.05314504984481699,
65 0.6365024991213986, 0.3103524510337845, 0.6365024991213986},
66 {0.05839313786318975, 0.05839313786318975, 0.05839313786318975,
67 0.0254224531851035, 0.0254224531851035, 0.0254224531851035,
68 0.04142553780918675, 0.04142553780918675, 0.04142553780918675,
69 0.04142553780918675, 0.04142553780918675, 0.04142553780918675}};
70
71} // namespace DNDS::Geom::Elem
double t_real
Definition Geometric.hpp:8