DAMASK with grid solvers  Revision: v2.0.3-2204-gdb1f2151
The Düsseldorf Advanced Material Simulation Kit with Grid Solvers
element.f90
Go to the documentation of this file.
1 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/element.f90"
2 # 1 "<built-in>"
3 # 1 "<command-line>"
4 # 1 "/home/damask_user/GitLabCI_Pipeline_4301/DAMASK/src/element.f90"
5 !--------------------------------------------------------------------------------------------------
8 !--------------------------------------------------------------------------------------------------
9 module element
10  use io
11 
12  implicit none
13  private
14 
15 !---------------------------------------------------------------------------------------------------
17 !---------------------------------------------------------------------------------------------------
18  type, public :: telement
19  integer :: &
20  elemtype, &
21  geomtype, & !< geometry type (same for same dimension and same number of integration points)
22  celltype, &
23  nnodes, &
24  ncellnodes, &
25  ncellnodespercell, &
26  nips, &
27  nipneighbors
28  character(len=:), allocatable :: &
29  vtktype
30  integer, dimension(:,:), allocatable :: &
31  cell, & !< intra-element (cell) nodes that constitute a cell
32  ipneighbor, &
33  cellface
34  integer, dimension(:,:), allocatable :: &
35  ! center of gravity of the weighted nodes gives the position of the cell node.
36  ! example: face-centered cell node with face nodes 1,2,5,6 to be used in,
37  ! e.g., an 8 node element, would be encoded: 1, 1, 0, 0, 1, 1, 0, 0
38  cellnodeparentnodeweights
39  contains
40  procedure :: init => telement_init
41  end type telement
42 
43 
44  integer, parameter :: &
45  nelemtype = 13
46 
47  integer, dimension(NELEMTYPE), parameter :: nnode = &
48  [ &
49  3, & ! 2D, 1 IP
50  6, & ! 2D, 3 IP
51  4, & ! 2D, 4 IP
52  8, & ! 2D, 9 IP
53  8, & ! 2D, 4 IP
54  !----------------------
55  4, & ! 3D, 1 IP
56  5, & ! 3D, 4 IP
57  10, & ! 3D, 4 IP
58  6, & ! 3D, 6 IP
59  8, & ! 3D, 1 IP
60  8, & ! 3D, 8 IP
61  20, & ! 3D, 8 IP
62  20 & ! 3D, 27 IP
63  ]
64 
65  integer, dimension(NELEMTYPE), parameter :: geomtype = &
66  [ &
67  1, & ! 1 triangle
68  2, & ! 3 quadrilaterals
69  3, & ! 4 quadrilaterals
70  4, & ! 9 quadrilaterals
71  3, & ! 4 quadrilaterals
72  !----------------------
73  5, & ! 1 tetrahedron
74  6, & ! 4 hexahedrons
75  6, & ! 4 hexahedrons
76  7, & ! 6 hexahedrons
77  8, & ! 1 hexahedron
78  9, & ! 8 hexahedrons
79  9, & ! 8 hexahedrons
80  10 & ! 27 hexahedrons
81  ]
82 
83  integer, dimension(maxval(GEOMTYPE)), parameter :: ncellnode = &
84  [ &
85  3, &
86  7, &
87  9, &
88  16, &
89  4, &
90  15, &
91  21, &
92  8, &
93  27, &
94  64 &
95  ]
96 
97  integer, dimension(maxval(GEOMTYPE)), parameter :: nip = &
98  [ &
99  1, &
100  3, &
101  4, &
102  9, &
103  1, &
104  4, &
105  6, &
106  1, &
107  8, &
108  27 &
109  ]
110 
111  integer, dimension(maxval(GEOMTYPE)), parameter :: celltype = &
112  [ &
113  1, & ! 2D, 3 node (Triangle)
114  2, & ! 2D, 4 node (Quadrilateral)
115  2, & ! - " -
116  2, & ! - " -
117  3, & ! 3D, 4 node (Tetrahedron)
118  4, & ! 3D, 4 node (Hexahedron)
119  4, & ! - " -
120  4, & ! - " -
121  4, & ! - " -
122  4 & ! - " -
123  ]
124 
125  integer, dimension(maxval(CELLTYPE)), parameter :: nipneighbor = &
126  [ &
127  3, &
128  4, &
129  4, &
130  6 &
131  ]
132 
133  integer, dimension(maxval(CELLTYPE)), parameter :: ncellnodepercellface = &
134  [ &
135  2, &
136  2, &
137  3, &
138  4 &
139  ]
140 
141  integer, dimension(maxval(CELLTYPE)), parameter :: ncellnodepercell = &
142  [ &
143  3, &
144  4, &
145  4, &
146  8 &
147  ]
148 
149  ! *** IPneighbor ***
150  ! list of the neighborhood of each IP.
151  ! It is sorted in (local) +x,-x, +y,-y, +z,-z direction.
152  ! Positive integers denote an intra-element IP identifier.
153  ! Negative integers denote the interface behind which the neighboring (extra-element) IP will be located.
154 
155  integer, dimension(NIPNEIGHBOR(CELLTYPE(1)),NIP(1)), parameter :: ipneighbor1 = &
156  reshape([&
157  -2,-3,-1 &
158 
159 
160 
161  ],[nipneighbor(celltype(1)),nip(1)])
162 
163 
164  integer, dimension(NIPNEIGHBOR(CELLTYPE(2)),NIP(2)), parameter :: ipneighbor2 = &
165  reshape([&
166  2,-3, 3,-1, &
167  -2, 1, 3,-1, &
168  2,-3,-2, 1 &
169 
170 
171 
172  ],[nipneighbor(celltype(2)),nip(2)])
173 
174 
175  integer, dimension(NIPNEIGHBOR(CELLTYPE(3)),NIP(3)), parameter :: ipneighbor3 = &
176  reshape([&
177  2,-4, 3,-1, &
178  -2, 1, 4,-1, &
179  4,-4,-3, 1, &
180  -2, 3,-3, 2 &
181 
182 
183 
184  ],[nipneighbor(celltype(3)),nip(3)])
185 
186 
187  integer, dimension(NIPNEIGHBOR(CELLTYPE(4)),NIP(4)), parameter :: ipneighbor4 = &
188  reshape([&
189  2,-4, 4,-1, &
190  3, 1, 5,-1, &
191  -2, 2, 6,-1, &
192  5,-4, 7, 1, &
193  6, 4, 8, 2, &
194  -2, 5, 9, 3, &
195  8,-4,-3, 4, &
196  9, 7,-3, 5, &
197  -2, 8,-3, 6 &
198 
199 
200 
201  ],[nipneighbor(celltype(4)),nip(4)])
202 
203 
204  integer, dimension(NIPNEIGHBOR(CELLTYPE(5)),NIP(5)), parameter :: ipneighbor5 = &
205  reshape([&
206  -1,-2,-3,-4 &
207 
208 
209 
210  ],[nipneighbor(celltype(5)),nip(5)])
211 
212 
213  integer, dimension(NIPNEIGHBOR(CELLTYPE(6)),NIP(6)), parameter :: ipneighbor6 = &
214  reshape([&
215  2,-4, 3,-2, 4,-1, &
216  -2, 1, 3,-2, 4,-1, &
217  2,-4,-3, 1, 4,-1, &
218  2,-4, 3,-2,-3, 1 &
219 
220 
221 
222  ],[nipneighbor(celltype(6)),nip(6)])
223 
224 
225  integer, dimension(NIPNEIGHBOR(CELLTYPE(7)),NIP(7)), parameter :: ipneighbor7 = &
226  reshape([&
227  2,-4, 3,-2, 4,-1, &
228  -3, 1, 3,-2, 5,-1, &
229  2,-4,-3, 1, 6,-1, &
230  5,-4, 6,-2,-5, 1, &
231  -3, 4, 6,-2,-5, 2, &
232  5,-4,-3, 4,-5, 3 &
233 
234 
235 
236  ],[nipneighbor(celltype(7)),nip(7)])
237 
238 
239  integer, dimension(NIPNEIGHBOR(CELLTYPE(8)),NIP(8)), parameter :: ipneighbor8 = &
240  reshape([&
241  -3,-5,-4,-2,-6,-1 &
242 
243 
244 
245  ],[nipneighbor(celltype(8)),nip(8)])
246 
247 
248  integer, dimension(NIPNEIGHBOR(CELLTYPE(9)),NIP(9)), parameter :: ipneighbor9 = &
249  reshape([&
250  2,-5, 3,-2, 5,-1, &
251  -3, 1, 4,-2, 6,-1, &
252  4,-5,-4, 1, 7,-1, &
253  -3, 3,-4, 2, 8,-1, &
254  6,-5, 7,-2,-6, 1, &
255  -3, 5, 8,-2,-6, 2, &
256  8,-5,-4, 5,-6, 3, &
257  -3, 7,-4, 6,-6, 4 &
258 
259 
260 
261  ],[nipneighbor(celltype(9)),nip(9)])
262 
263 
264  integer, dimension(NIPNEIGHBOR(CELLTYPE(10)),NIP(10)), parameter :: ipneighbor10 = &
265  reshape([&
266  2,-5, 4,-2,10,-1, &
267  3, 1, 5,-2,11,-1, &
268  -3, 2, 6,-2,12,-1, &
269  5,-5, 7, 1,13,-1, &
270  6, 4, 8, 2,14,-1, &
271  -3, 5, 9, 3,15,-1, &
272  8,-5,-4, 4,16,-1, &
273  9, 7,-4, 5,17,-1, &
274  -3, 8,-4, 6,18,-1, &
275  11,-5,13,-2,19, 1, &
276  12,10,14,-2,20, 2, &
277  -3,11,15,-2,21, 3, &
278  14,-5,16,10,22, 4, &
279  15,13,17,11,23, 5, &
280  -3,14,18,12,24, 6, &
281  17,-5,-4,13,25, 7, &
282  18,16,-4,14,26, 8, &
283  -3,17,-4,15,27, 9, &
284  20,-5,22,-2,-6,10, &
285  21,19,23,-2,-6,11, &
286  -3,20,24,-2,-6,12, &
287  23,-5,25,19,-6,13, &
288  24,22,26,20,-6,14, &
289  -3,23,27,21,-6,15, &
290  26,-5,-4,22,-6,16, &
291  27,25,-4,23,-6,17, &
292  -3,26,-4,24,-6,18 &
293 
294 
295 
296  ],[nipneighbor(celltype(10)),nip(10)])
297 
298 
299 
300  integer, dimension(NNODE(1),NCELLNODE(GEOMTYPE(1))), parameter :: cellnodeparentnodeweights1 = &
301  reshape([&
302  1, 0, 0, &
303  0, 1, 0, &
304  0, 0, 1 &
305 
306 
307 
308  ],[nnode(1),ncellnode(geomtype(1))])
309 
310 
311  integer, dimension(NNODE(2),NCELLNODE(GEOMTYPE(2))), parameter :: cellnodeparentnodeweights2 = &
312  reshape([&
313  1, 0, 0, 0, 0, 0, &
314  0, 1, 0, 0, 0, 0, &
315  0, 0, 1, 0, 0, 0, &
316  0, 0, 0, 1, 0, 0, &
317  0, 0, 0, 0, 1, 0, &
318  0, 0, 0, 0, 0, 1, &
319  1, 1, 1, 2, 2, 2 &
320 
321 
322 
323  ],[nnode(2),ncellnode(geomtype(2))])
324 
325 
326  integer, dimension(NNODE(3),NCELLNODE(GEOMTYPE(3))), parameter :: cellnodeparentnodeweights3 = &
327  reshape([&
328  1, 0, 0, 0, &
329  0, 1, 0, 0, &
330  0, 0, 1, 0, &
331  0, 0, 0, 1, &
332  1, 1, 0, 0, &
333  0, 1, 1, 0, &
334  0, 0, 1, 1, &
335  1, 0, 0, 1, &
336  1, 1, 1, 1 &
337 
338 
339 
340  ],[nnode(3),ncellnode(geomtype(3))])
341 
342 
343  integer, dimension(NNODE(4),NCELLNODE(GEOMTYPE(4))), parameter :: cellnodeparentnodeweights4 = &
344  reshape([&
345  1, 0, 0, 0, 0, 0, 0, 0, &
346  0, 1, 0, 0, 0, 0, 0, 0, &
347  0, 0, 1, 0, 0, 0, 0, 0, &
348  0, 0, 0, 1, 0, 0, 0, 0, &
349  1, 0, 0, 0, 2, 0, 0, 0, &
350  0, 1, 0, 0, 2, 0, 0, 0, &
351  0, 1, 0, 0, 0, 2, 0, 0, &
352  0, 0, 1, 0, 0, 2, 0, 0, &
353  0, 0, 1, 0, 0, 0, 2, 0, &
354  0, 0, 0, 1, 0, 0, 2, 0, &
355  0, 0, 0, 1, 0, 0, 0, 2, &
356  1, 0, 0, 0, 0, 0, 0, 2, &
357  4, 1, 1, 1, 8, 2, 2, 8, &
358  1, 4, 1, 1, 8, 8, 2, 2, &
359  1, 1, 4, 1, 2, 8, 8, 2, &
360  1, 1, 1, 4, 2, 2, 8, 8 &
361 
362 
363 
364  ],[nnode(4),ncellnode(geomtype(4))])
365 
366 
367  integer, dimension(NNODE(5),NCELLNODE(GEOMTYPE(5))), parameter :: cellnodeparentnodeweights5 = &
368  reshape([&
369  1, 0, 0, 0, 0, 0, 0, 0, &
370  0, 1, 0, 0, 0, 0, 0, 0, &
371  0, 0, 1, 0, 0, 0, 0, 0, &
372  0, 0, 0, 1, 0, 0, 0, 0, &
373  0, 0, 0, 0, 1, 0, 0, 0, &
374  0, 0, 0, 0, 0, 1, 0, 0, &
375  0, 0, 0, 0, 0, 0, 1, 0, &
376  0, 0, 0, 0, 0, 0, 0, 1, &
377  1, 1, 1, 1, 2, 2, 2, 2 &
378 
379 
380 
381  ],[nnode(5),ncellnode(geomtype(5))])
382 
383 
384  integer, dimension(NNODE(6),NcellNode(GEOMTYPE(6))), parameter :: cellnodeparentnodeweights6 = &
385  reshape([&
386  1, 0, 0, 0, &
387  0, 1, 0, 0, &
388  0, 0, 1, 0, &
389  0, 0, 0, 1 &
390 
391 
392 
393  ],[nnode(6),ncellnode(geomtype(6))])
394 
395 
396  integer, dimension(NNODE(7),NCELLNODE(GEOMTYPE(7))), parameter :: cellnodeparentnodeweights7 = &
397  reshape([&
398  1, 0, 0, 0, 0, &
399  0, 1, 0, 0, 0, &
400  0, 0, 1, 0, 0, &
401  0, 0, 0, 1, 0, &
402  1, 1, 0, 0, 0, &
403  0, 1, 1, 0, 0, &
404  1, 0, 1, 0, 0, &
405  1, 0, 0, 1, 0, &
406  0, 1, 0, 1, 0, &
407  0, 0, 1, 1, 0, &
408  1, 1, 1, 0, 0, &
409  1, 1, 0, 1, 0, &
410  0, 1, 1, 1, 0, &
411  1, 0, 1, 1, 0, &
412  0, 0, 0, 0, 1 &
413 
414 
415 
416  ],[nnode(7),ncellnode(geomtype(7))])
417 
418 
419  integer, dimension(NNODE(8),NCELLNODE(GEOMTYPE(8))), parameter :: cellnodeparentnodeweights8 = &
420  reshape([&
421  1, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
422  0, 1, 0, 0, 0, 0, 0, 0, 0, 0, &
423  0, 0, 1, 0, 0, 0, 0, 0, 0, 0, &
424  0, 0, 0, 1, 0, 0, 0, 0, 0, 0, &
425  0, 0, 0, 0, 1, 0, 0, 0, 0, 0, &
426  0, 0, 0, 0, 0, 1, 0, 0, 0, 0, &
427  0, 0, 0, 0, 0, 0, 1, 0, 0, 0, &
428  0, 0, 0, 0, 0, 0, 0, 1, 0, 0, &
429  0, 0, 0, 0, 0, 0, 0, 0, 1, 0, &
430  0, 0, 0, 0, 0, 0, 0, 0, 0, 1, &
431  1, 1, 1, 0, 2, 2, 2, 0, 0, 0, &
432  1, 1, 0, 1, 2, 0, 0, 2, 2, 0, &
433  0, 1, 1, 1, 0, 2, 0, 0, 2, 2, &
434  1, 0, 1, 1, 0, 0, 2, 2, 0, 2, &
435  3, 3, 3, 3, 4, 4, 4, 4, 4, 4 &
436 
437 
438 
439  ],[nnode(8),ncellnode(geomtype(8))])
440 
441 
442  integer, dimension(NNODE(9),NCELLNODE(GEOMTYPE(9))), parameter :: cellnodeparentnodeweights9 = &
443  reshape([&
444  1, 0, 0, 0, 0, 0, &
445  0, 1, 0, 0, 0, 0, &
446  0, 0, 1, 0, 0, 0, &
447  0, 0, 0, 1, 0, 0, &
448  0, 0, 0, 0, 1, 0, &
449  0, 0, 0, 0, 0, 1, &
450  1, 1, 0, 0, 0, 0, &
451  0, 1, 1, 0, 0, 0, &
452  1, 0, 1, 0, 0, 0, &
453  1, 0, 0, 1, 0, 0, &
454  0, 1, 0, 0, 1, 0, &
455  0, 0, 1, 0, 0, 1, &
456  0, 0, 0, 1, 1, 0, &
457  0, 0, 0, 0, 1, 1, &
458  0, 0, 0, 1, 0, 1, &
459  1, 1, 1, 0, 0, 0, &
460  1, 1, 0, 1, 1, 0, &
461  0, 1, 1, 0, 1, 1, &
462  1, 0, 1, 1, 0, 1, &
463  0, 0, 0, 1, 1, 1, &
464  1, 1, 1, 1, 1, 1 &
465 
466 
467 
468  ],[nnode(9),ncellnode(geomtype(9))])
469 
470 
471  integer, dimension(NNODE(10),NCELLNODE(GEOMTYPE(10))), parameter :: cellnodeparentnodeweights10 = &
472  reshape([&
473  1, 0, 0, 0, 0, 0, 0, 0, &
474  0, 1, 0, 0, 0, 0, 0, 0, &
475  0, 0, 1, 0, 0, 0, 0, 0, &
476  0, 0, 0, 1, 0, 0, 0, 0, &
477  0, 0, 0, 0, 1, 0, 0, 0, &
478  0, 0, 0, 0, 0, 1, 0, 0, &
479  0, 0, 0, 0, 0, 0, 1, 0, &
480  0, 0, 0, 0, 0, 0, 0, 1 &
481 
482 
483 
484  ],[nnode(10),ncellnode(geomtype(10))])
485 
486 
487  integer, dimension(NNODE(11),NCELLNODE(GEOMTYPE(11))), parameter :: cellnodeparentnodeweights11 = &
488  reshape([&
489  1, 0, 0, 0, 0, 0, 0, 0, & !
490  0, 1, 0, 0, 0, 0, 0, 0, & !
491  0, 0, 1, 0, 0, 0, 0, 0, & !
492  0, 0, 0, 1, 0, 0, 0, 0, & !
493  0, 0, 0, 0, 1, 0, 0, 0, & ! 5
494  0, 0, 0, 0, 0, 1, 0, 0, & !
495  0, 0, 0, 0, 0, 0, 1, 0, & !
496  0, 0, 0, 0, 0, 0, 0, 1, & !
497  1, 1, 0, 0, 0, 0, 0, 0, & !
498  0, 1, 1, 0, 0, 0, 0, 0, & ! 10
499  0, 0, 1, 1, 0, 0, 0, 0, & !
500  1, 0, 0, 1, 0, 0, 0, 0, & !
501  0, 0, 0, 0, 1, 1, 0, 0, & !
502  0, 0, 0, 0, 0, 1, 1, 0, & !
503  0, 0, 0, 0, 0, 0, 1, 1, & ! 15
504  0, 0, 0, 0, 1, 0, 0, 1, & !
505  1, 0, 0, 0, 1, 0, 0, 0, & !
506  0, 1, 0, 0, 0, 1, 0, 0, & !
507  0, 0, 1, 0, 0, 0, 1, 0, & !
508  0, 0, 0, 1, 0, 0, 0, 1, & ! 20
509  1, 1, 1, 1, 0, 0, 0, 0, & !
510  1, 1, 0, 0, 1, 1, 0, 0, & !
511  0, 1, 1, 0, 0, 1, 1, 0, & !
512  0, 0, 1, 1, 0, 0, 1, 1, & !
513  1, 0, 0, 1, 1, 0, 0, 1, & ! 25
514  0, 0, 0, 0, 1, 1, 1, 1, & !
515  1, 1, 1, 1, 1, 1, 1, 1 & !
516 
517 
518 
519  ],[nnode(11),ncellnode(geomtype(11))])
520 
521 
522  integer, dimension(NNODE(12),NCELLNODE(GEOMTYPE(12))), parameter :: cellnodeparentnodeweights12 = &
523  reshape([&
524  1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
525  0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
526  0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
527  0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
528  0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & ! 5
529  0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
530  0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
531  0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
532  0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
533  0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & ! 10
534  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
535  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, & !
536  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, & !
537  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, & !
538  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, & ! 15
539  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, & !
540  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, & !
541  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, & !
542  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, & !
543  0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, & ! 20
544  1, 1, 1, 1, 0, 0, 0, 0, 2, 2, 2, 2, 0, 0, 0, 0, 0, 0, 0, 0, & !
545  1, 1, 0, 0, 1, 1, 0, 0, 2, 0, 0, 0, 2, 0, 0, 0, 2, 2, 0, 0, & !
546  0, 1, 1, 0, 0, 1, 1, 0, 0, 2, 0, 0, 0, 2, 0, 0, 0, 2, 2, 0, & !
547  0, 0, 1, 1, 0, 0, 1, 1, 0, 0, 2, 0, 0, 0, 2, 0, 0, 0, 2, 2, & !
548  1, 0, 0, 1, 1, 0, 0, 1, 0, 0, 0, 2, 0, 0, 0, 2, 2, 0, 0, 2, & ! 25
549  0, 0, 0, 0, 1, 1, 1, 1, 0, 0, 0, 0, 2, 2, 2, 2, 0, 0, 0, 0, & !
550  3, 3, 3, 3, 3, 3, 3, 3, 4, 4, 4, 4, 4, 4, 4, 4, 4, 4, 4, 4 & !
551 
552 
553 
554  ],[nnode(12),ncellnode(geomtype(12))])
555 
556 
557  integer, dimension(NNODE(13),NCELLNODE(GEOMTYPE(13))), parameter :: cellnodeparentnodeweights13 = &
558  reshape([&
559  1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
560  0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
561  0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
562  0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
563  0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & ! 5
564  0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
565  0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
566  0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
567  1, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
568  0, 1, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & ! 10
569  0, 1, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
570  0, 0, 1, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
571  0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
572  0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, 0, & !
573  0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, & ! 15
574  1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, 0, & !
575  1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, & !
576  0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, & !
577  0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, & !
578  0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, & ! 20
579  0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, & !
580  0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, & !
581  0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, & !
582  0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, & !
583  0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, & ! 25
584  0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, 0, & !
585  0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, & !
586  0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, 0, & !
587  0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, & !
588  0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, 0, & ! 30
589  0, 0, 0, 0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, & !
590  0, 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2, 0, 0, 0, 0, & !
591  4, 1, 1, 1, 0, 0, 0, 0, 8, 2, 2, 8, 0, 0, 0, 0, 0, 0, 0, 0, & !
592  1, 4, 1, 1, 0, 0, 0, 0, 8, 8, 2, 2, 0, 0, 0, 0, 0, 0, 0, 0, & !
593  1, 1, 4, 1, 0, 0, 0, 0, 2, 8, 8, 2, 0, 0, 0, 0, 0, 0, 0, 0, & ! 35
594  1, 1, 1, 4, 0, 0, 0, 0, 2, 2, 8, 8, 0, 0, 0, 0, 0, 0, 0, 0, & !
595  4, 1, 0, 0, 1, 1, 0, 0, 8, 0, 0, 0, 2, 0, 0, 0, 8, 2, 0, 0, & !
596  1, 4, 0, 0, 1, 1, 0, 0, 8, 0, 0, 0, 2, 0, 0, 0, 2, 8, 0, 0, & !
597  0, 4, 1, 0, 0, 1, 1, 0, 0, 8, 0, 0, 0, 2, 0, 0, 0, 8, 2, 0, & !
598  0, 1, 4, 0, 0, 1, 1, 0, 0, 8, 0, 0, 0, 2, 0, 0, 0, 2, 8, 0, & ! 40
599  0, 0, 4, 1, 0, 0, 1, 1, 0, 0, 8, 0, 0, 0, 2, 0, 0, 0, 8, 2, & !
600  0, 0, 1, 4, 0, 0, 1, 1, 0, 0, 8, 0, 0, 0, 2, 0, 0, 0, 2, 8, & !
601  1, 0, 0, 4, 1, 0, 0, 1, 0, 0, 0, 8, 0, 0, 0, 2, 2, 0, 0, 8, & !
602  4, 0, 0, 1, 1, 0, 0, 1, 0, 0, 0, 8, 0, 0, 0, 2, 8, 0, 0, 2, & !
603  1, 1, 0, 0, 4, 1, 0, 0, 2, 0, 0, 0, 8, 0, 0, 0, 8, 2, 0, 0, & ! 45
604  1, 1, 0, 0, 1, 4, 0, 0, 2, 0, 0, 0, 8, 0, 0, 0, 2, 8, 0, 0, & !
605  0, 1, 1, 0, 0, 4, 1, 0, 0, 2, 0, 0, 0, 8, 0, 0, 0, 8, 2, 0, & !
606  0, 1, 1, 0, 0, 1, 4, 0, 0, 2, 0, 0, 0, 8, 0, 0, 0, 2, 8, 0, & !
607  0, 0, 1, 1, 0, 0, 4, 1, 0, 0, 2, 0, 0, 0, 8, 0, 0, 0, 8, 2, & !
608  0, 0, 1, 1, 0, 0, 1, 4, 0, 0, 2, 0, 0, 0, 8, 0, 0, 0, 2, 8, & ! 50
609  1, 0, 0, 1, 1, 0, 0, 4, 0, 0, 0, 2, 0, 0, 0, 8, 2, 0, 0, 8, & !
610  1, 0, 0, 1, 4, 0, 0, 1, 0, 0, 0, 2, 0, 0, 0, 8, 8, 0, 0, 2, & !
611  0, 0, 0, 0, 4, 1, 1, 1, 0, 0, 0, 0, 8, 2, 2, 8, 0, 0, 0, 0, & !
612  0, 0, 0, 0, 1, 4, 1, 1, 0, 0, 0, 0, 8, 8, 2, 2, 0, 0, 0, 0, & !
613  0, 0, 0, 0, 1, 1, 4, 1, 0, 0, 0, 0, 2, 8, 8, 2, 0, 0, 0, 0, & ! 55
614  0, 0, 0, 0, 1, 1, 1, 4, 0, 0, 0, 0, 2, 2, 8, 8, 0, 0, 0, 0, & !
615  24, 8, 4, 8, 8, 4, 3, 4, 32,12,12,32, 12, 4, 4,12, 32,12, 4,12, & !
616  8,24, 8, 4, 4, 8, 4, 3, 32,32,12,12, 12,12, 4, 4, 12,32,12, 4, & !
617  4, 8,24, 8, 3, 4, 8, 4, 12,32,32,12, 4,12,12, 4, 4,12,32,12, & !
618  8, 4, 8,24, 4, 3, 4, 8, 12,12,32,32, 4, 4,12,12, 12, 4,12,32, & ! 60
619  8, 4, 3, 4, 24, 8, 4, 8, 12, 4, 4,12, 32,12,12,32, 32,12, 4,12, & !
620  4, 8, 4, 3, 8,24, 8, 4, 12,12, 4, 4, 32,32,12,12, 12,32,12, 4, & !
621  3, 4, 8, 4, 4, 8,24, 8, 4,12,12, 4, 12,32,32,12, 4,12,32,12, & !
622  4, 3, 4, 8, 8, 4, 8,24, 4, 4,12,12, 12,12,32,32, 12, 4,12,32 & !
623 
624 
625 
626  ],[nnode(13),ncellnode(geomtype(13))])
627 
628 
629 
630  integer, dimension(NCELLNODEPERCELL(CELLTYPE(1)),NIP(1)), parameter :: cell1 = &
631  reshape([&
632  1,2,3 &
633 
634 
635 
636  ],[ncellnodepercell(celltype(1)),nip(1)])
637 
638 
639  integer, dimension(NCELLNODEPERCELL(CELLTYPE(2)),NIP(2)), parameter :: cell2 = &
640  reshape([&
641  1, 4, 7, 6, &
642  2, 5, 7, 4, &
643  3, 6, 7, 5 &
644 
645 
646 
647  ],[ncellnodepercell(celltype(2)),nip(2)])
648 
649 
650  integer, dimension(NCELLNODEPERCELL(CELLTYPE(3)),NIP(3)), parameter :: cell3 = &
651  reshape([&
652  1, 5, 9, 8, &
653  5, 2, 6, 9, &
654  8, 9, 7, 4, &
655  9, 6, 3, 7 &
656 
657 
658 
659  ],[ncellnodepercell(celltype(3)),nip(3)])
660 
661 
662  integer, dimension(NCELLNODEPERCELL(CELLTYPE(4)),NIP(4)), parameter :: cell4 = &
663  reshape([&
664  1, 5,13,12, &
665  5, 6,14,13, &
666  6, 2, 7,14, &
667  12,13,16,11, &
668  13,14,15,16, &
669  14, 7, 8,15, &
670  11,16,10, 4, &
671  16,15, 9,10, &
672  15, 8, 3, 9 &
673 
674 
675 
676  ],[ncellnodepercell(celltype(4)),nip(4)])
677 
678 
679  integer, dimension(NCELLNODEPERCELL(CELLTYPE(5)),NIP(5)), parameter :: cell5 = &
680  reshape([&
681  1, 2, 3, 4 &
682 
683 
684 
685  ],[ncellnodepercell(celltype(5)),nip(5)])
686 
687 
688  integer, dimension(NCELLNODEPERCELL(CELLTYPE(6)),NIP(6)), parameter :: cell6 = &
689  reshape([&
690  1, 5,11, 7, 8,12,15,14, &
691  5, 2, 6,11,12, 9,13,15, &
692  7,11, 6, 3,14,15,13,10, &
693  8,12,15, 4, 4, 9,13,10 &
694 
695 
696 
697  ],[ncellnodepercell(celltype(6)),nip(6)])
698 
699 
700  integer, dimension(NCELLNODEPERCELL(CELLTYPE(7)),NIP(7)), parameter :: cell7 = &
701  reshape([&
702  1, 7,16, 9,10,17,21,19, &
703  7, 2, 8,16,17,11,18,21, &
704  9,16, 8, 3,19,21,18,12, &
705  10,17,21,19, 4,13,20,15, &
706  17,11,18,21,13, 5,14,20, &
707  19,21,18,12,15,20,14, 6 &
708 
709 
710 
711  ],[ncellnodepercell(celltype(7)),nip(7)])
712 
713 
714  integer, dimension(NCELLNODEPERCELL(CELLTYPE(8)),NIP(8)), parameter :: cell8 = &
715  reshape([&
716  1, 2, 3, 4, 5, 6, 7, 8 &
717 
718 
719 
720  ],[ncellnodepercell(celltype(8)),nip(8)])
721 
722 
723  integer, dimension(NCELLNODEPERCELL(CELLTYPE(9)),NIP(9)), parameter :: cell9 = &
724  reshape([&
725  1, 9,21,12,17,22,27,25, &
726  9, 2,10,21,22,18,23,27, &
727  12,21,11, 4,25,27,24,20, &
728  21,10, 3,11,27,23,19,24, &
729  17,22,27,25, 5,13,26,16, &
730  22,18,23,27,13, 6,14,26, &
731  25,27,24,20,16,26,15, 8, &
732  27,23,19,24,26,14, 7,15 &
733 
734 
735 
736  ],[ncellnodepercell(celltype(9)),nip(9)])
737 
738 
739  integer, dimension(NCELLNODEPERCELL(CELLTYPE(10)),NIP(10)), parameter :: cell10 = &
740  reshape([&
741  1, 9,33,16,17,37,57,44, &
742  9,10,34,33,37,38,58,57, &
743  10, 2,11,34,38,18,39,58, &
744  16,33,36,15,44,57,60,43, &
745  33,34,35,36,57,58,59,60, &
746  34,11,12,35,58,39,40,59, &
747  15,36,14, 4,43,60,42,20, &
748  36,35,13,14,60,59,41,42, &
749  35,12, 3,13,59,40,19,41, &
750  17,37,57,44,21,45,61,52, &
751  37,38,58,57,45,46,62,61, &
752  38,18,39,58,46,22,47,62, &
753  44,57,60,43,52,61,64,51, &
754  57,58,59,60,61,62,63,64, &
755  58,39,40,59,62,47,48,63, &
756  43,60,42,20,51,64,50,24, &
757  60,59,41,42,64,63,49,50, &
758  59,40,19,41,63,48,23,49, &
759  21,45,61,52, 5,25,53,32, &
760  45,46,62,61,25,26,54,53, &
761  46,22,47,62,26, 6,27,54, &
762  52,61,64,51,32,53,56,31, &
763  61,62,63,64,53,54,55,56, &
764  62,47,48,63,54,27,28,55, &
765  51,64,50,24,31,56,30, 8, &
766  64,63,49,50,56,55,29,30, &
767  63,48,23,49,55,28, 7,29 &
768 
769 
770 
771  ],[ncellnodepercell(celltype(10)),nip(10)])
772 
773 
774 
775  integer, dimension(NCELLNODEPERCELLFACE(1),NIPNEIGHBOR(1)), parameter :: cellface1 = &
776  reshape([&
777  2,3, &
778  3,1, &
779  1,2 &
780 
781 
782 
784 
785 
786  integer, dimension(NCELLNODEPERCELLFACE(2),NIPNEIGHBOR(2)), parameter :: cellface2 = &
787  reshape([&
788  2,3, &
789  4,1, &
790  3,4, &
791  1,2 &
792 
793 
794 
796 
797 
798  integer, dimension(NCELLNODEPERCELLFACE(3),NIPNEIGHBOR(3)), parameter :: cellface3 = &
799  reshape([&
800  1,3,2, &
801  1,2,4, &
802  2,3,4, &
803  1,4,3 &
804 
805 
806 
808 
809 
810  integer, dimension(NCELLNODEPERCELLFACE(4),NIPNEIGHBOR(4)), parameter :: cellface4 = &
811  reshape([&
812  2,3,7,6, &
813  4,1,5,8, &
814  3,4,8,7, &
815  1,2,6,5, &
816  5,6,7,8, &
817  1,4,3,2 &
818 
819 
820 
822 
823 
824 
825 contains
826 
827 
828 !---------------------------------------------------------------------------------------------------
830 !---------------------------------------------------------------------------------------------------
831 subroutine telement_init(self,elemType)
832 
833  class(telement) :: self
834  integer, intent(in) :: elemType
835 
836  self%elemType = elemtype
837 
838  self%Nnodes = nnode(self%elemType)
839  self%geomType = geomtype(self%elemType)
840 
841  select case (self%elemType)
842  case(1)
843  self%cellNodeParentNodeWeights = cellnodeparentnodeweights1
844  case(2)
845  self%cellNodeParentNodeWeights = cellnodeparentnodeweights2
846  case(3)
847  self%cellNodeParentNodeWeights = cellnodeparentnodeweights3
848  case(4)
849  self%cellNodeParentNodeWeights = cellnodeparentnodeweights4
850  case(5)
851  self%cellNodeParentNodeWeights = cellnodeparentnodeweights5
852  case(6)
853  self%cellNodeParentNodeWeights = cellnodeparentnodeweights6
854  case(7)
855  self%cellNodeParentNodeWeights = cellnodeparentnodeweights7
856  case(8)
857  self%cellNodeParentNodeWeights = cellnodeparentnodeweights8
858  case(9)
859  self%cellNodeParentNodeWeights = cellnodeparentnodeweights9
860  case(10)
861  self%cellNodeParentNodeWeights = cellnodeparentnodeweights10
862  case(11)
863  self%cellNodeParentNodeWeights = cellnodeparentnodeweights11
864  case(12)
865  self%cellNodeParentNodeWeights = cellnodeparentnodeweights12
866  case(13)
867  self%cellNodeParentNodeWeights = cellnodeparentnodeweights13
868  case default
869  call io_error(0,ext_msg='invalid element type')
870  end select
871 
872 
873  self%NcellNodes = ncellnode(self%geomType)
874  self%nIPs = nip(self%geomType)
875  self%cellType = celltype(self%geomType)
876 
877  select case (self%geomType)
878  case(1)
879  self%IPneighbor = ipneighbor1
880  self%cell = cell1
881  case(2)
882  self%IPneighbor = ipneighbor2
883  self%cell = cell2
884  case(3)
885  self%IPneighbor = ipneighbor3
886  self%cell = cell3
887  case(4)
888  self%IPneighbor = ipneighbor4
889  self%cell = cell4
890  case(5)
891  self%IPneighbor = ipneighbor5
892  self%cell = cell5
893  case(6)
894  self%IPneighbor = ipneighbor6
895  self%cell = cell6
896  case(7)
897  self%IPneighbor = ipneighbor7
898  self%cell = cell7
899  case(8)
900  self%IPneighbor = ipneighbor8
901  self%cell = cell8
902  case(9)
903  self%IPneighbor = ipneighbor9
904  self%cell = cell9
905  case(10)
906  self%IPneighbor = ipneighbor10
907  self%cell = cell10
908  end select
909 
910  self%NcellnodesPerCell = ncellnodepercell(self%cellType)
911 
912  select case(self%cellType)
913  case(1)
914  self%cellFace = cellface1
915  self%vtkType = 'TRIANGLE'
916  case(2)
917  self%cellFace = cellface2
918  self%vtkType = 'QUAD'
919  case(3)
920  self%cellFace = cellface3
921  self%vtkType = 'TETRA'
922  case(4)
923  self%cellFace = cellface4
924  self%vtkType = 'HEXAHEDRON'
925  end select
926 
927  self%nIPneighbors = size(self%IPneighbor,1)
928 
929  write(6,'(/,a)') ' <<<+- element_init -+>>>'; flush(6)
930 
931  write(6,*) ' element type: ',self%elemType
932  write(6,*) ' geom type: ',self%geomType
933  write(6,*) ' cell type: ',self%cellType
934  write(6,*) ' # node: ',self%Nnodes
935  write(6,*) ' # IP: ',self%nIPs
936  write(6,*) ' # cellnode: ',self%Ncellnodes
937  write(6,*) ' # cellnode/cell: ',self%NcellnodesPerCell
938  write(6,*) ' # IP neighbor: ',self%nIPneighbors
939 
940 end subroutine telement_init
941 
942 end module element
element::cell8
integer, dimension(ncellnodepercell(celltype(8)), nip(8)), parameter cell8
Definition: element.f90:714
element::cellnodeparentnodeweights1
integer, dimension(nnode(1), ncellnode(geomtype(1))), parameter cellnodeparentnodeweights1
Definition: element.f90:300
element::cellnodeparentnodeweights9
integer, dimension(nnode(9), ncellnode(geomtype(9))), parameter cellnodeparentnodeweights9
Definition: element.f90:442
element::cellnodeparentnodeweights5
integer, dimension(nnode(5), ncellnode(geomtype(5))), parameter cellnodeparentnodeweights5
Definition: element.f90:367
element::ipneighbor1
integer, dimension(nipneighbor(celltype(1)), nip(1)), parameter ipneighbor1
Definition: element.f90:155
element::telement
Properties of a single element.
Definition: element.f90:18
element::cellnodeparentnodeweights7
integer, dimension(nnode(7), ncellnode(geomtype(7))), parameter cellnodeparentnodeweights7
Definition: element.f90:396
element::nnode
integer, dimension(nelemtype), parameter nnode
number of nodes that constitute a specific type of element
Definition: element.f90:47
io::io_error
subroutine, public io_error(error_ID, el, ip, g, instance, ext_msg)
write error statements to standard out and terminate the Marc/spectral run with exit #9xxx
Definition: IO.f90:305
element::cell1
integer, dimension(ncellnodepercell(celltype(1)), nip(1)), parameter cell1
Definition: element.f90:630
element::telement_init
subroutine telement_init(self, elemType)
define properties of an element
Definition: element.f90:832
element::nipneighbor
integer, dimension(maxval(celltype)), parameter nipneighbor
number of ip neighbors / cell faces
Definition: element.f90:125
element::cellnodeparentnodeweights4
integer, dimension(nnode(4), ncellnode(geomtype(4))), parameter cellnodeparentnodeweights4
Definition: element.f90:343
element::ncellnodepercellface
integer, dimension(maxval(celltype)), parameter ncellnodepercellface
number of cell nodes per face
Definition: element.f90:133
element::ipneighbor10
integer, dimension(nipneighbor(celltype(10)), nip(10)), parameter ipneighbor10
Definition: element.f90:264
element::ipneighbor8
integer, dimension(nipneighbor(celltype(8)), nip(8)), parameter ipneighbor8
Definition: element.f90:239
element::ipneighbor4
integer, dimension(nipneighbor(celltype(4)), nip(4)), parameter ipneighbor4
Definition: element.f90:187
element::cellnodeparentnodeweights13
integer, dimension(nnode(13), ncellnode(geomtype(13))), parameter cellnodeparentnodeweights13
Definition: element.f90:557
element::cell2
integer, dimension(ncellnodepercell(celltype(2)), nip(2)), parameter cell2
Definition: element.f90:639
element::ipneighbor9
integer, dimension(nipneighbor(celltype(9)), nip(9)), parameter ipneighbor9
Definition: element.f90:248
element::cellface2
integer, dimension(ncellnodepercellface(2), nipneighbor(2)), parameter cellface2
Definition: element.f90:786
element::ncellnodepercell
integer, dimension(maxval(celltype)), parameter ncellnodepercell
number of total cell nodes
Definition: element.f90:141
element::cellface1
integer, dimension(ncellnodepercellface(1), nipneighbor(1)), parameter cellface1
Definition: element.f90:775
element::cell3
integer, dimension(ncellnodepercell(celltype(3)), nip(3)), parameter cell3
Definition: element.f90:650
element::celltype
integer, dimension(maxval(geomtype)), parameter celltype
cell type
Definition: element.f90:111
io
input/output functions, partly depending on chosen solver
Definition: IO.f90:12
element::cellface4
integer, dimension(ncellnodepercellface(4), nipneighbor(4)), parameter cellface4
Definition: element.f90:810
element::cellnodeparentnodeweights3
integer, dimension(nnode(3), ncellnode(geomtype(3))), parameter cellnodeparentnodeweights3
Definition: element.f90:326
element::cell7
integer, dimension(ncellnodepercell(celltype(7)), nip(7)), parameter cell7
Definition: element.f90:700
element::cell4
integer, dimension(ncellnodepercell(celltype(4)), nip(4)), parameter cell4
Definition: element.f90:662
element::cell6
integer, dimension(ncellnodepercell(celltype(6)), nip(6)), parameter cell6
Definition: element.f90:688
element::cellnodeparentnodeweights11
integer, dimension(nnode(11), ncellnode(geomtype(11))), parameter cellnodeparentnodeweights11
Definition: element.f90:487
element::cellnodeparentnodeweights12
integer, dimension(nnode(12), ncellnode(geomtype(12))), parameter cellnodeparentnodeweights12
Definition: element.f90:522
element::ipneighbor5
integer, dimension(nipneighbor(celltype(5)), nip(5)), parameter ipneighbor5
Definition: element.f90:204
element::ncellnode
integer, dimension(maxval(geomtype)), parameter ncellnode
number of cell nodes
Definition: element.f90:83
element
Definition: element.f90:9
element::cellnodeparentnodeweights10
integer, dimension(nnode(10), ncellnode(geomtype(10))), parameter cellnodeparentnodeweights10
Definition: element.f90:471
element::cellnodeparentnodeweights2
integer, dimension(nnode(2), ncellnode(geomtype(2))), parameter cellnodeparentnodeweights2
Definition: element.f90:311
element::cellnodeparentnodeweights6
integer, dimension(nnode(6), ncellnode(geomtype(6))), parameter cellnodeparentnodeweights6
Definition: element.f90:384
element::ipneighbor3
integer, dimension(nipneighbor(celltype(3)), nip(3)), parameter ipneighbor3
Definition: element.f90:175
element::nelemtype
integer, parameter nelemtype
Definition: element.f90:44
element::geomtype
integer, dimension(nelemtype), parameter geomtype
geometry type (same number of cell nodes and IPs)
Definition: element.f90:65
element::cell10
integer, dimension(ncellnodepercell(celltype(10)), nip(10)), parameter cell10
Definition: element.f90:739
element::cell5
integer, dimension(ncellnodepercell(celltype(5)), nip(5)), parameter cell5
Definition: element.f90:679
element::cellnodeparentnodeweights8
integer, dimension(nnode(8), ncellnode(geomtype(8))), parameter cellnodeparentnodeweights8
Definition: element.f90:419
element::ipneighbor2
integer, dimension(nipneighbor(celltype(2)), nip(2)), parameter ipneighbor2
Definition: element.f90:164
element::ipneighbor6
integer, dimension(nipneighbor(celltype(6)), nip(6)), parameter ipneighbor6
Definition: element.f90:213
element::ipneighbor7
integer, dimension(nipneighbor(celltype(7)), nip(7)), parameter ipneighbor7
Definition: element.f90:225
element::cellface3
integer, dimension(ncellnodepercellface(3), nipneighbor(3)), parameter cellface3
Definition: element.f90:798
element::nip
integer, dimension(maxval(geomtype)), parameter nip
number of IPs
Definition: element.f90:97
element::cell9
integer, dimension(ncellnodepercell(celltype(9)), nip(9)), parameter cell9
Definition: element.f90:723