GemaCoreLib
The GeMA Core library
gmWedgeIntegrationRule.h
Go to the documentation of this file.
1 /************************************************************************
2 **
3 ** Copyright (C) 2014 by Carlos Augusto Teixera Mendes
4 ** All rights reserved.
5 **
6 ** This file is part of the "GeMA" software. It's use should respect
7 ** the terms in the license agreement that can be found together
8 ** with this source code.
9 ** It is provided AS IS, with NO WARRANTY OF ANY KIND,
10 ** INCLUDING THE WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR
11 ** A PARTICULAR PURPOSE.
12 **
13 ************************************************************************/
14 
26 #ifndef _GEMA_WEDGE_INT_RULE_H_
27 #define _GEMA_WEDGE_INT_RULE_H_
28 
29 #include "gmIntegrationRule.h"
31 #include "gmBarIntegrationRule.h" // Needed by edge rules
32 #include "gmQuadIntegrationRule.h" // Needed by face rules
33 #include "gmTriIntegrationRule.h" // Needed by face rules
34 
35 //--------------------------------------------------------
36 // Element rules
37 //--------------------------------------------------------
38 
39 
42 {
43 public:
44  GmWedgeGaussIntegrationRule(int triRule, int linearRule);
45 
47  GmWedgeGaussIntegrationRule(GmIntegrationRuleType irType, int rule1, int rule2, int rule3)
48  : GmWedgeGaussIntegrationRule(rule1, rule2)
49  {
50  assert(irType == GM_GAUSS_RULE_TYPE); Q_UNUSED(irType); Q_UNUSED(rule3);
51  }
52 
53  // See comments on the base class
54  bool isValid() const { return _tTable; }
55 
56  // See comments on the base class
57  virtual int numPoints() const { assert(isValid()); return _trule * _lrule; }
58 
59  virtual void integrationPoint(int index, GmVector& naturalCoord, double* weight) const;
60 
61  // See comments on the base class
62  virtual int numNaturalCoord() const { assert(isValid()); return 3; }
63 
64  // See comments on the base class.
65  virtual int degree() const { assert(isValid()); return _degree; }
66 
67  // See comments on the base class
69 
70  // See comments on the base class.
71  virtual QString ruleName() const { return makeRuleName(GM_GAUSS_RULE_TYPE, _trule, _lrule); }
72 
73  // See comments on the base class
74  virtual int rulePar(int num) const { return (num == 1) ? _trule : (num == 2 ? _lrule : -1); }
75 
76  // See comments on the base class
77  virtual int ruleParNumPoints(int num) const { return (num == 1) ? _tTable->size() : (num == 2 ? _lTable->size() : -1); }
78 
79  // See comments on the base class
80  virtual int cacheKey() const { return makeCacheKey(GM_WEDGE_GIR_KEY, false, _trule, _lrule); }
81 
82 private:
83  int _trule;
84  int _lrule;
85  int _degree;
88 };
89 
90 //--------------------------------------------------------
91 // Border rules
92 //--------------------------------------------------------
93 
96 {
97 public:
99  GmWedgeGaussFaceIntegrationRule(GmCellType type, GmCellFaceType faceType, int rule1, int rule2)
100  : GmProxyBorderIntegrationRule(type, 3, faceType == GM_QUAD_FACE ?
101  (GmIntegrationRule*) new GmQuadGaussIntegrationRule(rule1, rule2) :
102  (GmIntegrationRule*) new GmTriGaussIntegrationRule(rule1), true) {}
103 
106  : GmWedgeGaussFaceIntegrationRule(type, faceType, rule1, rule2) { Q_UNUSED(irType); }
107 };
108 
111 {
112 public:
115  : GmProxyBorderIntegrationRule(type, 3, new GmLineGaussIntegrationRule(rule), false) {}
116 
119  : GmWedgeGaussEdgeIntegrationRule(type, rule) { Q_UNUSED(irType); }
120 };
121 
122 
123 #endif
124 
Gauss Integration rule for triangle elements. Base on Felippa, IFEM, section 24.2.
Definition: gmTriIntegrationRule.h:42
GmIntegrationRuleType
The type of desired integration rule (Gauss quadrature, Lobatto quadrature, etc)
Definition: gmIntegrationRule.h:67
Declaration of the GmBorderIntegrationRule and GmProxyBorderIntegrationRule classes.
GmWedgeGaussFaceIntegrationRule(GmCellType type, GmCellFaceType faceType, int rule1, int rule2)
Constructor.
Definition: gmWedgeIntegrationRule.h:99
virtual int numNaturalCoord() const
Returns the number of natural coordinates used by this integration rule.
Definition: gmWedgeIntegrationRule.h:62
Gauss rules.
Definition: gmIntegrationRule.h:75
GmCellFaceType
The type of a 3D element's face.
Definition: gmCellType.h:135
virtual int cacheKey() const
Returns an unique integer that can be used to uniquelly represent this kind of integration rule (i....
Definition: gmWedgeIntegrationRule.h:80
int _trule
The triangular rule number. Options: 1, 3 or 6.
Definition: gmWedgeIntegrationRule.h:83
bool isValid() const
Returns true if this is a valid integration rule object, false if not.
Definition: gmWedgeIntegrationRule.h:54
const QVector< QPair< double, double > > * _lTable
Table with point & weight data for the given linear rule. NULL for an invalid object.
Definition: gmWedgeIntegrationRule.h:87
virtual int numPoints() const
Returns the number of integration points returned by this rule.
Definition: gmWedgeIntegrationRule.h:57
Declaration of the GmIntegrationRule class and its helper decendants.
2D Isotropic / anisotropic Gauss integration rules for Quad elements
Definition: gmQuadIntegrationRule.h:42
Integration rule base classe.
Definition: gmIntegrationRule.h:88
Gauss integration rule for Wedge element faces.
Definition: gmWedgeIntegrationRule.h:95
int makeCacheKey(GmIntegrationRuleCacheKeyBase base, bool closed, int rule1, int rule2=-1, int rule3=-1) const
Creates an unique cache key based on the rule type cache key and the rule parameters.
Definition: gmIntegrationRule.h:157
The element face is a quadrilateral.
Definition: gmCellType.h:137
virtual bool isValid() const =0
Returns true if this is a valid integration rule object, false if not.
GmWedgeGaussEdgeIntegrationRule(GmCellType type, GmIntegrationRuleType irType, int rule)
Standard constructor to enable the class use with GmCellIntegrationRuleSet.
Definition: gmWedgeIntegrationRule.h:118
Common code for border integration rules based on a given integration rule that is used to implement ...
Definition: gmBorderIntegrationRule.h:74
Base cache key for 3D Gauss wedge rules.
Definition: gmIntegrationRule.h:47
QString makeRuleName(GmIntegrationRuleType irType, int rule1, int rule2=-1, int rule3=-1) const
Returns a standard rule name based on rule type and rule parameters.
Definition: gmIntegrationRule.cpp:129
#define GMC_API_EXPORT
Macro for controling if the class is being exported (GEMA_CORE_LIB defined) or imported (GEMA_CORE_LI...
Definition: gmCoreConfig.h:35
const QVector< QPair< GmVector, double > > * _tTable
Table with point & weight data for the given triangle rule. NULL for an invalid object.
Definition: gmWedgeIntegrationRule.h:86
Declaration of the GmQuadXxxxIntegrationRule family of classes, including border rules....
GmWedgeGaussIntegrationRule(GmIntegrationRuleType irType, int rule1, int rule2, int rule3)
Standard constructor to enable the class use with GmCellIntegrationRuleSet.
Definition: gmWedgeIntegrationRule.h:47
virtual int ruleParNumPoints(int num) const
Definition: gmWedgeIntegrationRule.h:77
GmCellType
Mesh Cell types. Don't change type orders or add types without reading comments below.
Definition: gmCellType.h:30
virtual int rulePar(int num) const
Returns the value of the numbered rule parameter received in the rule constructor (-1 if unused)....
Definition: gmWedgeIntegrationRule.h:74
Gauss integration rule for Wedge element borders (edges)
Definition: gmWedgeIntegrationRule.h:110
3D Gauss integration rules for Wedge elements
Definition: gmWedgeIntegrationRule.h:41
GmWedgeGaussEdgeIntegrationRule(GmCellType type, int rule)
Constructor.
Definition: gmWedgeIntegrationRule.h:114
virtual int degree() const
Returns the integration degree for this rule (The polynomial degree for which this rule is exact)
Definition: gmWedgeIntegrationRule.h:65
Declaration of the GmTriXxxxIntegrationRule family of classes, including border rules....
arma::vec GmVector
The basic type for a GeMA vector object. Currently based on an Armadillo vector.
Definition: gmVector.h:34
int _degree
The rule degree (minimum of the triangular and linear rule degrees)
Definition: gmWedgeIntegrationRule.h:85
GmWedgeGaussFaceIntegrationRule(GmCellType type, GmCellFaceType faceType, GmIntegrationRuleType irType, int rule1, int rule2)
Standard constructor to enable the class use with GmCellIntegrationRuleSet.
Definition: gmWedgeIntegrationRule.h:105
1D Gauss Integration rules (for bar elements)
Definition: gmBarIntegrationRule.h:36
Declaration of the GmLineXxxxIntegrationRule family of classes Content previously in gmIntegrationRul...
virtual GmIntegrationRuleType ruleType() const
Returns the type of the integration rule.
Definition: gmWedgeIntegrationRule.h:68
virtual void integrationPoint(int index, GmVector &naturalCoord, double *weight) const =0
Given an index from 0 to numPoints() - 1, fills naturalCoord and weight with the position and weigth ...
int _lrule
The linear rule number. From 1 to 5.
Definition: gmWedgeIntegrationRule.h:84
virtual QString ruleName() const
Returns the current ruleName adopted by this integration rule.
Definition: gmWedgeIntegrationRule.h:71