31 #ifndef QUADRATUREBOUNDARY_H 32 #define QUADRATUREBOUNDARY_H 1
34 #include <lifev/core/LifeV.hpp> 36 #include <lifev/core/fem/QuadratureRule.hpp> 45 class QuadratureBoundary
49 QuadratureBoundary() {};
51 QuadratureBoundary (
const QuadratureBoundary& qrbd )
56 void setQuadrature (
const UInt i,
const QuadratureRule& qr)
58 while (i >= M_qrs.size() )
60 M_qrs.push_back (QuadratureRule() );
66 QuadratureRule qr (
const UInt i)
const 68 ASSERT (i < M_qrs.size(),
"Invalid quadrature number");
75 for (UInt i (0); i < M_qrs.size(); ++i)
77 std::cout <<
" ## QR " << i <<
" ## " << std::endl;
83 std::vector< QuadratureRule > M_qrs;
89 buildTetraBDQR (
const QuadratureRule& my_qr)
91 QuadratureBoundary qrbd;
94 qrbd.setQuadrature (0, my_qr);
97 QuadratureRule qf1 (
"none", TRIANGLE, 3, 0, 0);
98 for (UInt iq (0); iq < my_qr.nbQuadPt(); ++iq)
100 Real x (my_qr.quadPointCoor (iq, 0) );
101 Real y (my_qr.quadPointCoor (iq, 1) );
102 Real w (my_qr.weight (iq) );
104 qf1.addPoint (QuadraturePoint (x, 0, y, w) );
106 qrbd.setQuadrature (1, qf1);
109 QuadratureRule qf2 (
"none", TRIANGLE, 3, 0, 0);
110 for (UInt iq (0); iq < my_qr.nbQuadPt(); ++iq)
112 Real x (my_qr.quadPointCoor (iq, 0) );
113 Real y (my_qr.quadPointCoor (iq, 1) );
114 Real w (my_qr.weight (iq) );
116 qf2.addPoint (QuadraturePoint (x, y, 1 - x - y, w) );
118 qrbd.setQuadrature (2, qf2);
121 QuadratureRule qf3 (
"none", TRIANGLE, 3, 0, 0);
122 for (UInt iq (0); iq < my_qr.nbQuadPt(); ++iq)
124 Real x (my_qr.quadPointCoor (iq, 0) );
125 Real y (my_qr.quadPointCoor (iq, 1) );
126 Real w (my_qr.weight (iq) );
128 qf3.addPoint (QuadraturePoint (0, y, x, w) );
130 qrbd.setQuadrature (3, qf3);