FE-Project
Loading...
Searching...
No Matches
scale_atm_dyn_dgm_spongelayer.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
9!-------------------------------------------------------------------------------
10#include "scaleFElib.h"
12 !-----------------------------------------------------------------------------
13 !
14 !++ Used modules
15 !
16 use scale_precision
17 use scale_io
18
19 use scale_element_base, only: &
24
25 !-----------------------------------------------------------------------------
26 implicit none
27 private
28 !-----------------------------------------------------------------------------
29 !
30 !++ Public type & procedure
31 !
32
33 type, public :: atmdynspongelayer
34 real(rp) :: wdamp_tau
35 real(rp) :: wdamp_height
36 logical :: hveldamp_flag
37 contains
38 procedure :: init => atm_dyn_dgm_spongelayer_init
39 procedure :: final => atm_dyn_dgm_spongelayer_final
40 procedure :: addtend => atm_dyn_dgm_spongelayer_add_tend
41 end type atmdynspongelayer
42
43 !-----------------------------------------------------------------------------
44 !
45 !++ Public parameters & variables
46 !
47
48 !-----------------------------------------------------------------------------
49 !
50 !++ Private procedures & variables
51 !
52 private :: calc_wdampcoef
53
54contains
55
56!OCL SERIAL
57 subroutine atm_dyn_dgm_spongelayer_init( this, mesh3D, dtsec )
58 use scale_prc, only: prc_abort
61 implicit none
62
63 class(atmdynspongelayer), target, intent(inout) :: this
64 class(meshbase3d), target, intent(in) :: mesh3D
65 real(RP), intent(in) :: dtsec
66
67 real(RP) :: SL_WDAMP_TAU = -1.0_rp ! the maximum tau for Rayleigh damping of w [s]
68 real(RP) :: SL_WDAMP_HEIGHT = -1.0_rp ! the height to start apply Rayleigh damping [m]
69 integer :: SL_WDAMP_LAYER = -1 ! the vertical number of finite element to start apply Rayleigh damping [num]
70 logical :: SL_HORIVELDAMP_FLAG = .false. ! Is the horizontal velocity damped?
71
72 namelist /param_atmos_dyn_spongelayer/ &
73 sl_wdamp_tau, &
74 sl_wdamp_height, &
75 sl_wdamp_layer, &
76 sl_horiveldamp_flag
77
78 class(localmesh3d), pointer :: lcmesh3D
79 class(elementbase3d), pointer :: elem3D
80
81 integer :: NeGZ
82 integer :: ierr
83 !---------------------------------------------------------------
84
85 rewind(io_fid_conf)
86 read(io_fid_conf,nml=param_atmos_dyn_spongelayer,iostat=ierr)
87 if( ierr < 0 ) then !--- missing
88 log_info("ATMOS_DYN_setup_spongelayer",*) 'Not found namelist. Default used.'
89 else if( ierr > 0 ) then !--- fatal error
90 log_error("ATMOS_DYN_setup_spongelayer",*) 'Not appropriate names in namelist PARAM_ATMOS_DYN_SPONGELAYER. Check!'
91 call prc_abort
92 end if
93 log_nml(param_atmos_dyn_spongelayer)
94
95 this%wdamp_tau = sl_wdamp_tau
96 this%wdamp_height = sl_wdamp_height
97
98 lcmesh3d => mesh3d%lcmesh_list(1)
99 elem3d => lcmesh3d%refElem3D
100 select type(mesh3d)
101 type is (meshcubedom3d)
102 negz = mesh3d%NeGZ
103 type is (meshcubedspheredom3d)
104 negz = mesh3d%NeGZ
105 end select
106
107 if ( sl_wdamp_layer > negz ) then
108 log_error("ATMOS_DYN_setup_spongelayer",*) 'SL_wdamp_layer should be less than total of vertical elements (NeGZ). Check!'
109 call prc_abort
110 else if( sl_wdamp_layer > 0 ) then
111 this%wdamp_height = lcmesh3d%pos_en(1,1+(sl_wdamp_layer-1)*lcmesh3d%NeX*lcmesh3d%NeY,3)
112 end if
113 if ( this%wdamp_tau < 0.0_rp ) then
114 this%wdamp_tau = dtsec * 10.0_rp
115 else if ( this%wdamp_tau < dtsec ) then
116 log_error("ATMOS_DYN_setup_spongelayer",*) 'SL_wdamp_tau should be larger than TIME_DT (ATMOS_DYN). Check!'
117 call prc_abort
118 end if
119
120 this%hveldamp_flag = sl_horiveldamp_flag
121
122 return
123 end subroutine atm_dyn_dgm_spongelayer_init
124
125!OCL SERIAL
126 subroutine atm_dyn_dgm_spongelayer_final( this )
127 implicit none
128 class(atmdynspongelayer), target, intent(inout) :: this
129 !---------------------------------------------------------------
130
131 return
132 end subroutine atm_dyn_dgm_spongelayer_final
133
134 subroutine atm_dyn_dgm_spongelayer_add_tend( this, MOMX_dt, MOMY_dt, MOMZ_dt, &
135 MOMX_, MOMY_, MOMZ_, &
136 lmesh, elem )
137
138 implicit none
139
140 class(atmdynspongelayer), intent(in) :: this
141 class(localmesh3d), intent(in) :: lmesh
142 class(elementbase3d), intent(in) :: elem
143 real(RP), intent(inout) :: MOMX_dt(elem%Np,lmesh%NeA)
144 real(RP), intent(inout) :: MOMY_dt(elem%Np,lmesh%NeA)
145 real(RP), intent(inout) :: MOMZ_dt(elem%Np,lmesh%NeA)
146 real(RP), intent(in) :: MOMX_(elem%Np,lmesh%NeA)
147 real(RP), intent(in) :: MOMY_(elem%Np,lmesh%NeA)
148 real(RP), intent(in) :: MOMZ_(elem%Np,lmesh%NeA)
149
150 integer :: ke
151 integer :: ke_x, ke_y, ke_z
152 integer :: keZtop
153 real(RP) :: wdamp_coef(elem%Np)
154 real(RP) :: zTop(elem%Nnode_h1D**2)
155 real(RP) :: s
156 !-----------------------------------------------------------------
157
158 if ( this%hveldamp_flag ) then
159 s = 1.0_rp
160 else
161 s = 0.0_rp
162 end if
163
164 !$omp parallel do collapse(3) private(ke,keZtop,zTop,wdamp_coef)
165 do ke_z = 1, lmesh%NeZ
166 do ke_y = 1, lmesh%NeY
167 do ke_x = 1, lmesh%NeX
168 ke = ke_x + (ke_y-1)*lmesh%NeX + (ke_z-1)*lmesh%NeX*lmesh%NeY
169 keztop = ke_x + (ke_y-1)*lmesh%NeX + (lmesh%NeZ-1)*lmesh%NeX*lmesh%NeY
170 ztop(:) = lmesh%pos_en(elem%Hslice(:,elem%Nnode_v),keztop,3)
171
172 call calc_wdampcoef( &
173 this%wdamp_tau, this%wdamp_height, lmesh%pos_en(:,ke,3), ztop(:), &
174 elem%Nnode_h1D, elem%Nnode_v, &
175 wdamp_coef(:) )
176
177 momx_dt(:,ke) = momx_dt(:,ke) - s * wdamp_coef(:) * momx_(:,ke)
178 momy_dt(:,ke) = momy_dt(:,ke) - s * wdamp_coef(:) * momy_(:,ke)
179 momz_dt(:,ke) = momz_dt(:,ke) - wdamp_coef(:) * momz_(:,ke)
180 end do
181 end do
182 end do
183
184 return
185 end subroutine atm_dyn_dgm_spongelayer_add_tend
186
187!-- private ------------------------------
188
189!OCL SERIAL
190 subroutine calc_wdampcoef( &
191 wdamp_tau, wdamp_height, z, zTop, Nnode_h1D, Nnode_v, & ! (in)
192 wdamp_coef ) ! (out)
193
194 use scale_const, only: &
195 pi => const_pi
196 implicit none
197
198 integer, intent(in) :: Nnode_h1D
199 integer, intent(in) :: Nnode_v
200 real(RP), intent(out) :: wdamp_coef(Nnode_h1D**2,Nnode_v)
201 real(RP), intent(in) :: wdamp_tau
202 real(RP), intent(in) :: wdamp_height
203 real(RP), intent(in) :: z(Nnode_h1D**2,Nnode_v)
204 real(RP), intent(in) :: zTop(Nnode_h1D**2)
205
206 integer :: p_z
207 real(RP) :: sw(Nnode_h1D**2)
208 real(RP) :: r_wdamp_tau
209 !-----------------------------------------------------------------
210
211 r_wdamp_tau = 1.0_rp / wdamp_tau
212 do p_z=1, nnode_v
213 wdamp_coef(:,p_z) = 0.25_rp * r_wdamp_tau &
214 * ( 1.0_rp + sign( 1.0_rp, z(:,p_z) - wdamp_height ) ) &
215 * ( 1.0_rp - cos( pi * (z(:,p_z) - wdamp_height)/(ztop(:) - wdamp_height) ) )
216 end do
217
218 return
219 end subroutine calc_wdampcoef
220
221end module scale_atm_dyn_dgm_spongelayer
module FElib / Fluid dyn solver / Atmosphere / Common
subroutine atm_dyn_dgm_spongelayer_init(this, mesh3d, dtsec)
module FElib / Element / Base
module FElib / Mesh / Local 3D
module FElib / Mesh / Local, Base
module FElib / Mesh / Base 3D
module FElib / Mesh / Cubic 3D domain
module FElib / Mesh / Cubed-sphere 3D domain