FE-Project
Loading...
Searching...
No Matches
Functions/Subroutines
scale_meshutil_cubedsphere2d Module Reference

module FElib / Mesh / utility for 2D cubed-sphere mesh More...

Functions/Subroutines

subroutine, public meshutilcubedsphere2d_buildglobalmap (panelid_table, pi_table, pj_table, tileid_map, tilefaceid_map, tilepanelid_map, ntile)
 
subroutine, public meshutilcubedsphere2d_modifyconnectivity (tilepanelid_map, tileid_map, tilefaceid_map, panelid_table, pi_table, pj_table, nex, ney, ntile, nface)
 
subroutine, public meshutilcubedsphere2d_getpanelconnectivity (panel_connectivity, face_connectivity)
 
subroutine, public meshutilcubedsphere2d_getpanelid (panelid, lon, lat, np)
 

Detailed Description

module FElib / Mesh / utility for 2D cubed-sphere mesh

Description
A module useful for generating 2D cubed-sphere mesh
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ meshutilcubedsphere2d_buildglobalmap()

subroutine, public scale_meshutil_cubedsphere2d::meshutilcubedsphere2d_buildglobalmap ( integer, dimension(ntile), intent(out) panelid_table,
integer, dimension(ntile), intent(out) pi_table,
integer, dimension(ntile), intent(out) pj_table,
integer, dimension(4,ntile), intent(out) tileid_map,
integer, dimension(4,ntile), intent(out) tilefaceid_map,
integer, dimension(4,ntile), intent(out) tilepanelid_map,
integer, intent(in) ntile )

Definition at line 48 of file scale_meshutil_cubedsphere2d.F90.

52
53 use scale_meshutil_2d, only: &
55 implicit none
56
57 integer, intent(in) :: Ntile
58 integer, intent(out) :: panelID_table(Ntile)
59 integer, intent(out) :: pi_table(Ntile)
60 integer, intent(out) :: pj_table(Ntile)
61 integer, intent(out) :: tileID_map(4,Ntile)
62 integer, intent(out) :: tileFaceID_map(4,Ntile)
63 integer, intent(out) :: tilePanelID_map(4,Ntile)
64
65 integer :: NtilePerPanel
66 integer :: NeX, NeY, NvX, NvY
67 integer, allocatable :: nodesID_2d(:,:,:)
68 integer, allocatable :: EToV(:,:)
69 integer, allocatable :: EToE(:,:)
70 integer, allocatable :: EToF(:,:)
71 integer :: i, j, k, f
72 integer :: panelID
73 integer :: tileID, tileID_R
74 integer :: counter
75 !-----------------------------------------------------------------------------
76
77 ntileperpanel = ntile / 6
78 ney = int( sqrt(dble(ntileperpanel)) )
79 nex = ntileperpanel/ney
80 nvx = nex + 1
81 nvy = ney + 1
82 allocate( nodesid_2d(nvx,nvy,6) )
83 allocate( etov(ntile,4), etoe(ntile,4), etof(ntile,4) )
84
85 counter = 0
86 do panelid = 1, 6
87 do j = 1, nvy
88 do i = 1, nvx
89 counter = counter + 1
90 nodesid_2d(i,j,panelid) = counter
91 end do
92 end do
93 end do
94
95 !----
96
97 tileid = 0
98 do panelid = 1, 6
99 do j = 1, ney
100 do i = 1, nex
101 tileid = tileid + 1
102 panelid_table(tileid) = panelid
103 pi_table(tileid) = i; pj_table(tileid) = j;
104 etov(tileid,:) = (/ nodesid_2d(i,j ,panelid), nodesid_2d(i+1,j ,panelid), &
105 nodesid_2d(i,j+1,panelid), nodesid_2d(i+1,j+1,panelid) /)
106 end do
107 end do
108 end do
109
110 call meshutil2d_genconnectivity( etoe, etof, &
111 etov, ntile, 4 )
112 tileid_map(:,:) = transpose(etoe)
113 tilefaceid_map(:,:) = transpose(etof)
114
115 do tileid=1, ntile
116 do f=1, 4
117 tileid_r = tileid_map(f,tileid)
118 tilepanelid_map(f,tileid) = panelid_table(tileid_r)
119 end do
120 end do
121
122 call meshutilcubedsphere2d_modifyconnectivity( &
123 tilepanelid_map, tileid_map, tilefaceid_map, & ! (inout)
124 panelid_table, pi_table, pj_table, nex, ney, ntile, 4 ) ! (in)
125
126 return
module FElib / Mesh / utility for 2D mesh
subroutine, public meshutil2d_genconnectivity(etoe, etof, etov, ne, nfaces)

References scale_meshutil_2d::meshutil2d_genconnectivity(), and meshutilcubedsphere2d_modifyconnectivity().

Referenced by scale_mesh_cubedspheredom2d::meshcubedspheredom2d_setuplocaldom().

◆ meshutilcubedsphere2d_modifyconnectivity()

subroutine, public scale_meshutil_cubedsphere2d::meshutilcubedsphere2d_modifyconnectivity ( integer, dimension(nface,ntile), intent(out) tilepanelid_map,
integer, dimension(nface,ntile), intent(out) tileid_map,
integer, dimension(nface,ntile), intent(out) tilefaceid_map,
integer, dimension(ntile), intent(in) panelid_table,
integer, dimension(ntile), intent(in) pi_table,
integer, dimension(ntile), intent(in) pj_table,
integer, intent(in) nex,
integer, intent(in) ney,
integer, intent(in) ntile,
integer, intent(in) nface )

Definition at line 131 of file scale_meshutil_cubedsphere2d.F90.

133
134 integer, intent(in) :: Ntile
135 integer, intent(in) :: Nface
136 integer, intent(out) :: tileID_map(Nface,Ntile)
137 integer, intent(out) :: tileFaceID_map(Nface,Ntile)
138 integer, intent(out) :: tilePanelID_map(Nface,Ntile)
139 integer, intent(in) :: panelID_table(Ntile)
140 integer, intent(in) :: pi_table(Ntile)
141 integer, intent(in) :: pj_table(Ntile)
142 integer, intent(in) :: NeX, NeY
143
144 integer :: panel_connectivity(4,6)
145 integer :: face_connectivity (4,6)
146
147 integer :: tileID
148 integer :: panelID
149 integer :: f
150 integer :: pi_, pj_
151 !-----------------------------------------------------------------------------
152
153 call meshutilcubedsphere2d_getpanelconnectivity( &
154 panel_connectivity, face_connectivity )
155
156
157 do tileid=1, ntile
158 panelid = panelid_table(tileid)
159
160 do f=1, 4
161 if ( tilefaceid_map(f,tileid) /= f ) cycle ! Does the face correspond the boundary of panel of cubed sphere?
162
163 pi_ = pi_table(tileid)
164 pj_ = pj_table(tileid)
165
166 select case( panelid )
167 case ( 1 )
168 if ( mod(f,2) == 0 ) then ! West / East
169 if ( pi_table(tileid) == 1 ) pi_ = nex
170 if ( pi_table(tileid) == nex ) pi_ = 1
171 else ! North / South
172 if ( pj_table(tileid) == 1 ) pj_ = ney
173 if ( pj_table(tileid) == ney ) pj_ = 1
174 end if
175 case ( 2 )
176 if ( mod(f,2) == 0 ) then ! West / East
177 if ( pi_table(tileid) == 1 ) pi_ = nex
178 if ( pi_table(tileid) == nex ) pi_ = 1
179 else ! North / South
180 if ( pj_table(tileid) == 1 ) then
181 pi_ = nex; pj_ = ney - pi_table(tileid) + 1
182 end if
183 if ( pj_table(tileid) == ney ) then
184 pi_ = nex ; pj_ = pi_table(tileid)
185 end if
186 end if
187 case ( 3 )
188 if ( mod(f,2) == 0 ) then ! West / East
189 if ( pi_table(tileid) == 1 ) pi_ = nex
190 if ( pi_table(tileid) == nex ) pi_ = 1
191 else ! North / South
192 if ( pj_table(tileid) == 1 ) then
193 pi_ = nex - pi_table(tileid) + 1; pj_ = 1
194 end if
195 if ( pj_table(tileid) == ney ) then
196 pi_ = nex - pi_table(tileid) + 1; pj_ = ney
197 end if
198 end if
199 case ( 4 )
200 if ( mod(f,2) == 0 ) then ! West / East
201 if ( pi_table(tileid) == 1 ) pi_ = nex
202 if ( pi_table(tileid) == nex ) pi_ = 1
203 else ! North / South
204 if ( pj_table(tileid) == 1 ) then
205 pi_ = 1; pj_ = pi_table(tileid)
206 end if
207 if ( pj_table(tileid) == ney ) then
208 pi_ = 1; pj_ = ney - pi_table(tileid) + 1;
209 end if
210 end if
211 case ( 5 )
212 pj_ = ney
213 if ( mod(f,2) == 0 ) then
214 if ( pi_table(tileid) == 1 ) pi_ = ney - pj_table(tileid) + 1 ! West
215 if ( pi_table(tileid) == nex ) pi_ = pj_table(tileid) ! East
216 else
217 if ( pj_table(tileid) == 1 ) pi_ = pi_table(tileid) ! South
218 if ( pj_table(tileid) == ney ) pi_ = nex - pi_table(tileid) + 1 ! North
219 end if
220 case ( 6 )
221 pj_ = 1
222 if ( mod(f,2) == 0 ) then
223 if ( pi_table(tileid) == 1 ) pi_ = pj_table(tileid) ! West
224 if ( pi_table(tileid) == nex ) pi_ = ney - pj_table(tileid) + 1 ! East
225 else
226 if ( pj_table(tileid) == 1 ) pi_ = nex - pi_table(tileid) + 1 ! South
227 if ( pj_table(tileid) == ney ) pi_ = pi_table(tileid) ! North
228 end if
229 end select
230
231 tilepanelid_map(f,tileid) = panel_connectivity(f,panelid)
232 tileid_map(f,tileid) = pi_ + (pj_ - 1) * nex + (tilepanelid_map(f,tileid) - 1) * nex * ney
233 tilefaceid_map(f,tileid) = face_connectivity(f,panelid)
234 end do ! loop for f
235 end do ! loop for tile
236
237 ! if (PRC_myrank==0) then
238 ! write(*,*) " MeshUtilCubedSphere2D_modifyConnectivity"
239 ! do tileID=1, Ntile
240 ! write(*,'(a,i3,a,2i3)') "tileID=", tileID, ":", pi_table(tileID), pj_table(tileID)
241 ! write(*,'(a,6i4)') "map_tile:", tileID_map(:,tileID)
242 ! write(*,'(a,6i4)') "map_panel:", tilePanelID_map(:,tileID)
243 ! write(*,'(a,6i4)') "mac_face:", tileFaceID_map(:,tileID)
244 ! end do
245 ! end if
246
247 return

References meshutilcubedsphere2d_getpanelconnectivity().

Referenced by meshutilcubedsphere2d_buildglobalmap(), and scale_meshutil_cubedsphere3d::meshutilcubedsphere3d_buildglobalmap().

◆ meshutilcubedsphere2d_getpanelconnectivity()

subroutine, public scale_meshutil_cubedsphere2d::meshutilcubedsphere2d_getpanelconnectivity ( integer, dimension(4,6), intent(out) panel_connectivity,
integer, dimension(4,6), intent(out) face_connectivity )

Definition at line 251 of file scale_meshutil_cubedsphere2d.F90.

252
253 implicit none
254
255 integer, intent(out) :: panel_connectivity(4,6)
256 integer, intent(out) :: face_connectivity(4,6)
257
258 integer :: n
259 integer :: zonal_panelID_list(0:5)
260
261 zonal_panelid_list(:) = (/ 4, 1, 2, 3, 4, 1/)
262 do n=1, 4
263 panel_connectivity(1,n) = 6
264 if (zonal_panelid_list(4-n) > 2) then
265 face_connectivity(1,n) = +zonal_panelid_list(4-n)
266 else
267 face_connectivity(1,n) = -zonal_panelid_list(4-n)
268 end if
269
270 panel_connectivity(2,n) = zonal_panelid_list(n+1)
271 face_connectivity(2,n) = 4
272
273 panel_connectivity(3,n) = 5
274 face_connectivity(3,n) = n
275 if (n > 2) then
276 face_connectivity(3,n) = -n
277 else
278 face_connectivity(3,n) = n
279 end if
280
281 panel_connectivity(4,n) = zonal_panelid_list(n-1)
282 face_connectivity(4,n) = 2
283 end do
284
285 panel_connectivity(:,5) = (/ 1, 2, 3, 4 /)
286 face_connectivity(:,5) = (/ 3, 3, -3, -3 /)
287 panel_connectivity(:,6) = (/ 3, 2, 1, 4 /)
288 face_connectivity(:,6) = (/ -1, -1, 1, 1 /)
289

Referenced by meshutilcubedsphere2d_modifyconnectivity().

◆ meshutilcubedsphere2d_getpanelid()

subroutine, public scale_meshutil_cubedsphere2d::meshutilcubedsphere2d_getpanelid ( integer, dimension(np), intent(out) panelid,
real(rp), dimension(np), intent(in) lon,
real(rp), dimension(np), intent(in) lat,
integer, intent(in) np )

Definition at line 293 of file scale_meshutil_cubedsphere2d.F90.

294 use scale_cubedsphere_coord_cnv, only: &
296 implicit none
297
298 integer, intent(in) :: Np
299 integer, intent(out) :: panelID(Np)
300 real(RP), intent(in) :: lon(Np)
301 real(RP), intent(in) :: lat(Np)
302
303 integer :: p
304 integer :: pnl
305 real(RP) :: alph(Np), beta(Np)
306 real(RP) :: lon_(Np), lat_(Np)
307
308 real(RP), parameter :: EPS = 1.0e-64_rp
309 !------------------------------------------
310
311 !$omp parallel do
312 do p=1, np
313 if ( abs(cos(lon(p))) < eps ) then
314 lon_(p) = lon(p) + eps
315 else
316 lon_(p) = lon(p)
317 end if
318 if ( abs( lat(p) - 0.5_rp * pi ) < eps ) then
319 lat_(p) = lat(p) - sign(eps, lat(p))
320 else
321 lat_(p) = lat(p)
322 end if
323 end do
324
325 panelid(:) = -1
326 do pnl=1, 6
327 call cubedspherecoordcnv_lonlat2cspos( pnl, lon_, lat_, np, &
328 alph, beta )
329
330 select case(pnl)
331 case (5)
332 where ( lat(:) > 0.0_rp .and. abs(alph(:)) <= 0.25_rp * pi .and. abs(beta(:)) <= 0.25_rp * pi )
333 panelid(:) = pnl
334 end where
335 case (6)
336 where ( lat(:) < 0.0_rp .and. abs(alph(:)) <= 0.25_rp * pi .and. abs(beta(:)) <= 0.25_rp * pi )
337 panelid(:) = pnl
338 end where
339 case default
340 where ( abs(alph(:)) <= 0.25_rp * pi .and. abs(beta(:)) <= 0.25_rp * pi )
341 panelid(:) = pnl
342 end where
343 end select
344 end do
345
346 do p=1, np
347 if (panelid(p) < 0) then
348 log_error("MeshUtilCubedSphere2D_getPanelID",*) 'Fail to search a panel ID of cubed sphere grid!'
349 write(*,*) "p=", p, ": (lon,lat)=", lon(p), lat(p), "(alpha,beta)=", alph(p), beta(p)
350 call prc_abort
351 end if
352 end do
353
354 return
module common / Coordinate conversion with a cubed-sphere
subroutine, public cubedspherecoordcnv_lonlat2cspos(panelid, lon, lat, np, alpha, beta)

References scale_cubedsphere_coord_cnv::cubedspherecoordcnv_lonlat2cspos().