FrontISTR  5.7.0
Large-scale structural analysis program with finit element method
element.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
2 ! Copyright (c) 2019 FrontISTR Commons
3 ! This software is released under the MIT License, see LICENSE.txt
4 !-------------------------------------------------------------------------------
31 ! If you wish introduce new elements with new geometry or/and
32 ! new shape functions, you need do the following
33 !!
34 !!
35 !!- Introduce new element ID corresponding to your element in this module.
36 !!- Provide corresponding shape function and shape derivative in a new
37 !! module and include this new module here.
38 !!- Add the information of your new element into functions listed above.
39 !!- If new quadrature method needed, do some modification in MODULE
40 !! Quadrature.
41 !!- If you introduce new surface geometry and do contact calculation also,
42 !! you may need do some modification in MODULE mSurfElement.
44 
45  use shape_line2n
46  use shape_line3n
47  use shape_tri3n
48  use shape_tri6n
49  use shape_quad4n
50  use shape_quad8n
51  use shape_quad9n
52  use shape_hex8n
53  use shape_hex20n
54  use shape_tet4n
55  use shape_tet10n
56  use shape_prism6n
57  use shape_prism15n
58  implicit none
59 
60  integer, parameter, private :: kreal = kind(0.0d0)
61 
62  !-------------------------------------------
63  ! Fllowing ID of element types
64  !-------------------------------------------
65  integer, parameter :: fe_unknown = -1
66 
67  integer, parameter :: fe_line2n = 111
68  integer, parameter :: fe_line3n = 112
69  integer, parameter :: fe_tri3n = 231
70  integer, parameter :: fe_tri6n = 232
71  integer, parameter :: fe_tri6nc = 2322
72  integer, parameter :: fe_quad4n = 241
73  integer, parameter :: fe_quad8n = 242
74  integer, parameter :: fe_truss = 301
75  integer, parameter :: fe_tet4n = 341
76  integer, parameter :: fe_tet4n_pipi = 3414
77  integer, parameter :: fe_tet10n = 342
78  integer, parameter :: fe_tet10nc = 3422
79  integer, parameter :: fe_prism6n = 351
80  integer, parameter :: fe_prism15n = 352
81  integer, parameter :: fe_hex8n = 361
82  integer, parameter :: fe_hex20n = 362
83  integer, parameter :: fe_hex27n = 363
84 
85  integer, parameter :: fe_if_line2n = 511
86 
87  integer, parameter :: fe_beam2n = 611
88  integer, parameter :: fe_beam3n = 612
89  integer, parameter :: fe_beam341 = 641
90 
91  integer, parameter :: fe_tri6n_shell = 732
92  integer, parameter :: fe_dsg3_shell = 733
93  integer, parameter :: fe_mitc3_shell = 731
94  integer, parameter :: fe_mitc4_shell = 741
95  integer, parameter :: fe_mitc8_shell = 742
96  integer, parameter :: fe_mitc9_shell = 743
97 
98  integer, parameter :: fe_mitc3_shell361 = 761
99  integer, parameter :: fe_mitc4_shell361 = 781
100 
101  integer, parameter :: fe_nodesmooth_tet4n = 881
102  integer, parameter :: fe_edgesmooth_tet4n = 891
103 
104  integer, parameter :: fe_tri3n_patch = 1031
105  integer, parameter :: fe_tri6n_patch = 1032
106  integer, parameter :: fe_quad4n_patch = 1041
107  integer, parameter :: fe_quad8n_patch = 1042
108  ! ---------------------------------------------
109 
110 contains
111 
112  !************************************
113  ! Following geometric information
114  !************************************
116  integer(kind=kind(2)) function getspacedimension( etype )
117  integer, intent(in) :: etype
118 
119  select case( etype)
124  case default
126  end select
127  end function
128 
130  integer(kind=kind(2)) function getnumberofnodes( etype )
131  integer, intent(in) :: etype
132 
133  select case (etype)
135  getnumberofnodes = 2
136  case (fe_line3n, fe_beam3n)
137  getnumberofnodes = 3
139  getnumberofnodes = 3
141  getnumberofnodes = 6
143  getnumberofnodes = 4
145  getnumberofnodes = 8
146  case ( fe_mitc9_shell )
147  getnumberofnodes = 9
149  getnumberofnodes = 4
150  case ( fe_tet10n, fe_tet10nc )
151  getnumberofnodes = 10
152  case ( fe_prism6n )
153  getnumberofnodes = 6
154  case ( fe_prism15n )
155  getnumberofnodes = 15
156  case ( fe_hex8n )
157  getnumberofnodes = 8
158  case ( fe_hex20n )
159  getnumberofnodes = 20
160  case default
161  getnumberofnodes = -1
162  ! error message
163  end select
164  end function
165 
167  integer(kind=kind(2)) function getnumberofsubface( etype )
168  integer, intent(in) :: etype
169 
170  select case (etype)
179  case ( fe_prism6n, fe_prism15n )
181  case ( fe_hex8n, fe_hex20n)
185  case default
186  getnumberofsubface = -1
187  ! error message
188  end select
189  end function
190 
192  subroutine getsubface( intype, innumber, outtype, nodes )
193  integer, intent(in) :: intype
194  integer, intent(in) :: innumber
195  integer, intent(out) :: outtype
196  integer, intent(out) :: nodes(:)
197 
198  if( innumber>getnumberofsubface( intype ) ) stop "Error in getting subface"
199  select case ( intype )
201  outtype = fe_tri3n
202  select case ( innumber )
203  case (1)
204  nodes(1)=1; nodes(2)=2; nodes(3)=3
205  case (2)
206  nodes(1)=4; nodes(2)=2; nodes(3)=1
207  case (3)
208  nodes(1)=4; nodes(2)=3; nodes(3)=2
209  case (4)
210  nodes(1)=4; nodes(2)=1; nodes(3)=3
211  end select
212  case (fe_tet10n)
213  outtype = fe_tri6n
214  select case ( innumber )
215  case (1)
216  nodes(1)=1; nodes(2)=2; nodes(3)=3
217  nodes(4)=5; nodes(5)=6; nodes(6)=7
218  case (2)
219  nodes(1)=4; nodes(2)=2; nodes(3)=1
220  nodes(4)=9; nodes(5)=5; nodes(6)=8
221  case (3)
222  nodes(1)=4; nodes(2)=3; nodes(3)=2
223  nodes(4)=10; nodes(5)=6; nodes(6)=9
224  case (4)
225  nodes(1)=4; nodes(2)=1; nodes(3)=3
226  nodes(4)=8; nodes(5)=7; nodes(6)=10
227  end select
228  case (fe_tet10nc)
229  outtype = fe_tri6nc
230  select case ( innumber )
231  case (1)
232  nodes(1)=1; nodes(2)=2; nodes(3)=3
233  nodes(4)=5; nodes(5)=6; nodes(6)=7
234  case (2)
235  nodes(1)=4; nodes(2)=2; nodes(3)=1
236  nodes(4)=9; nodes(5)=5; nodes(6)=8
237  case (3)
238  nodes(1)=4; nodes(2)=3; nodes(3)=2
239  nodes(4)=10; nodes(5)=6; nodes(6)=9
240  case (4)
241  nodes(1)=4; nodes(2)=1; nodes(3)=3
242  nodes(4)=8; nodes(5)=7; nodes(6)=10
243  end select
244  case ( fe_hex8n )
245  outtype = fe_quad4n
246  select case ( innumber )
247  case (1)
248  nodes(1)=1; nodes(2)=2; nodes(3)=3; nodes(4)=4
249  case (2)
250  nodes(1)=8; nodes(2)=7; nodes(3)=6; nodes(4)=5
251  case (3)
252  nodes(1)=5; nodes(2)=6; nodes(3)=2; nodes(4)=1
253  case (4)
254  nodes(1)=6; nodes(2)=7; nodes(3)=3; nodes(4)=2
255  case (5)
256  nodes(1)=7; nodes(2)=8; nodes(3)=4; nodes(4)=3
257  case (6)
258  nodes(1)=8; nodes(2)=5; nodes(3)=1; nodes(4)=4
259  case default
260  ! error
261  end select
262  case (fe_hex20n)
263  outtype = fe_quad8n
264  select case ( innumber )
265  case (1)
266  nodes(1)=1; nodes(2)=2; nodes(3)=3; nodes(4)=4
267  nodes(5)=9; nodes(6)=10; nodes(7)=11; nodes(8)=12
268  case (2)
269  nodes(1)=8; nodes(2)=7; nodes(3)=6; nodes(4)=5
270  nodes(5)=15; nodes(6)=14; nodes(7)=13; nodes(8)=16
271  case (3)
272  nodes(1)=5; nodes(2)=6; nodes(3)=2; nodes(4)=1
273  nodes(5)=13; nodes(6)=18; nodes(7)=9; nodes(8)=17
274  case (4)
275  nodes(1)=6; nodes(2)=7; nodes(3)=3; nodes(4)=2
276  nodes(5)=14; nodes(6)=19; nodes(7)=10; nodes(8)=18
277  case (5)
278  nodes(1)=7; nodes(2)=8; nodes(3)=4; nodes(4)=3
279  nodes(5)=15; nodes(6)=20; nodes(7)=11; nodes(8)=19
280  case (6)
281  nodes(1)=8; nodes(2)=5; nodes(3)=1; nodes(4)=4
282  nodes(5)=16; nodes(6)=17; nodes(7)=12; nodes(8)=20
283  case default
284  ! error
285  end select
286  case (fe_prism6n)
287  select case ( innumber )
288  case (1)
289  outtype = fe_tri3n
290  nodes(1)=1; nodes(2)=2; nodes(3)=3
291  case (2)
292  outtype = fe_tri3n
293  nodes(1)=6; nodes(2)=5; nodes(3)=4
294  case (3)
295  outtype = fe_quad4n
296  nodes(1)=4; nodes(2)=5; nodes(3)=2; nodes(4)=1
297  case (4)
298  outtype = fe_quad4n
299  nodes(1)=5; nodes(2)=6; nodes(3)=3; nodes(4)=2
300  case (5)
301  outtype = fe_quad4n
302  nodes(1)=6; nodes(2)=4; nodes(3)=1; nodes(4)=3
303  end select
304  case (fe_prism15n)
305  select case ( innumber )
306  case (1)
307  outtype = fe_tri6n
308  nodes(1)=1; nodes(2)=2; nodes(3)=3
309  nodes(4)=7; nodes(5)=8; nodes(6)=9
310  case (2)
311  outtype = fe_tri6n
312  nodes(1)=6; nodes(2)=5; nodes(3)=4
313  nodes(4)=11; nodes(5)=10; nodes(6)=12
314  case (3)
315  outtype = fe_quad8n
316  nodes(1)=4; nodes(2)=5; nodes(3)=2; nodes(4)=1
317  nodes(5)=10; nodes(6)=14; nodes(7)=7; nodes(8)=13
318  case (4)
319  outtype = fe_quad8n
320  nodes(1)=5; nodes(2)=6; nodes(3)=3; nodes(4)=2
321  nodes(5)=11; nodes(6)=15; nodes(7)=8; nodes(8)=14
322  case (5)
323  outtype = fe_quad8n
324  nodes(1)=6; nodes(2)=4; nodes(3)=1; nodes(4)=3
325  nodes(5)=12; nodes(6)=13; nodes(7)=9; nodes(8)=15
326  end select
327  case ( fe_tri3n, fe_mitc3_shell )
328  outtype = fe_line2n
329  select case (innumber )
330  case (1)
331  nodes(1) = 1; nodes(2)=2
332  case (2)
333  nodes(1) = 2; nodes(2)=3
334  case (3)
335  nodes(1) = 3; nodes(2)=1
336  end select
338  outtype = fe_line3n
339  select case (innumber )
340  case (1)
341  nodes(1) = 1; nodes(2)=2; nodes(3)=4
342  case (2)
343  nodes(1) = 2; nodes(2)=3; nodes(3)=5
344  case (3)
345  nodes(1) = 3; nodes(2)=1; nodes(3)=6
346  end select
347  case ( fe_quad4n, fe_mitc4_shell )
348  outtype = fe_line2n
349  select case (innumber )
350  case (1)
351  nodes(1) = 1; nodes(2)=2
352  case (2)
353  nodes(1) = 2; nodes(2)=3
354  case (3)
355  nodes(1) = 3; nodes(2)=4
356  case (4)
357  nodes(1) = 4; nodes(2)=1
358  end select
360  outtype = fe_line3n
361  select case (innumber )
362  case (1)
363  nodes(1) = 1; nodes(2)=2; nodes(3)=5
364  case (2)
365  nodes(1) = 2; nodes(2)=3; nodes(3)=6
366  case (3)
367  nodes(1) = 3; nodes(2)=4; nodes(3)=7
368  case (4)
369  nodes(1) = 4; nodes(2)=1; nodes(3)=8
370  end select
371  case (fe_mitc3_shell361)
372  select case ( innumber )
373  case (1)
374  outtype = fe_tri3n
375  nodes(1)=1; nodes(2)=2; nodes(3)=3
376  case (2)
377  outtype = fe_tri3n
378  nodes(1)=6; nodes(2)=5; nodes(3)=4
379  case (3)
380  outtype = fe_quad4n
381  nodes(1)=4; nodes(2)=5; nodes(3)=2; nodes(4)=1
382  case (4)
383  outtype = fe_quad4n
384  nodes(1)=5; nodes(2)=6; nodes(3)=3; nodes(4)=2
385  case (5)
386  outtype = fe_quad4n
387  nodes(1)=6; nodes(2)=4; nodes(3)=1; nodes(4)=3
388  end select
389  case ( fe_mitc4_shell361 )
390  outtype = fe_quad4n
391  select case ( innumber )
392  case (1)
393  nodes(1)=1; nodes(2)=2; nodes(3)=3; nodes(4)=4
394  case (2)
395  nodes(1)=8; nodes(2)=7; nodes(3)=6; nodes(4)=5
396  case (3)
397  nodes(1)=5; nodes(2)=6; nodes(3)=2; nodes(4)=1
398  case (4)
399  nodes(1)=6; nodes(2)=7; nodes(3)=3; nodes(4)=2
400  case (5)
401  nodes(1)=7; nodes(2)=8; nodes(3)=4; nodes(4)=3
402  case (6)
403  nodes(1)=8; nodes(2)=5; nodes(3)=1; nodes(4)=4
404  case default
405  ! error
406  end select
407  case ( fe_tri3n_patch )
408  outtype = fe_tri3n
409  select case ( innumber )
410  case (1)
411  nodes(1)=1; nodes(2)=2; nodes(3)=3
412  case default
413  !error
414  end select
415  case ( fe_tri6n_patch )
416  outtype = fe_tri6n
417  select case ( innumber )
418  case (1)
419  nodes(1)=1; nodes(2)=2; nodes(3)=3
420  nodes(4)=4; nodes(5)=5; nodes(6)=6
421  case default
422  !error
423  end select
424  case ( fe_quad4n_patch )
425  outtype = fe_quad4n
426  select case ( innumber )
427  case (1)
428  nodes(1)=1; nodes(2)=2; nodes(3)=3; nodes(4)=4
429  case default
430  !error
431  end select
432  case ( fe_quad8n_patch )
433  outtype = fe_quad8n
434  select case ( innumber )
435  case (1)
436  nodes(1)=1; nodes(2)=2; nodes(3)=3; nodes(4)=4
437  nodes(5)=5; nodes(6)=6; nodes(7)=7; nodes(8)=8
438  case default
439  !error
440  end select
441  case default
442  outtype = fe_unknown
443  stop "element type not defined-sbs"
444  ! error message
445  end select
446  end subroutine
447 
449  integer function numofquadpoints( fetype )
450  integer, intent(in) :: fetype
451  select case (fetype)
453  numofquadpoints = 1
454  case ( fe_tri6n )
455  numofquadpoints = 3
456  case ( fe_tri6nc )
457  numofquadpoints = 4
458  case (fe_line3n )
459  numofquadpoints = 2
461  numofquadpoints = 4
462  case ( fe_quad8n, fe_mitc9_shell )
463  numofquadpoints = 9
464  case ( fe_hex8n )
465  numofquadpoints = 8
466  case ( fe_hex20n, fe_mitc8_shell )
467  numofquadpoints = 27
468  case ( fe_prism6n )
469  numofquadpoints = 2
471  numofquadpoints = 3
472  case ( fe_prism15n, fe_tri6n_shell )
473  numofquadpoints = 9
474  case ( fe_tet10n, fe_tet4n_pipi )
475  numofquadpoints = 4
476  case ( fe_tet10nc )
477  numofquadpoints = 12
479  numofquadpoints = 1
480  case default
481  numofquadpoints = -1
482  ! error message
483  stop "element type not defined-np"
484  end select
485  end function
486 
488  subroutine getquadpoint( fetype, np, pos )
490  integer, intent(in) :: fetype
491  integer, intent(in) :: np
492  real(kind=kreal), intent(out) :: pos(:)
493 
494  if( np<1 .or. np>numofquadpoints(fetype) ) then
495  ! error
496  endif
497 
498  select case (fetype)
499  case (fe_tri3n)
500  pos(1:2)=gauss2d4(:,np)
501  case ( fe_tri6n, fe_mitc3_shell )
502  pos(1:2)=gauss2d5(:,np)
503  case (fe_tri6nc )
504  pos(1:2)=gauss2d6(:,np)
505  case ( fe_quad4n, fe_mitc4_shell )
506  pos(1:2)=gauss2d2(:,np)
507  case ( fe_quad8n, fe_mitc9_shell )
508  pos(1:2)=gauss2d3(:,np)
509  case ( fe_hex8n, fe_mitc4_shell361 )
510  pos(1:3)=gauss3d2(:,np)
511  case ( fe_hex20n, fe_mitc8_shell )
512  pos(1:3)=gauss3d3(:,np)
513  case ( fe_prism6n, fe_mitc3_shell361 )
514  pos(1:3)=gauss3d7(:,np)
515  case ( fe_prism15n, fe_tri6n_shell )
516  pos(1:3)=gauss3d8(:,np)
517  case ( fe_tet4n, fe_beam341 )
518  pos(1:3)=gauss3d4(:,np)
519  case ( fe_tet10n, fe_tet4n_pipi )
520  pos(1:3)=gauss3d5(:,np)
521  case ( fe_tet10nc )
522  pos(1:3)=np
523  case ( fe_line2n, fe_if_line2n )
524  pos(1:1)=gauss1d1(:,np)
525  case ( fe_line3n )
526  pos(1:1)=gauss1d2(:,np)
527  case default
528  ! error message
529  stop "element type not defined-qp"
530  end select
531  end subroutine
532 
534  real(kind=kreal) function getweight( fetype, np )
536  integer, intent(in) :: fetype
537  integer, intent(in) :: np
538  if( np<1 .or. np>numofquadpoints(fetype) ) then
539  ! error
540  endif
541 
542  select case (fetype)
543  case (fe_tri3n)
544  getweight = weight2d4(1)
545  case ( fe_tri6n, fe_mitc3_shell )
546  getweight = weight2d5(np)
547  case ( fe_quad4n, fe_mitc4_shell )
548  getweight = weight2d2(np)
549  case ( fe_quad8n, fe_mitc9_shell )
550  getweight = weight2d3(np)
551  case ( fe_hex8n, fe_mitc4_shell361 )
552  getweight = weight3d2(np)
553  case ( fe_hex20n)
554  getweight = weight3d3(np)
555  case ( fe_prism6n, fe_mitc3_shell361 )
556  getweight = weight3d7(np)
557  case ( fe_prism15n )
558  getweight = weight3d8(np)
559  case ( fe_tet4n, fe_beam341 )
560  getweight = weight3d4(1)
561  case ( fe_tet10n, fe_tet4n_pipi )
562  getweight = weight3d5(np)
563  case ( fe_line2n, fe_if_line2n )
564  getweight = weight1d1(1)
565  case ( fe_line3n )
566  getweight = weight1d2(np)
567  case default
568  getweight = 0.d0
569  ! error message
570  end select
571  end function
572 
573  !************************************
574  ! Following shape function information
575  !************************************
577  subroutine getshapederiv( fetype, localcoord, shapederiv )
578  integer, intent(in) :: fetype
579  real(kind=kreal), intent(in) :: localcoord(:)
580  real(kind=kreal), intent(out) :: shapederiv(:,:)
581 
582  select case (fetype)
583  case ( fe_tri3n, fe_mitc3_shell )
584  !error check
585  call shapederiv_tri3n(shapederiv(1:3,1:2))
586  case (fe_tri6n)
587  !error check
588  call shapederiv_tri6n(localcoord,shapederiv(1:6,1:2) )
589  case ( fe_quad4n, fe_mitc4_shell )
590  !error check
591  call shapederiv_quad4n(localcoord,shapederiv(1:4,1:2))
592  case (fe_quad8n)
593  !error check
594  call shapederiv_quad8n(localcoord,shapederiv(1:8,1:2))
595  case ( fe_mitc9_shell )
596  !error check
597  call shapederiv_quad9n(localcoord,shapederiv(1:9,1:2))
599  ! error check
600  call shapederiv_hex8n(localcoord,shapederiv(1:8,1:3))
601  case (fe_hex20n)
602  ! error check
603  call shapederiv_hex20n(localcoord, shapederiv(1:20,1:3))
605  call shapederiv_prism6n(localcoord,shapederiv(1:6,1:3))
606  case (fe_prism15n)
607  call shapederiv_prism15n(localcoord,shapederiv(1:15,1:3))
609  ! error check
610  call shapederiv_tet4n(shapederiv(1:4,1:3))
611  case (fe_tet10n)
612  ! error check
613  call shapederiv_tet10n(localcoord,shapederiv(1:10,1:3))
614  case default
615  ! error message
616  stop "Element type not defined-sde"
617  end select
618  end subroutine
619 
621  subroutine getshape2ndderiv( fetype, localcoord, shapederiv )
622  integer, intent(in) :: fetype
623  real(kind=kreal), intent(in) :: localcoord(:)
624  real(kind=kreal), intent(out) :: shapederiv(:,:,:)
625 
626  select case (fetype)
627  case ( fe_tri3n, fe_mitc3_shell )
628  !error check
629  call shape2ndderiv_tri3n(shapederiv(1:3,1:2,1:2))
630  case (fe_tri6n)
631  !error check
632  call shape2ndderiv_tri6n(shapederiv(1:6,1:2,1:2))
633  case ( fe_quad4n, fe_mitc4_shell )
634  !error check
635  call shape2ndderiv_quad4n(shapederiv(1:4,1:2,1:2))
636  case (fe_quad8n)
637  !error check
638  call shape2ndderiv_quad8n(localcoord,shapederiv(1:8,1:2,1:2))
639  case default
640  ! error message
641  stop "Cannot calculate second derivatives of shape function"
642  end select
643  end subroutine
644 
646  subroutine getshapefunc( fetype, localcoord, func )
647  integer, intent(in) :: fetype
648  real(kind=kreal), intent(in) :: localcoord(:)
649  real(kind=kreal), intent(out) :: func(:)
650 
651  select case (fetype)
652  case ( fe_tri3n, fe_mitc3_shell )
653  !error check
654  call shapefunc_tri3n(localcoord,func(1:3))
655  case (fe_tri6n)
656  !error check
657  call shapefunc_tri6n(localcoord,func(1:6))
658  case ( fe_quad4n, fe_mitc4_shell )
659  !error check
660  call shapefunc_quad4n(localcoord,func(1:4))
661  case (fe_quad8n)
662  !error check
663  call shapefunc_quad8n(localcoord,func(1:8))
665  ! error check
666  call shapefunc_hex8n(localcoord,func(1:8))
667  case ( fe_mitc9_shell )
668  !error check
669  call shapefunc_quad9n(localcoord,func(1:9))
670  case (fe_hex20n)
671  ! error check
672  call shapefunc_hex20n(localcoord,func(1:20))
674  call shapefunc_prism6n(localcoord,func(1:6))
675  case (fe_prism15n)
676  call shapefunc_prism15n(localcoord,func(1:15))
678  ! error check
679  call shapefunc_tet4n(localcoord,func(1:4))
680  case (fe_tet10n)
681  ! error check
682  call shapefunc_tet10n(localcoord,func(1:10))
683  case (fe_line2n, fe_if_line2n)
684  !error check
685  call shapefunc_line2n(localcoord,func(1:2))
686  case (fe_line3n)
687  !error check
688  call shapefunc_line3n(localcoord,func(1:3))
689  case default
690  stop "Element type not defined-sf"
691  ! error message
692  end select
693  end subroutine
694 
695 
696  ! (Gaku Hashimoto, The University of Tokyo, 2012/11/15) <
697  !####################################################################
698  subroutine getnodalnaturalcoord(fetype, nncoord)
699  !####################################################################
700 
701  integer, intent(in) :: fetype
702  real(kind = kreal), intent(out) :: nncoord(:, :)
703 
704  !--------------------------------------------------------------------
705 
706  select case( fetype )
708 
709  !error check
710  call nodalnaturalcoord_tri3n( nncoord(1:3, 1:2) )
711 
713 
714  !error check
715  call nodalnaturalcoord_quad4n( nncoord(1:4, 1:2) )
716 
717  case( fe_mitc9_shell )
718 
719  !error check
720  call nodalnaturalcoord_quad9n( nncoord(1:9, 1:2) )
721 
722  case default
723 
724  ! error message
725  stop "Element type not defined-sde"
726 
727  end select
728 
729  !--------------------------------------------------------------------
730 
731  return
732 
733  !####################################################################
734  end subroutine getnodalnaturalcoord
735  !####################################################################
736  ! > (Gaku Hashimoto, The University of Tokyo, 2012/11/15)
737 
738 
740  subroutine getglobalderiv( fetype, nn, localcoord, elecoord, det, gderiv )
741  integer, intent(in) :: fetype
742  integer, intent(in) :: nn
743  real(kind=kreal), intent(in) :: localcoord(:)
744  real(kind=kreal), intent(in) :: elecoord(:,:)
745  real(kind=kreal), intent(out) :: det
746  real(kind=kreal), intent(out) :: gderiv(:,:)
747 
748  real(kind=kreal) :: dum, xj(3,3), xji(3,3), deriv(nn,3)
749  integer :: nspace
750 
751  nspace = getspacedimension( fetype )
752  call getshapederiv( fetype, localcoord(:), deriv(1:nn,:) )
753 
754  if( nspace==2 ) then
755  xj(1:2,1:2)=matmul( elecoord(1:2,1:nn), deriv(1:nn,1:2) )
756  det=xj(1,1)*xj(2,2)-xj(2,1)*xj(1,2)
757  if( det==0.d0 ) stop "Math error in GetGlobalDeriv! Determinant==0.0"
758  dum=1.d0/det
759  xji(1,1)= xj(2,2)*dum
760  xji(1,2)=-xj(1,2)*dum
761  xji(2,1)=-xj(2,1)*dum
762  xji(2,2)= xj(1,1)*dum
763  else
764  ! JACOBI MATRIX
765  xj(1:3,1:3)= matmul( elecoord(1:3,1:nn), deriv(1:nn,1:3) )
766  !DETERMINANT OF JACOBIAN
767  det=xj(1,1)*xj(2,2)*xj(3,3) &
768  +xj(2,1)*xj(3,2)*xj(1,3) &
769  +xj(3,1)*xj(1,2)*xj(2,3) &
770  -xj(3,1)*xj(2,2)*xj(1,3) &
771  -xj(2,1)*xj(1,2)*xj(3,3) &
772  -xj(1,1)*xj(3,2)*xj(2,3)
773  if( det==0.d0 ) stop "Math error in GetGlobalDeriv! Determinant==0.0"
774  ! INVERSION OF JACOBIAN
775  dum=1.d0/det
776  xji(1,1)=dum*( xj(2,2)*xj(3,3)-xj(3,2)*xj(2,3) )
777  xji(1,2)=dum*(-xj(1,2)*xj(3,3)+xj(3,2)*xj(1,3) )
778  xji(1,3)=dum*( xj(1,2)*xj(2,3)-xj(2,2)*xj(1,3) )
779  xji(2,1)=dum*(-xj(2,1)*xj(3,3)+xj(3,1)*xj(2,3) )
780  xji(2,2)=dum*( xj(1,1)*xj(3,3)-xj(3,1)*xj(1,3) )
781  xji(2,3)=dum*(-xj(1,1)*xj(2,3)+xj(2,1)*xj(1,3) )
782  xji(3,1)=dum*( xj(2,1)*xj(3,2)-xj(3,1)*xj(2,2) )
783  xji(3,2)=dum*(-xj(1,1)*xj(3,2)+xj(3,1)*xj(1,2) )
784  xji(3,3)=dum*( xj(1,1)*xj(2,2)-xj(2,1)*xj(1,2) )
785  endif
786 
787  gderiv(1:nn,1:nspace)=matmul( deriv(1:nn,1:nspace), xji(1:nspace,1:nspace) )
788  end subroutine
789 
791  real(kind=kreal) function getdeterminant( fetype, nn, localcoord, elecoord )
792  integer, intent(in) :: fetype
793  integer, intent(in) :: nn
794  real(kind=kreal), intent(in) :: localcoord(:)
795  real(kind=kreal), intent(in) :: elecoord(:,:)
796 
797  real(kind=kreal) :: xj(3,3), deriv(nn,3)
798  integer :: nspace
799 
800  nspace = getspacedimension( fetype )
801  call getshapederiv( fetype, localcoord(:), deriv(1:nn,:) )
802 
803  if( nspace==2 ) then
804  xj(1:2,1:2)=matmul( elecoord(1:2,1:nn), deriv(1:nn,1:2) )
805  getdeterminant=xj(1,1)*xj(2,2)-xj(2,1)*xj(1,2)
806  else
807  xj(1:3,1:3)= matmul( elecoord(1:3,1:nn), deriv(1:nn,1:3) )
808  getdeterminant=xj(1,1)*xj(2,2)*xj(3,3) &
809  +xj(2,1)*xj(3,2)*xj(1,3) &
810  +xj(3,1)*xj(1,2)*xj(2,3) &
811  -xj(3,1)*xj(2,2)*xj(1,3) &
812  -xj(2,1)*xj(1,2)*xj(3,3) &
813  -xj(1,1)*xj(3,2)*xj(2,3)
814  endif
815 
816  end function
817 
819  subroutine getjacobian( fetype, nn, localcoord, elecoord, det, jacobian, inverse )
820  integer, intent(in) :: fetype
821  integer, intent(in) :: nn
822  real(kind=kreal), intent(in) :: localcoord(:)
823  real(kind=kreal), intent(in) :: elecoord(:,:)
824  real(kind=kreal), intent(out) :: det
825  real(kind=kreal), intent(out) :: jacobian(:,:)
826  real(kind=kreal), intent(out) :: inverse(:,:)
827 
828  real(kind=kreal) :: dum, deriv(nn,3)
829  integer :: nspace
830 
831  nspace = getspacedimension( fetype )
832  call getshapederiv( fetype, localcoord(:), deriv(1:nn,:) )
833 
834  if( nspace==2 ) then
835  jacobian(1:2,1:2)=matmul( elecoord(1:2,1:nn), deriv(1:nn,1:2) )
836  det=jacobian(1,1)*jacobian(2,2)-jacobian(2,1)*jacobian(1,2)
837  if( det==0.d0 ) stop "Math error in getJacobain! Determinant==0.0"
838  dum=1.0/det
839  inverse(1,1)= jacobian(2,2)*dum
840  inverse(1,2)=-jacobian(1,2)*dum
841  inverse(2,1)=-jacobian(2,1)*dum
842  inverse(2,2)= jacobian(1,1)*dum
843  else
844  ! JACOBI MATRIX
845  jacobian(1:3,1:3)= matmul( elecoord(1:3,1:nn), deriv(1:nn,1:3) )
846  !DETERMINANT OF JACOBIAN
847  det=jacobian(1,1)*jacobian(2,2)*jacobian(3,3) &
848  +jacobian(2,1)*jacobian(3,2)*jacobian(1,3) &
849  +jacobian(3,1)*jacobian(1,2)*jacobian(2,3) &
850  -jacobian(3,1)*jacobian(2,2)*jacobian(1,3) &
851  -jacobian(2,1)*jacobian(1,2)*jacobian(3,3) &
852  -jacobian(1,1)*jacobian(3,2)*jacobian(2,3)
853  if( det==0.d0 ) stop "Math error in getJacobain! Determinant==0.0"
854  ! INVERSION OF JACOBIAN
855  dum=1.d0/det
856  inverse(1,1)=dum*( jacobian(2,2)*jacobian(3,3)-jacobian(3,2)*jacobian(2,3) )
857  inverse(1,2)=dum*(-jacobian(1,2)*jacobian(3,3)+jacobian(3,2)*jacobian(1,3) )
858  inverse(1,3)=dum*( jacobian(1,2)*jacobian(2,3)-jacobian(2,2)*jacobian(1,3) )
859  inverse(2,1)=dum*(-jacobian(2,1)*jacobian(3,3)+jacobian(3,1)*jacobian(2,3) )
860  inverse(2,2)=dum*( jacobian(1,1)*jacobian(3,3)-jacobian(3,1)*jacobian(1,3) )
861  inverse(2,3)=dum*(-jacobian(1,1)*jacobian(2,3)+jacobian(2,1)*jacobian(1,3) )
862  inverse(3,1)=dum*( jacobian(2,1)*jacobian(3,2)-jacobian(3,1)*jacobian(2,2) )
863  inverse(3,2)=dum*(-jacobian(1,1)*jacobian(3,2)+jacobian(3,1)*jacobian(1,2) )
864  inverse(3,3)=dum*( jacobian(1,1)*jacobian(2,2)-jacobian(2,1)*jacobian(1,2) )
865  endif
866  end subroutine
867 
869  function surfacenormal( fetype, nn, localcoord, elecoord ) result( normal )
870  integer, intent(in) :: fetype
871  integer, intent(in) :: nn
872  real(kind=kreal), intent(in) :: localcoord(2)
873  real(kind=kreal), intent(in) :: elecoord(3,nn)
874  real(kind=kreal) :: normal(3)
875  real(kind=kreal) :: deriv(nn,2), gderiv(3,2)
876 
877  select case (fetype)
878  case (fe_tri3n)
879  !error check
880  call shapederiv_tri3n(deriv(1:3,1:2))
881  case (fe_tri6n)
882  !error check
883  call shapederiv_tri6n(localcoord,deriv(1:6,1:2))
884  case (fe_quad4n)
885  !error check
886  call shapederiv_quad4n(localcoord,deriv(1:4,1:2))
887  case (fe_quad8n)
888  !error check
889  call shapederiv_quad8n(localcoord,deriv(1:8,1:2))
890  case default
891  ! error message
892  normal =0.d0
893  return
894  end select
895 
896  gderiv = matmul( elecoord, deriv )
897  normal(1) = gderiv(2,1)*gderiv(3,2) - gderiv(3,1)*gderiv(2,2)
898  normal(2) = gderiv(3,1)*gderiv(1,2) - gderiv(1,1)*gderiv(3,2)
899  normal(3) = gderiv(1,1)*gderiv(2,2) - gderiv(2,1)*gderiv(1,2)
900  ! normal = normal/dsqrt(dot_product(normal, normal))
901  end function
902 
904  function edgenormal( fetype, nn, localcoord, elecoord ) result( normal )
905  integer, intent(in) :: fetype
906  integer, intent(in) :: nn
907  real(kind=kreal), intent(in) :: localcoord(1)
908  real(kind=kreal), intent(in) :: elecoord(2,nn)
909  real(kind=kreal) :: normal(2)
910  real(kind=kreal) :: deriv(nn,1), gderiv(2,1)
911 
912  select case (fetype)
913  case (fe_line2n, fe_if_line2n)
914  !error check
915  call shapederiv_line2n(deriv(1:nn,:))
916  case (fe_line3n)
917  !error check
918  call shapederiv_line3n(localcoord,deriv(1:nn,:))
919  case default
920  ! error message
921  normal =0.d0
922  return
923  end select
924 
925  gderiv = matmul( elecoord, deriv )
926  normal(1) = -gderiv(2,1)
927  normal(2) = gderiv(1,1)
928  ! normal = normal/dsqrt(dot_product(normal, normal))
929  end function
930 
932  subroutine tangentbase( fetype, nn, localcoord, elecoord, tangent )
933  integer, intent(in) :: fetype
934  integer, intent(in) :: nn
935  real(kind=kreal), intent(in) :: localcoord(2)
936  real(kind=kreal), intent(in) :: elecoord(3,nn)
937  real(kind=kreal), intent(out) :: tangent(3,2)
938  real(kind=kreal) :: deriv(nn,2)
939 
940  select case (fetype)
941  case (fe_tri3n)
942  !error check
943  call shapederiv_tri3n(deriv(1:3,1:2))
944  case (fe_tri6n)
945  !error check
946  call shapederiv_tri6n(localcoord,deriv(1:6,1:2))
947  case (fe_tri6nc)
948  !error check
949  call shapederiv_tri6n(localcoord,deriv(1:6,1:2))
950  case (fe_quad4n)
951  !error check
952  call shapederiv_quad4n(localcoord,deriv(1:4,1:2))
953  case (fe_quad8n)
954  !error check
955  call shapederiv_quad8n(localcoord,deriv(1:8,1:2))
956  case default
957  ! error message
958  tangent =0.d0
959  return
960  end select
961 
962  tangent = matmul( elecoord, deriv )
963  end subroutine tangentbase
964 
966  subroutine curvature( fetype, nn, localcoord, elecoord, l2ndderiv, normal, curv )
967  integer, intent(in) :: fetype
968  integer, intent(in) :: nn
969  real(kind=kreal), intent(in) :: localcoord(2)
970  real(kind=kreal), intent(in) :: elecoord(3,nn)
971  real(kind=kreal), intent(out) :: l2ndderiv(3,2,2)
972  real(kind=kreal), intent(in), optional :: normal(3)
973  real(kind=kreal), intent(out), optional :: curv(2,2)
974  real(kind=kreal) :: deriv2(nn,2,2)
975 
976  select case (fetype)
977  case (fe_tri3n)
978  !error check
979  call shape2ndderiv_tri3n(deriv2(1:3,1:2,1:2))
980  case (fe_tri6n)
981  !error check
982  call shape2ndderiv_tri6n(deriv2(1:6,1:2,1:2))
983  case (fe_tri6nc)
984  !error check
985  call shape2ndderiv_tri6n(deriv2(1:6,1:2,1:2))
986  ! deriv2=0.d0
987  case (fe_quad4n)
988  !error check
989  call shape2ndderiv_quad4n(deriv2(1:4,1:2,1:2))
990  case (fe_quad8n)
991  !error check
992  call shape2ndderiv_quad8n(localcoord,deriv2(1:8,1:2,1:2))
993  case default
994  ! error message
995  stop "Cannot calculate second derivatives of shape function"
996  end select
997 
998  l2ndderiv(1:3,1,1) = matmul( elecoord(1:3,1:nn), deriv2(1:nn,1,1) )
999  l2ndderiv(1:3,1,2) = matmul( elecoord(1:3,1:nn), deriv2(1:nn,1,2) )
1000  l2ndderiv(1:3,2,1) = matmul( elecoord(1:3,1:nn), deriv2(1:nn,2,1) )
1001  l2ndderiv(1:3,2,2) = matmul( elecoord(1:3,1:nn), deriv2(1:nn,2,2) )
1002  if( present(curv) ) then
1003  curv(1,1) = dot_product( l2ndderiv(:,1,1), normal(:) )
1004  curv(1,2) = dot_product( l2ndderiv(:,1,2), normal(:) )
1005  curv(2,1) = dot_product( l2ndderiv(:,2,1), normal(:) )
1006  curv(2,2) = dot_product( l2ndderiv(:,2,2), normal(:) )
1007  endif
1008  end subroutine curvature
1009 
1011  subroutine getelementcenter( fetype, localcoord )
1012  integer, intent(in) :: fetype
1013  real(kind=kreal), intent(out) :: localcoord(2)
1014 
1015  select case (fetype)
1016  case (fe_tri3n, fe_tri6n, fe_tri6nc)
1017  localcoord(:) = 1.d0/3.d0
1018  case (fe_quad4n, fe_quad8n)
1019  localcoord(:) = 0.d0
1020  case default
1021  ! error message
1022  localcoord(:) = 0.d0
1023  end select
1024  end subroutine getelementcenter
1025 
1028  integer function isinsideelement( fetype, localcoord, clearance )
1029  integer, intent(in) :: fetype
1030  real(kind=kreal), intent(inout) :: localcoord(2)
1031  real(kind=kreal), optional :: clearance
1032  real(kind=kreal) :: clr, coord3
1033 
1034  clr = 1.d-6
1035  if( present(clearance) ) clr = clearance
1036  if( dabs(localcoord(1))<clr ) localcoord(1)=0.d0
1037  if( dabs(localcoord(2))<clr ) localcoord(2)=0.d0
1038  if( dabs(dabs(localcoord(1))-1.d0)<clr ) &
1039  localcoord(1)=sign(1.d0,localcoord(1))
1040  if( dabs(dabs(localcoord(2))-1.d0)<clr ) &
1041  localcoord(2)=sign(1.d0,localcoord(2))
1042  isinsideelement = -1
1043  select case (fetype)
1044  case (fe_tri3n, fe_tri6n, fe_tri6nc)
1045  !error check
1046  coord3 = 1.d0-(localcoord(1)+localcoord(2))
1047  if( dabs(coord3)<clr ) coord3=0.d0
1048  if( localcoord(1)>=0.d0 .and. localcoord(1)<=1.d0 .and. &
1049  localcoord(2)>=0.d0 .and. localcoord(2)<=1.d0 .and. &
1050  coord3>=0.d0 .and. coord3<=1.d0 ) then
1051  isinsideelement = 0
1052  if( localcoord(1)==1.d0 ) then
1053  isinsideelement = 1
1054  elseif( localcoord(2)==1.d0 ) then
1055  isinsideelement = 2
1056  elseif( coord3==1.d0 ) then
1057  isinsideelement = 3
1058  elseif( coord3==0.d0 ) then
1059  isinsideelement = 12
1060  elseif( localcoord(1)==0.d0 ) then
1061  isinsideelement = 23
1062  elseif( localcoord(2)==0.d0 ) then
1063  isinsideelement = 31
1064  endif
1065  endif
1066  case (fe_quad4n, fe_quad8n)
1067  !error check
1068  if( all(dabs(localcoord)<=1.d0) ) then
1069  isinsideelement = 0
1070  if( localcoord(1)==-1.d0 .and. localcoord(2)==-1.d0 ) then
1071  isinsideelement = 1
1072  elseif( localcoord(1)==1.d0 .and. localcoord(2)==-1.d0 ) then
1073  isinsideelement = 2
1074  elseif( localcoord(1)==1.d0 .and. localcoord(2)==1.d0 ) then
1075  isinsideelement = 3
1076  elseif( localcoord(1)==-1.d0 .and. localcoord(2)==1.d0 ) then
1077  isinsideelement = 4
1078  elseif( localcoord(2)==-1.d0 ) then
1079  isinsideelement = 12
1080  elseif( localcoord(1)==1.d0 ) then
1081  isinsideelement = 23
1082  elseif( localcoord(2)==1.d0 ) then
1083  isinsideelement = 34
1084  elseif( localcoord(1)==-1.d0 ) then
1085  isinsideelement = 41
1086  endif
1087  endif
1088  end select
1089  end function isinsideelement
1090 
1093  integer function isinside3delement( fetype, localcoord, clearance )
1094  integer, intent(in) :: fetype
1095  real(kind=kreal), intent(inout) :: localcoord(3)
1096  real(kind=kreal), optional :: clearance
1097  real(kind=kreal) :: clr, coord4
1098 
1099  integer :: idof
1100 
1101  clr = 1.d-6
1102  if( present(clearance) ) clr = clearance
1103  do idof=1,3
1104  if( dabs(localcoord(idof))<clr ) localcoord(idof)=0.d0
1105  if( dabs(dabs(localcoord(idof))-1.d0)<clr ) &
1106  & localcoord(idof)=sign(1.d0,localcoord(idof))
1107  enddo
1108 
1109  isinside3delement = -1
1110  select case (fetype)
1112  !error check
1113  coord4 = 1.d0-(localcoord(1)+localcoord(2)+localcoord(3))
1114  if( dabs(coord4)<clr ) coord4=0.d0
1115  isinside3delement = 0
1116  do idof=1,3
1117  if( localcoord(idof) < 0.d0 .or. localcoord(idof) > 1.d0 ) isinside3delement = -1
1118  enddo
1119  if( coord4 < 0.d0 .or. coord4 > 1.d0 ) isinside3delement = -1
1120  case (fe_prism6n, fe_prism15n)
1121  !error check
1122  coord4 = 1.d0-(localcoord(1)+localcoord(2))
1123  isinside3delement = 0
1124  do idof=1,2
1125  if( localcoord(idof) < 0.d0 .or. localcoord(idof) > 1.d0 ) isinside3delement = -1
1126  enddo
1127  if( localcoord(3) < -1.d0 .or. localcoord(3) > 1.d0 ) isinside3delement = -1
1128  if( coord4 < 0.d0 .or. coord4 > 1.d0 ) isinside3delement = -1
1129  case (fe_hex8n, fe_hex20n, fe_hex27n)
1130  if( all(dabs(localcoord)<=1.d0) ) isinside3delement = 0
1131  end select
1132  end function
1133 
1135  subroutine getvertexcoord( fetype, cnode, localcoord )
1136  integer, intent(in) :: fetype
1137  integer, intent(in) :: cnode
1138  real(kind=kreal), intent(out) :: localcoord(2)
1139 
1140  select case (fetype)
1141  case (fe_tri3n, fe_tri6n, fe_tri6nc)
1142  if( cnode==1 ) then
1143  localcoord(1) =1.d0
1144  localcoord(2) =0.d0
1145  elseif( cnode==2 ) then
1146  localcoord(1) =0.d0
1147  localcoord(2) =1.d0
1148  else
1149  localcoord(1) =0.d0
1150  localcoord(2) =0.d0
1151  endif
1152  case (fe_quad4n, fe_quad8n)
1153  if( cnode==1 ) then
1154  localcoord(1) =-1.d0
1155  localcoord(2) =-1.d0
1156  elseif( cnode==2 ) then
1157  localcoord(1) =1.d0
1158  localcoord(2) =-1.d0
1159  elseif( cnode==3 ) then
1160  localcoord(1) =1.d0
1161  localcoord(2) =1.d0
1162  else
1163  localcoord(1) =-1.d0
1164  localcoord(2) =1.d0
1165  endif
1166  end select
1167  end subroutine
1168 
1170  subroutine extrapolatevalue( lpos, fetype, nnode, pvalue, ndvalue )
1171  real(kind=kreal), intent(in) :: lpos(:)
1172  integer, intent(in) :: fetype
1173  integer, intent(in) :: nnode
1174  real(kind=kreal), intent(in) :: pvalue(:)
1175  real(kind=kreal), intent(out) :: ndvalue(:,:)
1176 
1177  integer :: i
1178  real(kind=kreal) :: shapefunc(nnode)
1179  call getshapefunc( fetype, lpos, shapefunc )
1180  do i=1,nnode
1181  ndvalue(i,:) = shapefunc(i)*pvalue(:)
1182  enddo
1183  end subroutine
1184 
1186  subroutine interapolatevalue( lpos, fetype, nnode, pvalue, ndvalue )
1187  real(kind=kreal), intent(in) :: lpos(:)
1188  integer, intent(in) :: fetype
1189  integer, intent(in) :: nnode
1190  real(kind=kreal), intent(out) :: pvalue(:)
1191  real(kind=kreal), intent(in) :: ndvalue(:,:)
1192 
1193  integer :: i
1194  real(kind=kreal) :: shapefunc(nnode)
1195  call getshapefunc( fetype, lpos, shapefunc )
1196  pvalue(:) = 0
1197  do i=1,nnode
1198  pvalue(:) = pvalue(:)+ shapefunc(i)*ndvalue(i,:)
1199  enddo
1200  end subroutine
1201 
1203  subroutine gauss2node( fetype, gaussv, nodev )
1204  integer, intent(in) :: fetype
1205  real(kind=kreal), intent(in) :: gaussv(:,:)
1206  real(kind=kreal), intent(out) :: nodev(:,:)
1207 
1208  integer :: i, ngauss, nnode
1209  real(kind=kreal) :: localcoord(3), func(100)
1210  ngauss = numofquadpoints( fetype )
1211  nnode = getnumberofnodes( fetype )
1212  ! error checking
1213  select case (fetype)
1214  case (fe_tri3n)
1215  !error check
1216  do i=1,nnode
1217  nodev(i,:) = gaussv(1,:)
1218  enddo
1219  case (fe_tri6n)
1220  !error check
1221  ! func(1:6) = ShapeFunc_tri6n(localcoord)
1222  case (fe_quad4n)
1223  !error check
1224  ! nodev(:,:) = gaussv(1,:)
1225  case (fe_quad8n)
1226  !error check
1227  call shapefunc_quad8n(localcoord,func(1:8))
1228  case (fe_hex8n, fe_mitc4_shell361)
1229  ! error check
1230  call shapefunc_hex8n(localcoord,func(1:8))
1231  case (fe_hex20n)
1232  ! error check
1233  call shapefunc_hex20n(localcoord,func(1:20))
1235  do i=1,3
1236  nodev(i,:) = gaussv(1,:)
1237  enddo
1238  do i=1,3
1239  nodev(i+3,:) = gaussv(2,:)
1240  enddo
1241  case (fe_prism15n)
1242  call shapefunc_prism15n(localcoord,func(1:15))
1244  ! error check
1245  do i=1,nnode
1246  nodev(i,:) = gaussv(1,:)
1247  enddo
1248  case (fe_tet10n)
1249  ! error check
1250  call shapefunc_tet10n(localcoord,func(1:10))
1251  case default
1252  stop "Element type not defined"
1253  ! error message
1254  end select
1255  end subroutine
1256 
1258  real(kind=kreal) function getreferencelength( fetype, nn, localcoord, elecoord )
1259  integer, intent(in) :: fetype
1260  integer, intent(in) :: nn
1261  real(kind=kreal),intent(in) :: localcoord(2)
1262  real(kind=kreal),intent(in) :: elecoord(3,nn)
1263  real(kind=kreal) :: detjxy, detjyz, detjxz, detj
1264  detjxy = getdeterminant( fetype, nn, localcoord, elecoord(1:2,1:nn) )
1265  detjyz = getdeterminant( fetype, nn, localcoord, elecoord(2:3,1:nn) )
1266  detjxz = getdeterminant( fetype, nn, localcoord, elecoord(1:3:2,1:nn) )
1267  detj = dsqrt( detjxy **2 + detjyz **2 + detjxz **2 )
1268  getreferencelength = dsqrt( detj )
1269  end function getreferencelength
1270 
1271 
1272 end module
elementinfo::fe_edgesmooth_tet4n
integer, parameter fe_edgesmooth_tet4n
Definition: element.f90:102
shape_prism15n::shapederiv_prism15n
subroutine shapederiv_prism15n(ncoord, func)
Definition: prism15n.f90:37
quadrature::gauss3d4
real(kind=kreal), dimension(3, 1) gauss3d4
Definition: quadrature.f90:32
elementinfo::getshapederiv
subroutine getshapederiv(fetype, localcoord, shapederiv)
Calculate derivatives of shape function in natural coordinate system.
Definition: element.f90:578
elementinfo::fe_hex8n
integer, parameter fe_hex8n
Definition: element.f90:81
quadrature::gauss3d5
real(kind=kreal), dimension(3, 4) gauss3d5
Definition: quadrature.f90:32
shape_hex20n::shapefunc_hex20n
subroutine shapefunc_hex20n(localcoord, func)
Definition: hex20n.f90:12
elementinfo::fe_quad4n
integer, parameter fe_quad4n
Definition: element.f90:72
elementinfo::fe_tet10nc
integer, parameter fe_tet10nc
Definition: element.f90:78
elementinfo::fe_tet4n_pipi
integer, parameter fe_tet4n_pipi
Definition: element.f90:76
quadrature::weight1d1
real(kind=kreal), dimension(1) weight1d1
Definition: quadrature.f90:32
elementinfo::getweight
real(kind=kreal) function getweight(fetype, np)
Fetch the weight value in given gauss point.
Definition: element.f90:535
elementinfo::curvature
subroutine curvature(fetype, nn, localcoord, elecoord, l2ndderiv, normal, curv)
Calculate curvature tensor at a point along 3d surface.
Definition: element.f90:967
elementinfo::numofquadpoints
integer function numofquadpoints(fetype)
Obtains the number of quadrature points of the element.
Definition: element.f90:450
shape_quad9n::nodalnaturalcoord_quad9n
subroutine nodalnaturalcoord_quad9n(nncoord)
Definition: quad9n.f90:174
quadrature::gauss1d1
real(kind=kreal), dimension(1, 1) gauss1d1
Definition: quadrature.f90:32
elementinfo::gauss2node
subroutine gauss2node(fetype, gaussv, nodev)
This subroutine extroplate value in quadrature point to element nodes.
Definition: element.f90:1204
elementinfo::fe_nodesmooth_tet4n
integer, parameter fe_nodesmooth_tet4n
Definition: element.f90:101
quadrature::weight2d4
real(kind=kreal), dimension(1) weight2d4
Definition: quadrature.f90:32
quadrature::gauss1d2
real(kind=kreal), dimension(1, 2) gauss1d2
Definition: quadrature.f90:32
shape_tri6n::shapederiv_tri6n
subroutine shapederiv_tri6n(areacoord, func)
Definition: tri6n.f90:26
elementinfo::getjacobian
subroutine getjacobian(fetype, nn, localcoord, elecoord, det, jacobian, inverse)
calculate Jacobian matrix, its determinant and inverse
Definition: element.f90:820
quadrature::weight2d2
real(kind=kreal), dimension(4) weight2d2
Definition: quadrature.f90:32
elementinfo::fe_tet4n
integer, parameter fe_tet4n
Definition: element.f90:75
shape_line2n
This module contains functions for interpolation in 2 node line element (Langrange interpolation)
Definition: line2n.f90:7
elementinfo::fe_line2n
integer, parameter fe_line2n
Definition: element.f90:67
shape_hex8n
This module contains functions for interpolation in 8 node hexahedral element (Langrange interpolatio...
Definition: hex8n.f90:7
elementinfo::getnumberofsubface
integer(kind=kind(2)) function getnumberofsubface(etype)
Obtain number of sub-surface.
Definition: element.f90:168
elementinfo::fe_prism15n
integer, parameter fe_prism15n
Definition: element.f90:80
elementinfo::fe_tri3n
integer, parameter fe_tri3n
Definition: element.f90:69
shape_tet10n
This module contains functions for interpolation in 10 node tetrahedron element (Langrange interpolat...
Definition: tet10n.f90:7
shape_tet10n::shapederiv_tet10n
subroutine shapederiv_tet10n(volcoord, shp)
Definition: tet10n.f90:30
shape_line3n::shapederiv_line3n
subroutine shapederiv_line3n(lcoord, func)
Definition: line3n.f90:20
elementinfo::fe_truss
integer, parameter fe_truss
Definition: element.f90:74
elementinfo::fe_quad8n
integer, parameter fe_quad8n
Definition: element.f90:73
elementinfo::surfacenormal
real(kind=kreal) function, dimension(3) surfacenormal(fetype, nn, localcoord, elecoord)
Calculate normal of 3d-surface.
Definition: element.f90:870
quadrature::gauss3d7
real(kind=kreal), dimension(3, 2) gauss3d7
Definition: quadrature.f90:32
shape_quad9n
This module contains functions for interpolation in 9 node quadrilateral element.
Definition: quad9n.f90:7
quadrature::gauss2d6
real(kind=kreal), dimension(2, 4) gauss2d6
Definition: quadrature.f90:32
shape_tri3n::shape2ndderiv_tri3n
subroutine shape2ndderiv_tri3n(func)
Definition: tri3n.f90:30
shape_tri6n::shape2ndderiv_tri6n
subroutine shape2ndderiv_tri6n(func)
Definition: tri6n.f90:48
shape_quad8n::shape2ndderiv_quad8n
subroutine shape2ndderiv_quad8n(lcoord, func)
Definition: quad8n.f90:56
shape_prism15n
This module contains functions for interpolation in 15 node prism element (Langrange interpolation)
Definition: prism15n.f90:7
quadrature::weight3d8
real(kind=kreal), dimension(9) weight3d8
Definition: quadrature.f90:32
shape_hex20n
This module contains functions for interpolation in 20 node hexahedral element (Serendipity interpola...
Definition: hex20n.f90:7
elementinfo::fe_tet10n
integer, parameter fe_tet10n
Definition: element.f90:77
elementinfo::fe_mitc3_shell
integer, parameter fe_mitc3_shell
Definition: element.f90:93
elementinfo::extrapolatevalue
subroutine extrapolatevalue(lpos, fetype, nnode, pvalue, ndvalue)
This subroutine extrapolate a point value into elemental nodes.
Definition: element.f90:1171
shape_quad8n::shapederiv_quad8n
subroutine shapederiv_quad8n(lcoord, func)
Definition: quad8n.f90:33
shape_hex8n::shapefunc_hex8n
subroutine shapefunc_hex8n(localcoord, func)
Definition: hex8n.f90:12
elementinfo::fe_tri6nc
integer, parameter fe_tri6nc
Definition: element.f90:71
shape_quad4n::shape2ndderiv_quad4n
subroutine shape2ndderiv_quad4n(func)
Definition: quad4n.f90:35
quadrature::gauss3d8
real(kind=kreal), dimension(3, 9) gauss3d8
Definition: quadrature.f90:32
elementinfo::getshape2ndderiv
subroutine getshape2ndderiv(fetype, localcoord, shapederiv)
Calculate the 2nd derivative of shape function in natural coordinate system.
Definition: element.f90:622
quadrature::weight3d2
real(kind=kreal), dimension(8) weight3d2
Definition: quadrature.f90:32
quadrature::weight3d4
real(kind=kreal), dimension(1) weight3d4
Definition: quadrature.f90:32
elementinfo::isinside3delement
integer function isinside3delement(fetype, localcoord, clearance)
if a point is inside a surface element -1: No; 0: Yes; >0: Node's (vertex) number
Definition: element.f90:1094
shape_prism6n
This module contains functions for interpolation in 6 node prism element (Langrange interpolation)
Definition: prism6n.f90:7
elementinfo::fe_beam2n
integer, parameter fe_beam2n
Definition: element.f90:87
elementinfo::fe_tri3n_patch
integer, parameter fe_tri3n_patch
Definition: element.f90:104
quadrature::gauss2d4
real(kind=kreal), dimension(2, 1) gauss2d4
Definition: quadrature.f90:32
shape_tet10n::shapefunc_tet10n
subroutine shapefunc_tet10n(volcoord, shp)
Definition: tet10n.f90:12
shape_line2n::shapefunc_line2n
subroutine shapefunc_line2n(lcoord, func)
Definition: line2n.f90:12
elementinfo::fe_unknown
integer, parameter fe_unknown
Definition: element.f90:65
elementinfo::tangentbase
subroutine tangentbase(fetype, nn, localcoord, elecoord, tangent)
Calculate base vector of tangent space of 3d surface.
Definition: element.f90:933
elementinfo::edgenormal
real(kind=kreal) function, dimension(2) edgenormal(fetype, nn, localcoord, elecoord)
Calculate normal of 2d-edge.
Definition: element.f90:905
shape_tri6n::shapefunc_tri6n
subroutine shapefunc_tri6n(areacoord, func)
Definition: tri6n.f90:12
elementinfo::fe_beam341
integer, parameter fe_beam341
Definition: element.f90:89
shape_prism6n::shapederiv_prism6n
subroutine shapederiv_prism6n(ncoord, func)
Definition: prism6n.f90:27
quadrature::gauss2d2
real(kind=kreal), dimension(2, 4) gauss2d2
Definition: quadrature.f90:32
elementinfo::fe_hex20n
integer, parameter fe_hex20n
Definition: element.f90:82
elementinfo::getglobalderiv
subroutine getglobalderiv(fetype, nn, localcoord, elecoord, det, gderiv)
Calculate shape derivative in global coordinate system.
Definition: element.f90:741
quadrature::weight2d5
real(kind=kreal), dimension(3) weight2d5
Definition: quadrature.f90:32
elementinfo::fe_mitc9_shell
integer, parameter fe_mitc9_shell
Definition: element.f90:96
elementinfo::getnodalnaturalcoord
subroutine getnodalnaturalcoord(fetype, nncoord)
Definition: element.f90:699
shape_tri3n::shapederiv_tri3n
subroutine shapederiv_tri3n(func)
Definition: tri3n.f90:19
quadrature::weight3d7
real(kind=kreal), dimension(2) weight3d7
Definition: quadrature.f90:32
shape_hex8n::shapederiv_hex8n
subroutine shapederiv_hex8n(localcoord, func)
Definition: hex8n.f90:25
shape_tet4n
This module contains functions for interpolation in 4 node tetrahedron element (Langrange interpolati...
Definition: tet4n.f90:7
quadrature::gauss3d2
real(kind=kreal), dimension(3, 8) gauss3d2
Definition: quadrature.f90:32
elementinfo::fe_prism6n
integer, parameter fe_prism6n
Definition: element.f90:79
elementinfo::fe_tri6n_patch
integer, parameter fe_tri6n_patch
Definition: element.f90:105
elementinfo::fe_mitc4_shell
integer, parameter fe_mitc4_shell
Definition: element.f90:94
elementinfo::getreferencelength
real(kind=kreal) function getreferencelength(fetype, nn, localcoord, elecoord)
This function calculates reference length at a point in surface.
Definition: element.f90:1259
elementinfo::fe_quad4n_patch
integer, parameter fe_quad4n_patch
Definition: element.f90:106
shape_line3n::shapefunc_line3n
subroutine shapefunc_line3n(lcoord, func)
Definition: line3n.f90:12
elementinfo::getdeterminant
real(kind=kreal) function getdeterminant(fetype, nn, localcoord, elecoord)
Calculate shape derivative in global coordinate system.
Definition: element.f90:792
quadrature::weight2d3
real(kind=kreal), dimension(9) weight2d3
Definition: quadrature.f90:32
elementinfo
This module encapsulate the basic functions of all elements provide by this software.
Definition: element.f90:43
shape_quad8n::shapefunc_quad8n
subroutine shapefunc_quad8n(lcoord, func)
Definition: quad8n.f90:14
elementinfo::isinsideelement
integer function isinsideelement(fetype, localcoord, clearance)
if a point is inside a surface element -1: No; 0: Yes; >0: Node's (vertex) number
Definition: element.f90:1029
shape_quad9n::shapederiv_quad9n
subroutine shapederiv_quad9n(lcoord, func)
Definition: quad9n.f90:91
elementinfo::getquadpoint
subroutine getquadpoint(fetype, np, pos)
Fetch the coordinate of gauss point.
Definition: element.f90:489
elementinfo::fe_hex27n
integer, parameter fe_hex27n
Definition: element.f90:83
elementinfo::getshapefunc
subroutine getshapefunc(fetype, localcoord, func)
Calculate the shape function in natural coordinate system.
Definition: element.f90:647
shape_quad4n::shapefunc_quad4n
subroutine shapefunc_quad4n(lcoord, func)
Definition: quad4n.f90:12
elementinfo::fe_line3n
integer, parameter fe_line3n
Definition: element.f90:68
quadrature::weight1d2
real(kind=kreal), dimension(2) weight1d2
Definition: quadrature.f90:32
elementinfo::fe_if_line2n
integer, parameter fe_if_line2n
Definition: element.f90:85
elementinfo::fe_beam3n
integer, parameter fe_beam3n
Definition: element.f90:88
shape_quad8n
This module contains functions for interpolation in 8 node quadrilateral element (Serendipity interpo...
Definition: quad8n.f90:7
shape_tri3n
This module contains functions for interpolation in 3 node trianglar element (Langrange interpolation...
Definition: tri3n.f90:7
elementinfo::fe_tri6n_shell
integer, parameter fe_tri6n_shell
Definition: element.f90:91
elementinfo::getsubface
subroutine getsubface(intype, innumber, outtype, nodes)
Find the definition of surface of the element.
Definition: element.f90:193
elementinfo::interapolatevalue
subroutine interapolatevalue(lpos, fetype, nnode, pvalue, ndvalue)
This subroutine interapolate element nodes value into a point value.
Definition: element.f90:1187
shape_tri3n::nodalnaturalcoord_tri3n
subroutine nodalnaturalcoord_tri3n(nncoord)
Definition: tri3n.f90:38
shape_line3n
This module contains functions for interpolation in 3 nodes line element (Langrange interpolation)
Definition: line3n.f90:7
quadrature::gauss2d5
real(kind=kreal), dimension(2, 3) gauss2d5
Definition: quadrature.f90:32
elementinfo::getspacedimension
integer(kind=kind(2)) function getspacedimension(etype)
Obtain the space dimension of the element.
Definition: element.f90:117
elementinfo::getvertexcoord
subroutine getvertexcoord(fetype, cnode, localcoord)
Get the natural coord of a vertex node.
Definition: element.f90:1136
shape_hex20n::shapederiv_hex20n
subroutine shapederiv_hex20n(localcoord, func)
Definition: hex20n.f90:41
shape_tet4n::shapederiv_tet4n
subroutine shapederiv_tet4n(func)
Definition: tet4n.f90:19
elementinfo::fe_mitc8_shell
integer, parameter fe_mitc8_shell
Definition: element.f90:95
shape_prism15n::shapefunc_prism15n
subroutine shapefunc_prism15n(ncoord, shp)
Definition: prism15n.f90:13
quadrature::weight3d5
real(kind=kreal), dimension(4) weight3d5
Definition: quadrature.f90:32
shape_tri3n::shapefunc_tri3n
subroutine shapefunc_tri3n(areacoord, func)
Definition: tri3n.f90:12
elementinfo::fe_quad8n_patch
integer, parameter fe_quad8n_patch
Definition: element.f90:107
shape_line2n::shapederiv_line2n
subroutine shapederiv_line2n(func)
Definition: line2n.f90:19
quadrature::gauss2d3
real(kind=kreal), dimension(2, 9) gauss2d3
Definition: quadrature.f90:32
shape_quad4n::shapederiv_quad4n
subroutine shapederiv_quad4n(lcoord, func)
Definition: quad4n.f90:21
shape_quad4n
This module contains functions for interpolation in 4 node qudrilateral element (Langrange interpolat...
Definition: quad4n.f90:7
shape_quad4n::nodalnaturalcoord_quad4n
subroutine nodalnaturalcoord_quad4n(nncoord)
Definition: quad4n.f90:53
shape_tet4n::shapefunc_tet4n
subroutine shapefunc_tet4n(volcoord, func)
Definition: tet4n.f90:12
quadrature
This module contains Gauss point information.
Definition: quadrature.f90:28
shape_prism6n::shapefunc_prism6n
subroutine shapefunc_prism6n(ncoord, func)
Definition: prism6n.f90:12
quadrature::weight3d3
real(kind=kreal), dimension(27) weight3d3
Definition: quadrature.f90:32
elementinfo::fe_tri6n
integer, parameter fe_tri6n
Definition: element.f90:70
elementinfo::getelementcenter
subroutine getelementcenter(fetype, localcoord)
Return natural coordinate of the center of surface element.
Definition: element.f90:1012
elementinfo::getnumberofnodes
integer(kind=kind(2)) function getnumberofnodes(etype)
Obtain number of nodes of the element.
Definition: element.f90:131
elementinfo::fe_mitc3_shell361
integer, parameter fe_mitc3_shell361
Definition: element.f90:98
shape_quad9n::shapefunc_quad9n
subroutine shapefunc_quad9n(lcoord, func)
Definition: quad9n.f90:23
elementinfo::fe_dsg3_shell
integer, parameter fe_dsg3_shell
Definition: element.f90:92
quadrature::gauss3d3
real(kind=kreal), dimension(3, 27) gauss3d3
Definition: quadrature.f90:32
shape_tri6n
This module contains functions for interpolation in 6 node trianglar element (Langrange interpolation...
Definition: tri6n.f90:7
elementinfo::fe_mitc4_shell361
integer, parameter fe_mitc4_shell361
Definition: element.f90:99