FE-Project
Loading...
Searching...
No Matches
mod_atmos_phy_tb.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
10!-------------------------------------------------------------------------------
11#include "scaleFElib.h"
13 !-----------------------------------------------------------------------------
14 !
15 !++ used modules
16 !
17 use scale_precision
18 use scale_prc
19 use scale_io
20 use scale_prof
21 use scale_const, only: &
22 undef8 => const_undef8
23 use scale_tracer, only: qa
24
25 use scale_mesh_base, only: meshbase
28
32 use scale_element_base, only: &
34
36 use scale_localmeshfield_base, only: &
38
39 use scale_model_mesh_manager, only: modelmeshbase
42
44
48
49 !-----------------------------------------------------------------------------
50 implicit none
51 private
52 !-----------------------------------------------------------------------------
53 !
54 !++ Public type & procedure
55 !
56
59 type, extends(modelcomponentproc), public :: atmosphytb
60 integer :: tb_typeid
61 type(atmphytbdgmdriver) :: tb_driver
62
63 type(atmosphytbvars) :: vars
64 type(atmdynbnd), pointer :: dyn_bnd
65 contains
66 procedure, public :: setup => atmosphytb_setup
67 procedure, public :: calc_tendency => atmosphytb_calc_tendency
68 procedure, public :: update => atmosphytb_update
69 procedure, public :: finalize => atmosphytb_finalize
70 procedure, public :: setdynbc => atmosphytb_setdynbc
71 end type atmosphytb
72
73 !-----------------------------------------------------------------------------
74 !++ Public parameters & variables
75 !-----------------------------------------------------------------------------
76
77 !-----------------------------------------------------------------------------
78 !
79 !++ Private procedure
80 !
81 !-----------------------------------------------------------------------------
82 !
83 !++ Private parameters & variables
84 !
85 !-----------------------------------------------------------------------------
86
87contains
93!OCL SERIAL
94 subroutine atmosphytb_setup( this, model_mesh, tm_parent_comp )
95 use mod_atmos_mesh, only: atmosmesh
97
98 implicit none
99 class(atmosphytb), intent(inout) :: this
100 class(modelmeshbase), target, intent(in) :: model_mesh
101 class(time_manager_component), intent(inout) :: tm_parent_comp
102
103 real(DP) :: TIME_DT = undef8
104 character(len=H_SHORT) :: TIME_DT_UNIT = 'SEC'
105
106 character(len=H_MID) :: TB_TYPE = 'SMAGORINSKY'
107 namelist /param_atmos_phy_tb/ &
108 time_dt, &
109 time_dt_unit, &
110 tb_type
111
112 class(atmosmesh), pointer :: atm_mesh
113 class(meshbase), pointer :: ptr_mesh
114 real(DP) :: dtsec
115
116 integer :: ierr
117 !--------------------------------------------------
118
119 if (.not. this%IsActivated()) return
120
121 log_newline
122 log_info("ATMOS_PHY_TB_setup",*) 'Setup'
123
124 !--- read namelist
125 rewind(io_fid_conf)
126 read(io_fid_conf,nml=param_atmos_phy_tb,iostat=ierr)
127 if( ierr < 0 ) then !--- missing
128 log_info("ATMOS_PHY_TB_setup",*) 'Not found namelist. Default used.'
129 elseif( ierr > 0 ) then !--- fatal error
130 log_error("ATMOS_PHY_TB_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_TB. Check!'
131 call prc_abort
132 endif
133 log_nml(param_atmos_phy_tb)
134
135 !- get mesh --------------------------------------------------
136
137 call model_mesh%GetModelMesh( ptr_mesh )
138 select type(model_mesh)
139 class is (atmosmesh)
140 atm_mesh => model_mesh
141 end select
142
143 !--- Regist this compoent in the time manager
144
145 call tm_parent_comp%Regist_process( 'ATMOS_PHY_TB', time_dt, time_dt_unit, & ! (in)
146 this%tm_process_id ) ! (out)
147
148 dtsec = tm_parent_comp%process_list(this%tm_process_id)%dtsec
149
150 !- initialize the variables
151 call this%vars%Init( model_mesh )
152
153 !- Initialize a module for turbulence model
154 call this%tb_driver%Init( tb_type, dtsec, atm_mesh )
155
156 !--
157 this%dyn_bnd => null()
158
159 return
160 end subroutine atmosphytb_setup
161
171!OCL SERIAL
172 subroutine atmosphytb_calc_tendency( &
173 this, model_mesh, prgvars_list, trcvars_list, &
174 auxvars_list, forcing_list, is_update )
175
176 use scale_tracer, only: &
177 tracer_advc
180 use mod_atmos_vars, only: &
186 use mod_atmos_phy_tb_vars, only: &
188
189 implicit none
190
191 class(atmosphytb), intent(inout) :: this
192 class(modelmeshbase), intent(in) :: model_mesh
193 class(modelvarmanager), intent(inout) :: prgvars_list
194 class(modelvarmanager), intent(inout) :: trcvars_list
195 class(modelvarmanager), intent(inout) :: auxvars_list
196 class(modelvarmanager), intent(inout) :: forcing_list
197 logical, intent(in) :: is_update
198
199 class(meshbase), pointer :: mesh
200 class(meshbase3d), pointer :: mesh3D
201 class(localmesh3d), pointer :: lcmesh
202 integer :: n
203 integer :: ke
204 integer :: iq
205
206 class(localmeshfieldbase), pointer :: DENS_tp, MOMX_tp, MOMY_tp, MOMZ_tp, RHOT_tp, RHOH_P
207 type(localmeshfieldbaselist) :: RHOQ_tp(QA)
208 class(localmeshfieldbase), pointer :: tb_MOMX_t, tb_MOMY_t, tb_MOMZ_t, tb_RHOT_t
209 class(localmeshfieldbase), pointer :: tb_RHOQ_t
210 type(localmeshfieldbaselist) :: tb_RHOQ_t_list(QA)
211 !--------------------------------------------------
212
213 if (.not. this%IsActivated()) return
214
215 !LOG_INFO('AtmosDyn_tendency',*)
216
217 log_progress(*) 'atmosphere / physics / turbulence'
218
219 call model_mesh%GetModelMesh( mesh )
220 select type(mesh)
221 class is (meshbase3d)
222 mesh3d => mesh
223 end select
224
225 if (is_update) then
226 call this%tb_driver%Tendency( this%vars%tends_manager, &
227 prgvars_list, trcvars_list, auxvars_list, &
228 this%vars%auxvars_manager, this%vars%diagvars_manager, &
229 this%dyn_bnd, &
230 model_mesh%DOptrMat(1), model_mesh%DOptrMat(2), model_mesh%DOptrMat(3), &
231 model_mesh%SOptrMat(1), model_mesh%SOptrMat(2), model_mesh%SOptrMat(3), &
232 model_mesh%LiftOptrMat, mesh3d )
233 end if
234
235 call prof_rapstart('ATM_PHY_TB_add_tend', 2)
236 do n=1, mesh%LOCAL_MESH_NUM
238 mesh, forcing_list, &
239 dens_tp, momx_tp, momy_tp, momz_tp, rhot_tp, &
240 rhoh_p, rhoq_tp )
241
243 mesh, this%vars%tends_manager, &
244 tb_momx_t, tb_momy_t, tb_momz_t, tb_rhot_t, tb_rhoq_t_list, &
245 lcmesh )
246
247 !$omp parallel private(ke, iq)
248 !$omp do
249 do ke=lcmesh%NeS, lcmesh%NeE
250 momx_tp%val(:,ke) = momx_tp%val(:,ke) + tb_momx_t%val(:,ke)
251 momy_tp%val(:,ke) = momy_tp%val(:,ke) + tb_momy_t%val(:,ke)
252 momz_tp%val(:,ke) = momz_tp%val(:,ke) + tb_momz_t%val(:,ke)
253 rhot_tp%val(:,ke) = rhot_tp%val(:,ke) + tb_rhot_t%val(:,ke)
254 end do
255 !$omp end do
256 do iq = 1, qa
257 if ( .not. tracer_advc(iq) ) cycle
258 !$omp do
259 do ke = lcmesh%NeS, lcmesh%NeE
260 rhoq_tp(iq)%ptr%val(:,ke) = rhoq_tp(iq)%ptr%val(:,ke) &
261 + tb_rhoq_t_list(iq)%ptr%val(:,ke)
262 end do
263 !$omp end do
264 end do
265 !$omp end parallel
266 end do
267 call prof_rapend('ATM_PHY_TB_add_tend', 2)
268
269 return
270 end subroutine atmosphytb_calc_tendency
271
280!OCL SERIAL
281 subroutine atmosphytb_update( this, model_mesh, prgvars_list, trcvars_list, auxvars_list, forcing_list, is_update )
282 implicit none
283
284 class(atmosphytb), intent(inout) :: this
285 class(modelmeshbase), intent(in) :: model_mesh
286 class(modelvarmanager), intent(inout) :: prgvars_list
287 class(modelvarmanager), intent(inout) :: trcvars_list
288 class(modelvarmanager), intent(inout) :: auxvars_list
289 class(modelvarmanager), intent(inout) :: forcing_list
290 logical, intent(in) :: is_update
291 !--------------------------------------------------
292
293 return
294 end subroutine atmosphytb_update
295
298!OCL SERIAL
299 subroutine atmosphytb_finalize( this )
300 implicit none
301 class(atmosphytb), intent(inout) :: this
302
303 !--------------------------------------------------
304 if (.not. this%IsActivated()) return
305
306 call this%tb_driver%Final()
307 call this%vars%Final()
308
309 return
310 end subroutine atmosphytb_finalize
311
316!OCL SERIAL
317 subroutine atmosphytb_setdynbc( this, dyn_bnd )
318 implicit none
319 class(atmosphytb), intent(inout) :: this
320 type(atmdynbnd), intent(in), target :: dyn_bnd
321 !--------------------------------------------------
322
323 this%dyn_bnd => dyn_bnd
324
325 return
326 end subroutine atmosphytb_setdynbc
327
328!- private ------------------------------------------------
329
330end module mod_atmos_phy_tb
module Atmosphere / Mesh
module ATMOSPHERE physics / sub-grid scale turbulence
subroutine, public atmosphytbvars_getlocalmeshfields_tend(domid, mesh, tb_tends_list, tb_momx_t, tb_momy_t, tb_momz_t, tb_rhot_t, tb_rhoq_t, lcmesh3d)
module ATMOSPHERE physics / sub-grid scale turbulence process
subroutine atmosphytb_setup(this, model_mesh, tm_parent_comp)
Setup a component of SGS turbulence process.
module ATMOSPHERE / Variables
subroutine, public atmosvars_getlocalmeshprgvar(domid, mesh, prgvars_list, auxvars_list, varid, var, dens_hyd, pres_hyd, lcmesh3d)
subroutine, public atmosvars_getlocalmeshprgvars(domid, mesh, prgvars_list, auxvars_list, ddens, momx, momy, momz, therm, dens_hyd, pres_hyd, rtot, cvtot, cptot, lcmesh3d)
subroutine, public atmosvars_getlocalmeshphytends(domid, mesh, phytends_list, dens_tp, momx_tp, momy_tp, momz_tp, rhot_tp, rhoh_p, rhoq_tp, lcmesh3d)
subroutine, public atmosvars_getlocalmeshqtrcvar(domid, mesh, trcvars_list, varid, var, lcmesh3d)
subroutine, public atmosvars_getlocalmeshphyauxvars(domid, mesh, phyauxvars_list, pres, pt, lcmesh3d)
module FElib / Fluid dyn solver / Atmosphere / Boundary
module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / Common
module FElib / Physics turbulence / Atmosphere / driver
module FElib / Element / Base
module FElib / Mesh / Local 2D
module FElib / Mesh / Local 3D
module FElib / Mesh / Local, Base
module FElib / Mesh / Base 2D
module FElib / Mesh / Base 3D
module FElib / Mesh / Base
module FElib / Data / base
FElib / model framework / physics process.
FElib / model framework / mesh manager.
FElib / model framework / variable manager.
module common / time
Derived type to manage a component of sub-grid scale turbulent process.