Class documentation of Concepts

Loading...
Searching...
No Matches
quad.hh
Go to the documentation of this file.
1
6#ifndef quad_hh
7#define quad_hh
8
9#include <iostream>
10#include <memory>
11#include "basics/typedefs.hh"
14#include "geometry/cell2D.hh"
15#include "geometry/integral.hh"
17#include "space/hpMethod.hh"
18#include "element.hh"
21
22namespace hp2D {
23
24 using concepts::Real;
25
26 // ******************************************************** IntegrableQuad **
27
32 public:
34
40 inline concepts::Real2d chi(const Real x, const Real y) const {
41 return cell_.elemMap( concepts::Real2d(x,y) );
42 }
43
45 inline concepts::MapReal2d jacobian(const Real x, const Real y) const
46 {
47 concepts::Quad2d* cell = dynamic_cast<concepts::Quad2d*>(&cell_);
49 concepts::MapReal2d jac = cell->jacobian(x, y);
50
52 if(intRule_.get()->domain())
53 jac *= 0.5; // factor comes due to integration over [-1,1]^2
54
55 return jac;
56 }
57
60 const Real y) const
61 {
62 return jacobian(x, y).inverse();
63 }
64
65 // Computes the second partial derivatives of the inverse of the Mapping.
66 inline concepts::MapReal2d inverseLaplace(const Real x, const Real y) const
67 {
68 // TODO: move up the inverseLaplace to QuadNd and then
69 // remove this assert and dynamic cast again
70 const concepts::Quad2d *cell = dynamic_cast<concepts::Quad2d*>(&cell_);
72
73 concepts::MapReal2d invLap = cell->inverseLaplace(x,y);
74
76 if(intRule_.get()->domain())
77 invLap *= 2;
78
79 return invLap;
80 }
81
82
83// /// Computes the determinant of the Jacobian
84// inline Real jacobianDeterminant(const Real x, const Real y) const
85// {
86// return jacobian(x, y).determinant();
87// }
88
91 inline Real gramDeterminantRoot(const Real x, const Real y) const
92 {
93 Real gramian = cell_.gramDeterminantRoot(x,y);
94
96 //if integration is not on [0,1]^2, we need to scale
97 if(intRule_.get()->domain())
98 gramian *= 0.25; // factor due to the integration over [-1,1]^2
99 return gramian;
100 }
101
109 void setStrategy(const concepts::Quad2dSubdivision* strategy = 0)
110 {
111 // TODO: move up the subdivision strategy to QuadNd and then
112 // remove this assert and dynamic cast again
113 concepts::Quad2d *cell = dynamic_cast<concepts::Quad2d*>(&cell_);
115 cell->setStrategy(strategy);
116 }
117
122 {
123 // TODO: move up the subdivision strategy to QuadNd and then
124 // remove this assert and dynamic cast again
125 const concepts::Quad2d *cell = dynamic_cast<concepts::Quad2d*>(&cell_);
127 return cell->getStrategy();
128 }
129
130 virtual bool quadraturePoint(uint i, intPoint& p, intFormType form = ZERO,
131 bool localCoord = false) const;
132
133
134 // Returns the integration rule
135 const concepts::QuadratureRule2d* const intRule() const { return intRule_.get(); }
136
144 static std::unique_ptr<concepts::QuadRuleFactoryTensor2d>& factory() { return factory_; }
145
148 static concepts::QuadRuleFactoryTensor2d* factory_rp() { return factory_.get(); }
149
150 protected:
152 //std::unique_ptr<concepts::QuadratureRule1d> intX_;
153 // std::unique_ptr<concepts::QuadratureRule1d> intY_;
154 // Factory for the integration rule
155 //static concepts::QuadRuleFactory rule_;
156
157 // integration rule
158 std::unique_ptr<concepts::QuadratureRule2d> intRule_;
159
162
163
164 // factory for the integration rule, i.e tensored structural
165 static std::unique_ptr<concepts::QuadRuleFactoryTensor2d> factory_;
166 };
167
168
174 void setQuadrature(enum concepts::intRule rule, uint noP);
175
176 // ************************************************************** BaseQuad **
177
180 template<class F = Real>
181 class BaseQuad : public Element<F>, public IntegrableQuad {
182 public:
190 BaseQuad(concepts::QuadNd& cell, const ushort* p,
192 virtual ~BaseQuad();
193
194 virtual const concepts::Quad& support() const { return cell_.connector(); }
195 virtual concepts::Real3d vertex(uint i) const { return cell_.vertex(i); }
196 virtual const concepts::QuadNd& cell() const { return cell_; }
197
199 virtual const concepts::ElementGraphics<F>* graphics() const = 0;
200
201 protected:
202 virtual std::ostream& info(std::ostream& os) const;
203 };
204
205 // **************************************************** QuadShapeFunctions **
206
213 public:
218 QuadShapeFunctions(const ushort p, const concepts::QuadratureRule2d *intRule);
224 QuadShapeFunctions(const ushort* p, const concepts::QuadratureRule2d *intRule);
225 virtual ~QuadShapeFunctions();
226
229 inline const ushort* p() const { return p_; }
230
232 inline const concepts::Karniadakis<1,0>* shpfctX() const {
233 return shpfctX_.get();
234 }
236 inline const concepts::Karniadakis<1,0>* shpfctY() const {
237 return shpfctY_.get(); }
238
240 inline const concepts::Karniadakis<1,1>* shpfctDX() const {
241 return shpfctDX_.get(); }
243 inline const concepts::Karniadakis<1,1>* shpfctDY() const {
244 return shpfctDY_.get(); }
245 protected:
247 // void computeShapefunctions_(const concepts::QuadratureRule1d *intX,
248 // const concepts::QuadratureRule1d *intY);
249
252
253 private:
255 ushort p_[2];
257 std::unique_ptr<concepts::Karniadakis<1,0> > shpfctX_, shpfctY_;
259 std::unique_ptr<concepts::Karniadakis<1,1> > shpfctDX_, shpfctDY_;
260 };
261
262 // ****************************************************************** Quad **
263
270 template<class F = Real>
271 class Quad : public BaseQuad<F>, public QuadShapeFunctions {
272 public:
279 Quad(concepts::QuadNd& cell, const ushort* p,
281
282 virtual ~Quad();
283
289 void recomputeShapefunctions(const uint nq[2]);
290
292 void edgeP(const uint i, const uint& p) {
293 edges_[i] = &p;
294 }
295 const uint edgeP(const uint i) const {
298 return *edges_[i];
299 }
300 protected:
301 virtual std::ostream& info(std::ostream& os) const;
302 private:
304 static std::unique_ptr<concepts::ElementGraphics<F> > graphics_;
306 const uint* edges_[4];
307 };
308
309
310 // ********************************************************** InfiniteQuad **
311
316 class InfiniteQuad : public Element<Real> {
317 public:
326
327 virtual ~InfiniteQuad();
328
329 const concepts::Connector2& support() const { return cell_.connector(); }
330 concepts::Real3d vertex(uint i) const { return cell_.vertex(i); }
331 const concepts::InfiniteQuad2d& cell() const { return cell_; }
332
335 inline const ushort* p() const { return p_; }
336
339 { return intX_.get(); }
340
343 { return intY_.get(); }
344
351 static void setIntegrationRuleX(enum concepts::intRule rule, bool constant,
352 uint points);
353
360 static void setIntegrationRuleY(bool constant, uint points,
361 const Real border);
362
364 static void resetIntegrationRule();
365
367 static std::string integrationRule();
368
370 inline const concepts::Karniadakis<1,0>* shpfctX() const {
371 return shpfctX_.get();
372 }
373
375 inline const concepts::Karniadakis<1,1>* shpfctDX() const {
376 return shpfctDX_.get(); }
377
379 virtual const concepts::LaguerreBasis<0>* shpfctY() const = 0;
380
382 virtual const concepts::LaguerreBasis<1>* shpfctDY() const = 0;
383
384 virtual const concepts::ElementGraphics<Real>* graphics() const = 0;
385
387 static Real& borderY() { return borderY_; }
388 protected:
389 void computeIntegrationRuleFromP_(const ushort* p);
390
391 void computeIntegrationRule_(const uint nq[2]);
392
395
396 virtual std::ostream& info(std::ostream& os) const;
397
399 ushort p_[2];
400
402 std::unique_ptr<concepts::QuadratureRule1d> intX_;
403 std::unique_ptr<concepts::QuadratureRule1d> intY_;
404 private:
407
413 static enum concepts::intRule integrationTypeX_;
416 static enum concepts::intRule integrationTypeY_;
420 static uint constNumerOfPointsX_, constNumerOfPointsY_;
424 static uint addNumberOfPointsX_, addNumberOfPointsY_;
428 static bool useConstantNumberOfPointsX_, useConstantNumberOfPointsY_;
430 static Real borderY_;
431
433 std::unique_ptr<concepts::Karniadakis<1,0> > shpfctX_;
435 std::unique_ptr<concepts::Karniadakis<1,1> > shpfctDX_;
436 };
437
438 // ************************************************** InfiniteLaguerreQuad **
439
452 public:
462 virtual ~InfiniteLaguerreQuad();
463
465 inline const concepts::LaguerreBasis<0>* shpfctY() const {
466 return shpfctY_.get(); }
467
469 inline const concepts::LaguerreBasis<1>* shpfctDY() const {
470 return shpfctDY_.get(); }
471
472 void recomputeShapefunctions();
473
474 virtual const concepts::ElementGraphics<Real>* graphics() const;
475 protected:
476 virtual std::ostream& info(std::ostream& os) const;
477 private:
479 static std::unique_ptr<concepts::ElementGraphics<Real> > graphics_;
480
482 std::unique_ptr<concepts::LaguerreBasis<0> > shpfctY_;
484 std::unique_ptr<concepts::LaguerreBasis<1> > shpfctDY_;
485
487 void computeShapefunctionsY_();
488 };
489
490 // ****************************************************** ArrayQuadWeights **
491
494 class ArrayQuadWeights : public concepts::Array<Real> {
495 public:
501 };
502
503 // ***************************************************** KarniadakisDeriv2 **
504
517 public:
525 KarniadakisDeriv2(const int P, const Real* xP, const int NxP,
526 const bool cache = true);
527
529 // commented, because there is not an assignement operator= . See Stroustrup "The Meaning of Copy" and "Explicit Defaults" chapters.
530 //KarniadakisDeriv2(const KarniadakisDeriv2& arg);
531
533 protected:
534 virtual std::ostream& info(std::ostream& os) const;
535 private:
538 };
539
540 namespace l2 {
541
542 // **************************************************** QuadShapeFunctions **
543
549 class QuadShapeFunctions {
550 public:
555 QuadShapeFunctions(const ushort p, const concepts::QuadratureRule2d *intRule);
561 QuadShapeFunctions(const ushort* p, const concepts::QuadratureRule2d *intRule);
562 virtual ~QuadShapeFunctions();
563
566 inline const ushort* p() const { return p_; }
567
569 inline const KarniadakisDeriv2* shpfctX() const {
570 return shpfctX_.get();
571 }
573 inline const KarniadakisDeriv2* shpfctY() const {
574 return shpfctY_.get(); }
575 protected:
577 void computeShapefunctions_(const concepts::QuadratureRule2d *intRule);
578 private:
580 ushort p_[2];
582 std::unique_ptr<KarniadakisDeriv2 > shpfctX_, shpfctY_;
583 };
584
585 // ****************************************************************** Quad **
586
593 template<class F = Real>
594 class Quad : public BaseQuad<F>, public QuadShapeFunctions {
595 public:
604 Quad(concepts::QuadNd& cell, const ushort* p,
606
607 virtual ~Quad();
608
609 virtual const concepts::ElementGraphics<F>* graphics() const;
613 void recomputeShapefunctions();
614 void recomputeShapefunctions(const uint nq[2]);
615 protected:
616 virtual std::ostream& info(std::ostream& os) const;
617 private:
619 static std::unique_ptr<concepts::ElementGraphics<F> > graphics_;
620 };
621
622 } // namespace hp2D::l2
623
624} // namespace hp2D
625
626#endif // quad_hh
virtual Real3d vertex(uint i) const =0
Returns the coordinates of the ith vertex.
virtual Real3d elemMap(const Real2d &coord_local) const =0
Element map from point local coordinates in 2D.
Real3d vertex(uint i) const
Returns the coordinates of the ith vertex.
InfiniteQuad & connector() const
Returns the quadrilateral connector (topology)
Definition cell2D.hh:628
Mapping< F, DimX, DimY > inverse() const
Returns the inverse of the matrix.
MapReal2d inverseLaplace(const Real xi, const Real eta) const
MapReal2d jacobian(const Real xi, const Real eta) const
const Quad2dSubdivision * getStrategy() const
Definition cell2D.hh:558
virtual void setStrategy(const Quad2dSubdivision *strategy=0)
Quad & connector() const
Returns the quadrilateral connector (topology)
Definition cell2D.hh:234
virtual Real gramDeterminantRoot(const Real xi, const Real eta) const =0
ArrayQuadWeights(const Quad<> &quad)
virtual const concepts::ElementGraphics< F > * graphics() const =0
Returns element graphics class.
virtual concepts::Real3d vertex(uint i) const
Definition quad.hh:195
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
virtual const concepts::Quad & support() const
Definition quad.hh:194
BaseQuad(concepts::QuadNd &cell, const ushort *p, concepts::TColumn< F > *T0, concepts::TColumn< F > *T1)
virtual const concepts::QuadNd & cell() const
Returns the cell on which the element is built.
Definition quad.hh:196
const concepts::LaguerreBasis< 1 > * shpfctDY() const
Returns the derivatives of the shape functions in x direction.
Definition quad.hh:469
InfiniteLaguerreQuad(concepts::InfiniteQuad2d &cell, const ushort *p, concepts::TColumn< Real > *T0, concepts::TColumn< Real > *T1)
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
const concepts::LaguerreBasis< 0 > * shpfctY() const
Returns the shape functions in x direction.
Definition quad.hh:465
std::unique_ptr< concepts::QuadratureRule1d > intX_
The integration rules.
Definition quad.hh:402
virtual const concepts::LaguerreBasis< 1 > * shpfctDY() const =0
Returns the derivatives of the shape functions in x direction.
void computeShapefunctionsX_()
gets the shapefunctions in x direction
ushort p_[2]
Polynomial degree.
Definition quad.hh:399
virtual const concepts::LaguerreBasis< 0 > * shpfctY() const =0
Returns the shape functions in x direction.
static Real & borderY()
Largest value in local y coordinates. Default is 1.0.
Definition quad.hh:387
const concepts::QuadratureRule1d * integrationX() const
Returns the integration rule in x direction.
Definition quad.hh:338
const concepts::Karniadakis< 1, 0 > * shpfctX() const
Returns the shape functions in x direction.
Definition quad.hh:370
const concepts::Karniadakis< 1, 1 > * shpfctDX() const
Returns the derivatives of the shape functions in x direction.
Definition quad.hh:375
static void setIntegrationRuleX(enum concepts::intRule rule, bool constant, uint points)
InfiniteQuad(concepts::InfiniteQuad2d &cell, const ushort *p, concepts::TColumn< Real > *T0, concepts::TColumn< Real > *T1)
static void setIntegrationRuleY(bool constant, uint points, const Real border)
const ushort * p() const
Definition quad.hh:335
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
const concepts::QuadratureRule1d * integrationY() const
Returns the integration rule in y direction.
Definition quad.hh:342
concepts::Real3d vertex(uint i) const
Definition quad.hh:330
const concepts::Connector2 & support() const
Definition quad.hh:329
static void resetIntegrationRule()
Set the standard type of integration.
static std::string integrationRule()
Returns information on the settings of the quadrature rule.
const concepts::InfiniteQuad2d & cell() const
Returns the cell on which the element is built.
Definition quad.hh:331
virtual bool quadraturePoint(uint i, intPoint &p, intFormType form=ZERO, bool localCoord=false) const
concepts::MapReal2d jacobianInverse(const Real x, const Real y) const
Computes the inverse of the Jacobian.
Definition quad.hh:59
static concepts::QuadRuleFactoryTensor2d * factory_rp()
Definition quad.hh:148
concepts::QuadNd & cell_
The cell.
Definition quad.hh:161
concepts::MapReal2d jacobian(const Real x, const Real y) const
Computes the Jacobian.
Definition quad.hh:45
std::unique_ptr< concepts::QuadratureRule2d > intRule_
The integration rules.
Definition quad.hh:158
Real gramDeterminantRoot(const Real x, const Real y) const
Definition quad.hh:91
concepts::Real2d chi(const Real x, const Real y) const
Definition quad.hh:40
static std::unique_ptr< concepts::QuadRuleFactoryTensor2d > & factory()
Definition quad.hh:144
const concepts::Quad2dSubdivision * getStrategy() const
Definition quad.hh:121
void setStrategy(const concepts::Quad2dSubdivision *strategy=0)
Definition quad.hh:109
~KarniadakisDeriv2()
Copy constructor.
KarniadakisDeriv2(const int P, const Real *xP, const int NxP, const bool cache=true)
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
void computeShapefunctions_(const concepts::QuadratureRule2d *intRule)
gets the shapefunctions, used in both constructors
const concepts::Karniadakis< 1, 1 > * shpfctDY() const
Returns the shape functions in y direction.
Definition quad.hh:243
QuadShapeFunctions(const ushort *p, const concepts::QuadratureRule2d *intRule)
QuadShapeFunctions(const ushort p, const concepts::QuadratureRule2d *intRule)
const concepts::Karniadakis< 1, 0 > * shpfctX() const
Returns the shape functions in x direction.
Definition quad.hh:232
const concepts::Karniadakis< 1, 0 > * shpfctY() const
Returns the shape functions in y direction.
Definition quad.hh:236
const concepts::Karniadakis< 1, 1 > * shpfctDX() const
Returns the derivatives of the shape functions in x direction.
Definition quad.hh:240
const ushort * p() const
Definition quad.hh:229
void edgeP(const uint i, const uint &p)
Set polynomial degree of edge i to p.
Definition quad.hh:292
virtual const concepts::ElementGraphics< F > * graphics() const
Returns element graphics class.
virtual std::ostream & info(std::ostream &os) const
Returns information in an output stream.
void recomputeShapefunctions()
Quad(concepts::QuadNd &cell, const ushort *p, concepts::TColumn< F > *T0, concepts::TColumn< F > *T1)
#define conceptsAssert(cond, exc)
double Real
Definition typedefs.hh:39
unsigned short ushort
Abbreviation for unsigned short.
Definition typedefs.hh:51
intRule
Types of integration rules to choose from.
Definition defines.hh:13
void setQuadrature(enum concepts::intRule rule, uint noP)