SPEED
MAKE_SPX_CON_LOC_BOUND.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
37
38
39 subroutine make_spx_con_loc_bound(cs_nnz,cs,ne_bc,cm_bc,nm,tm,sd,&
40 Ennz,Ebin,cs_nnz_bc,cs_bc,loc_el_num, nel_loc, &
41 nf_loc)
42
43 implicit none
44
45 integer*4 :: cs_nnz,ne_bc,nm,Ennz,cs_nnz_bc,nel_loc,nf_loc
46 integer*4 :: im,ibc,ielem,nn, icur, iloc
47 integer*4 :: i,j,k,an,bn,cn,ia,ib,ic,ja,jb,jc,ka,kb,kc
48 integer*4 :: i1,i2,j1,j2,k1,k2,is,js
49
50 integer*4, dimension(0:cs_nnz) :: cs
51 integer*4, dimension(nm) :: tm
52 integer*4, dimension(nm) :: sd
53 integer*4, dimension(0:Ennz) :: Ebin
54 integer*4, dimension(0:cs_nnz_bc) :: cs_bc
55 integer*4, dimension(nel_loc) :: loc_el_num
56
57 integer*4, dimension(ne_bc,5) :: cm_bc
58
59
60 iloc = 0
61
62 do ibc = 1,ne_bc
63 an = cm_bc(ibc,1 +1)
64 bn = cm_bc(ibc,2 +1)
65 cn = cm_bc(ibc,4 +1)
66
67 call get_elem_from_face(ennz,ebin,an,bn,cn,ielem)
68 call get_indloc_from_indglo(loc_el_num, nel_loc, ielem, icur)
69
70
71 if (icur .ne. 0) then
72 iloc = iloc + 1
73 do im = 1,nm
74 if (tm(im).eq.cs(cs(icur -1))) then
75 nn = sd(im) +1
76 cs_bc(iloc) = cs_bc(iloc -1) +nn*nn +1
77 cs_bc(cs_bc(iloc -1)) = cm_bc(ibc,1)
78 endif
79 enddo
80 endif
81 enddo
82
83
84
85 iloc = 0
86 do ibc = 1,ne_bc
87 an = cm_bc(ibc,1 +1)
88 bn = cm_bc(ibc,2 +1)
89 cn = cm_bc(ibc,4 +1)
90
91
92 call get_elem_from_face(ennz,ebin,an,bn,cn,ielem)
93 call get_indloc_from_indglo(loc_el_num, nel_loc, ielem, icur)
94
95
96 if(icur.ne.0) then
97 iloc = iloc + 1
98 do im = 1,nm
99 if (tm(im).eq.cs(cs(icur -1) +0)) nn = sd(im) +1
100 enddo
101
102 ia = 0
103 ja = 0
104 ka = 0
105 ib = 0
106 jb = 0
107 kb = 0
108 ic = 0
109 jc = 0
110 kc = 0
111
112 do k = 1,nn
113 do j = 1,nn
114 do i = 1,nn
115 js = nn*nn*(k -1) +nn*(j -1) +i
116
117 if (cs(cs(icur -1) +js).eq.an) then
118 ia = i
119 ja = j
120 ka = k
121 endif
122
123 if (cs(cs(icur -1) +js).eq.bn) then
124 ib = i
125 jb = j
126 kb = k
127 endif
128
129 if (cs(cs(icur -1) +js).eq.cn) then
130 ic = i
131 jc = j
132 kc = k
133 endif
134 enddo
135 enddo
136 enddo
137
138 i1 = min(ia,ib,ic)
139 i2 = max(ia,ib,ic)
140 j1 = min(ja,jb,jc)
141 j2 = max(ja,jb,jc)
142 k1 = min(ka,kb,kc)
143 k2 = max(ka,kb,kc)
144
145 is = 1
146 do k = k1,k2
147 do j = j1,j2
148 do i = i1,i2
149 js = nn*nn*(k -1) +nn*(j -1) +i
150
151 cs_bc(cs_bc(iloc -1) +is) = cs(cs(icur -1) +js)
152 is = is +1
153 enddo
154 enddo
155 enddo
156 endif
157
158 enddo
159
160 return
161
162 end subroutine make_spx_con_loc_bound
163
subroutine get_elem_from_face(nb_nz_el, list_el, v1, v2, v3, ie)
Find element index from verteces number.
subroutine get_indloc_from_indglo(local_el, nel_loc, ie, ic)
Returns local id from global id.
subroutine make_spx_con_loc_bound(cs_nnz, cs, ne_bc, cm_bc, nm, tm, sd, ennz, ebin, cs_nnz_bc, cs_bc, loc_el
Makes spectral connectivity vector for the boundary.