1 |
|
2 |
/******************************************************* |
3 |
* |
4 |
* Copyright (c) 2003-2009 by University of Queensland |
5 |
* Earth Systems Science Computational Center (ESSCC) |
6 |
* http://www.uq.edu.au/esscc |
7 |
* |
8 |
* Primary Business: Queensland, Australia |
9 |
* Licensed under the Open Software License version 3.0 |
10 |
* http://www.opensource.org/licenses/osl-3.0.php |
11 |
* |
12 |
*******************************************************/ |
13 |
|
14 |
|
15 |
/**************************************************************/ |
16 |
|
17 |
/* Finley: Reference elements */ |
18 |
|
19 |
/**************************************************************/ |
20 |
|
21 |
#ifndef INC_FINLEY_REFERENCEELEMENTS |
22 |
#define INC_FINLEY_REFERENCEELEMENTS |
23 |
|
24 |
|
25 |
/**************************************************************/ |
26 |
|
27 |
#include "Finley.h" |
28 |
#include "ShapeFunctions.h" |
29 |
#include "Quadrature.h" |
30 |
|
31 |
/**************************************************************/ |
32 |
|
33 |
/* The ids of the allowed reference ellements: */ |
34 |
|
35 |
#define MAX_numNodes 64 |
36 |
#define MAX_numSubElements 8 |
37 |
#define MAX_numSides 2 |
38 |
|
39 |
typedef enum { |
40 |
Point1, |
41 |
Line2, |
42 |
Line3, |
43 |
Line4, |
44 |
Tri3, |
45 |
Tri6, |
46 |
Tri9, |
47 |
Tri10, |
48 |
Rec4, |
49 |
Rec8, |
50 |
Rec9, |
51 |
Rec12, |
52 |
Rec16, |
53 |
Tet4, |
54 |
Tet10, |
55 |
Tet16, |
56 |
Hex8, |
57 |
Hex20, |
58 |
Hex27, |
59 |
Hex32, |
60 |
Line2Face, |
61 |
Line3Face, |
62 |
Line4Face, |
63 |
Tri3Face, |
64 |
Tri6Face, |
65 |
Tri9Face, |
66 |
Tri10Face, |
67 |
Rec4Face, |
68 |
Rec8Face, |
69 |
Rec9Face, |
70 |
Rec12Face, |
71 |
Rec16Face, |
72 |
Tet4Face, |
73 |
Tet10Face, |
74 |
Tet16Face, |
75 |
Hex8Face, |
76 |
Hex20Face, |
77 |
Hex27Face, |
78 |
Hex32Face, |
79 |
Point1_Contact, |
80 |
Line2_Contact, |
81 |
Line3_Contact, |
82 |
Line4_Contact, |
83 |
Tri3_Contact, |
84 |
Tri6_Contact, |
85 |
Tri9_Contact, |
86 |
Tri10_Contact, |
87 |
Rec4_Contact, |
88 |
Rec8_Contact, |
89 |
Rec9_Contact, |
90 |
Rec12_Contact, |
91 |
Rec16_Contact, |
92 |
Line2Face_Contact, |
93 |
Line3Face_Contact, |
94 |
Line4Face_Contact, |
95 |
Tri3Face_Contact, |
96 |
Tri6Face_Contact, |
97 |
Tri9Face_Contact, |
98 |
Tri10Face_Contact, |
99 |
Rec4Face_Contact, |
100 |
Rec8Face_Contact, |
101 |
Rec9Face_Contact, |
102 |
Rec12Face_Contact, |
103 |
Rec16Face_Contact, |
104 |
Tet4Face_Contact, |
105 |
Tet10Face_Contact, |
106 |
Tet16Face_Contact, |
107 |
Hex8Face_Contact, |
108 |
Hex20Face_Contact, |
109 |
Hex27Face_Contact, |
110 |
Hex32Face_Contact, |
111 |
Line3Macro, |
112 |
Tri6Macro, |
113 |
Rec9Macro, |
114 |
Tet10Macro, |
115 |
Hex27Macro, |
116 |
NoRef /* marks end of list */ |
117 |
} ElementTypeId; |
118 |
|
119 |
/**************************************************************/ |
120 |
|
121 |
/* this struct holds the definition of the reference element: */ |
122 |
|
123 |
typedef struct Finley_ReferenceElementInfo { |
124 |
ElementTypeId TypeId; /* the id */ |
125 |
char* Name; /* the name in text form e.g. Line1,Rec12,... */ |
126 |
dim_t numNodes; /* number of nodes defining the element*/ |
127 |
dim_t numSubElements; /* number of subelements. >1 is macro elements are used. */ |
128 |
dim_t numSides; /* specifies the number of sides the element supports. This =2 if contatact elements are used |
129 |
otherwise =1. */ |
130 |
|
131 |
|
132 |
index_t offsets[MAX_numSides+1]; /* offset to the side nodes: offsets[s]...offset[s+1]-1 referes to the nodes to be used for side s*/ |
133 |
|
134 |
|
135 |
ElementTypeId LinearTypeId; /* id of the linear version of the element */ |
136 |
|
137 |
index_t linearNodes[MAX_numNodes*MAX_numSides]; /* gives the list of nodes defining the linear or macro element */ |
138 |
|
139 |
Finley_QuadTypeId Quadrature; /* quadrature scheme */ |
140 |
Finley_ShapeFunctionTypeId Parametrization; /* shape function for parametrization of the element */ |
141 |
Finley_ShapeFunctionTypeId BasisFunctions; /* shape function for the basis functions */ |
142 |
|
143 |
index_t subElementNodes[MAX_numNodes*MAX_numSides*MAX_numSubElements]; /* gives the list of nodes defining the subelements: |
144 |
subElementNodes[INDEX2(i,s,BasisFunctions->numShape*numSides)] is the i-th node in the s-th subelement.*/ |
145 |
/*********************************************************************************************************************************** */ |
146 |
dim_t numRelevantGeoNodes; /* number of nodes used to describe the geometry of the geometrically relevant part of the element |
147 |
typically this is numNodes but for 'Face' elements where the quadrature points are defined on face of the element |
148 |
this is the number of nodes on the particular face. */ |
149 |
index_t relevantGeoNodes[MAX_numNodes]; /* list to gather the geometrically relevant nodes (length used is numRelevantGeoNodes) |
150 |
this list is used for VTK interface */ |
151 |
|
152 |
dim_t numNodesOnFace; /* if the element is allowed as a face element, numNodesOnFace defines the number of nodes defining the face */ |
153 |
/* the following lists are only used for face elements defined by numNodesOnFace>0 */ |
154 |
index_t faceNodes[MAX_numNodes]; /* list of the nodes defining the face */ |
155 |
index_t shiftNodes[MAX_numNodes]; /* defines a permutation of the nodes which rotates the nodes on the face */ |
156 |
index_t reverseNodes[MAX_numNodes]; /* reverses the order of the nodes on a face. the permutation has keep 0 fixed. */ |
157 |
/* shiftNodes={-1} or reverseNodes={-1} are ignored. */ |
158 |
} Finley_ReferenceElementInfo; |
159 |
|
160 |
|
161 |
/**************************************************************/ |
162 |
|
163 |
/* this struct holds the realization of a reference element */ |
164 |
|
165 |
typedef struct Finley_ReferenceElement { |
166 |
Finley_ReferenceElementInfo* Type; /* type of the reference element */ |
167 |
Finley_ReferenceElementInfo* LinearType; /* type of the linear reference element */ |
168 |
index_t reference_counter; /* reference counter */ |
169 |
dim_t numNodes; |
170 |
dim_t numLocalDim; |
171 |
dim_t numLinearNodes; |
172 |
Finley_ShapeFunction* Parametrization; |
173 |
Finley_ShapeFunction* BasisFunctions; |
174 |
Finley_ShapeFunction* LinearBasisFunctions; |
175 |
double* DBasisFunctionDv; /* pointer to derivatives to basis function corresponding to the Parametrization quad points */ |
176 |
bool_t DBasisFunctionDvShared; /* TRUE to indicate that DBasisFunctionDv is shared with another object which is managing it */ |
177 |
|
178 |
} Finley_ReferenceElement; |
179 |
|
180 |
/**************************************************************/ |
181 |
|
182 |
/* interfaces: */ |
183 |
|
184 |
Finley_ReferenceElement* Finley_ReferenceElement_alloc(ElementTypeId,int); |
185 |
void Finley_ReferenceElement_dealloc(Finley_ReferenceElement*); |
186 |
ElementTypeId Finley_ReferenceElement_getTypeId(char*); |
187 |
Finley_ReferenceElement* Finley_ReferenceElement_reference(Finley_ReferenceElement* in); |
188 |
Finley_ReferenceElementInfo* Finley_ReferenceElement_getInfo(ElementTypeId id); |
189 |
|
190 |
|
191 |
#define Finley_ReferenceElement_getNumNodes(__in__) (__in__)->Type->numNodes |
192 |
|
193 |
#endif /* #ifndef INC_FINLEY_REFERENCEELEMENTS */ |