FV3-JEDI
fv3jedi_tlm_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 ! iso
9 use iso_c_binding
10 
11 ! fckit
12 use fckit_configuration_module, only: fckit_configuration
13 
14 ! oops
15 use duration_mod
16 use oops_variables_mod
17 
18 ! fv3-jedi-linearmodel
19 use fv3jedi_lm_mod, only: fv3jedi_lm_type
20 
21 ! fv3-jedi
26 use fv3jedi_traj_mod, only: fv3jedi_traj
27 
28 implicit none
29 private
30 public :: fv3jedi_tlm
31 
32 ! --------------------------------------------------------------------------------------------------
33 
34 !> Fortran derived type to hold tlm definition
35 type:: fv3jedi_tlm
36  type(fv3jedi_lm_type) :: fv3jedi_lm !<Linearized model object
37  contains
38  procedure :: create
39  procedure :: delete
40  procedure :: initialize_tl
41  procedure :: initialize_ad
42  procedure :: step_tl
43  procedure :: step_ad
44  procedure :: finalize_tl
45  procedure :: finalize_ad
46 end type fv3jedi_tlm
47 
48 ! --------------------------------------------------------------------------------------------------
49 
50 contains
51 
52 ! --------------------------------------------------------------------------------------------------
53 
54 subroutine create(self, geom, conf)
55 
56 implicit none
57 class(fv3jedi_tlm), intent(inout) :: self
58 type(fv3jedi_geom), intent(in) :: geom
59 type(fckit_configuration), intent(in) :: conf
60 
61 !Locals
62 character(len=20) :: ststep
63 type(duration) :: dtstep
64 real(kind=kind_real) :: dt
65 character(len=:), allocatable :: str
66 
67 ! Model time step
68 ! ---------------
69 call conf%get_or_die("tstep",str)
70 ststep = str
71 deallocate(str)
72 dtstep = trim(ststep)
73 dt = real(duration_seconds(dtstep),kind_real)
74 
75 ! Model configuration and creation
76 ! --------------------------------
77 call conf%get_or_die("lm_do_dyn",self%fv3jedi_lm%conf%do_dyn)
78 call conf%get_or_die("lm_do_trb",self%fv3jedi_lm%conf%do_phy_trb)
79 call conf%get_or_die("lm_do_mst",self%fv3jedi_lm%conf%do_phy_mst)
80 
81 call self%fv3jedi_lm%create(dt,geom%npx,geom%npy,geom%npz,geom%ptop,geom%ak,geom%bk)
82 
83 end subroutine create
84 
85 ! --------------------------------------------------------------------------------------------------
86 
87 subroutine delete(self)
88 
89 implicit none
90 class(fv3jedi_tlm), intent(inout) :: self
91 
92 !Delete the model
93 !----------------
94 call self%fv3jedi_lm%delete()
95 
96 end subroutine delete
97 
98 ! --------------------------------------------------------------------------------------------------
99 
100 subroutine initialize_ad(self, inc)
101 
102 implicit none
103 class(fv3jedi_tlm), intent(inout) :: self
104 type(fv3jedi_increment), intent(inout) :: inc
105 
106 call inc_to_lm(inc,self%fv3jedi_lm)
107 call self%fv3jedi_lm%init_ad()
108 call lm_to_inc(self%fv3jedi_lm,inc)
109 
110 end subroutine initialize_ad
111 
112 ! --------------------------------------------------------------------------------------------------
113 
114 subroutine initialize_tl(self, inc)
115 
116 implicit none
117 class(fv3jedi_tlm), intent(inout) :: self
118 type(fv3jedi_increment), intent(inout) :: inc
119 
120 call inc_to_lm(inc,self%fv3jedi_lm)
121 call self%fv3jedi_lm%init_tl()
122 call lm_to_inc(self%fv3jedi_lm,inc)
123 
124 end subroutine initialize_tl
125 
126 ! --------------------------------------------------------------------------------------------------
127 
128 subroutine step_ad(self, inc, traj)
129 
130 implicit none
131 class(fv3jedi_tlm), intent(inout) :: self
132 type(fv3jedi_increment), intent(inout) :: inc
133 type(fv3jedi_traj), intent(in) :: traj
134 
135 call traj_to_traj(traj,self%fv3jedi_lm)
136 
137 call inc_to_lm(inc,self%fv3jedi_lm)
138 call self%fv3jedi_lm%step_ad()
139 call lm_to_inc(self%fv3jedi_lm,inc)
140 
141 end subroutine step_ad
142 
143 ! --------------------------------------------------------------------------------------------------
144 
145 subroutine step_tl(self, inc, traj)
146 
147 implicit none
148 class(fv3jedi_tlm), intent(inout) :: self
149 type(fv3jedi_increment), intent(inout) :: inc
150 type(fv3jedi_traj), intent(in) :: traj
151 
152 call traj_to_traj(traj,self%fv3jedi_lm)
153 
154 call inc_to_lm(inc,self%fv3jedi_lm)
155 call self%fv3jedi_lm%step_tl()
156 call lm_to_inc(self%fv3jedi_lm,inc)
157 
158 end subroutine step_tl
159 
160 ! --------------------------------------------------------------------------------------------------
161 
162 subroutine finalize_ad(self, inc)
163 
164 implicit none
165 class(fv3jedi_tlm), intent(inout) :: self
166 type(fv3jedi_increment), intent(inout) :: inc
167 
168 call inc_to_lm(inc,self%fv3jedi_lm)
169 call self%fv3jedi_lm%final_ad()
170 call lm_to_inc(self%fv3jedi_lm,inc)
171 
172 end subroutine finalize_ad
173 
174 ! --------------------------------------------------------------------------------------------------
175 
176 subroutine finalize_tl(self, inc)
177 
178 implicit none
179 class(fv3jedi_tlm), intent(inout) :: self
180 type(fv3jedi_increment), intent(inout) :: inc
181 
182 call inc_to_lm(inc,self%fv3jedi_lm)
183 call self%fv3jedi_lm%final_tl()
184 call lm_to_inc(self%fv3jedi_lm,inc)
185 
186 end subroutine finalize_tl
187 
188 ! --------------------------------------------------------------------------------------------------
189 
190 subroutine inc_to_lm(inc, lm)
191 
192 implicit none
193 type(fv3jedi_increment), intent(in) :: inc
194 type(fv3jedi_lm_type), intent(inout) :: lm
195 
196 real(kind=kind_real), pointer, dimension(:,:,:) :: ud
197 real(kind=kind_real), pointer, dimension(:,:,:) :: vd
198 
199 ! Bounds mismatch for D-Grid winds so first get pointer
200 call inc%get_field('ud', ud )
201 call inc%get_field('vd', vd )
202 lm%pert%u = ud(inc%isc:inc%iec,inc%jsc:inc%jec,1:inc%npz)
203 lm%pert%v = vd(inc%isc:inc%iec,inc%jsc:inc%jec,1:inc%npz)
204 
205 call inc%get_field('t' , lm%pert%t )
206 call inc%get_field('delp' , lm%pert%delp)
207 call inc%get_field('sphum' , lm%pert%qv )
208 call inc%get_field('ice_wat', lm%pert%qi )
209 call inc%get_field('liq_wat', lm%pert%ql )
210 
211 ! Optional fields
212 if (inc%has_field('o3mr' )) call inc%get_field('o3mr' , lm%pert%o3 )
213 if (inc%has_field('o3ppmv' )) call inc%get_field('o3ppmv' , lm%pert%o3 )
214 if (inc%has_field('w' )) call inc%get_field('w' , lm%pert%w )
215 if (inc%has_field('delz' )) call inc%get_field('delz' , lm%pert%delz)
216 
217 end subroutine inc_to_lm
218 
219 ! --------------------------------------------------------------------------------------------------
220 
221 subroutine lm_to_inc(lm, inc)
222 
223 implicit none
224 type(fv3jedi_lm_type), intent(in) :: lm
225 type(fv3jedi_increment), intent(inout) :: inc
226 
227 real(kind=kind_real), pointer, dimension(:,:,:) :: ud
228 real(kind=kind_real), pointer, dimension(:,:,:) :: vd
229 
230 ! Bounds mismatch for D-Grid winds so first get pointer
231 call inc%get_field('ud' , ud )
232 call inc%get_field('vd' , vd )
233 ud(inc%isc:inc%iec,inc%jsc:inc%jec,1:inc%npz) = lm%pert%u
234 vd(inc%isc:inc%iec,inc%jsc:inc%jec,1:inc%npz) = lm%pert%v
235 
236 call inc%put_field('t' , lm%pert%t )
237 call inc%put_field('delp' , lm%pert%delp)
238 call inc%put_field('sphum' , lm%pert%qv )
239 call inc%put_field('ice_wat', lm%pert%qi )
240 call inc%put_field('liq_wat', lm%pert%ql )
241 
242 ! Optional fields
243 if (inc%has_field('o3mr' )) call inc%put_field('o3mr' , lm%pert%o3 )
244 if (inc%has_field('o3ppmv' )) call inc%put_field('o3ppmv' , lm%pert%o3 )
245 if (inc%has_field('w' )) call inc%put_field('w' , lm%pert%w )
246 if (inc%has_field('delz' )) call inc%put_field('delz' , lm%pert%delz)
247 
248 end subroutine lm_to_inc
249 
250 ! --------------------------------------------------------------------------------------------------
251 
252 subroutine traj_to_traj( traj, lm )
253 
254 implicit none
255 type(fv3jedi_traj), intent(in) :: traj
256 type(fv3jedi_lm_type), intent(inout) :: lm
257 
258 lm%traj%u = traj%u
259 lm%traj%v = traj%v
260 lm%traj%ua = traj%ua
261 lm%traj%va = traj%va
262 lm%traj%t = traj%t
263 lm%traj%delp = traj%delp
264 lm%traj%qv = traj%qv
265 lm%traj%ql = traj%ql
266 lm%traj%qi = traj%qi
267 lm%traj%o3 = traj%o3
268 
269 if (.not. lm%conf%hydrostatic) then
270  lm%traj%w = traj%w
271  lm%traj%delz = traj%delz
272 endif
273 
274 if (lm%conf%do_phy_mst .ne. 0) then
275  lm%traj%qls = traj%qls
276  lm%traj%qcn = traj%qcn
277  lm%traj%cfcn = traj%cfcn
278 endif
279 
280 !> Rank two
281 lm%traj%phis = traj%phis
282 lm%traj%frocean = traj%frocean
283 lm%traj%frland = traj%frland
284 lm%traj%varflt = traj%varflt
285 lm%traj%ustar = traj%ustar
286 lm%traj%bstar = traj%bstar
287 lm%traj%zpbl = traj%zpbl
288 lm%traj%cm = traj%cm
289 lm%traj%ct = traj%ct
290 lm%traj%cq = traj%cq
291 lm%traj%kcbl = traj%kcbl
292 lm%traj%ts = traj%ts
293 lm%traj%khl = traj%khl
294 lm%traj%khu = traj%khu
295 
296 end subroutine traj_to_traj
297 
298 ! --------------------------------------------------------------------------------------------------
299 
300 end module fv3jedi_tlm_mod
fv3jedi_state_mod::fv3jedi_state
Fortran derived type to hold FV3JEDI state.
Definition: fv3jedi_state_mod.F90:30
fv3jedi_tlm_mod
Definition: fv3jedi_tlm_mod.f90:6
fv3jedi_tlm_mod::step_tl
subroutine step_tl(self, inc, traj)
Definition: fv3jedi_tlm_mod.f90:146
fv3jedi_tlm_mod::initialize_ad
subroutine initialize_ad(self, inc)
Definition: fv3jedi_tlm_mod.f90:101
fv3jedi_state_mod
Definition: fv3jedi_state_mod.F90:6
fv3jedi_traj_mod
Definition: fv3jedi_traj_mod.f90:6
fv3jedi_tlm_mod::fv3jedi_tlm
Fortran derived type to hold tlm definition.
Definition: fv3jedi_tlm_mod.f90:35
fv3jedi_tlm_mod::traj_to_traj
subroutine traj_to_traj(traj, lm)
Definition: fv3jedi_tlm_mod.f90:253
fv3jedi_geom_mod
Fortran module handling geometry for the FV3 model.
Definition: fv3jedi_geom_mod.f90:8
fv3jedi_tlm_mod::create
subroutine create(self, geom, conf)
Definition: fv3jedi_tlm_mod.f90:55
fv3jedi_increment_mod
Definition: fv3jedi_increment_mod.F90:6
fv3jedi_tlm_mod::lm_to_inc
subroutine lm_to_inc(lm, inc)
Definition: fv3jedi_tlm_mod.f90:222
fv3jedi_geom_mod::fv3jedi_geom
Fortran derived type to hold geometry data for the FV3JEDI model.
Definition: fv3jedi_geom_mod.f90:46
fv3jedi_tlm_mod::inc_to_lm
subroutine inc_to_lm(inc, lm)
Definition: fv3jedi_tlm_mod.f90:191
fv3jedi_tlm_mod::initialize_tl
subroutine initialize_tl(self, inc)
Definition: fv3jedi_tlm_mod.f90:115
fv3jedi_tlm_mod::step_ad
subroutine step_ad(self, inc, traj)
Definition: fv3jedi_tlm_mod.f90:129
fv3jedi_tlm_mod::finalize_ad
subroutine finalize_ad(self, inc)
Definition: fv3jedi_tlm_mod.f90:163
fv3jedi_tlm_mod::delete
subroutine delete(self)
Definition: fv3jedi_tlm_mod.f90:88
fv3jedi_kinds_mod::kind_real
integer, parameter, public kind_real
Definition: fv3jedi_kinds_mod.f90:14
fv3jedi_tlm_mod::finalize_tl
subroutine finalize_tl(self, inc)
Definition: fv3jedi_tlm_mod.f90:177
fv3jedi_kinds_mod
Definition: fv3jedi_kinds_mod.f90:6
fv3jedi_increment_mod::fv3jedi_increment
Definition: fv3jedi_increment_mod.F90:34