SPEED
GET_HIGHEST_NODE.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine get_highest_node (nn_loc, ne_loc, zz_loc, loc_n_num, nnz
 Computes the highest node of the mesh (z-dir)
 

Function/Subroutine Documentation

◆ get_highest_node()

subroutine get_highest_node ( integer*4  nn_loc,
integer*4  ne_loc,
real*8, dimension(nn_loc)  zz_loc,
integer*4, dimension(nn_loc)  loc_n_num,
  nnz 
)

Computes the highest node of the mesh (z-dir)

Author
Ilario Mazzieri
Date
September, 2013
Version
1.0
Parameters
[in]nn_locnumber of local nodes
[in]ne_locnumber of local element
[in]zz_locz-coordinate for the nodes
[in]local_n_numvector for local node numbering
[in]nnz_loclength of cs_loc
[in]cs_loclocal spectral connectivity
[in]nmnumber of materials
[in]sdpolynomial degree vector
[in]tmlabels for material vector

Definition at line 34 of file GET_HIGHEST_NODE.f90.

35 nm, tm, sd, highest)
36
37
38 implicit none
39
40 integer*4 :: nn_loc, nnz_loc, ne_loc
41 integer*4 :: ie, nn, nm, im
42 integer*4 :: n1,n2,n3,n4,n5,n6,n7,n8
43 integer*4 :: ic1,ic2,ic3,ic4,ic5,ic6,ic7,ic8
44
45 integer*4, dimension(nn_loc) :: loc_n_num
46 integer*4, dimension(nm) :: tm
47 integer*4, dimension(nm) :: sd
48
49 integer*4, dimension(0:nnz_loc) :: cs_loc
50
51 real*8 :: zz1,zz2,zz3,zz4,zz5,zz6,zz7,zz8
52
53 real*8, dimension(:), allocatable :: ct,ww
54 real*8, dimension(nn_loc) :: zz_loc
55 real*8, dimension(ne_loc) :: highest
56
57 real*8, dimension(:,:), allocatable :: dd
58
59
60 nn = 2
61 allocate(ct(nn),ww(nn),dd(nn,nn))
62 call make_lgl_nw(nn,ct,ww,dd)
63
64 do im = 1,nm
65 if ((sd(im) +1).ne.nn) then
66 deallocate(ct,ww,dd)
67 nn = sd(im) +1
68 allocate(ct(nn),ww(nn),dd(nn,nn))
69 call make_lgl_nw(nn,ct,ww,dd)
70 endif
71
72 do ie = 1,ne_loc
73 if (cs_loc(cs_loc(ie -1) +0).eq.tm(im)) then
74
75 n1 = nn*nn*(1 -1) +nn*(1 -1) +1
76 n2 = nn*nn*(1 -1) +nn*(1 -1) +nn
77 n3 = nn*nn*(1 -1) +nn*(nn -1) +nn
78 n4 = nn*nn*(1 -1) +nn*(nn -1) +1
79 n5 = nn*nn*(nn -1) +nn*(1 -1) +1
80 n6 = nn*nn*(nn -1) +nn*(1 -1) +nn
81 n7 = nn*nn*(nn -1) +nn*(nn -1) +nn
82 n8 = nn*nn*(nn -1) +nn*(nn -1) +1
83
84 ic1 = cs_loc(cs_loc(ie -1) +n1)
85 ic2 = cs_loc(cs_loc(ie -1) +n2)
86 ic3 = cs_loc(cs_loc(ie -1) +n3)
87 ic4 = cs_loc(cs_loc(ie -1) +n4)
88 ic5 = cs_loc(cs_loc(ie -1) +n5)
89 ic6 = cs_loc(cs_loc(ie -1) +n6)
90 ic7 = cs_loc(cs_loc(ie -1) +n7)
91 ic8 = cs_loc(cs_loc(ie -1) +n8)
92
93 zz1 = zz_loc(ic1)
94 zz2 = zz_loc(ic2)
95 zz3 = zz_loc(ic3)
96 zz4 = zz_loc(ic4)
97 zz5 = zz_loc(ic5)
98 zz6 = zz_loc(ic6)
99 zz7 = zz_loc(ic7)
100 zz8 = zz_loc(ic8)
101
102 highest(ie) = max(zz1,zz2,zz3,zz4,zz5,zz6,zz7,zz8)
103
104 endif
105 enddo
106 enddo
107
108
109 deallocate(ct,ww,dd)
110
111 return
112
subroutine make_lgl_nw(nb_pnt, xq, wq, dd)
Makes Gauss-Legendre-Lobatto nodes, weigths and spectral derivatives.

References make_lgl_nw().

Referenced by find_monitor_position(), and read_system_position().

Here is the call graph for this function:
Here is the caller graph for this function: