FV3 Bundle
fv3jedi_state_mod.f90
Go to the documentation of this file.
1 ! (C) Copyright 2017-2018 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 
6 !> Handle state for the FV3JEDI odel
7 
9 
10 use iso_c_binding
11 use config_mod
12 use datetime_mod
13 
14 use ioda_locs_mod
15 use ufo_vars_mod
17 
26 
27 implicit none
28 
29 private
30 public :: create, delete, zeros, copy, axpy, add_incr, &
33 public :: fv3jedi_state
34 
35 ! ------------------------------------------------------------------------------
36 
37 contains
38 
39 ! ------------------------------------------------------------------------------
40 
41 subroutine create(self, geom, vars)
42 
43 implicit none
44 type(fv3jedi_state), intent(inout) :: self
45 type(fv3jedi_geom), target, intent(in) :: geom
46 type(fv3jedi_vars), intent(in) :: vars
47 
48 integer :: var
49 
50 ! Copy the variable names
51 self%vars%nv = vars%nv
52 allocate(self%vars%fldnames(self%vars%nv))
53 self%vars%fldnames = vars%fldnames
54 
55 ! Allocate variables based on names
56 do var = 1, self%vars%nv
57 
58  select case (trim(self%vars%fldnames(var)))
59 
60  case("ud")
61  if (.not.allocated( self%ud)) allocate ( self%ud(geom%isc:geom%iec, geom%jsc:geom%jec+1, geom%npz))
62  case("vd")
63  if (.not.allocated( self%vd)) allocate ( self%vd(geom%isc:geom%iec+1,geom%jsc:geom%jec , geom%npz))
64  case("ua")
65  if (.not.allocated( self%ua)) allocate ( self%ua(geom%isc:geom%iec, geom%jsc:geom%jec , geom%npz))
66  case("va")
67  if (.not.allocated( self%va)) allocate ( self%va(geom%isc:geom%iec, geom%jsc:geom%jec , geom%npz))
68  case("t")
69  if (.not.allocated( self%t)) allocate ( self%t(geom%isc:geom%iec, geom%jsc:geom%jec , geom%npz))
70  case("delp")
71  if (.not.allocated(self%delp)) allocate (self%delp(geom%isc:geom%iec, geom%jsc:geom%jec , geom%npz))
72  case("q")
73  if (.not.allocated( self%q)) allocate ( self%q(geom%isc:geom%iec, geom%jsc:geom%jec , geom%npz))
74  case("qi")
75  if (.not.allocated( self%qi)) allocate ( self%qi(geom%isc:geom%iec, geom%jsc:geom%jec , geom%npz))
76  case("ql")
77  if (.not.allocated( self%ql)) allocate ( self%ql(geom%isc:geom%iec, geom%jsc:geom%jec , geom%npz))
78  case("o3")
79  if (.not.allocated( self%o3)) allocate ( self%o3(geom%isc:geom%iec, geom%jsc:geom%jec , geom%npz))
80  case("w")
81  if (.not.allocated( self%w)) allocate ( self%w(geom%isc:geom%iec, geom%jsc:geom%jec , geom%npz))
82  case("delz")
83  if (.not.allocated(self%delz)) allocate (self%delz(geom%isc:geom%iec, geom%jsc:geom%jec , geom%npz))
84  case("ps")
85  case default
86  call abor1_ftn("Create: unknown variable "//trim(self%vars%fldnames(var)))
87 
88  end select
89 
90 enddo
91 
92 self%hydrostatic = .true.
93 if (allocated(self%w).and.allocated(self%delz)) self%hydrostatic = .false.
94 
95 if (.not.allocated(self%phis)) allocate(self%phis(geom%isc:geom%iec,geom%jsc:geom%jec))
96 
97 !CRTM surface variables
98 if (.not.allocated(self%slmsk )) allocate(self%slmsk (geom%isc:geom%iec,geom%jsc:geom%jec))
99 if (.not.allocated(self%sheleg)) allocate(self%sheleg(geom%isc:geom%iec,geom%jsc:geom%jec))
100 if (.not.allocated(self%tsea )) allocate(self%tsea (geom%isc:geom%iec,geom%jsc:geom%jec))
101 if (.not.allocated(self%vtype )) allocate(self%vtype (geom%isc:geom%iec,geom%jsc:geom%jec))
102 if (.not.allocated(self%stype )) allocate(self%stype (geom%isc:geom%iec,geom%jsc:geom%jec))
103 if (.not.allocated(self%vfrac )) allocate(self%vfrac (geom%isc:geom%iec,geom%jsc:geom%jec))
104 if (.not.allocated(self%stc )) allocate(self%stc (geom%isc:geom%iec,geom%jsc:geom%jec,4))
105 if (.not.allocated(self%smc )) allocate(self%smc (geom%isc:geom%iec,geom%jsc:geom%jec,4))
106 if (.not.allocated(self%snwdph)) allocate(self%snwdph(geom%isc:geom%iec,geom%jsc:geom%jec))
107 if (.not.allocated(self%u_srf )) allocate(self%u_srf (geom%isc:geom%iec,geom%jsc:geom%jec))
108 if (.not.allocated(self%v_srf )) allocate(self%v_srf (geom%isc:geom%iec,geom%jsc:geom%jec))
109 if (.not.allocated(self%f10m )) allocate(self%f10m (geom%isc:geom%iec,geom%jsc:geom%jec))
110 
111 !Linearized model trajectory
112 if (.not.allocated(self%qls )) allocate(self%qls (geom%isc:geom%iec,geom%jsc:geom%jec,geom%npz))
113 if (.not.allocated(self%qcn )) allocate(self%qcn (geom%isc:geom%iec,geom%jsc:geom%jec,geom%npz))
114 if (.not.allocated(self%cfcn )) allocate(self%cfcn (geom%isc:geom%iec,geom%jsc:geom%jec,geom%npz))
115 if (.not.allocated(self%frocean)) allocate(self%frocean(geom%isc:geom%iec,geom%jsc:geom%jec))
116 if (.not.allocated(self%frland )) allocate(self%frland (geom%isc:geom%iec,geom%jsc:geom%jec))
117 if (.not.allocated(self%varflt )) allocate(self%varflt (geom%isc:geom%iec,geom%jsc:geom%jec))
118 if (.not.allocated(self%ustar )) allocate(self%ustar (geom%isc:geom%iec,geom%jsc:geom%jec))
119 if (.not.allocated(self%bstar )) allocate(self%bstar (geom%isc:geom%iec,geom%jsc:geom%jec))
120 if (.not.allocated(self%zpbl )) allocate(self%zpbl (geom%isc:geom%iec,geom%jsc:geom%jec))
121 if (.not.allocated(self%cm )) allocate(self%cm (geom%isc:geom%iec,geom%jsc:geom%jec))
122 if (.not.allocated(self%ct )) allocate(self%ct (geom%isc:geom%iec,geom%jsc:geom%jec))
123 if (.not.allocated(self%cq )) allocate(self%cq (geom%isc:geom%iec,geom%jsc:geom%jec))
124 if (.not.allocated(self%kcbl )) allocate(self%kcbl (geom%isc:geom%iec,geom%jsc:geom%jec))
125 if (.not.allocated(self%ts )) allocate(self%ts (geom%isc:geom%iec,geom%jsc:geom%jec))
126 if (.not.allocated(self%khl )) allocate(self%khl (geom%isc:geom%iec,geom%jsc:geom%jec))
127 if (.not.allocated(self%khu )) allocate(self%khu (geom%isc:geom%iec,geom%jsc:geom%jec))
128 
129 ! Initialize all domain arrays to zero
130 call zeros(self)
131 self%phis = 0.0_kind_real
132 
133 ! For convenience
134 self%isc = geom%isc
135 self%iec = geom%iec
136 self%jsc = geom%jsc
137 self%jec = geom%jec
138 self%isd = geom%isd
139 self%ied = geom%ied
140 self%jsd = geom%jsd
141 self%jed = geom%jed
142 self%npx = geom%npx
143 self%npy = geom%npy
144 self%npz = geom%npz
145 
146 self%ntile = geom%ntile
147 self%ntiles = geom%ntiles
148 
149 end subroutine create
150 
151 ! ------------------------------------------------------------------------------
152 
153 subroutine delete(self)
154 implicit none
155 type(fv3jedi_state), intent(inout) :: self
156 
157 if (allocated(self%ud )) deallocate (self%ud )
158 if (allocated(self%vd )) deallocate (self%vd )
159 if (allocated(self%ua )) deallocate (self%ua )
160 if (allocated(self%va )) deallocate (self%va )
161 if (allocated(self%t )) deallocate (self%t )
162 if (allocated(self%delp)) deallocate (self%delp)
163 if (allocated(self%q )) deallocate (self%q )
164 if (allocated(self%qi )) deallocate (self%qi )
165 if (allocated(self%ql )) deallocate (self%ql )
166 if (allocated(self%o3 )) deallocate (self%o3 )
167 if (allocated(self%phis)) deallocate (self%phis)
168 if (allocated(self%w )) deallocate (self%w )
169 if (allocated(self%delz)) deallocate (self%delz)
170 
171 if (allocated(self%slmsk )) deallocate(self%slmsk )
172 if (allocated(self%sheleg)) deallocate(self%sheleg)
173 if (allocated(self%tsea )) deallocate(self%tsea )
174 if (allocated(self%vtype )) deallocate(self%vtype )
175 if (allocated(self%stype )) deallocate(self%stype )
176 if (allocated(self%vfrac )) deallocate(self%vfrac )
177 if (allocated(self%stc )) deallocate(self%stc )
178 if (allocated(self%smc )) deallocate(self%smc )
179 if (allocated(self%snwdph)) deallocate(self%snwdph)
180 if (allocated(self%u_srf )) deallocate(self%u_srf )
181 if (allocated(self%v_srf )) deallocate(self%v_srf )
182 if (allocated(self%f10m )) deallocate(self%f10m )
183 
184 if (allocated(self%qls )) deallocate(self%qls )
185 if (allocated(self%qcn )) deallocate(self%qcn )
186 if (allocated(self%cfcn )) deallocate(self%cfcn )
187 if (allocated(self%frocean)) deallocate(self%frocean)
188 if (allocated(self%frland )) deallocate(self%frland )
189 if (allocated(self%varflt )) deallocate(self%varflt )
190 if (allocated(self%ustar )) deallocate(self%ustar )
191 if (allocated(self%bstar )) deallocate(self%bstar )
192 if (allocated(self%zpbl )) deallocate(self%zpbl )
193 if (allocated(self%cm )) deallocate(self%cm )
194 if (allocated(self%ct )) deallocate(self%ct )
195 if (allocated(self%cq )) deallocate(self%cq )
196 if (allocated(self%kcbl )) deallocate(self%kcbl )
197 if (allocated(self%ts )) deallocate(self%ts )
198 if (allocated(self%khl )) deallocate(self%khl )
199 if (allocated(self%khu )) deallocate(self%khu )
200 
201 end subroutine delete
202 
203 ! ------------------------------------------------------------------------------
204 
205 subroutine zeros(self)
206 implicit none
207 type(fv3jedi_state), intent(inout) :: self
208 
209 !Zero out the entire domain
210 
211 !Model
212 if(allocated(self%ud )) self%ud = 0.0_kind_real
213 if(allocated(self%vd )) self%vd = 0.0_kind_real
214 if(allocated(self%ua )) self%ua = 0.0_kind_real
215 if(allocated(self%va )) self%va = 0.0_kind_real
216 if(allocated(self%t )) self%t = 0.0_kind_real
217 if(allocated(self%delp)) self%delp = 0.0_kind_real
218 if(allocated(self%q )) self%q = 0.0_kind_real
219 if(allocated(self%qi )) self%qi = 0.0_kind_real
220 if(allocated(self%ql )) self%ql = 0.0_kind_real
221 if(allocated(self%o3 )) self%o3 = 0.0_kind_real
222 if(allocated(self%w )) self%w = 0.0_kind_real
223 if(allocated(self%delz)) self%delz = 0.0_kind_real
224 
225 end subroutine zeros
226 
227 ! ------------------------------------------------------------------------------
228 
229 subroutine copy(self,rhs)
230 implicit none
231 type(fv3jedi_state), intent(inout) :: self
232 type(fv3jedi_state), intent(in) :: rhs
233 
234 self%isc = rhs%isc
235 self%iec = rhs%iec
236 self%jsc = rhs%jsc
237 self%jec = rhs%jec
238 self%isd = rhs%isd
239 self%ied = rhs%ied
240 self%jsd = rhs%jsd
241 self%jed = rhs%jed
242 self%npx = rhs%npx
243 self%npy = rhs%npy
244 self%npz = rhs%npz
245 self%havecrtmfields = rhs%havecrtmfields
246 self%hydrostatic = rhs%hydrostatic
247 self%calendar_type = rhs%calendar_type
248 self%date = rhs%date
249 self%date_init = rhs%date_init
250 self%ntile = rhs%ntile
251 self%ntiles = rhs%ntiles
252 
253 if(allocated(self%ud ).and.allocated(rhs%ud )) self%ud = rhs%ud
254 if(allocated(self%vd ).and.allocated(rhs%vd )) self%vd = rhs%vd
255 if(allocated(self%ua ).and.allocated(rhs%ua )) self%ua = rhs%ua
256 if(allocated(self%va ).and.allocated(rhs%va )) self%va = rhs%va
257 if(allocated(self%t ).and.allocated(rhs%t )) self%t = rhs%t
258 if(allocated(self%delp).and.allocated(rhs%delp)) self%delp = rhs%delp
259 if(allocated(self%q ).and.allocated(rhs%q )) self%q = rhs%q
260 if(allocated(self%qi ).and.allocated(rhs%qi )) self%qi = rhs%qi
261 if(allocated(self%ql ).and.allocated(rhs%ql )) self%ql = rhs%ql
262 if(allocated(self%o3 ).and.allocated(rhs%o3 )) self%o3 = rhs%o3
263 if(allocated(self%w ).and.allocated(rhs%w )) self%w = rhs%w
264 if(allocated(self%delz).and.allocated(rhs%delz)) self%delz = rhs%delz
265 
266 self%phis = rhs%phis
267 self%slmsk = rhs%slmsk
268 self%sheleg = rhs%sheleg
269 self%tsea = rhs%tsea
270 self%vtype = rhs%vtype
271 self%stype = rhs%stype
272 self%vfrac = rhs%vfrac
273 self%stc = rhs%stc
274 self%smc = rhs%smc
275 self%snwdph = rhs%snwdph
276 self%u_srf = rhs%u_srf
277 self%v_srf = rhs%v_srf
278 self%f10m = rhs%f10m
279 
280 self%qls = rhs%qls
281 self%qcn = rhs%qcn
282 self%cfcn = rhs%cfcn
283 self%frocean = rhs%frocean
284 self%frland = rhs%frland
285 self%varflt = rhs%varflt
286 self%ustar = rhs%ustar
287 self%bstar = rhs%bstar
288 self%zpbl = rhs%zpbl
289 self%cm = rhs%cm
290 self%ct = rhs%ct
291 self%cq = rhs%cq
292 self%kcbl = rhs%kcbl
293 self%ts = rhs%ts
294 self%khl = rhs%khl
295 self%khu = rhs%khu
296 
297 return
298 end subroutine copy
299 
300 ! ------------------------------------------------------------------------------
301 
302 subroutine axpy(self,zz,rhs)
303 implicit none
304 type(fv3jedi_state), intent(inout) :: self
305 real(kind=kind_real), intent(in) :: zz
306 type(fv3jedi_state), intent(in) :: rhs
307 
308 if(allocated(self%ud ).and.allocated(rhs%ud )) self%ud = self%ud + zz * rhs%ud
309 if(allocated(self%vd ).and.allocated(rhs%vd )) self%vd = self%vd + zz * rhs%vd
310 if(allocated(self%ua ).and.allocated(rhs%ua )) self%ua = self%ua + zz * rhs%ua
311 if(allocated(self%va ).and.allocated(rhs%va )) self%va = self%va + zz * rhs%va
312 if(allocated(self%t ).and.allocated(rhs%t )) self%t = self%t + zz * rhs%t
313 if(allocated(self%delp).and.allocated(rhs%delp)) self%delp = self%delp + zz * rhs%delp
314 if(allocated(self%q ).and.allocated(rhs%q )) self%q = self%q + zz * rhs%q
315 if(allocated(self%qi ).and.allocated(rhs%qi )) self%qi = self%qi + zz * rhs%qi
316 if(allocated(self%ql ).and.allocated(rhs%ql )) self%ql = self%ql + zz * rhs%ql
317 if(allocated(self%o3 ).and.allocated(rhs%o3 )) self%o3 = self%o3 + zz * rhs%o3
318 if(allocated(self%w ).and.allocated(rhs%w )) self%w = self%w + zz * rhs%w
319 if(allocated(self%delz).and.allocated(rhs%delz)) self%delz = self%delz + zz * rhs%delz
320 
321 return
322 end subroutine axpy
323 
324 ! ------------------------------------------------------------------------------
325 
326 subroutine add_incr(geom,self,rhs)
328 use wind_vt_mod, only: a2d
329 
330 implicit none
331 type(fv3jedi_geom), intent(inout) :: geom
332 type(fv3jedi_state), intent(inout) :: self
333 type(fv3jedi_increment), intent(in) :: rhs
334 
335 integer :: isc,iec,jsc,jec,isd,ied,jsd,jed,npz,k
336 
337 real(kind=kind_real), allocatable, dimension(:,:,:) :: ud, vd
338 
339 !Check for matching resolution between state and increment
340 if ((rhs%iec-rhs%isc+1)-(self%iec-self%isc+1)==0) then
341 
342  isc = rhs%isc
343  iec = rhs%iec
344  jsc = rhs%jsc
345  jec = rhs%jec
346  isd = rhs%isd
347  ied = rhs%ied
348  jsd = rhs%jsd
349  jed = rhs%jed
350  npz = rhs%npz
351 
352  !Convert A-Grid increment to D-Grid
353  allocate(ud(isc:iec ,jsc:jec+1,1:npz))
354  allocate(vd(isc:iec+1,jsc:jec ,1:npz))
355  ud = 0.0_kind_real
356  vd = 0.0_kind_real
357 
358  call a2d(geom, rhs%ua(isc:iec,jsc:jec,1:npz), rhs%va(isc:iec,jsc:jec,1:npz), ud, vd)
359 
360  if(allocated(self%ud )) self%ud(isc:iec ,jsc:jec+1,:) = self%ud(isc:iec ,jsc:jec+1,:) + ud(isc:iec ,jsc:jec+1,:)
361  if(allocated(self%vd )) self%vd(isc:iec+1,jsc:jec ,:) = self%vd(isc:iec+1,jsc:jec ,:) + vd(isc:iec+1,jsc:jec ,:)
362 
363  deallocate(ud,vd)
364 
365  if(allocated(self%ua )) self%ua = self%ua + rhs%ua
366  if(allocated(self%va )) self%va = self%va + rhs%va
367  if(allocated(self%t )) self%t = self%t + rhs%t
368  if(allocated(self%delp)) then
369  if (allocated(rhs%ps)) then
370  do k = 1,geom%npz
371  self%delp(:,:,k) = self%delp(:,:,k) + (geom%bk(k+1)-geom%bk(k))*rhs%ps
372  enddo
373  elseif (allocated(rhs%delp)) then
374  self%delp = self%delp + rhs%delp
375  endif
376  endif
377  if(allocated(self%q )) self%q = self%q + rhs%q
378  if(allocated(self%qi )) self%qi = self%qi + rhs%qi
379  if(allocated(self%ql )) self%ql = self%ql + rhs%ql
380  if(allocated(self%o3 )) self%o3 = self%o3 + rhs%o3
381  if(allocated(self%w )) self%w = self%w + rhs%w
382  if(allocated(self%delz)) self%delz = self%delz + rhs%delz
383 else
384  call abor1_ftn("fv3jedi state: add_incr not implemented for low res increment yet")
385 endif
386 
387 return
388 end subroutine add_incr
389 
390 ! ------------------------------------------------------------------------------
391 
392 subroutine change_resol(state,rhs)
393 implicit none
394 type(fv3jedi_state), intent(inout) :: state
395 type(fv3jedi_state), intent(in) :: rhs
396 
397 integer :: check
398 check = (rhs%iec-rhs%isc+1) - (state%iec-state%isc+1)
399 
400 if (check==0) then
401  call copy(state, rhs)
402 else
403  call abor1_ftn("fv3jedi_state: change_resol not implmeneted yet")
404 endif
405 
406 return
407 end subroutine change_resol
408 
409 ! ------------------------------------------------------------------------------
410 !> Analytic Initialization for the FV3 Model
411 !!
412 !! \details **analytic_IC()** initializes the FV3JEDI state and State objects using one of
413 !! several alternative idealized analytic models. This is intended to facilitate testing by
414 !! eliminating the need to read in the initial state from a file and by providing exact expressions
415 !! to test interpolations. This function is activated by setting the "analytic_init" state in the
416 !! "initial" or "StateFile" section of the configuration file.
417 !!
418 !! Initialization options that begin with "dcmip" refer to tests defined by the multi-institutional
419 !! 2012 [Dynamical Core Intercomparison Project](https://earthsystealcmcog.org/projects/dcmip-2012)
420 !! and the associated Summer School, sponsored by NOAA, NSF, DOE, NCAR, and the University of Michigan.
421 !!
422 !! Currently implemented options for analytic_init include:
423 !! * invent-state: Backward compatibility with original analytic init option
424 !! * dcmip-test-1-1: 3D deformational flow
425 !! * dcmip-test-1-2: 3D Hadley-like meridional circulation
426 !! * dcmip-test-3-1: Non-hydrostatic gravity wave
427 !! * dcmip-test-4-0: Baroclinic instability
428 !!
429 !! \author M. Miesch (adapted from a pre-existing call to invent_state)
430 !! \date March, 2018: Created
431 !! \date May, 2018: Added dcmip-test-3-1
432 !! \date June, 2018: Added dcmip-test-4-0
433 !!
434 !! \warning This routine initializes the fv3jedi_state object. However, since the fv_atmos_type
435 !! component of fv3jedi_state is a subset of the corresponding object in the fv3 model,
436 !! this initialization routine is not sufficient to comprehensively define the full fv3 state.
437 !! So, this intitialization can be used for interpolation and other tests within JEDI but it is
438 !! cannot currently be used to initiate a forecast with gfs.
439 !!
440 !! \warning This routine does not initialize the fv3jedi_interp member of the fv3jedi_state object
441 !!
442 !! \warning Though an input state file is not required for these analytic initialization routines,
443 !! some grid information (in particular the hybrid vertical grid coefficients ak and bk)
444 !! is still read in from an input file when creating the geometry object that is a required
445 !! member of fv3jedi_state; see c_fv3jedi_geo_setup() in fv3jedi_geom_mod.F90.
446 !!
447 !! \warning It's unclear whether the pt member of the fv_atmos_type structure is potential temperature
448 !! or temperature. This routine assumes the latter. If this is not correct, then we will need to
449 !! implement a conversion
450 !!
451 subroutine analytic_ic(state, geom, c_conf, vdate)
454  use iso_c_binding
455  use datetime_mod
456  use fckit_log_module, only : log
457  use constants_mod, only: pi=>pi_8
461 
462  !FV3 Test Cases
466 
467  implicit none
468 
469  type(fv3jedi_state), intent(inout) :: state !< State
470  type(fv3jedi_geom), intent(inout) :: geom !< Geometry
471  type(c_ptr), intent(in) :: c_conf !< Configuration
472  type(datetime), intent(inout) :: vdate !< DateTime
473 
474  character(len=30) :: ic
475  character(len=20) :: sdate
476  character(len=1024) :: buf
477  Integer :: i,j,k
478  real(kind=kind_real) :: rlat, rlon, z
479  real(kind=kind_real) :: pk,pe1,pe2,ps
480  real(kind=kind_real) :: u0,v0,w0,t0,phis0,ps0,rho0,hum0,q1,q2,q3,q4
481 
482  type(fv_atmos_type), allocatable :: fv_atmic(:)
483  real(kind=kind_real) :: dtdummy = 900.0
484  logical, allocatable :: grids_on_this_pe(:)
485  integer :: p_split = 1
486 
487  If (config_element_exists(c_conf,"analytic_init")) Then
488  ic = trim(config_get_string(c_conf,len(ic),"analytic_init"))
489  Else
490  ! This default value is for backward compatibility
491  ic = "invent-state"
492  EndIf
493 
494  call log%warning("fv3jedi_state:analytic_init: "//ic)
495  sdate = config_get_string(c_conf,len(sdate),"date")
496  WRITE(buf,*) 'validity date is: '//sdate
497  call log%info(buf)
498  call datetime_set(sdate, vdate)
499 
500  !===================================================================
501  int_option: Select Case (ic)
502 
503  Case("invent-state")
504 
505  call invent_state(state,c_conf,geom)
506 
507  Case("fv3_init_case")
508 
509  !Initialize temporary FV_Atm fv3 construct
510  call fv_init(fv_atmic, dtdummy, grids_on_this_pe, p_split)
511  deallocate(pelist_all)
512 
513  !Test case to run, see fv3: /tools/test_cases.F90 for possibilities
514  test_case = config_get_int(c_conf,"fv3_test_case")
515 
516  call init_case( fv_atmic(1)%u,fv_atmic(1)%v,fv_atmic(1)%w,fv_atmic(1)%pt,fv_atmic(1)%delp,fv_atmic(1)%q, &
517  fv_atmic(1)%phis, fv_atmic(1)%ps,fv_atmic(1)%pe, fv_atmic(1)%peln,fv_atmic(1)%pk,fv_atmic(1)%pkz, &
518  fv_atmic(1)%uc,fv_atmic(1)%vc, fv_atmic(1)%ua,fv_atmic(1)%va, &
519  fv_atmic(1)%ak, fv_atmic(1)%bk, fv_atmic(1)%gridstruct, fv_atmic(1)%flagstruct,&
520  fv_atmic(1)%npx, fv_atmic(1)%npy, fv_atmic(1)%npz, fv_atmic(1)%ng, &
521  fv_atmic(1)%flagstruct%ncnst, fv_atmic(1)%flagstruct%nwat, &
522  fv_atmic(1)%flagstruct%ndims, fv_atmic(1)%flagstruct%ntiles, &
523  fv_atmic(1)%flagstruct%dry_mass, &
524  fv_atmic(1)%flagstruct%mountain, &
525  fv_atmic(1)%flagstruct%moist_phys, fv_atmic(1)%flagstruct%hydrostatic, &
526  fv_atmic(1)%flagstruct%hybrid_z, fv_atmic(1)%delz, fv_atmic(1)%ze0, &
527  fv_atmic(1)%flagstruct%adiabatic, fv_atmic(1)%ks, fv_atmic(1)%neststruct%npx_global, &
528  fv_atmic(1)%ptop, fv_atmic(1)%domain, fv_atmic(1)%tile, fv_atmic(1)%bd )
529 
530  !Copy from temporary structure into state
531  state%ud = fv_atmic(1)%u
532  state%vd = fv_atmic(1)%v
533  state%t = fv_atmic(1)%pt
534  state%delp = fv_atmic(1)%delp
535  state%q = fv_atmic(1)%q(:,:,:,1)
536  state%phis = fv_atmic(1)%phis
537  geom%ak = fv_atmic(1)%ak
538  geom%ak = fv_atmic(1)%ak
539  geom%ptop = fv_atmic(1)%ptop
540  if (.not. state%hydrostatic) then
541  state%w = fv_atmic(1)%w
542  state%delz = fv_atmic(1)%delz
543  endif
544 
545  !Deallocate temporary FV_Atm fv3 structure
546  call deallocate_fv_atmos_type(fv_atmic(1))
547  deallocate(fv_atmic)
548  deallocate(grids_on_this_pe)
549 
550  Case ("dcmip-test-1-1")
551 
552  do i = geom%isc,geom%iec
553  do j = geom%jsc,geom%jec
554  rlat = geom%grid_lat(i,j)
555  rlon = geom%grid_lon(i,j)
556 
557  ! Call the routine first just to get the surface pressure
558  Call test1_advection_deformation(rlon,rlat,pk,0.d0,1,u0,v0,w0,t0,&
559  phis0,ps,rho0,hum0,q1,q2,q3,q4)
560 
561  state%phis(i,j) = phis0
562 
563  ! Now loop over all levels
564  do k = 1, geom%npz
565 
566  pe1 = geom%ak(k) + geom%bk(k)*ps
567  pe2 = geom%ak(k+1) + geom%bk(k+1)*ps
568  pk = 0.5_kind_real * (pe1+pe2)
569  Call test1_advection_deformation(rlon,rlat,pk,0.d0,0,u0,v0,w0,t0,&
570  phis0,ps0,rho0,hum0,q1,q2,q3,q4)
571 
572  state%ud(i,j,k) = u0 !ATTN Not going to necessary keep a-grid winds, u can be either a or d grid
573  state%vd(i,j,k) = v0 !so this needs to be generic. You cannot drive the model with A grid winds
574  If (.not.state%hydrostatic) state%w(i,j,k) = w0
575  state%t(i,j,k) = t0
576  state%delp(i,j,k) = pe2-pe1
577  state%q (i,j,k) = hum0
578  state%qi(i,j,k) = q1
579  state%ql(i,j,k) = q2
580  state%o3(i,j,k) = q3
581 
582  enddo
583  enddo
584  enddo
585 
586  Case ("dcmip-test-1-2")
587 
588  do i = geom%isc,geom%iec
589  do j = geom%jsc,geom%jec
590  rlat = geom%grid_lat(i,j)
591  rlon = geom%grid_lon(i,j)
592 
593  ! Call the routine first just to get the surface pressure
594  Call test1_advection_hadley(rlon,rlat,pk,0.d0,1,u0,v0,w0,&
595  t0,phis0,ps,rho0,hum0,q1)
596 
597  state%phis(i,j) = phis0
598 
599  ! Now loop over all levels
600  do k = 1, geom%npz
601 
602  pe1 = geom%ak(k) + geom%bk(k)*ps
603  pe2 = geom%ak(k+1) + geom%bk(k+1)*ps
604  pk = 0.5_kind_real * (pe1+pe2)
605  Call test1_advection_hadley(rlon,rlat,pk,0.d0,0,u0,v0,w0,&
606  t0,phis0,ps,rho0,hum0,q1)
607 
608  state%ud(i,j,k) = u0 !ATTN comment above
609  state%vd(i,j,k) = v0
610  If (.not.state%hydrostatic) state%w(i,j,k) = w0
611  state%t(i,j,k) = t0
612  state%delp(i,j,k) = pe2-pe1
613  state%q (i,j,k) = hum0
614  state%qi(i,j,k) = q1
615 
616  enddo
617  enddo
618  enddo
619 
620  Case ("dcmip-test-3-1")
621 
622  do i = geom%isc,geom%iec
623  do j = geom%jsc,geom%jec
624  rlat = geom%grid_lat(i,j)
625  rlon = geom%grid_lon(i,j)
626 
627  ! Call the routine first just to get the surface pressure
628  Call test3_gravity_wave(rlon,rlat,pk,0.d0,1,u0,v0,w0,&
629  t0,phis0,ps,rho0,hum0)
630 
631  state%phis(i,j) = phis0
632 
633  ! Now loop over all levels
634  do k = 1, geom%npz
635 
636  pe1 = geom%ak(k) + geom%bk(k)*ps
637  pe2 = geom%ak(k+1) + geom%bk(k+1)*ps
638  pk = 0.5_kind_real * (pe1+pe2)
639  Call test3_gravity_wave(rlon,rlat,pk,0.d0,0,u0,v0,w0,&
640  t0,phis0,ps,rho0,hum0)
641 
642  state%ud(i,j,k) = u0 !ATTN comment above
643  state%vd(i,j,k) = v0
644  If (.not.state%hydrostatic) state%w(i,j,k) = w0
645  state%t(i,j,k) = t0
646  state%delp(i,j,k) = pe2-pe1
647  state%q(i,j,k) = hum0
648 
649  enddo
650  enddo
651  enddo
652 
653  Case ("dcmip-test-4-0")
654 
655  do i = geom%isc,geom%iec
656  do j = geom%jsc,geom%jec
657  rlat = geom%grid_lat(i,j)
658  rlon = geom%grid_lon(i,j)
659 
660  ! Call the routine first just to get the surface pressure
661  Call test4_baroclinic_wave(0,1.0_kind_real,rlon,rlat,pk,0.d0,1,u0,v0,w0,&
662  t0,phis0,ps,rho0,hum0,q1,q2)
663 
664  state%phis(i,j) = phis0
665 
666  ! Now loop over all levels
667  do k = 1, geom%npz
668 
669  pe1 = geom%ak(k) + geom%bk(k)*ps
670  pe2 = geom%ak(k+1) + geom%bk(k+1)*ps
671  pk = 0.5_kind_real * (pe1+pe2)
672  Call test4_baroclinic_wave(0,1.0_kind_real,rlon,rlat,pk,0.d0,0,u0,v0,w0,&
673  t0,phis0,ps,rho0,hum0,q1,q2)
674 
675  state%ud(i,j,k) = u0 !ATTN comment above
676  state%vd(i,j,k) = v0
677  If (.not.state%hydrostatic) state%w(i,j,k) = w0
678  state%t(i,j,k) = t0
679  state%delp(i,j,k) = pe2-pe1
680  state%q(i,j,k) = hum0
681 
682  enddo
683  enddo
684  enddo
685 
686  Case Default
687 
688  call invent_state(state,c_conf,geom)
689 
690  End Select int_option
691 
692 end subroutine analytic_ic
693 
694 ! ------------------------------------------------------------------------------
695 subroutine invent_state(state,config,geom)
698 
699 implicit none
700 
701 type(fv3jedi_state), intent(inout) :: state !< Model state
702 type(c_ptr), intent(in) :: config !< Configuration structure
703 type(fv3jedi_geom), intent(in) :: geom
704 
705 integer :: i,j,k
706 
707 !ud
708 do k = 1,geom%npz
709  do j = geom%jsc,geom%jec
710  do i = geom%isc,geom%iec
711  state%ud(i,j,k) = cos(0.25*geom%grid_lon(i,j)) + cos(0.25*geom%grid_lat(i,j))
712  enddo
713  enddo
714 enddo
715 
716 !vd
717 do k = 1,geom%npz
718  do j = geom%jsc,geom%jec
719  do i = geom%isc,geom%iec
720  state%vd(i,j,k) = 1.0_kind_real
721  enddo
722  enddo
723 enddo
724 
725 !ua
726 do k = 1,geom%npz
727  do j = geom%jsc,geom%jec
728  do i = geom%isc,geom%iec
729  state%ua(i,j,k) = cos(0.25*geom%grid_lon(i,j)) + cos(0.25*geom%grid_lat(i,j))
730  enddo
731  enddo
732 enddo
733 
734 !va
735 do k = 1,geom%npz
736  do j = geom%jsc,geom%jec
737  do i = geom%isc,geom%iec
738  state%va(i,j,k) = 1.0_kind_real
739  enddo
740  enddo
741 enddo
742 
743 !t
744 do k = 1,geom%npz
745  do j = geom%jsc,geom%jec
746  do i = geom%isc,geom%iec
747  state%t(i,j,k) = cos(0.25*geom%grid_lon(i,j)) + cos(0.25*geom%grid_lat(i,j))
748  enddo
749  enddo
750 enddo
751 
752 !delp
753 do k = 1,geom%npz
754  do j = geom%jsc,geom%jec
755  do i = geom%isc,geom%iec
756  state%delp(i,j,k) = k
757  enddo
758  enddo
759 enddo
760 
761 !q
762 do k = 1,geom%npz
763  do j = geom%jsc,geom%jec
764  do i = geom%isc,geom%iec
765  state%q(i,j,k) = 0.0
766  enddo
767  enddo
768 enddo
769 
770 !qi
771 do k = 1,geom%npz
772  do j = geom%jsc,geom%jec
773  do i = geom%isc,geom%iec
774  state%q(i,j,k) = 0.0
775  enddo
776  enddo
777 enddo
778 
779 !ql
780 do k = 1,geom%npz
781  do j = geom%jsc,geom%jec
782  do i = geom%isc,geom%iec
783  state%q(i,j,k) = 0.0
784  enddo
785  enddo
786 enddo
787 
788 !o3
789 do k = 1,geom%npz
790  do j = geom%jsc,geom%jec
791  do i = geom%isc,geom%iec
792  state%q(i,j,k) = 0.0
793  enddo
794  enddo
795 enddo
796 
797 return
798 end subroutine invent_state
799 
800 ! ------------------------------------------------------------------------------
801 
802 subroutine read_file(geom, state, c_conf, vdate)
804  implicit none
805 
806  type(fv3jedi_geom), intent(inout) :: geom !< Geometry
807  type(fv3jedi_state), intent(inout) :: state !< State
808  type(c_ptr), intent(in) :: c_conf !< Configuration
809  type(datetime), intent(inout) :: vdate !< DateTime
810 
811  character(len=10) :: restart_type
812 
813  restart_type = config_get_string(c_conf,len(restart_type),"restart_type")
814 
815  if (trim(restart_type) == 'gfs') then
816  call read_fms_restart(geom, state, c_conf, vdate)
817  elseif (trim(restart_type) == 'geos') then
818  call read_geos_restart(state, c_conf, vdate)
819  else
820  call abor1_ftn("fv3jedi_state read: restart type not supported")
821  endif
822 
823  return
824 
825 end subroutine read_file
826 
827 ! ------------------------------------------------------------------------------
828 
829 subroutine write_file(geom, state, c_conf, vdate)
831  implicit none
832 
833  type(fv3jedi_geom), intent(inout) :: geom !< Geometry
834  type(fv3jedi_state), intent(in) :: state !< State
835  type(c_ptr), intent(in) :: c_conf !< Configuration
836  type(datetime), intent(inout) :: vdate !< DateTime
837 
838  character(len=10) :: restart_type
839 
840  restart_type = config_get_string(c_conf,len(restart_type),"restart_type")
841 
842  if (trim(restart_type) == 'gfs') then
843  call write_fms_restart(geom, state, c_conf, vdate)
844  elseif (trim(restart_type) == 'geos') then
845  call write_geos_restart(geom, state, c_conf, vdate)
846  else
847  call abor1_ftn("fv3jedi_state write: restart type not supported")
848  endif
849 
850  return
851 
852 end subroutine write_file
853 
854 ! ------------------------------------------------------------------------------
855 
856 subroutine gpnorm(state, nf, pstat)
857 implicit none
858 type(fv3jedi_state), intent(in) :: state
859 integer, intent(in) :: nf
860 real(kind=kind_real), intent(inout) :: pstat(3, nf)
861 
862 integer :: isc, iec, jsc, jec, gs
863 
864 !1. Min
865 !2. Max
866 !3. RMS
867 
868 isc = state%isc
869 iec = state%iec
870 jsc = state%jsc
871 jec = state%jec
872 
873 gs = (iec-isc+1)*(jec-jsc+1)*state%npz
874 
875 pstat = 0.0_kind_real
876 
877 !ud
878 if (allocated(state%ud)) then
879  pstat(1,1) = minval(state%ud(isc:iec,jsc:jec,:))
880  pstat(2,1) = maxval(state%ud(isc:iec,jsc:jec,:))
881  pstat(3,1) = sqrt((sum(state%ud(isc:iec,jsc:jec,:))/gs)**2)
882 endif
883 
884 !vd
885 if (allocated(state%vd)) then
886  pstat(1,2) = minval(state%vd(isc:iec,jsc:jec,:))
887  pstat(2,2) = maxval(state%vd(isc:iec,jsc:jec,:))
888  pstat(3,2) = sqrt((sum(state%vd(isc:iec,jsc:jec,:))/gs)**2)
889 endif
890 
891 !ua
892 if (allocated(state%ua)) then
893  pstat(1,3) = minval(state%ua(isc:iec,jsc:jec,:))
894  pstat(2,3) = maxval(state%ua(isc:iec,jsc:jec,:))
895  pstat(3,3) = sqrt((sum(state%ua(isc:iec,jsc:jec,:))/gs)**2)
896 endif
897 
898 !va
899 if (allocated(state%va)) then
900  pstat(1,4) = minval(state%va(isc:iec,jsc:jec,:))
901  pstat(2,4) = maxval(state%va(isc:iec,jsc:jec,:))
902  pstat(3,4) = sqrt((sum(state%va(isc:iec,jsc:jec,:))/gs)**2)
903 endif
904 
905 !t
906 if (allocated(state%t)) then
907  pstat(1,5) = minval(state%t(isc:iec,jsc:jec,:))
908  pstat(2,5) = maxval(state%t(isc:iec,jsc:jec,:))
909  pstat(3,5) = sqrt((sum(state%t(isc:iec,jsc:jec,:))/gs)**2)
910 endif
911 
912 !delp
913 if (allocated(state%delp)) then
914  pstat(1,6) = minval(state%delp(isc:iec,jsc:jec,:))
915  pstat(2,6) = maxval(state%delp(isc:iec,jsc:jec,:))
916  pstat(3,6) = sqrt((sum(state%delp(isc:iec,jsc:jec,:))/gs)**2)
917 endif
918 
919 !q
920 if (allocated(state%q)) then
921  pstat(1,7) = minval(state%q(isc:iec,jsc:jec,:))
922  pstat(2,7) = maxval(state%q(isc:iec,jsc:jec,:))
923  pstat(3,7) = sqrt((sum(state%q(isc:iec,jsc:jec,:))/gs)**2)
924 endif
925 
926 !qi
927 if (allocated(state%qi)) then
928  pstat(1,8) = minval(state%qi(isc:iec,jsc:jec,:))
929  pstat(2,8) = maxval(state%qi(isc:iec,jsc:jec,:))
930  pstat(3,8) = sqrt((sum(state%qi(isc:iec,jsc:jec,:))/gs)**2)
931 endif
932 
933 !ql
934 if (allocated(state%ql)) then
935  pstat(1,9) = minval(state%ql(isc:iec,jsc:jec,:))
936  pstat(2,9) = maxval(state%ql(isc:iec,jsc:jec,:))
937  pstat(3,9) = sqrt((sum(state%ql(isc:iec,jsc:jec,:))/gs)**2)
938 endif
939 
940 !o3
941 if (allocated(state%o3)) then
942  pstat(1,10) = minval(state%o3(isc:iec,jsc:jec,:))
943  pstat(2,10) = maxval(state%o3(isc:iec,jsc:jec,:))
944  pstat(3,10) = sqrt((sum(state%o3(isc:iec,jsc:jec,:))/gs)**2)
945 endif
946 
947 return
948 
949 end subroutine gpnorm
950 
951 ! ------------------------------------------------------------------------------
952 
953 subroutine staterms(state, prms)
954 use fckit_mpi_module, only : fckit_mpi_comm, fckit_mpi_sum
955 implicit none
956 type(fv3jedi_state), intent(in) :: state
957 real(kind=kind_real), intent(out) :: prms
958 
959 real(kind=kind_real) :: zz
960 integer i,j,k,ii,nt,ierr,npes,iisum
961 integer :: isc,iec,jsc,jec,npz
962 type(fckit_mpi_comm) :: f_comm
963 
964 isc = state%isc
965 iec = state%iec
966 jsc = state%jsc
967 jec = state%jec
968 npz = state%npz
969 
970 f_comm = fckit_mpi_comm()
971 
972 zz = 0.0_kind_real
973 prms = 0.0_kind_real
974 ii = 0
975 
976 !ud
977 if (allocated(state%ud)) then
978  do k = 1,npz
979  do j = jsc,jec
980  do i = isc,iec
981  zz = zz + state%ud(i,j,k)**2
982  ii = ii + 1
983  enddo
984  enddo
985  enddo
986 endif
987 
988 !vd
989 if (allocated(state%vd)) then
990  do k = 1,npz
991  do j = jsc,jec
992  do i = isc,iec
993  zz = zz + state%vd(i,j,k)**2
994  ii = ii + 1
995  enddo
996  enddo
997  enddo
998 endif
999 
1000 !ua
1001 if (allocated(state%ua)) then
1002  do k = 1,npz
1003  do j = jsc,jec
1004  do i = isc,iec
1005  zz = zz + state%ua(i,j,k)**2
1006  ii = ii + 1
1007  enddo
1008  enddo
1009  enddo
1010 endif
1011 
1012 !va
1013 if (allocated(state%va)) then
1014  do k = 1,npz
1015  do j = jsc,jec
1016  do i = isc,iec
1017  zz = zz + state%va(i,j,k)**2
1018  ii = ii + 1
1019  enddo
1020  enddo
1021  enddo
1022 endif
1023 
1024 !t
1025 if (allocated(state%t)) then
1026  do k = 1,npz
1027  do j = jsc,jec
1028  do i = isc,iec
1029  zz = zz + state%t(i,j,k)**2
1030  ii = ii + 1
1031  enddo
1032  enddo
1033  enddo
1034 endif
1035 
1036 !delp
1037 if (allocated(state%delp)) then
1038  do k = 1,npz
1039  do j = jsc,jec
1040  do i = isc,iec
1041  zz = zz + state%delp(i,j,k)**2
1042  ii = ii + 1
1043  enddo
1044  enddo
1045  enddo
1046 endif
1047 
1048 !q
1049 if (allocated(state%q)) then
1050  do k = 1,npz
1051  do j = jsc,jec
1052  do i = isc,iec
1053  zz = zz + state%q(i,j,k)**2
1054  ii = ii + 1
1055  enddo
1056  enddo
1057  enddo
1058 endif
1059 
1060 !qi
1061 if (allocated(state%qi)) then
1062  do k = 1,npz
1063  do j = jsc,jec
1064  do i = isc,iec
1065  zz = zz + state%qi(i,j,k)**2
1066  ii = ii + 1
1067  enddo
1068  enddo
1069  enddo
1070 endif
1071 
1072 !ql
1073 if (allocated(state%ql)) then
1074  do k = 1,npz
1075  do j = jsc,jec
1076  do i = isc,iec
1077  zz = zz + state%ql(i,j,k)**2
1078  ii = ii + 1
1079  enddo
1080  enddo
1081  enddo
1082 endif
1083 
1084 !o3
1085 if (allocated(state%o3)) then
1086  do k = 1,npz
1087  do j = jsc,jec
1088  do i = isc,iec
1089  zz = zz + state%o3(i,j,k)**2
1090  ii = ii + 1
1091  enddo
1092  enddo
1093  enddo
1094 endif
1095 
1096 !Get global values
1097 call f_comm%allreduce(zz,prms,fckit_mpi_sum())
1098 call f_comm%allreduce(ii,iisum,fckit_mpi_sum())
1099 
1100 !if (ierr .ne. 0) then
1101 ! print *,'error in staterms/mpi_allreduce, error code=',ierr
1102 !endif
1103 prms = sqrt(prms/real(iisum,kind_real))
1104 
1105 end subroutine staterms
1106 
1107 ! ------------------------------------------------------------------------------
1108 
1109 end module fv3jedi_state_mod
subroutine, public gpnorm(inc, nf, pstat)
subroutine, public staterms(state, prms)
Fortran derived type to hold FV3JEDI increment.
subroutine, public change_resol(inc, rhs)
real(kind=8), parameter, public pi_8
Ratio of circle circumference to diameter [N/A].
Definition: constants.F90:73
real(kind=kind_real), parameter, public rad2deg
subroutine, public write_file(geom, inc, c_conf, vdate)
subroutine, public init_case(u, v, w, pt, delp, q, phis, ps, pe, peln, pk, pkz, uc, vc, ua, va, ak, bk, gridstruct, flagstruct, npx, npy, npz, ng, ncnst, nwat, ndims, nregions, dry_mass, mountain, moist_phys, hydrostatic, hybrid_z, delz, ze0, adiabatic, ks, npx_global, ptop, domain_in, tile_in, bd)
subroutine, public copy(self, rhs)
Fortran derived type to hold FV3JEDI state.
subroutine, public getvalues(geom, state, locs, vars, gom, traj)
subroutine, public read_fms_restart(geom, state, c_conf, vdate)
subroutine test1_advection_deformation(lon, lat, p, z, zcoords, u, v, w, t, phis, ps, rho, q, q1, q2, q3, q4)
integer, dimension(:), allocatable, public pelist_all
subroutine test4_baroclinic_wave(moist, X, lon, lat, p, z, zcoords, u, v, w, t, phis, ps, rho, q, q1, q2)
Fortran derived type to hold geometry data for the FV3JEDI model.
subroutine, public read_geos_restart(state, c_conf, vdate)
Fortran module to handle variables for the FV3JEDI model.
subroutine, public delete(self)
subroutine, public fv_init(Atm, dt_atmos, grids_on_this_pe, p_split)
Handle state for the FV3JEDI odel.
subroutine, public write_fms_restart(geom, state, c_conf, vdate)
subroutine invent_state(state, config, geom)
subroutine, public a2d(geom, ua, va, ud, vd)
subroutine, public add_incr(self, rhs)
Utilities for increment for the FV3JEDI model.
Variable transforms on wind variables for fv3-jedi Daniel Holdaway, NASA/JCSDA.
Utilities for state for the FV3JEDI model.
subroutine test1_advection_hadley(lon, lat, p, z, zcoords, u, v, w, t, phis, ps, rho, q, q1)
subroutine, public read_file(geom, inc, c_conf, vdate)
subroutine, public analytic_ic(state, geom, c_conf, vdate)
Analytic Initialization for the FV3 Model.
Fortran module handling geometry for the FV3 model.
subroutine, public zeros(self)
subroutine deallocate_fv_atmos_type(Atm)
subroutine test3_gravity_wave(lon, lat, p, z, zcoords, u, v, w, t, phis, ps, rho, q)
real(kind=kind_real), parameter, public constoz
integer, public test_case
Fortran derived type to represent fv3jedi model variables.
Fortran module handling observation locations.
integer, parameter, public kind_real
subroutine, public axpy(self, zz, rhs)
subroutine, public create(self, geom, vars)
subroutine, public write_geos_restart(geom, state, c_conf, vdate)
real(fp), parameter, public pi