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

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

Functions/Subroutines

subroutine, public meshutilcubedsphere3d_buildglobalmap (panelid_table, pi_table, pj_table, pk_table, tileid_map, tilefaceid_map, tilepanelid_map, ntile, nez)
 

Detailed Description

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

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

Function/Subroutine Documentation

◆ meshutilcubedsphere3d_buildglobalmap()

subroutine, public scale_meshutil_cubedsphere3d::meshutilcubedsphere3d_buildglobalmap ( integer, dimension(ntile), intent(out) panelid_table,
integer, dimension(ntile), intent(out) pi_table,
integer, dimension(ntile), intent(out) pj_table,
integer, dimension(ntile), intent(out) pk_table,
integer, dimension(6,ntile), intent(out) tileid_map,
integer, dimension(6,ntile), intent(out) tilefaceid_map,
integer, dimension(6,ntile), intent(out) tilepanelid_map,
integer, intent(in) ntile,
integer, intent(in) nez )

Definition at line 45 of file scale_meshutil_cubedsphere3d.F90.

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

References scale_meshutil_3d::meshutil3d_genconnectivity(), and scale_meshutil_cubedsphere2d::meshutilcubedsphere2d_modifyconnectivity().

Referenced by scale_mesh_cubedspheredom3d::meshcubedspheredom3d_init().