SPEED
MAKE_SPX_NODES.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine make_spx_nodes (nb_node, spx_con_nnz, spx_con, node_wgt, nnz
 Makes a pointer for spectral connectivity vector.
 

Function/Subroutine Documentation

◆ make_spx_nodes()

subroutine make_spx_nodes ( integer*4  nb_node,
integer*4  spx_con_nnz,
integer*4, dimension(0:spx_con_nnz)  spx_con,
integer*4, dimension(nb_node)  node_wgt,
integer*4  nnz 
)

Makes a pointer for spectral connectivity vector.

Author
Ilario Mazzieri
Date
September, 2013
Version
1.0
Parameters
[in]nb_nodeGLL nodes
[in]spx_conspectral connectivity vector
[in]node_wgtnode_wgt(i) contains the multiplicity of the GLL node i (e.g. if the GLL node i is shared by 4 elements node_wgt(i) = 4
[in]nnznumber of GLL nodes + number of GLL nodes including repetitions + 1
[out]node_pointeras explained in MAKE_GRID_NODES.f90

Definition at line 30 of file MAKE_SPX_NODES.f90.

31
32
33
34 implicit none
35
36 integer*4 :: nb_node,spx_con_nnz,nnz
37 integer*4 :: i,j,ie,nb_elem,nn3
38
39 integer*4, dimension(:), allocatable :: i4count
40 integer*4, dimension(0:spx_con_nnz) :: spx_con
41 integer*4, dimension(nb_node) :: node_wgt
42 integer*4, dimension(0:nnz) :: node_pointer
43
44 allocate(i4count(nb_node))
45
46 node_pointer(0) = nb_node + 1
47 do j = 1, nb_node
48 node_pointer(j) = node_pointer(j -1) + node_wgt(j)
49 enddo
50
51 do j = 1, nb_node
52 i4count(j) = node_pointer(j -1)
53 enddo
54
55 nb_elem = spx_con(0) - 1
56 do ie = 1, nb_elem
57
58 nn3 = spx_con(ie) - spx_con(ie -1) -1
59 do i = 1, nn3
60
61 j = spx_con(spx_con(ie -1) +i)
62 node_pointer(i4count(j)) = ie
63 i4count(j) = i4count(j) + 1
64 enddo
65 enddo
66
67 deallocate(i4count)
68
69 return
70