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

module FElib / Fluid dyn solver / Atmosphere / Global shallow water More...

Functions/Subroutines

subroutine, public atm_dyn_dgm_globalsw_init (mesh)
 
subroutine, public atm_dyn_dgm_globalsw_final ()
 
subroutine, public atm_dyn_dgm_globalsw_cal_tend (h_dt, u_dt, v_dt, h_, u_, v_, hs_, u1_, u2_, coriolis, dx, dy, sx, sy, lift, lmesh, elem)
 

Detailed Description

module FElib / Fluid dyn solver / Atmosphere / Global shallow water

Description
DGM scheme for Atmospheric dynamical process.
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ atm_dyn_dgm_globalsw_init()

subroutine, public scale_atm_dyn_dgm_globalsw::atm_dyn_dgm_globalsw_init ( class(meshbase2d), intent(in) mesh)

Definition at line 63 of file scale_atm_dyn_dgm_globalsw.F90.

64 implicit none
65 class(MeshBase2D), intent(in) :: mesh
66 !--------------------------------------------
67 return

◆ atm_dyn_dgm_globalsw_final()

subroutine, public scale_atm_dyn_dgm_globalsw::atm_dyn_dgm_globalsw_final

Definition at line 71 of file scale_atm_dyn_dgm_globalsw.F90.

72 implicit none
73 !--------------------------------------------
74 return

◆ atm_dyn_dgm_globalsw_cal_tend()

subroutine, public scale_atm_dyn_dgm_globalsw::atm_dyn_dgm_globalsw_cal_tend ( real(rp), dimension(elem%np,lmesh%nea), intent(out) h_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) u_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(out) v_dt,
real(rp), dimension(elem%np,lmesh%nea), intent(in) h_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) u_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) v_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) hs_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) u1_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) u2_,
real(rp), dimension(elem%np,lmesh%nea), intent(in) coriolis,
type(sparsemat), intent(in) dx,
type(sparsemat), intent(in) dy,
type(sparsemat), intent(in) sx,
type(sparsemat), intent(in) sy,
type(sparsemat), intent(in) lift,
class(localmesh2d), intent(in) lmesh,
class(elementbase2d), intent(in) elem )

Definition at line 80 of file scale_atm_dyn_dgm_globalsw.F90.

84
85 implicit none
86
87 class(LocalMesh2D), intent(in) :: lmesh
88 class(ElementBase2D), intent(in) :: elem
89 type(SparseMat), intent(in) :: Dx, Dy, Sx, Sy, Lift
90 real(RP), intent(out) :: h_dt(elem%Np,lmesh%NeA)
91 real(RP), intent(out) :: U_dt(elem%Np,lmesh%NeA)
92 real(RP), intent(out) :: V_dt(elem%Np,lmesh%NeA)
93 real(RP), intent(in) :: h_(elem%Np,lmesh%NeA)
94 real(RP), intent(in) :: U_(elem%Np,lmesh%NeA)
95 real(RP), intent(in) :: V_(elem%Np,lmesh%NeA)
96 real(RP), intent(in) :: hs_(elem%Np,lmesh%NeA)
97 real(RP), intent(in) :: u1_(elem%Np,lmesh%NeA)
98 real(RP), intent(in) :: u2_(elem%Np,lmesh%NeA)
99 real(RP), intent(in) :: CORIOLIS(elem%Np,lmesh%NeA)
100
101 real(RP) :: Fx(elem%Np), Fy(elem%Np), LiftDelFlx(elem%Np)
102 real(RP) :: del_flux(elem%NfpTot,lmesh%Ne,PROG_VARS_NUM)
103 real(RP) :: del_flux_aux(elem%NfpTot,lmesh%Ne,1)
104 real(RP) :: u1_dt(elem%Np), u2_dt(elem%Np)
105 real(RP) :: VOR(elem%Np), E(elem%Np)
106 integer :: ke
107 !------------------------------------------------------------------------
108
109 call prof_rapstart('cal_dyn_tend_bndflux', 3)
110 call cal_del_flux_dyn( del_flux, del_flux_aux, & ! (out)
111 h_, u_, v_, hs_, u1_, u2_, & ! (in)
112 lmesh%Gsqrt, lmesh%GIJ(:,:,1,1), lmesh%GIJ(:,:,2,2), & ! (in)
113 lmesh%normal_fn(:,:,1), lmesh%normal_fn(:,:,2), & ! (in)
114 lmesh%vmapM, lmesh%vmapP, & ! (in)
115 lmesh, elem ) ! (in)
116 call prof_rapend('cal_dyn_tend_bndflux', 3)
117
118 !-----
119 call prof_rapstart('cal_dyn_tend_interior', 3)
120
121 !$omp parallel do private( &
122 !$omp E, VOR, &
123 !$omp u1_dt, u2_dt, Fx, Fy, LiftDelFlx )
124 do ke = lmesh%NeS, lmesh%NeE
125 !--
126 call sparsemat_matmul(dx, u2_(:,ke), fx)
127 call sparsemat_matmul(dy, u1_(:,ke), fy)
128 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux_aux(:,ke,1), liftdelflx)
129
130 vor(:) = ( &
131 lmesh%Escale(:,ke,1,1) * fx(:) &
132 - lmesh%Escale(:,ke,2,2) * fy(:) &
133 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
134
135 e(:) = grav * ( hs_(:,ke) + h_(:,ke) ) &
136 + 0.5_rp * ( u1_(:,ke) * u_(:,ke) + u2_(:,ke) * v_(:,ke) )
137
138 !-- h
139 call sparsemat_matmul(dx, lmesh%Gsqrt(:,ke) * h_(:,ke) * u_(:,ke), fx)
140 call sparsemat_matmul(dy, lmesh%Gsqrt(:,ke) * h_(:,ke) * v_(:,ke), fy)
141 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,vars_h_id), liftdelflx)
142
143 h_dt(:,ke) = - ( &
144 lmesh%Escale(:,ke,1,1) * fx(:) &
145 + lmesh%Escale(:,ke,2,2) * fy(:) &
146 + liftdelflx(:) ) / lmesh%Gsqrt(:,ke)
147
148 !-- u1
149 call sparsemat_matmul(dx, e(:), fx)
150 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,vars_u1_id), liftdelflx)
151
152 u1_dt(:) = - ( &
153 lmesh%Escale(:,ke,1,1) * fx(:) &
154 + liftdelflx(:) &
155 - lmesh%Gsqrt(:,ke) * ( coriolis(:,ke) + vor(:) ) * v_(:,ke) )
156
157 !-- u2
158 call sparsemat_matmul(dy, e(:), fy)
159 call sparsemat_matmul(lift, lmesh%Fscale(:,ke) * del_flux(:,ke,vars_u2_id), liftdelflx)
160
161 u2_dt(:) = - ( &
162 lmesh%Escale(:,ke,2,2) * fy(:) &
163 + liftdelflx(:) &
164 + lmesh%Gsqrt(:,ke) * ( coriolis(:,ke) + vor(:) ) * u_(:,ke) )
165
166 !--
167 u_dt(:,ke) = lmesh%GIJ(:,ke,1,1) * u1_dt(:) + lmesh%GIJ(:,ke,1,2) * u2_dt(:)
168 v_dt(:,ke) = lmesh%GIJ(:,ke,2,1) * u1_dt(:) + lmesh%GIJ(:,ke,2,2) * u2_dt(:)
169 end do
170
171 call prof_rapend('cal_dyn_tend_interior', 3)
172
173 return