FV3-JEDI
fv3jedi_fv3lm_mod.f90
Go to the documentation of this file.
1 ! (C) Copyright 2017-2020 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
7 
8 ! oops uses
9 use datetime_mod
10 use duration_mod
11 
12 ! fckit
13 use fckit_configuration_module, only: fckit_configuration
14 
15 ! Linear model
16 use fv3jedi_lm_mod, only: fv3jedi_lm_type
17 
18 ! fv3-jedi uses
22 
23 implicit none
24 private
25 public :: fv3lm_model
26 
27 ! --------------------------------------------------------------------------------------------------
28 
29 type :: fv3lm_model
30  type(fv3jedi_lm_type) :: fv3jedi_lm !<Linearized model object
31  contains
32  procedure, public :: create
33  procedure, public :: delete
34  procedure, public :: initialize
35  procedure, public :: step
36  procedure, public :: finalize
37 end type fv3lm_model
38 
39 ! --------------------------------------------------------------------------------------------------
40 
41 contains
42 
43 ! --------------------------------------------------------------------------------------------------
44 
45 subroutine create(self, geom, conf)
46 
47 implicit none
48 class(fv3lm_model), intent(inout) :: self
49 type(fv3jedi_geom), intent(in) :: geom
50 type(fckit_configuration), intent(in) :: conf
51 
52 !Locals
53 character(len=20) :: ststep
54 type(duration) :: dtstep
55 real(kind=kind_real) :: dt
56 character(len=:), allocatable :: str
57 
58 
59 ! Model time step
60 ! ---------------
61 call conf%get_or_die("tstep",str)
62 ststep = str
63 deallocate(str)
64 
65 dtstep = trim(ststep)
66 dt = real(duration_seconds(dtstep),kind_real)
67 
68 
69 ! Model configuration and creation
70 ! --------------------------------
71 call conf%get_or_die("lm_do_dyn",self%fv3jedi_lm%conf%do_dyn)
72 call conf%get_or_die("lm_do_trb",self%fv3jedi_lm%conf%do_phy_trb)
73 call conf%get_or_die("lm_do_mst",self%fv3jedi_lm%conf%do_phy_mst)
74 
75 call self%fv3jedi_lm%create(dt,geom%npx,geom%npy,geom%npz,geom%ptop,geom%ak,geom%bk)
76 
77 
78 ! Safety checks
79 ! -------------
80 
81 !The full trajecotory of the tlm/adm is not output by this simplified model
82 !so if being used to generate the trajectry with physics the traj must be read
83 !from file or obtained by running GEOS or GFS.
84 if ((self%fv3jedi_lm%conf%do_phy_trb .ne. 0) .or. &
85  (self%fv3jedi_lm%conf%do_phy_mst .ne. 0) ) then
86  call abor1_ftn("fv3lm_model | FV3LM : unless reading the trajecotory physics should be off")
87 endif
88 
89 end subroutine create
90 
91 ! --------------------------------------------------------------------------------------------------
92 
93 subroutine delete(self)
94 
95 implicit none
96 class(fv3lm_model), intent(inout) :: self
97 
98 !Delete the model
99 !----------------
100 call self%fv3jedi_lm%delete()
101 
102 end subroutine delete
103 
104 ! --------------------------------------------------------------------------------------------------
105 
106 subroutine initialize(self, state)
107 
108 implicit none
109 class(fv3lm_model), intent(inout) :: self
110 type(fv3jedi_state), intent(in) :: state
111 
112 call self%fv3jedi_lm%init_nl()
113 
114 end subroutine initialize
115 
116 ! --------------------------------------------------------------------------------------------------
117 
118 subroutine step(self, state, geom, vdate)
119 
120 implicit none
121 class(fv3lm_model), intent(inout) :: self
122 type(fv3jedi_state), intent(inout) :: state
123 type(fv3jedi_geom), intent(inout) :: geom
124 type(datetime), intent(inout) :: vdate !< Valid datetime after step
125 
126 call state_to_lm(state,self%fv3jedi_lm)
127 call self%fv3jedi_lm%step_nl()
128 call lm_to_state(self%fv3jedi_lm,state)
129 
130 end subroutine step
131 
132 ! --------------------------------------------------------------------------------------------------
133 
134 subroutine finalize(self, state)
135 
136 implicit none
137 class(fv3lm_model), intent(inout) :: self
138 type(fv3jedi_state), intent(inout) :: state
139 
140 call self%fv3jedi_lm%final_nl()
141 
142 end subroutine finalize
143 
144 ! --------------------------------------------------------------------------------------------------
145 
146 subroutine state_to_lm( state, lm )
147 
148 implicit none
149 type(fv3jedi_state), intent(in) :: state
150 type(fv3jedi_lm_type), intent(inout) :: lm
151 
152 real(kind=kind_real), pointer, dimension(:,:,:) :: ud
153 real(kind=kind_real), pointer, dimension(:,:,:) :: vd
154 real(kind=kind_real), pointer, dimension(:,:,:) :: ua
155 real(kind=kind_real), pointer, dimension(:,:,:) :: va
156 real(kind=kind_real), pointer, dimension(:,:,:) :: t
157 real(kind=kind_real), pointer, dimension(:,:,:) :: delp
158 real(kind=kind_real), pointer, dimension(:,:,:) :: q
159 real(kind=kind_real), pointer, dimension(:,:,:) :: qi
160 real(kind=kind_real), pointer, dimension(:,:,:) :: ql
161 real(kind=kind_real), pointer, dimension(:,:,:) :: o3
162 real(kind=kind_real), pointer, dimension(:,:,:) :: w
163 real(kind=kind_real), pointer, dimension(:,:,:) :: delz
164 real(kind=kind_real), pointer, dimension(:,:,:) :: phis
165 
166 call state%get_field('ud' , ud )
167 call state%get_field('vd' , vd )
168 call state%get_field('t' , t )
169 call state%get_field('delp' , delp)
170 call state%get_field('sphum' , q )
171 call state%get_field('ice_wat', qi )
172 call state%get_field('liq_wat', ql )
173 if (state%has_field('o3mr' )) call state%get_field('o3mr' , o3)
174 if (state%has_field('o3ppmv')) call state%get_field('o3ppmv', o3)
175 lm%traj%ua = 0.0_kind_real
176 lm%traj%va = 0.0_kind_real
177 
178 lm%traj%u = ud(state%isc:state%iec,state%jsc:state%jec,:)
179 lm%traj%v = vd(state%isc:state%iec,state%jsc:state%jec,:)
180 lm%traj%t = t
181 lm%traj%delp = delp
182 lm%traj%qv = q
183 lm%traj%qi = qi
184 lm%traj%ql = ql
185 lm%traj%o3 = o3
186 
187 if (state%has_field('ua')) then
188  call state%get_field('ua', ua )
189  call state%get_field('va', va )
190  lm%traj%ua = ua
191  lm%traj%va = va
192 endif
193 
194 if (.not. lm%conf%hydrostatic) then
195  call state%get_field('w ', w )
196  call state%get_field('delz', delz)
197  lm%traj%w = w
198  lm%traj%delz = delz
199 endif
200 
201 call state%get_field('phis', phis )
202 lm%traj%phis = phis(:,:,1)
203 
204 end subroutine state_to_lm
205 
206 ! --------------------------------------------------------------------------------------------------
207 
208 subroutine lm_to_state( lm, state )
209 
210 implicit none
211 type(fv3jedi_lm_type), intent(in) :: lm
212 type(fv3jedi_state), intent(inout) :: state
213 
214 real(kind=kind_real), pointer, dimension(:,:,:) :: ud
215 real(kind=kind_real), pointer, dimension(:,:,:) :: vd
216 real(kind=kind_real), pointer, dimension(:,:,:) :: ua
217 real(kind=kind_real), pointer, dimension(:,:,:) :: va
218 real(kind=kind_real), pointer, dimension(:,:,:) :: t
219 real(kind=kind_real), pointer, dimension(:,:,:) :: delp
220 real(kind=kind_real), pointer, dimension(:,:,:) :: q
221 real(kind=kind_real), pointer, dimension(:,:,:) :: qi
222 real(kind=kind_real), pointer, dimension(:,:,:) :: ql
223 real(kind=kind_real), pointer, dimension(:,:,:) :: o3
224 real(kind=kind_real), pointer, dimension(:,:,:) :: w
225 real(kind=kind_real), pointer, dimension(:,:,:) :: delz
226 real(kind=kind_real), pointer, dimension(:,:,:) :: phis
227 
228 call state%get_field('ud' , ud )
229 call state%get_field('vd' , vd )
230 call state%get_field('t' , t )
231 call state%get_field('delp' , delp)
232 call state%get_field('sphum' , q )
233 call state%get_field('ice_wat', qi )
234 call state%get_field('liq_wat', ql )
235 if ( state%has_field('o3mr' ) )call state%get_field('o3mr' , o3 )
236 if ( state%has_field('o3ppmv' ) )call state%get_field('o3ppmv' , o3 )
237 ud(state%isc:state%iec,state%jsc:state%jec,:) = lm%traj%u
238 vd(state%isc:state%iec,state%jsc:state%jec,:) = lm%traj%v
239 t = lm%traj%t
240 delp = lm%traj%delp
241 q = lm%traj%qv
242 qi = lm%traj%qi
243 ql = lm%traj%ql
244 o3 = lm%traj%o3
245 
246 if (state%has_field('ua')) then
247  call state%get_field('ua', ua )
248  call state%get_field('va', va )
249  ua = lm%traj%ua
250  va = lm%traj%va
251 endif
252 
253 if (.not. lm%conf%hydrostatic) then
254  call state%get_field('w ', w )
255  call state%get_field('delz', delz)
256  w = lm%traj%w
257  delz = lm%traj%delz
258 endif
259 
260 call state%get_field('phis', phis )
261 phis(:,:,1) = lm%traj%phis
262 
263 end subroutine lm_to_state
264 
265 ! --------------------------------------------------------------------------------------------------
266 
267 end module fv3jedi_fv3lm_mod
fv3jedi_fv3lm_mod::fv3lm_model
Definition: fv3jedi_fv3lm_mod.f90:29
fv3jedi_fv3lm_mod::state_to_lm
subroutine state_to_lm(state, lm)
Definition: fv3jedi_fv3lm_mod.f90:147
fv3jedi_state_mod::fv3jedi_state
Fortran derived type to hold FV3JEDI state.
Definition: fv3jedi_state_mod.F90:30
fv3jedi_fv3lm_mod::finalize
subroutine finalize(self, state)
Definition: fv3jedi_fv3lm_mod.f90:135
fv3jedi_fv3lm_mod::initialize
subroutine initialize(self, state)
Definition: fv3jedi_fv3lm_mod.f90:107
fv3jedi_state_mod
Definition: fv3jedi_state_mod.F90:6
fv3jedi_geom_mod
Fortran module handling geometry for the FV3 model.
Definition: fv3jedi_geom_mod.f90:8
fv3jedi_fv3lm_mod::step
subroutine step(self, state, geom, vdate)
Definition: fv3jedi_fv3lm_mod.f90:119
fv3jedi_fv3lm_mod::create
subroutine create(self, geom, conf)
Definition: fv3jedi_fv3lm_mod.f90:46
fv3jedi_geom_mod::fv3jedi_geom
Fortran derived type to hold geometry data for the FV3JEDI model.
Definition: fv3jedi_geom_mod.f90:46
fv3jedi_fv3lm_mod
Definition: fv3jedi_fv3lm_mod.f90:6
fv3jedi_kinds_mod::kind_real
integer, parameter, public kind_real
Definition: fv3jedi_kinds_mod.f90:14
fv3jedi_kinds_mod
Definition: fv3jedi_kinds_mod.f90:6
fv3jedi_fv3lm_mod::delete
subroutine delete(self)
Definition: fv3jedi_fv3lm_mod.f90:94
fv3jedi_fv3lm_mod::lm_to_state
subroutine lm_to_state(lm, state)
Definition: fv3jedi_fv3lm_mod.f90:209