![]() |
GemaCoreLib
The GeMA Core library
|
3D Gauss integration rules for Wedge elements More...
#include <gmWedgeIntegrationRule.h>
Public Member Functions | |
GmWedgeGaussIntegrationRule (int triRule, int linearRule) | |
Constructor. Receives the rule number for the triangle "base" (1, 3 or 6) plus the rule for the "linear" part (1 to 5). More... | |
GmWedgeGaussIntegrationRule (GmIntegrationRuleType irType, int rule1, int rule2, int rule3) | |
Standard constructor to enable the class use with GmCellIntegrationRuleSet. | |
bool | isValid () const |
Returns true if this is a valid integration rule object, false if not. More... | |
virtual int | numPoints () const |
Returns the number of integration points returned by this rule. | |
virtual void | integrationPoint (int index, GmVector &naturalCoord, double *weight) const |
This function returns weights and points to be sampled for 3D Gauss integration over a wedge with triangular axis (zeta1, zeta2) and linear axis (Xi). More... | |
virtual int | numNaturalCoord () const |
Returns the number of natural coordinates used by this integration rule. | |
virtual int | degree () const |
Returns the integration degree for this rule (The polynomial degree for which this rule is exact) | |
virtual GmIntegrationRuleType | ruleType () const |
Returns the type of the integration rule. | |
virtual QString | ruleName () const |
Returns the current ruleName adopted by this integration rule. | |
virtual int | rulePar (int num) const |
Returns the value of the numbered rule parameter received in the rule constructor (-1 if unused). Num should be a value between 1 and 3. | |
virtual int | ruleParNumPoints (int num) const |
virtual int | cacheKey () const |
Returns an unique integer that can be used to uniquelly represent this kind of integration rule (i.e, two objects with the same cache key should return exactly the same set of integration points). More... | |
![]() | |
virtual | ~GmIntegrationRule () |
Virtual destructor. | |
virtual const QVector< int > & | closestIntegrationIndex (GmCellType type, int P, int Q, int nodeIndex) const |
Returns a list with the index of the closest(s) integraion points from nodeIndex. More... | |
Private Attributes | |
int | _trule |
The triangular rule number. Options: 1, 3 or 6. | |
int | _lrule |
The linear rule number. From 1 to 5. | |
int | _degree |
The rule degree (minimum of the triangular and linear rule degrees) | |
const QVector< QPair< GmVector, double > > * | _tTable |
Table with point & weight data for the given triangle rule. NULL for an invalid object. | |
const QVector< QPair< double, double > > * | _lTable |
Table with point & weight data for the given linear rule. NULL for an invalid object. | |
Additional Inherited Members | |
![]() | |
static GmIntegrationRule * | instance (GmCellType type, int P, int Q, GmIntegrationRuleType irType, int rule1=-1, int rule2=-1, int rule3=-1) |
Instanciates an integration rule for the specified element type using the provided parameters. Returns NULL if the parameter set results in an invalid integration rule. More... | |
![]() | |
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. | |
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. | |
3D Gauss integration rules for Wedge elements
GmWedgeGaussIntegrationRule::GmWedgeGaussIntegrationRule | ( | int | tRule, |
int | lRule | ||
) |
Constructor. Receives the rule number for the triangle "base" (1, 3 or 6) plus the rule for the "linear" part (1 to 5).
Creates an invalid object for rules outside the valid range. See comments on GmIntegrationRule::isValid()
|
inlinevirtual |
Returns an unique integer that can be used to uniquelly represent this kind of integration rule (i.e, two objects with the same cache key should return exactly the same set of integration points).
The returned value is used together with an element type as a cache key in closestIntegrationIndex() and in GmShape::gaussExtrapolationMatrix().
Implements GmIntegrationRule.
|
virtual |
This function returns weights and points to be sampled for 3D Gauss integration over a wedge with triangular axis (zeta1, zeta2) and linear axis (Xi).
See Felippa AFM_Apendix I
index | Index (i, from 0 to numPoints()-1) of the point to be sampled |
naturalCoord | Returns the point where the function F should be sampled (Zeta1, Zeta2, Xi) |
weight | Returns the weight of this point w = w_zeta * w_xi) |
Implements GmIntegrationRule.
|
inlinevirtual |
Returns true if this is a valid integration rule object, false if not.
An invalid object can be created when an invalid rule number is passed to the constructor of an integration rule. The result of calling a method other then isValid() over an invalid integration rule is undefined.
Implements GmIntegrationRule.
|
inlinevirtual |
Returns the number of integration points associated with the given rule parameter. This is usefull for compound rules to access the number of points in each "dimension". For non compound rules, returns the number of points. Keep in mind that following the rulePar way, it returns -1 for a unknown parameter. In a 2 x 3 quad rule, for example, it should return 2 for parameter 0 and 3 for parameter 1. This method is necessary since for Newton Coates rules for tri/tet elements the number of points is NOT equal to the parameter and when combined in a Wedge this information is necessary.
Implements GmIntegrationRule.