FE-Project
Loading...
Searching...
No Matches
mod_atmos_dyn_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
25 use scale_mesh_base, only: meshbase
28
31
33
34 use scale_model_var_manager, only: &
35 modelvarmanager, variableinfo
37
38 use mod_atmos_mesh, only: atmosmesh
39
40 !-----------------------------------------------------------------------------
41 implicit none
42 private
43
44 !-----------------------------------------------------------------------------
45 !
46 !++ Public type & procedures
47 !
48 type, public :: atmosdynvars
49 type(meshfield2d), allocatable :: aux_vars2d(:)
50 type(modelvarmanager) :: auxvars2d_manager
51 contains
52 procedure :: init => atmosdynvars_init
53 procedure :: final => atmosdynvars_final
54 procedure :: history => atmosdynvars_history
55 end type atmosdynvars
56
58 !public :: AtmosDynVars_GetLocalMeshFields_analysis
59
60 !-----------------------------------------------------------------------------
61 !
62 !++ Public parameters & variables
63 !
64
65 !-
66 integer, public, parameter :: atmos_dyn_auxvars2d_num = 1
67 integer, public, parameter :: atmos_dyn_auxvars2d_coriolis_id = 1
68
71 variableinfo( atmos_dyn_auxvars2d_coriolis_id, 'CORIOLIS', 'coriolis parameter', &
72 's-1', 2, 'XY', '' ) /
73
74
75 ! integer, public, parameter :: ATMOS_DYN_ANALYSISVARS_NUM = 6
76 ! integer, public, parameter :: ATMOS_DYN_ANALYSISVAR_MOMZ_t = 1
77 ! integer, public, parameter :: ATMOS_DYN_ANALYSISVAR_MOMZ_t_advX = 2
78 ! integer, public, parameter :: ATMOS_DYN_ANALYSISVAR_MOMZ_t_advY = 3
79 ! integer, public, parameter :: ATMOS_DYN_ANALYSISVAR_MOMZ_t_advZ = 4
80 ! integer, public, parameter :: ATMOS_DYN_ANALYSISVAR_MOMZ_t_lift = 5
81 ! integer, public, parameter :: ATMOS_DYN_ANALYSISVAR_MOMZ_t_buoy = 6
82
83 ! type(VariableInfo), public :: ATMOS_DYN_ANALYSISVAR_VINFO(ATMOS_DYN_ANALYSISVARS_NUM)
84 ! DATA ATMOS_DYN_ANALYSISVAR_VINFO / &
85 ! VariableInfo( ATMOS_DYN_ANALYSISVAR_MOMZ_t, 'MOMZ_t', 'MOMZ_t', &
86 ! 'kg.m-2.s-1', 3, 'XYZ', '' ), &
87 ! VariableInfo( ATMOS_DYN_ANALYSISVAR_MOMZ_t_advX, 'MOMZ_t_advx', 'MOMZ_t_advX', &
88 ! 'kg.m-2.s-1', 3, 'XYZ', '' ), &
89 ! VariableInfo( ATMOS_DYN_ANALYSISVAR_MOMZ_t_advY, 'MOMZ_t_advy', 'MOMZ_t_advY', &
90 ! 'kg.m-2.s-1', 3, 'XYZ', '' ), &
91 ! VariableInfo( ATMOS_DYN_ANALYSISVAR_MOMZ_t_advZ, 'MOMZ_t_advz', 'MOMZ_t_advZ', &
92 ! 'kg.m-2.s-1', 3, 'XYZ', '' ), &
93 ! VariableInfo( ATMOS_DYN_ANALYSISVAR_MOMZ_t_lift, 'MOMZ_t_lift', 'MOMZ_t_lift', &
94 ! 'kg.m-2.s-1', 3, 'XYZ', '' ), &
95 ! VariableInfo( ATMOS_DYN_ANALYSISVAR_MOMZ_t_buoy, 'MOMZ_t_buoy', 'MOMZ_t_buo', &
96 ! 'kg.m-2.s-1', 3, 'XYZ', '' ) &
97 ! /
98
99 !-----------------------------------------------------------------------------
100 !
101 !++ Private procedures
102 !
103 !-------------------
104
105contains
106!OCL SERIAL
107 subroutine atmosdynvars_init( this, model_mesh )
109 implicit none
110 class(atmosdynvars), 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 !--------------------------------------------------
122
123 log_info('AtmosDynVars_Init',*)
124
125 nullify( atm_mesh )
126 select type(model_mesh)
127 class is (atmosmesh)
128 atm_mesh => model_mesh
129 end select
130 mesh3d => atm_mesh%ptr_mesh
131
132 call mesh3d%GetMesh2D( mesh2d )
133
134 !- Initialize 2D auxiliary variables
135
136 call this%AUXVARS2D_manager%Init()
137 allocate( this%AUX_VARS2D(atmos_dyn_auxvars2d_num) )
138
139 reg_file_hist = .false.
141 call this%AUXVARS2D_manager%Regist( &
142 atmos_dyn_auxvars2d_vinfo(v), mesh2d, & ! (in)
143 this%AUX_VARS2D(v), reg_file_hist ) ! (out)
144
145 do n = 1, mesh3d%LOCAL_MESH_NUM
146 this%AUX_VARS2D(v)%local(n)%val(:,:) = 0.0_rp
147 end do
148 end do
149
150 !-
151 ! call this%ANALYSISVARS_manager%Init()
152 ! allocate( this%ANALYSIS_VARS3D(ATMOS_DYN_ANALYSISVARS_NUM) )
153
154 ! reg_file_hist = .true.
155 ! do v = 1, ATMOS_DYN_ANALYSISVARS_NUM
156 ! call this%ANALYSISVARS_manager%Regist( &
157 ! ATMOS_DYN_ANALYSISVAR_VINFO(v), atm_mesh%mesh, &
158 ! this%ANALYSIS_VARS3D(v), reg_file_hist )
159 ! end do
160
161 return
162 end subroutine atmosdynvars_init
163
164 subroutine atmosdynvars_final( this )
165 implicit none
166 class(atmosdynvars), intent(inout) :: this
167
168 !--------------------------------------------------
169
170 log_info('AtmosDynVars_Final',*)
171
172 call this%AUXVARS2D_manager%Final()
173 deallocate( this%AUX_VARS2D )
174
175 !call this%ANALYSISVARS_manager%Final()
176
177 return
178 end subroutine atmosdynvars_final
179
180 subroutine atmosdynvars_history( this )
182 implicit none
183 class(atmosdynvars), intent(in) :: this
184
185 ! integer :: v
186 ! integer :: hst_id
187 !-------------------------------------------------------------------------
188
189 ! do v = 1, ATMOS_DYN_ANALYSISVARS_NUM
190 ! hst_id = this%ANALYSIS_VARS3D(v)%hist_id
191 ! if ( hst_id > 0 ) call FILE_HISTORY_meshfield_put( hst_id, this%ANALYSIS_VARS3D(v) )
192 ! end do
193
194 return
195 end subroutine atmosdynvars_history
196
197 subroutine atmosdynauxvars_getlocalmeshfields( domID, mesh, auxvars_list, &
198 Coriolis, &
199 lcmesh3D &
200 )
201
202 use scale_mesh_base, only: meshbase
204 implicit none
205
206 integer, intent(in) :: domid
207 class(meshbase), intent(in) :: mesh
208 class(modelvarmanager), intent(inout) :: auxvars_list
209 class(localmeshfieldbase), pointer, intent(out) :: coriolis
210 class(localmesh3d), pointer, intent(out), optional :: lcmesh3d
211
212 class(meshfieldbase), pointer :: field
213 class(localmeshbase), pointer :: lcmesh
214 !-------------------------------------------------------
215
216 !--
217 call auxvars_list%Get(atmos_dyn_auxvars2d_coriolis_id, field)
218 call field%GetLocalMeshField(domid, coriolis)
219 !---
220
221 if (present(lcmesh3d)) then
222 call mesh%GetLocalMesh( domid, lcmesh )
223 nullify( lcmesh3d )
224
225 select type(lcmesh)
226 type is (localmesh3d)
227 if (present(lcmesh3d)) lcmesh3d => lcmesh
228 end select
229 end if
230
231 return
233
234 ! subroutine AtmosDynVars_GetLocalMeshFields_analysis( domID, mesh, analysis_list, &
235 ! MOMZ_t, MOMZ_t_advx, MOMZ_t_advY, MOMZ_t_advZ, MOMZ_t_lift, MOMZ_t_buoy, &
236 ! lcmesh3D &
237 ! )
238
239 ! use scale_mesh_base, only: MeshBase
240 ! use scale_meshfield_base, only: MeshFieldBase
241 ! implicit none
242
243 ! integer, intent(in) :: domID
244 ! class(MeshBase), intent(in) :: mesh
245 ! class(ModelVarManager), intent(inout) :: analysis_list
246 ! class(LocalMeshFieldBase), pointer, intent(out) :: &
247 ! MOMZ_t, MOMZ_t_advx, MOMZ_t_advY, MOMZ_t_advZ, MOMZ_t_lift, MOMZ_t_buoy
248 ! class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
249
250 ! class(MeshFieldBase), pointer :: field
251 ! class(LocalMeshBase), pointer :: lcmesh
252 ! !-------------------------------------------------------
253
254 ! !--
255 ! call analysis_list%Get(ATMOS_DYN_ANALYSISVAR_MOMZ_t, field)
256 ! call field%GetLocalMeshField(domID, MOMZ_t)
257 ! !---
258 ! call analysis_list%Get(ATMOS_DYN_ANALYSISVAR_MOMZ_t_advx, field)
259 ! call field%GetLocalMeshField(domID, MOMZ_t_advx)
260 ! !---
261 ! call analysis_list%Get(ATMOS_DYN_ANALYSISVAR_MOMZ_t_advy, field)
262 ! call field%GetLocalMeshField(domID, MOMZ_t_advy)
263 ! !---
264 ! call analysis_list%Get(ATMOS_DYN_ANALYSISVAR_MOMZ_t_advz, field)
265 ! call field%GetLocalMeshField(domID, MOMZ_t_advz)
266 ! !---
267 ! call analysis_list%Get(ATMOS_DYN_ANALYSISVAR_MOMZ_t_lift, field)
268 ! call field%GetLocalMeshField(domID, MOMZ_t_lift)
269 ! !---
270 ! call analysis_list%Get(ATMOS_DYN_ANALYSISVAR_MOMZ_t_buoy, field)
271 ! call field%GetLocalMeshField(domID, MOMZ_t_buoy)
272
273
274 ! if (present(lcmesh3D)) then
275 ! call mesh%GetLocalMesh( domID, lcmesh )
276 ! nullify( lcmesh3D )
277
278 ! select type(lcmesh)
279 ! type is (LocalMesh3D)
280 ! if (present(lcmesh3D)) lcmesh3D => lcmesh
281 ! end select
282 ! end if
283
284 ! return
285 ! end subroutine AtmosDynVars_GetLocalMeshFields_analysis
286
287
288end module mod_atmos_dyn_vars
module ATMOSPHERE dynamics
type(variableinfo), dimension(atmos_dyn_auxvars2d_num), public atmos_dyn_auxvars2d_vinfo
integer, parameter, public atmos_dyn_auxvars2d_coriolis_id
integer, parameter, public atmos_dyn_auxvars2d_num
subroutine, public atmosdynauxvars_getlocalmeshfields(domid, mesh, auxvars_list, coriolis, lcmesh3d)
module Atmosphere / Mesh
module FElib / Element / Base
module FElib / Mesh / Local 2D
module FElib / Mesh / Local 3D
module FElib / Mesh / Local, Base
integer, parameter, public local_meshfield_type_nodes_faceval
module FElib / Mesh / Base 2D
module FElib / Mesh / Base 3D
module FElib / Mesh / Base
module FElib / Data / base
module FElib / Data / Communication base
FElib / model framework / mesh manager (base)
FElib / model framework / variable manager.