FE-Project
Loading...
Searching...
No Matches
mod_atmos_phy_sfc_vars.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_io
19 use scale_prc
20
22 use scale_mesh_base, only: meshbase
24 use scale_mesh_base3d, only: &
25 meshbase3d, &
26 dimtype_xyz => meshbase3d_dimtypeid_xyz
30 use scale_meshfield_base, only: &
32
35
36 use scale_model_var_manager, only: &
37 modelvarmanager, variableinfo
38 use scale_model_mesh_manager, only: modelmeshbase
39
40 use mod_atmos_mesh, only: atmosmesh
41
42 !-----------------------------------------------------------------------------
43 implicit none
44 private
45
46
47 !-----------------------------------------------------------------------------
48 !
49 !++ Public type & procedures
50 !
51 type, public :: atmosphysfcvars
52 type(meshfield2d), allocatable :: sfc_vars(:)
53 type(modelvarmanager) :: sfcvars_manager
54
55 type(meshfield2d), allocatable :: sfc_flx(:)
56 type(modelvarmanager) :: sfcflx_manager
57
58 integer :: sfcflx_num_tot
59 contains
60 procedure :: init => atmosphysfcvars_init
61 procedure :: final => atmosphysfcvars_final
62 procedure :: history => atmosphysfcvars_history
63 end type atmosphysfcvars
64
65 integer, public, parameter :: atmos_phy_sf_svar_temp_id = 1
66 integer, public, parameter :: atmos_phy_sf_svar_num = 1
67 type(variableinfo), public :: atmos_phy_sf_svar_vinfo(atmos_phy_sf_svar_num)
69 variableinfo( atmos_phy_sf_svar_num, 'SFC_TEMP', 'surface skin temperature', &
70 'K', 2, 'XY', '' ) /
71
72 integer, public, parameter :: atmos_phy_sf_sflx_mu_id = 1
73 integer, public, parameter :: atmos_phy_sf_sflx_mv_id = 2
74 integer, public, parameter :: atmos_phy_sf_sflx_mw_id = 3
75 integer, public, parameter :: atmos_phy_sf_sflx_sh_id = 4
76 integer, public, parameter :: atmos_phy_sf_sflx_lh_id = 5
77 integer, public, parameter :: atmos_phy_sf_sflx_qv_id = 6
78 integer, public, parameter :: atmos_phy_sf_sflx_num1 = 6
79
80 type(variableinfo), public :: atmos_phy_sf_sflx_vinfo(atmos_phy_sf_sflx_num1)
82 variableinfo( atmos_phy_sf_sflx_mu_id, 'SFLX_MU', 'x-momentum flux', &
83 'm/s*kg/m2/s', 2, 'XY', '' ), &
84 variableinfo( atmos_phy_sf_sflx_mv_id, 'SFLX_MV', 'y-momentum flux', &
85 'm/s*kg/m2/s', 2, 'XY', '' ), &
86 variableinfo( atmos_phy_sf_sflx_mw_id, 'SFLX_MW', 'z-momentum flux', &
87 'm/s*kg/m2/s', 2, 'XY', '' ), &
88 variableinfo( atmos_phy_sf_sflx_sh_id, 'SFLX_SH', 'sensible heat flux', &
89 'J/m2/s', 2, 'XY', '' ), &
90 variableinfo( atmos_phy_sf_sflx_lh_id, 'SFLX_LH', 'latent heat flux', &
91 'J/m2/s', 2, 'XY', '' ), &
92 variableinfo( atmos_phy_sf_sflx_qv_id, 'SFLX_QV', 'water vapor flux', &
93 'kg/m2/s', 2, 'XY', '' ) /
94
96
97 !-----------------------------------------------------------------------------
98 !
99 !++ Private procedures
100 !
101 !-------------------
102 !
103 !++ Private parameters & variables
104 !
105 real(rp), private :: atmos_phy_sfc_default_sfc_temp = 300.0_rp
106
107contains
108 subroutine atmosphysfcvars_init( this, model_mesh )
109 implicit none
110 class(atmosphysfcvars), target, intent(inout) :: this
111 class(modelmeshbase), target, intent(in) :: model_mesh
112
113 integer :: v
114 integer :: n
115 logical :: reg_file_hist
116
117 class(atmosmesh), pointer :: atm_mesh
118 class(meshbase2d), pointer :: mesh2d
119 class(meshbase3d), pointer :: mesh3d
120
121 namelist / param_atmos_phy_sfc_vars / &
122 atmos_phy_sfc_default_sfc_temp
123
124 integer :: ierr
125 !--------------------------------------------------
126
127 log_info('AtmosPhySfcVars_Init',*)
128
129 this%SFCFLX_NUM_TOT = atmos_phy_sf_sflx_num1
130
131 !--- read namelist
132 rewind(io_fid_conf)
133 read(io_fid_conf,nml=param_atmos_phy_sfc_vars,iostat=ierr)
134 if( ierr < 0 ) then !--- missing
135 log_info("ATMOS_PHY_SFC_vars_setup",*) 'Not found namelist. Default used.'
136 elseif( ierr > 0 ) then !--- fatal error
137 log_error("ATMOS_PHY_SFC_vars_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_SFC_VARS. Check!'
138 call prc_abort
139 endif
140 log_nml(param_atmos_phy_sfc_vars)
141
142 !- Initialize auxiliary and diagnostic variables
143
144 nullify( atm_mesh )
145 select type(model_mesh)
146 class is (atmosmesh)
147 atm_mesh => model_mesh
148 end select
149 mesh3d => atm_mesh%ptr_mesh
150
151 call mesh3d%GetMesh2D( mesh2d )
152
153 !----
154 call this%SFCVARS_manager%Init()
155 allocate( this%SFC_VARS(atmos_phy_sf_svar_num) )
156
157 reg_file_hist = .true.
158 do v = 1, atmos_phy_sf_svar_num
159 call this%SFCVARS_manager%Regist( &
160 atmos_phy_sf_svar_vinfo(v), mesh2d, &
161 this%SFC_VARS(v), reg_file_hist, fill_zero=.true. )
162 end do
163
164 do n=1, mesh2d%LOCAL_MESH_NUM
165 this%SFC_VARS(atmos_phy_sf_svar_temp_id)%local(n)%val(:,:) = atmos_phy_sfc_default_sfc_temp
166 end do
167
168 !----
169 call this%SFCFLX_manager%Init()
170 allocate( this%SFC_FLX(this%SFCFLX_NUM_TOT) )
171
172 reg_file_hist = .true.
173 do v = 1, this%SFCFLX_NUM_TOT
174 call this%SFCFLX_manager%Regist( &
175 atmos_phy_sf_sflx_vinfo(v), mesh2d, &
176 this%SFC_FLX(v), reg_file_hist, fill_zero=.true. )
177 end do
178
179 return
180 end subroutine atmosphysfcvars_init
181
182 subroutine atmosphysfcvars_final( this )
183 implicit none
184 class(atmosphysfcvars), intent(inout) :: this
185
186 !--------------------------------------------------
187
188 log_info('AtmosPhySfcVars_Final',*)
189
190 call this%SFCVARS_manager%Final()
191 deallocate( this%SFC_VARS )
192
193 call this%SFCFLX_manager%Final()
194 deallocate( this%SFC_FLX )
195
196 return
197 end subroutine atmosphysfcvars_final
198
199!OCL SERIAL
200 subroutine atmosphysfcvars_history( this )
202 implicit none
203 class(atmosphysfcvars), intent(inout) :: this
204
205 integer :: v
206 integer :: iq
207 integer :: hst_id
208 !-------------------------------------------------------------------------
209
210 do v = 1, atmos_phy_sf_svar_num
211 hst_id = this%SFC_VARS(v)%hist_id
212 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%SFC_VARS(v) )
213 end do
214
215 do v = 1, this%SFCFLX_NUM_TOT
216 hst_id = this%SFC_FLX(v)%hist_id
217 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%SFC_FLX(v) )
218 end do
219
220 return
221 end subroutine atmosphysfcvars_history
222
223 subroutine atmosphysfcvars_getlocalmeshfields( domID, mesh, svars_list, sflx_list, &
224 SFC_TEMP, SFLX_MU, SFLX_MV, SFLX_MW, SFLX_SH, SFLX_LH, SFLX_QV, &
225 lcmesh3D &
226 )
227
228 use scale_mesh_base, only: meshbase
230 implicit none
231
232 integer, intent(in) :: domid
233 class(meshbase), intent(in) :: mesh
234 class(modelvarmanager), intent(inout) :: svars_list
235 class(modelvarmanager), intent(inout) :: sflx_list
236 class(localmeshfieldbase), pointer, intent(out) :: sfc_temp
237 class(localmeshfieldbase), pointer, intent(out) :: sflx_mu
238 class(localmeshfieldbase), pointer, intent(out) :: sflx_mv
239 class(localmeshfieldbase), pointer, intent(out) :: sflx_mw
240 class(localmeshfieldbase), pointer, intent(out) :: sflx_sh
241 class(localmeshfieldbase), pointer, intent(out) :: sflx_lh
242 class(localmeshfieldbase), pointer, intent(out) :: sflx_qv
243 class(localmesh3d), pointer, intent(out), optional :: lcmesh3d
244
245 class(meshfieldbase), pointer :: field
246 class(localmeshbase), pointer :: lcmesh
247 !-------------------------------------------------------
248
249 !--
250 call svars_list%Get(atmos_phy_sf_svar_temp_id, field)
251 call field%GetLocalMeshField(domid, sfc_temp)
252
253 call sflx_list%Get(atmos_phy_sf_sflx_mu_id, field)
254 call field%GetLocalMeshField(domid, sflx_mu)
255
256 call sflx_list%Get(atmos_phy_sf_sflx_mv_id, field)
257 call field%GetLocalMeshField(domid, sflx_mv)
258
259 call sflx_list%Get(atmos_phy_sf_sflx_mw_id, field)
260 call field%GetLocalMeshField(domid, sflx_mw)
261
262 call sflx_list%Get(atmos_phy_sf_sflx_sh_id, field)
263 call field%GetLocalMeshField(domid, sflx_sh)
264
265 call sflx_list%Get(atmos_phy_sf_sflx_lh_id, field)
266 call field%GetLocalMeshField(domid, sflx_lh)
267
268 call sflx_list%Get(atmos_phy_sf_sflx_qv_id, field)
269 call field%GetLocalMeshField(domid, sflx_qv)
270 !---
271
272 if (present(lcmesh3d)) then
273 call mesh%GetLocalMesh( domid, lcmesh )
274 nullify( lcmesh3d )
275
276 select type(lcmesh)
277 type is (localmesh3d)
278 if (present(lcmesh3d)) lcmesh3d => lcmesh
279 end select
280 end if
281
282 return
284
285end module mod_atmos_phy_sfc_vars
module Atmosphere / Mesh
module ATMOSPHERE physics / surface process
type(variableinfo), dimension(atmos_phy_sf_sflx_num1), public atmos_phy_sf_sflx_vinfo
subroutine, public atmosphysfcvars_getlocalmeshfields(domid, mesh, svars_list, sflx_list, sfc_temp, sflx_mu, sflx_mv, sflx_mw, sflx_sh, sflx_lh, sflx_qv, lcmesh3d)
integer, parameter, public atmos_phy_sf_sflx_mv_id
type(variableinfo), dimension(atmos_phy_sf_svar_num), public atmos_phy_sf_svar_vinfo
integer, parameter, public atmos_phy_sf_sflx_qv_id
integer, parameter, public atmos_phy_sf_sflx_mu_id
integer, parameter, public atmos_phy_sf_sflx_num1
integer, parameter, public atmos_phy_sf_sflx_sh_id
integer, parameter, public atmos_phy_sf_svar_num
integer, parameter, public atmos_phy_sf_sflx_lh_id
integer, parameter, public atmos_phy_sf_svar_temp_id
integer, parameter, public atmos_phy_sf_sflx_mw_id
module FElib / Element / Base
module FElib / Mesh / Local 3D
module FElib / Mesh / Local, Base
module FElib / Mesh / Base 2D
module FElib / Mesh / Base 3D
integer, public meshbase3d_dimtypeid_xyz
module FElib / Mesh / Base
module FElib / Data / base
FElib / model framework / mesh manager.
FElib / model framework / variable manager.