FE-Project
Loading...
Searching...
No Matches
Functions/Subroutines | Variables
scale_meshutil_vcoord Module Reference

module FElib / Mesh / utility for general vertical coordinate More...

Functions/Subroutines

integer function, public meshutil_get_vcoord_typeid (vcoord_type)
 
subroutine, public meshutil_vcoord_getmetric (g13, g23, zlev, gsqrtv, topo, ztop, vcoord_id, lcmesh, elem, lcmesh2d, elem2d, dx2d, dy2d, lift2d)
 

Variables

character(*), parameter, public mesh_vcoord_terrain_following_name = "TERRAIN_FOLLOWING"
 
integer, parameter, public mesh_vcoord_terrain_following_id = 1
 

Detailed Description

module FElib / Mesh / utility for general vertical coordinate

Description
A module useful for general vertical coordinate
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ meshutil_get_vcoord_typeid()

integer function, public scale_meshutil_vcoord::meshutil_get_vcoord_typeid ( character(len=*), intent(in) vcoord_type)

Definition at line 61 of file scale_meshutil_vcoord.F90.

62 implicit none
63
64 character(len=*), intent(in) :: vcoord_type
65 integer :: vcoord_id
66
67 select case( vcoord_type )
68 case( mesh_vcoord_terrain_following_name )
69 vcoord_id = mesh_vcoord_terrain_following_id
70 case default
71 log_error("MeshUtil_VCoord_TypeID",*) "vcoord_type is inappropriate. Check!", vcoord_type
72 call prc_abort
73 end select
74
75 return

References mesh_vcoord_terrain_following_id, and mesh_vcoord_terrain_following_name.

Referenced by mod_atmos_mesh_gm::atmosmeshgm_init(), and mod_atmos_mesh_rm::atmosmeshrm_init().

◆ meshutil_vcoord_getmetric()

subroutine, public scale_meshutil_vcoord::meshutil_vcoord_getmetric ( real(rp), dimension(elem%np,lcmesh%nea), intent(out) g13,
real(rp), dimension(elem%np,lcmesh%nea), intent(out) g23,
real(rp), dimension(elem%np,lcmesh%nea), intent(out) zlev,
real(rp), dimension(elem%np,lcmesh%nea), intent(inout) gsqrtv,
real(rp), dimension(elem2d%np,lcmesh2d%nea), intent(in) topo,
real(rp), intent(in) ztop,
integer, intent(in) vcoord_id,
class(localmesh3d), intent(in) lcmesh,
class(elementbase3d), intent(in) elem,
class(localmesh2d), intent(in) lcmesh2d,
class(elementbase2d), intent(in) elem2d,
type(sparsemat), intent(in) dx2d,
type(sparsemat), intent(in) dy2d,
type(sparsemat), intent(in) lift2d )

Definition at line 79 of file scale_meshutil_vcoord.F90.

82 implicit none
83
84 class(LocalMesh3D), intent(in) :: lcmesh
85 class(ElementBase3D), intent(in) :: elem
86 class(LocalMesh2D), intent(in) :: lcmesh2D
87 class(ElementBase2D), intent(in) :: elem2D
88 real(RP), intent(out) :: G13(elem%Np,lcmesh%NeA)
89 real(RP), intent(out) :: G23(elem%Np,lcmesh%NeA)
90 real(RP), intent(out) :: zlev(elem%Np,lcmesh%NeA)
91 real(RP), intent(inout) :: GsqrtV(elem%Np,lcmesh%NeA)
92 real(RP), intent(in) :: topo(elem2D%Np,lcmesh2D%NeA)
93 integer, intent(in) :: vcoord_id
94 real(RP), intent(in) :: zTop
95 type(SparseMat), intent(in) :: Dx2D
96 type(SparseMat), intent(in) :: Dy2D
97 type(SparseMat), intent(in) :: Lift2D
98
99 integer :: ke, ke2D
100 real(RP) :: del_flux(elem2D%NfpTot,lcmesh2D%Ne,2)
101 real(RP) :: Fx2D(elem2D%Np), Fy2D(elem2D%Np), LiftDelFlux2D(elem2D%Np)
102 real(RP) :: GradZs(elem2D%Np,lcmesh2D%Ne,2)
103 real(RP) :: coef3D(elem%NP)
104 !------------------------------------------------
105
106 if ( vcoord_id == mesh_vcoord_terrain_following_id ) then
107
108 ! * z = topo + (1 - topo / zTop ) * zeta
109 ! * zeta = zTop * (z - topo)/(zTop - topo)
110 ! * Gi3 = (dzeta(x1,x2,z)/dxi)_z = d (zeta,z) / d (xi,z) = - d(z,zeta)/ d(xi,zeta) * d(xi,zeta)/d(xi,z)
111 ! = - (dz/dxi)_zeta * dzeta/dz
112 ! = (GsqrtV)^-1 * [ - 1 + zeta / zTop ] * d topo /dxi (i=1, 2)
113
114 call cal_del_flux( del_flux, &
115 topo, lcmesh2d%normal_fn(:,:,1), lcmesh2d%normal_fn(:,:,2), &
116 lcmesh2d%VMapM, lcmesh2d%VMapP, lcmesh2d, elem2d )
117
118 !$omp parallel private(ke2D, ke, &
119 !$omp Fx2D, Fy2D, LiftDelFlux2D, coef3D )
120
121 !$omp do
122 do ke2d=1, lcmesh2d%Ne
123 call sparsemat_matmul( dx2d, topo(:,ke2d), fx2d )
124 call sparsemat_matmul( lift2d, lcmesh2d%Fscale(:,ke2d) * del_flux(:,ke2d,1), liftdelflux2d)
125 gradzs(:,ke2d,1) = lcmesh2d%Escale(:,ke2d,1,1) * fx2d(:) + liftdelflux2d(:)
126
127 call sparsemat_matmul( dy2d, topo(:,ke2d), fy2d )
128 call sparsemat_matmul( lift2d, lcmesh2d%Fscale(:,ke2d) * del_flux(:,ke2d,2), liftdelflux2d)
129 gradzs(:,ke2d,2) = lcmesh2d%Escale(:,ke2d,2,2) * fy2d(:) + liftdelflux2d(:)
130 end do
131 !$omp end do
132
133 !$omp do
134 do ke=1, lcmesh%Ne
135 ke2d = lcmesh%EMap3Dto2D(ke)
136 coef3d(:) = 1.0_rp - lcmesh%pos_en(:,ke,3) / ztop
137 zlev(:,ke) = lcmesh%pos_en(:,ke,3) &
138 + coef3d(:) * topo(elem%IndexH2Dto3D,ke2d)
139
140 gsqrtv(:,ke) = 1.0_rp - topo(elem%IndexH2Dto3D,ke2d) / ztop ! dz/dzeta
141 coef3d(:) = - coef3d(:) / gsqrtv(:,ke)
142 g13(:,ke) = coef3d(:) * gradzs(elem%IndexH2Dto3D(:),ke2d,1)
143 g23(:,ke) = coef3d(:) * gradzs(elem%IndexH2Dto3D(:),ke2d,2)
144 end do
145 !$omp end do
146 !$omp end parallel
147 else
148 log_error("Mesh_VCoord_GetMetric",*) "vcoord_id is inappropriate. Check!", vcoord_id
149 call prc_abort
150 end if
151
152 return

References mesh_vcoord_terrain_following_id.

Referenced by scale_mesh_topography::meshtopography_init().

Variable Documentation

◆ mesh_vcoord_terrain_following_name

character(*), parameter, public scale_meshutil_vcoord::mesh_vcoord_terrain_following_name = "TERRAIN_FOLLOWING"

Definition at line 44 of file scale_meshutil_vcoord.F90.

44 character(*), public, parameter :: MESH_VCOORD_TERRAIN_FOLLOWING_NAME = "TERRAIN_FOLLOWING"

Referenced by meshutil_get_vcoord_typeid().

◆ mesh_vcoord_terrain_following_id

integer, parameter, public scale_meshutil_vcoord::mesh_vcoord_terrain_following_id = 1

Definition at line 45 of file scale_meshutil_vcoord.F90.

45 integer, public, parameter :: MESH_VCOORD_TERRAIN_FOLLOWING_ID = 1

Referenced by meshutil_get_vcoord_typeid(), and meshutil_vcoord_getmetric().