shape_function.F90 coverage: 100.00 %func 38.01 %block
1) module Shape_Function_module
2)
3) use Gauss_module
4) use Grid_Unstructured_Cell_module
5) use PFLOTRAN_Constants_module
6)
7) implicit none
8)
9) #include "petsc/finclude/petscsys.h"
10)
11) private
12)
13) type, public :: shapefunction_type
14) PetscInt :: EleType ! element type
15) PetscReal, pointer :: zeta(:) ! coordinates of point in reference element
16) PetscReal, pointer :: N(:) ! shape function for all nodes evaluated at zeta (N is a row vector)
17) PetscReal, pointer :: DN(:,:) ! derivatives of shape function with respect to zeta
18) PetscReal, pointer :: coord(:,:) ! local coordinates of the vertices in the reference element
19) end type shapefunction_type
20)
21) public :: ShapeFunctionInitialize, ShapeFunctionCalculate, &
22) ShapeFunctionDestroy
23)
24) contains
25)
26) ! ************************************************************************** !
27)
28) subroutine ShapeFunctionInitialize(shapefunction)
29) !
30) ! Allocate memory for shapefunction type
31) !
32) ! Author: Satish Karra, LANL
33) ! Date: 5/17/2013
34) !
35)
36) type(shapefunction_type) :: shapefunction
37) PetscReal, pointer :: coord(:,:)
38)
39) select case(shapefunction%EleType)
40) case(LINE_TYPE)
41) allocate(shapefunction%N(TWO_INTEGER))
42) allocate(shapefunction%DN(TWO_INTEGER,ONE_INTEGER))
43) allocate(shapefunction%zeta(ONE_INTEGER))
44) allocate(shapefunction%coord(TWO_INTEGER,ONE_INTEGER))
45) case(QUAD_TYPE)
46) allocate(shapefunction%N(FOUR_INTEGER))
47) allocate(shapefunction%DN(FOUR_INTEGER,TWO_INTEGER))
48) allocate(shapefunction%zeta(TWO_INTEGER))
49) allocate(shapefunction%coord(FOUR_INTEGER,TWO_INTEGER))
50) case(WEDGE_TYPE)
51) allocate(shapefunction%N(SIX_INTEGER))
52) allocate(shapefunction%DN(SIX_INTEGER,THREE_INTEGER))
53) allocate(shapefunction%zeta(THREE_INTEGER))
54) allocate(shapefunction%coord(SIX_INTEGER,THREE_INTEGER))
55) case(TET_TYPE)
56) allocate(shapefunction%N(FOUR_INTEGER))
57) allocate(shapefunction%DN(FOUR_INTEGER,THREE_INTEGER))
58) allocate(shapefunction%zeta(THREE_INTEGER))
59) allocate(shapefunction%coord(FOUR_INTEGER,THREE_INTEGER))
60) case(PYR_TYPE)
61) allocate(shapefunction%N(FIVE_INTEGER))
62) allocate(shapefunction%DN(FIVE_INTEGER,THREE_INTEGER))
63) allocate(shapefunction%zeta(THREE_INTEGER))
64) allocate(shapefunction%coord(FIVE_INTEGER,THREE_INTEGER))
65) case(HEX_TYPE)
66) allocate(shapefunction%N(EIGHT_INTEGER))
67) allocate(shapefunction%DN(EIGHT_INTEGER,THREE_INTEGER))
68) allocate(shapefunction%zeta(THREE_INTEGER))
69) allocate(shapefunction%coord(EIGHT_INTEGER,THREE_INTEGER))
70) case default
71) print *, 'Error: Invalid EleType. Enter an EleType from L2,' // &
72) ' T3, Q4, B8 only.'
73) end select
74)
75) shapefunction%zeta = 0.d0
76) shapefunction%N = 0.d0
77) shapefunction%DN = 0.d0
78) shapefunction%coord = 0.d0
79)
80) coord => shapefunction%coord
81)
82) select case(shapefunction%EleType)
83) case(LINE_TYPE)
84) coord(1,1) = -1.d0
85) coord(2,1) = 1.d0
86) case(QUAD_TYPE)
87) coord(1,1) = -1.d0
88) coord(1,2) = -1.d0
89) coord(2,1) = 1.d0
90) coord(2,2) = -1.d0
91) coord(3,1) = 1.d0
92) coord(3,2) = 1.d0
93) coord(4,1) = -1.d0
94) coord(4,2) = 1.d0
95) case(WEDGE_TYPE)
96) coord(1,:) = -(/0.d0,1.d0,1.d0/)
97) coord(2,:) = -(/-1.d0,-1.d0,1.d0/)
98) coord(3,:) = -(/1.d0,-1.d0,1.d0/)
99) coord(4,:) = -(/0.d0,1.d0,-1.d0/)
100) coord(5,:) = -(/-1.d0,-1.d0,-1.d0/)
101) coord(6,:) = -(/1.d0,-1.d0,-1.d0/)
102) case(TET_TYPE)
103) coord(1,:) = -(/0.d0,1.d0,1.d0/)
104) coord(2,:) = -(/-1.d0,-1.d0,1.d0/)
105) coord(3,:) = -(/1.d0,-1.d0,1.d0/)
106) coord(4,:) = -(/0.d0,0.d0,-1.d0/)
107) case(PYR_TYPE)
108) coord(1,:) = -(/1.d0,1.d0,1.d0/)
109) coord(2,:) = -(/-1.d0,1.d0,1.d0/)
110) coord(3,:) = -(/-1.d0,-1.d0,1.d0/)
111) coord(4,:) = -(/1.d0,-1.d0,1.d0/)
112) coord(5,:) = -(/0.d0,0.d0,-1.d0/)
113) case(HEX_TYPE)
114) coord(1,:) = -(/1.d0,1.d0,1.d0/)
115) coord(2,:) = -(/-1.d0,1.d0,1.d0/)
116) coord(3,:) = -(/-1.d0,-1.d0,1.d0/)
117) coord(4,:) = -(/1.d0,-1.d0,1.d0/)
118) coord(5,:) = -(/1.d0,1.d0,-1.d0/)
119) coord(6,:) = -(/-1.d0,1.d0,-1.d0/)
120) coord(7,:) = -(/-1.d0,-1.d0,-1.d0/)
121) coord(8,:) = -(/1.d0,-1.d0,-1.d0/)
122) case default
123) print *, 'Error: Invalid EleType. Enter an EleType from L2,' // &
124) ' T3, Q4, B8 only.'
125) end select
126)
127) end subroutine ShapeFunctionInitialize
128)
129) ! ************************************************************************** !
130)
131) subroutine ShapeFunctionCalculate(shapefunction)
132) !
133) ! Subroutine provides shape functions and its derivatives
134) ! at a given spatial point (in the reference element) 'zeta' for various finite
135) ! elements.
136) ! Input variables
137) ! EleType: element type
138) ! L2: one-dimensional element -1 <= x <= +1
139) ! Q4: four node quadrilateral element -1 <= x, y <= +1
140) ! Q9: nine node quadrilateral element
141) ! T3: three node right-angled triangle with h = b = 1
142) ! T6: six node right-angled triangle
143) ! B8: eight node right-angled tetrahedron element
144) ! zeta: coordinates of a point in the reference element
145) ! Output variables
146) ! N: shape functions for all nodes evaluated at zeta (N is a row
147) ! vector!)
148) ! DN: derivatives of shape functions with respect to zeta
149) !
150) ! Author: Satish Karra, LANL
151) ! Date: 5/17/2013
152) !
153)
154) type(shapefunction_type) :: shapefunction
155) PetscReal, pointer :: zeta(:)
156) PetscReal, pointer :: N(:)
157) PetscReal, pointer :: DN(:,:)
158) PetscInt :: i
159)
160) N => shapefunction%N
161) DN => shapefunction%DN
162) zeta => shapefunction%zeta
163)
164) do i = 1, size(zeta)
165) if (zeta(i) < -1.d0 .or. zeta(i) > 1.d0) then
166) print *, 'Error: Enter values between -1 and 1 for zeta'
167) stop
168) endif
169) enddo
170)
171) select case(shapefunction%EleType)
172) case(LINE_TYPE)
173) N(1) = 0.5d0*(1.d0 - zeta(1))
174) N(2) = 0.5d0*(1.d0 + zeta(1))
175) DN(1,1) = -0.5d0
176) DN(2,1) = 0.5d0
177) case(QUAD_TYPE)
178) N(1) = 1.d0/4.d0*(1.d0 - zeta(1))*(1.d0 - zeta(2))
179) N(2) = 1.d0/4.d0*(1.d0 + zeta(1))*(1.d0 - zeta(2))
180) N(3) = 1.d0/4.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))
181) N(4) = 1.d0/4.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))
182) DN(1,1) = -1.d0/4.d0*(1.d0 - zeta(2))
183) DN(1,2) = -1.d0/4.d0*(1.d0 - zeta(1))
184) DN(2,1) = 1.d0/4.d0*(1.d0 - zeta(2))
185) DN(2,2) = -1.d0/4.d0*(1.d0 + zeta(1))
186) DN(3,1) = 1.d0/4.d0*(1.d0 + zeta(2))
187) DN(3,2) = 1.d0/4.d0*(1.d0 + zeta(1))
188) DN(4,1) = -1.d0/4.d0*(1.d0 + zeta(2))
189) DN(4,2) = 1.d0/4.d0*(1.d0 - zeta(1))
190) case(WEDGE_TYPE)
191) N(1) = 1.d0/4.d0*(1.d0 - zeta(2))*(1.d0 - zeta(3))
192) N(2) = 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(1.d0 - zeta(3))
193) N(3) = 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(1.d0 - zeta(3))
194) N(4) = 1.d0/4.d0*(1.d0 - zeta(2))*(1.d0 + zeta(3))
195) N(5) = 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(1.d0 + zeta(3))
196) N(6) = 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(1.d0 + zeta(3))
197)
198) DN(1,:) = (/0.d0, &
199) 1.d0/4.d0*(-1.d0)*(1.d0 - zeta(3)), &
200) 1.d0/4.d0*(1.d0 - zeta(2))*(-1.d0)/)
201) DN(2,:) = (/1.d0/8.d0*(+1.d0)*(1.d0 + zeta(2))*(1.d0 - zeta(3)), &
202) 1.d0/8.d0*(1.d0 + zeta(1))*(+1.d0)*(1.d0 - zeta(3)), &
203) 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(-1.d0)/)
204) DN(3,:) = (/1.d0/8.d0*(-1.d0)*(1.d0 + zeta(2))*(1.d0 - zeta(3)), &
205) 1.d0/8.d0*(1.d0 - zeta(1))*(+1.d0)*(1.d0 - zeta(3)), &
206) 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(-1.d0)/)
207) DN(4,:) = (/0.d0, &
208) 1.d0/4.d0*(-1.d0)*(1.d0 + zeta(3)), &
209) 1.d0/4.d0*(1.d0 - zeta(2))*(+1.d0)/)
210) DN(5,:) = (/1.d0/8.d0*(+1.d0)*(1.d0 + zeta(2))*(1.d0 + zeta(3)), &
211) 1.d0/8.d0*(1.d0 + zeta(1))*(+1.d0)*(1.d0 + zeta(3)), &
212) 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(+1.d0)/)
213) DN(6,:) = (/1.d0/8.d0*(-1.d0)*(1.d0 + zeta(2))*(1.d0 + zeta(3)), &
214) 1.d0/8.d0*(1.d0 - zeta(1))*(+1.d0)*(1.d0 + zeta(3)), &
215) 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(+1.d0)/)
216)
217) case(TET_TYPE)
218) N(1) = 1.d0/4.d0*(1.d0 - zeta(2))*(1.d0 - zeta(3))
219) N(2) = 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(1.d0 - zeta(3))
220) N(3) = 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(1.d0 - zeta(3))
221) N(4) = 1.d0/2.d0*(1.d0 + zeta(3))
222)
223) DN(1,:) = (/0.d0, &
224) 1.d0/4.d0*(-1.d0)*(1.d0 - zeta(3)), &
225) 1.d0/4.d0*(1.d0 - zeta(2))*(-1.d0)/)
226) DN(2,:) = (/1.d0/8.d0*(+1.d0)*(1.d0 + zeta(2))*(1.d0 - zeta(3)), &
227) 1.d0/8.d0*(1.d0 + zeta(1))*(+1.d0)*(1.d0 - zeta(3)), &
228) 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(-1.d0)/)
229) DN(3,:) = (/1.d0/8.d0*(-1.d0)*(1.d0 + zeta(2))*(1.d0 - zeta(3)), &
230) 1.d0/8.d0*(1.d0 - zeta(1))*(+1.d0)*(1.d0 - zeta(3)), &
231) 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(-1.d0)/)
232) DN(4,:) = (/0.d0,0.d0,1.d0/2.d0/)
233)
234) case(PYR_TYPE)
235) N(1) = 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 - zeta(2))*(1.d0 - zeta(3))
236) N(2) = 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 - zeta(2))*(1.d0 - zeta(3))
237) N(3) = 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(1.d0 - zeta(3))
238) N(4) = 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(1.d0 - zeta(3))
239) N(5) = 1.d0/2.d0*(1.d0 + zeta(3))
240)
241) DN(1,:) = (/1.d0/8.d0*(-1.d0)*(1.d0 - zeta(2))*(1.d0 - zeta(3)), &
242) 1.d0/8.d0*(1.d0 - zeta(1))*(-1.d0)*(1.d0 - zeta(3)), &
243) 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 - zeta(2))*(-1.d0)/)
244) DN(2,:) = (/1.d0/8.d0*(+1.d0)*(1.d0 - zeta(2))*(1.d0 - zeta(3)), &
245) 1.d0/8.d0*(1.d0 + zeta(1))*(-1.d0)*(1.d0 - zeta(3)), &
246) 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 - zeta(2))*(-1.d0)/)
247) DN(3,:) = (/1.d0/8.d0*(+1.d0)*(1.d0 + zeta(2))*(1.d0 - zeta(3)), &
248) 1.d0/8.d0*(1.d0 + zeta(1))*(+1.d0)*(1.d0 - zeta(3)), &
249) 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(-1.d0)/)
250) DN(4,:) = (/1.d0/8.d0*(-1.d0)*(1.d0 + zeta(2))*(1.d0 - zeta(3)), &
251) 1.d0/8.d0*(1.d0 - zeta(1))*(+1.d0)*(1.d0 - zeta(3)), &
252) 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(-1.d0)/)
253) DN(5,:) = (/0.d0,0.d0,1.d0/2.d0/)
254)
255) case(HEX_TYPE)
256) N(1) = 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 - zeta(2))*(1.d0 - zeta(3))
257) N(2) = 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 - zeta(2))*(1.d0 - zeta(3))
258) N(3) = 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(1.d0 - zeta(3))
259) N(4) = 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(1.d0 - zeta(3))
260) N(5) = 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 - zeta(2))*(1.d0 + zeta(3))
261) N(6) = 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 - zeta(2))*(1.d0 + zeta(3))
262) N(7) = 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(1.d0 + zeta(3))
263) N(8) = 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(1.d0 + zeta(3))
264)
265) DN(1,:) = (/1.d0/8.d0*(-1.d0)*(1.d0 - zeta(2))*(1.d0 - zeta(3)), &
266) 1.d0/8.d0*(1.d0 - zeta(1))*(-1.d0)*(1.d0 - zeta(3)), &
267) 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 - zeta(2))*(-1.d0)/)
268) DN(2,:) = (/1.d0/8.d0*(+1.d0)*(1.d0 - zeta(2))*(1.d0 - zeta(3)), &
269) 1.d0/8.d0*(1.d0 + zeta(1))*(-1.d0)*(1.d0 - zeta(3)), &
270) 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 - zeta(2))*(-1.d0)/)
271) DN(3,:) = (/1.d0/8.d0*(+1.d0)*(1.d0 + zeta(2))*(1.d0 - zeta(3)), &
272) 1.d0/8.d0*(1.d0 + zeta(1))*(+1.d0)*(1.d0 - zeta(3)), &
273) 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(-1.d0)/)
274) DN(4,:) = (/1.d0/8.d0*(-1.d0)*(1.d0 + zeta(2))*(1.d0 - zeta(3)), &
275) 1.d0/8.d0*(1.d0 - zeta(1))*(+1.d0)*(1.d0 - zeta(3)), &
276) 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(-1.d0)/)
277) DN(5,:) = (/1.d0/8.d0*(-1.d0)*(1.d0 - zeta(2))*(1.d0 + zeta(3)), &
278) 1.d0/8.d0*(1.d0 - zeta(1))*(-1.d0)*(1.d0 + zeta(3)), &
279) 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 - zeta(2))*(+1.d0)/)
280) DN(6,:) = (/1.d0/8.d0*(+1.d0)*(1.d0 - zeta(2))*(1.d0 + zeta(3)), &
281) 1.d0/8.d0*(1.d0 + zeta(1))*(-1.d0)*(1.d0 + zeta(3)), &
282) 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 - zeta(2))*(+1.d0)/)
283) DN(7,:) = (/1.d0/8.d0*(+1.d0)*(1.d0 + zeta(2))*(1.d0 + zeta(3)), &
284) 1.d0/8.d0*(1.d0 + zeta(1))*(+1.d0)*(1.d0 + zeta(3)), &
285) 1.d0/8.d0*(1.d0 + zeta(1))*(1.d0 + zeta(2))*(+1.d0)/)
286) DN(8,:) = (/1.d0/8.d0*(-1.d0)*(1.d0 + zeta(2))*(1.d0 + zeta(3)), &
287) 1.d0/8.d0*(1.d0 - zeta(1))*(+1.d0)*(1.d0 + zeta(3)), &
288) 1.d0/8.d0*(1.d0 - zeta(1))*(1.d0 + zeta(2))*(+1.d0)/)
289) case default
290) print *, 'Error: Invalid EleType. Enter an EleType from L2,' // &
291) ' T3, Q4, B8 only.'
292) end select
293)
294) end subroutine ShapeFunctionCalculate
295)
296) ! ************************************************************************** !
297)
298) subroutine ShapeFunctionDestroy(shapefunction)
299) !
300) ! Dellocate memory for shapefunction type
301) !
302) ! Author: Satish Karra, LANL
303) ! Date: 5/17/2013
304) !
305)
306) type(shapefunction_type) :: shapefunction
307)
308) deallocate(shapefunction%N)
309) nullify(shapefunction%N)
310) deallocate(shapefunction%DN)
311) nullify(shapefunction%DN)
312) deallocate(shapefunction%zeta)
313) nullify(shapefunction%zeta)
314) deallocate(shapefunction%coord)
315) nullify(shapefunction%coord)
316)
317) end subroutine ShapeFunctionDestroy
318)
319) end module Shape_Function_module
320)