SPEED
MAKE_SPX_LOC_NUMERATION.f90
Go to the documentation of this file.
1! Copyright (C) 2012 The SPEED FOUNDATION
2! Author: Ilario Mazzieri
3!
4! This file is part of SPEED.
5!
6! SPEED is free software; you can redistribute it and/or modify it
7! under the terms of the GNU Affero General Public License as
8! published by the Free Software Foundation, either version 3 of the
9! License, or (at your option) any later version.
10!
11! SPEED is distributed in the hope that it will be useful, but
12! WITHOUT ANY WARRANTY; without even the implied warranty of
13! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14! Affero General Public License for more details.
15!
16! You should have received a copy of the GNU Affero General Public License
17! along with SPEED. If not, see <http://www.gnu.org/licenses/>.
18
19
31
32
33 subroutine make_spx_loc_numeration(nnloc, loc_n_num, &
34 cs_nnz, cs_loc, &
35 cs_nnz_bc, cs_bc_loc)
36
37
38 implicit none
39
40 integer*4 :: cs_nnz, cs_nnz_bc, nnloc
41 integer*4, dimension(0:cs_nnz), intent(inout) :: cs_loc
42 integer*4, dimension(0:cs_nnz_bc), intent(inout) :: cs_bc_loc
43 integer*4, dimension(nnloc), intent(inout) :: loc_n_num
44
45 integer*4 :: ne_loc, ne_loc_bc
46 integer*4 :: ie, j, iglo, iloc
47
48
49 ne_loc = cs_loc(0) - 1
50
51 do ie = 1, ne_loc
52 do j = cs_loc(ie -1) + 1, cs_loc(ie) - 1
53
54 iglo = cs_loc(j)
55 call get_indloc_from_indglo(loc_n_num, nnloc, iglo, iloc)
56 cs_loc(j) = iloc
57
58 if(iloc .eq. 0) write(*,*) '1 Error in MAKE_SPX_LOC_NUMERATION'
59 enddo
60 enddo
61
62 ne_loc_bc = cs_bc_loc(0) - 1
63
64 do ie = 1, ne_loc_bc
65 do j = cs_bc_loc(ie -1) + 1, cs_bc_loc(ie) - 1
66
67 iglo = cs_bc_loc(j)
68 call get_indloc_from_indglo(loc_n_num, nnloc, iglo, iloc)
69 cs_bc_loc(j) = iloc
70
71 if(iloc .eq. 0) write(*,*) '2 Error in MAKE_SPX_LOC_NUMERATION'
72 enddo
73 enddo
74
75
76
77 end subroutine make_spx_loc_numeration
78
subroutine get_indloc_from_indglo(local_el, nel_loc, ie, ic)
Returns local id from global id.
subroutine make_spx_loc_numeration(nnloc, loc_n_num, cs_nnz, cs_loc, cs_nnz_bc, cs_bc_loc)
Makes local numeration of spectral nodes.