FE-Project
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines
scale_mesh_cubedspheredom3d Module Reference

module FElib / Mesh / Cubed-sphere 3D domain More...

Data Types

type  meshcubedspheredom3d
 

Functions/Subroutines

subroutine meshcubedspheredom3d_init (this, negx, negy, negz, rplanet, dom_zmin, dom_zmax, refelem, nlocalmeshperprc, nproc, myrank, fz, shallow_approx)
 

Detailed Description

module FElib / Mesh / Cubed-sphere 3D domain

Description
Manage mesh data of cubed-sphere 3D domain for element-based methods
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ meshcubedspheredom3d_init()

subroutine scale_mesh_cubedspheredom3d::meshcubedspheredom3d_init ( class(meshcubedspheredom3d), intent(inout) this,
integer, intent(in) negx,
integer, intent(in) negy,
integer, intent(in) negz,
real(rp), intent(in) rplanet,
real(rp), intent(in) dom_zmin,
real(rp), intent(in) dom_zmax,
type(hexahedralelement), intent(in), target refelem,
integer, intent(in) nlocalmeshperprc,
integer, intent(in), optional nproc,
integer, intent(in), optional myrank,
real(rp), dimension(negz+1), intent(in), optional fz,
logical, intent(in), optional shallow_approx )

Definition at line 96 of file scale_mesh_cubedspheredom3d.F90.

102
103 use scale_const, only: &
104 pi => const_pi
105 implicit none
106
107 class(MeshCubedSphereDom3D), intent(inout) :: this
108 integer, intent(in) :: NeGX
109 integer, intent(in) :: NeGY
110 integer, intent(in) :: NeGZ
111 real(RP), intent(in) :: RPlanet
112 real(RP), intent(in) :: dom_Zmin
113 real(RP), intent(in) :: dom_zmax
114 type(HexahedralElement), intent(in), target :: refElem
115 integer, intent(in) :: NLocalMeshPerPrc
116 integer, intent(in), optional :: nproc
117 integer, intent(in), optional :: myrank
118 real(RP), intent(in), optional :: FZ(NeGZ+1)
119 logical, intent(in), optional :: shallow_approx
120
121 integer :: k
122 real(RP) :: dz
123 !-----------------------------------------------------------------------------
124
125 this%NeGX = negx
126 this%NeGY = negy
127 this%NeGZ = negz
128
129 this%xmin_gl = - 0.25_rp * pi
130 this%xmax_gl = + 0.25_rp * pi
131 this%ymin_gl = - 0.25_rp * pi
132 this%ymax_gl = + 0.25_rp * pi
133 this%zmin_gl = dom_zmin
134 this%zmax_gl = dom_zmax
135 this%RPlanet = rplanet
136 this%dom_vol = 4.0_rp / 3.0_rp * pi * ( ( dom_zmax + rplanet )**3 - ( dom_zmin + rplanet )**3 )
137
138
139 !- Fz
140 allocate( this%FZ(this%NeGZ+1) )
141 if ( present(fz) ) then
142 this%FZ(:) = fz(:)
143 else
144 this%FZ(1 ) = dom_zmin
145 this%FZ(this%NeGZ+1) = dom_zmax
146 dz = (dom_zmax - dom_zmin) / dble(this%NeGZ)
147 do k=2, this%NeGZ
148 this%FZ(k) = this%FZ(k-1) + dz
149 end do
150 end if
151
152 !-
153 if ( present(shallow_approx) ) then
154 this%shallow_approx = shallow_approx
155 else
156 this%shallow_approx = .true.
157 end if
158
159 !--
160 call meshbase3d_init( this, refelem, nlocalmeshperprc, 6, &
161 nproc, myrank )
162
163 !---
164 call this%refElem2D%Init( this%refElem3D%PolyOrder_h, refelem%IsLumpedMatrix() )
165 call this%mesh2D%Init( negx, negy, rplanet, this%refElem2D, nlocalmeshperprc, &
166 nproc, myrank )
167
168 !-- Modify the information of dimension for the cubed sphere mesh
169 call this%SetDimInfo( meshbase3d_dimtypeid_x, "x", "1", "X-coordinate" )
170 call this%SetDimInfo( meshbase3d_dimtypeid_y, "y", "1", "Y-coordinate" )
171 call this%SetDimInfo( meshbase3d_dimtypeid_z, "z", "m", "Z-coordinate" )
172 call this%SetDimInfo( meshbase3d_dimtypeid_xyz, "xyz", "1", "XYZ-coordinate" )
173 call this%SetDimInfo( meshbase3d_dimtypeid_zt, "zt", "1", "XYZ-coordinate" )
174 call this%SetDimInfo( meshbase3d_dimtypeid_xyzt, "xyzt", "1", "XYZ-coordinate" )
175
176 return

References scale_localmesh_base::bctype_interior, scale_cubedsphere_coord_cnv::cubedspherecoordcnv_cs2lonlatpos(), scale_cubedsphere_coord_cnv::cubedspherecoordcnv_getmetric(), scale_mesh_base3d::meshbase3d_dimtypeid_x, scale_mesh_base3d::meshbase3d_dimtypeid_xyz, scale_mesh_base3d::meshbase3d_dimtypeid_xyzt, scale_mesh_base3d::meshbase3d_dimtypeid_y, scale_mesh_base3d::meshbase3d_dimtypeid_z, scale_mesh_base3d::meshbase3d_dimtypeid_zt, scale_mesh_base3d::meshbase3d_final(), scale_mesh_base3d::meshbase3d_init(), scale_mesh_base3d::meshbase3d_setgeometricinfo(), scale_mesh_cubedspheredom2d::meshcubedspheredom2d_check_division_params(), scale_mesh_cubedspheredom2d::meshcubedspheredom2d_setuplocaldom(), and scale_meshutil_cubedsphere3d::meshutilcubedsphere3d_buildglobalmap().