FV3-JEDI
fv3jedi_varcha_c2a_mod.f90
Go to the documentation of this file.
1 ! (C) Copyright 2018-2019 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 datetime_mod
16 
17 ! femps
18 use femps_grid_mod
19 use femps_operators_mod
20 use femps_testgrid_mod
21 use femps_solve_mod
22 use femps_fv3_mod
23 
24 ! fv3jedi
30 
34 use wind_vt_mod
35 
36 implicit none
37 private
38 
39 public :: fv3jedi_varcha_c2a
40 public :: create
41 public :: delete
42 public :: changevar
43 public :: changevarinverse
44 
46  type(fempsgrid) :: grid
47  type(fempsoprs) :: oprs
48  integer :: lprocs
49  integer, allocatable :: lev_start(:), lev_final(:)
50 end type fv3jedi_varcha_c2a
51 
52 ! ------------------------------------------------------------------------------
53 
54 contains
55 
56 ! ------------------------------------------------------------------------------
57 
58 subroutine create(self, geom, conf)
59 
60 implicit none
61 type(fv3jedi_varcha_c2a), intent(inout) :: self
62 type(fv3jedi_geom), intent(inout) :: geom
63 type(fckit_configuration), intent(in) :: conf
64 
65 integer :: ngrids, niter, lprocs, lstart
66 logical :: check_convergence
67 character(len=:), allocatable :: str
68 character(len=2055) :: path2fv3gridfiles
69 
70 integer :: n, levs_per_proc
71 
72 ! Grid and operators for the femps Poisson solver
73 ! -----------------------------------------------
74 
75 ! Configuration
76 call conf%get_or_die("femps_iterations",niter)
77 call conf%get_or_die("femps_ngrids",ngrids)
78 call conf%get_or_die("femps_path2fv3gridfiles",str); path2fv3gridfiles = str
79 if( .not. conf%get('femps_levelprocs',lprocs) ) then
80  lprocs = -1
81 endif
82 if( .not. conf%get('femps_checkconvergence',check_convergence) ) then
83  check_convergence = .false.
84 endif
85 
86 ! Processors that will do the work
87 ! --------------------------------
88 lprocs = min(lprocs,geom%f_comm%size())
89 lprocs = min(lprocs,geom%npz)
90 
91 if (lprocs == -1) then
92  self%lprocs = min(geom%npz,geom%f_comm%size())
93 else
94  self%lprocs = lprocs
95 endif
96 
97 if (geom%f_comm%rank() == 0 ) print*, "Running femps with ", self%lprocs, " processors."
98 
99 allocate(self%lev_start(self%lprocs))
100 allocate(self%lev_final(self%lprocs))
101 
102 if (self%lprocs == geom%npz) then
103  do n = 1,self%lprocs
104  self%lev_start(n) = n
105  self%lev_final(n) = n
106  enddo
107 else
108  levs_per_proc = floor(real(geom%npz,kind_real)/real(self%lprocs,kind_real))
109  lstart = 0
110  do n = 1,self%lprocs
111  self%lev_start(n) = lstart+1
112  self%lev_final(n) = self%lev_start(n) + levs_per_proc - 1
113  if (n .le. mod(geom%npz, self%lprocs)) self%lev_final(n) = self%lev_final(n) + 1
114  lstart = self%lev_final(n)
115  enddo
116 endif
117 
118 if (self%lev_final(self%lprocs) .ne. geom%npz) &
119  call abor1_ftn("fv3jedi_varcha_c2a_mod.create: last level not equal to number of levels.")
120 
121 ! Processors doing the work need grid and operators
122 if (geom%f_comm%rank() < self%lprocs ) then
123 
124  if (geom%f_comm%rank() == 0 ) print*, 'Creating FEMPS grid object'
125  call self%grid%setup('cs',ngrids=ngrids,cube=geom%npx-1,niter=niter,&
126  comm = geom%f_comm%communicator(), &
127  rank = geom%f_comm%rank(), &
128  csize = geom%f_comm%size(), &
129  check_convergence = check_convergence )
130 
131  if (geom%f_comm%rank() == 0 ) print*, 'Creating FEMPS grid hierarchy from files'
132  call fv3grid_to_ugrid(self%grid,path2fv3gridfiles)
133 
134  ! Build the connectivity and extra geom
135  if (geom%f_comm%rank() == 0 ) print*, 'Creating FEMPS cubed-sphere connectivity'
136  call self%grid%build_cs(1,1)
137 
138  ! Perform all the setup
139  if (geom%f_comm%rank() == 0 ) print*, 'Creating FEMPS static operators'
140  call preliminary(self%grid,self%oprs)
141 
142  ! Partial delete of operators not needed
143  if (geom%f_comm%rank() == 0 ) print*, 'FEMPS partial deallocate'
144  call self%oprs%pdelete()
145 
146 endif
147 
148 end subroutine create
149 
150 ! ------------------------------------------------------------------------------
151 
152 subroutine delete(self)
153 
154 implicit none
155 type(fv3jedi_varcha_c2a), intent(inout) :: self
156 
157 call self%oprs%delete()
158 call self%grid%delete()
159 
160 deallocate(self%lev_start,self%lev_final)
161 
162 end subroutine delete
163 
164 ! ------------------------------------------------------------------------------
165 
166 subroutine changevar(self,geom,xctl,xana)
167 
168 implicit none
169 type(fv3jedi_varcha_c2a), intent(inout) :: self
170 type(fv3jedi_geom), intent(inout) :: geom
171 type(fv3jedi_state), intent(in) :: xctl
172 type(fv3jedi_state), intent(inout) :: xana
173 
174 integer :: f
175 character(len=field_clen), allocatable :: fields_to_do(:)
176 real(kind=kind_real), pointer :: field_ptr(:,:,:)
177 
178 ! Stream function/velocity potential
179 logical :: have_udvd, have_vodi
180 real(kind=kind_real), pointer :: psi(:,:,:) ! Stream function
181 real(kind=kind_real), pointer :: chi(:,:,:) ! Velocity potentail
182 real(kind=kind_real), allocatable :: ud(:,:,:) ! D-grid u wind
183 real(kind=kind_real), allocatable :: vd(:,:,:) ! D-grid v wind
184 real(kind=kind_real), allocatable :: vort(:,:,:) ! Vorticity
185 real(kind=kind_real), allocatable :: divg(:,:,:) ! Divergence
186 
187 ! Pressure
188 logical :: have_pres
189 real(kind=kind_real), pointer :: ps(:,:,:) ! Surface pressure
190 real(kind=kind_real), allocatable :: delp(:,:,:) ! Pressure thickness
191 real(kind=kind_real), allocatable :: pe(:,:,:) ! Pressure edges
192 real(kind=kind_real), allocatable :: p(:,:,:) ! Pressure mid
193 real(kind=kind_real), allocatable :: pkz(:,:,:) ! Pressure ^ kappa
194 
195 ! Temperaure
196 logical :: have_t, have_pt
197 real(kind=kind_real), pointer :: t(:,:,:) ! Temperature
198 real(kind=kind_real), allocatable :: pt(:,:,:) ! Potential temperature
199 
200 ! Clouds
201 logical :: have_cld4
202 real(kind=kind_real), pointer :: qi(:,:,:) ! Cloud liquid ice
203 real(kind=kind_real), pointer :: ql(:,:,:) ! Cloud liquid water
204 real(kind=kind_real), pointer :: qilsf(:,:,:) ! Fraction ice large scale
205 real(kind=kind_real), pointer :: qicnf(:,:,:) ! Fraction ice convective
206 real(kind=kind_real), allocatable :: qils(:,:,:) ! Cloud liquid ice large scale
207 real(kind=kind_real), allocatable :: qicn(:,:,:) ! Cloud liquid ice convective
208 real(kind=kind_real), allocatable :: qlls(:,:,:) ! Cloud liquid water large scale
209 real(kind=kind_real), allocatable :: qlcn(:,:,:) ! Cloud liquid water convective
210 
211 ! Copy fields that are the same in both
212 ! -------------------------------------
213 call copy_subset(xctl%fields, xana%fields, fields_to_do)
214 
215 ! If variable change is the identity early exit
216 ! ---------------------------------------------
217 if (.not.allocated(fields_to_do)) return
218 
219 ! Wind variables
220 ! --------------
221 have_udvd = .false.
222 have_vodi = .false.
223 if (xctl%has_field('psi') .and. xctl%has_field('chi')) then
224  call xctl%get_field('psi', psi)
225  call xctl%get_field('chi', chi)
226  allocate(ud(geom%isc:geom%iec ,geom%jsc:geom%jec+1,1:geom%npz))
227  allocate(vd(geom%isc:geom%iec+1,geom%jsc:geom%jec ,1:geom%npz))
228  call psichi_to_udvd(geom, psi, chi, ud, vd)
229  have_udvd = .true.
230  allocate(vort(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
231  allocate(divg(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
232  call psichi_to_vortdivg(geom, self%grid, self%oprs, psi, chi, self%lprocs, self%lev_start, &
233  self%lev_final, vort, divg)
234  have_vodi = .true.
235 endif
236 
237 ! Pressure
238 ! --------
239 have_pres = .false.
240 if (xctl%has_field('ps')) then
241  call xctl%get_field('ps', ps)
242  allocate(delp(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
243  allocate( pe(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz+1))
244  allocate( p(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
245  allocate( pkz(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
246  call ps_to_delp(geom, ps, delp)
247  call delp_to_pe_p_logp(geom, delp, pe, p)
248  call pe_to_pkz(geom, pe, pkz)
249  have_pres = .true.
250 endif
251 
252 ! Temperature
253 ! -----------
254 have_t = .false.
255 have_pt = .false.
256 if (xctl%has_field('t')) then
257  call xctl%get_field('t', t)
258  have_t = .true.
259  if (have_pres) then
260  allocate(pt(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
261  call t_to_pt(geom, pkz, t, pt)
262  have_pt = .true.
263  endif
264 endif
265 
266 ! Clouds
267 ! ------
268 have_cld4 = .false.
269 if ( xctl%has_field('ice_wat') .and. xctl%has_field('liq_wat') .and. &
270  xctl%has_field('qilsf') .and. xctl%has_field('qicnf')) then
271  call xctl%get_field('ice_wat', qi)
272  call xctl%get_field('liq_wat', ql)
273  call xctl%get_field('qilsf', qilsf)
274  call xctl%get_field('qicnf', qicnf)
275  allocate(qils(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
276  allocate(qicn(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
277  allocate(qlls(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
278  allocate(qlcn(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
279  call q2_to_q4(geom, qi, ql, qilsf, qicnf, qils, qicn, qlls, qlcn)
280  have_cld4 = .true.
281 endif
282 
283 ! Loop over the fields not found in the input state and work through cases
284 ! ------------------------------------------------------------------------
285 do f = 1, size(fields_to_do)
286 
287  call xana%get_field(trim(fields_to_do(f)), field_ptr)
288 
289  select case (trim(fields_to_do(f)))
290 
291  case ("ud")
292 
293  if (.not. have_udvd) call field_fail(fields_to_do(f))
294  field_ptr = ud
295 
296  case ("vd")
297 
298  if (.not. have_udvd) call field_fail(fields_to_do(f))
299  field_ptr = vd
300 
301  case ("vort")
302 
303  if (.not. have_vodi) call field_fail(fields_to_do(f))
304  field_ptr = vort
305 
306  case ("divg")
307 
308  if (.not. have_vodi) call field_fail(fields_to_do(f))
309  field_ptr = divg
310 
311  case ("ps")
312 
313  if (.not. have_pres) call field_fail(fields_to_do(f))
314  field_ptr = ps
315 
316  case ("delp")
317 
318  if (.not. have_pres) call field_fail(fields_to_do(f))
319  field_ptr = delp
320 
321  case ("pkz")
322 
323  if (.not. have_pres) call field_fail(fields_to_do(f))
324  field_ptr = pkz
325 
326  case ("qils")
327 
328  if (.not. have_pres) call field_fail(fields_to_do(f))
329  field_ptr = qils
330 
331  case ("qicn")
332 
333  if (.not. have_pres) call field_fail(fields_to_do(f))
334  field_ptr = qicn
335 
336  case ("qlls")
337 
338  if (.not. have_pres) call field_fail(fields_to_do(f))
339  field_ptr = qlls
340 
341  case ("qlcn")
342 
343  if (.not. have_pres) call field_fail(fields_to_do(f))
344  field_ptr = qlcn
345 
346  end select
347 
348 enddo
349 
350 ! Copy calendar infomation
351 xana%calendar_type = xctl%calendar_type
352 xana%date_init = xctl%date_init
353 
354 end subroutine changevar
355 
356 ! ------------------------------------------------------------------------------
357 
358 subroutine changevarinverse(self,geom,xana,xctl)
359 
360 implicit none
361 type(fv3jedi_varcha_c2a), intent(inout) :: self
362 type(fv3jedi_geom), intent(inout) :: geom
363 type(fv3jedi_state), intent(in) :: xana
364 type(fv3jedi_state), intent(inout) :: xctl
365 
366 integer :: f
367 character(len=field_clen), allocatable :: fields_to_do(:)
368 real(kind=kind_real), pointer :: field_ptr(:,:,:)
369 
370 ! Stream function/velocity potential
371 logical :: have_pcvd
372 real(kind=kind_real), pointer :: ud(:,:,:) ! D-grid u wind
373 real(kind=kind_real), pointer :: vd(:,:,:) ! D-grid v wind
374 real(kind=kind_real), allocatable :: psi(:,:,:) ! Stream function
375 real(kind=kind_real), allocatable :: chi(:,:,:) ! Velocity potentail
376 real(kind=kind_real), allocatable :: vort(:,:,:) ! Vorticity
377 real(kind=kind_real), allocatable :: divg(:,:,:) ! Divergence
378 
379 ! Pressure
380 logical :: have_pres
381 real(kind=kind_real), allocatable :: delp(:,:,:) ! Pressure thickness
382 real(kind=kind_real), pointer :: ps(:,:,:) ! Pressure edges
383 real(kind=kind_real), pointer :: pe(:,:,:) ! Pressure mid
384 real(kind=kind_real), allocatable :: p(:,:,:) ! Pressure mid
385 
386 ! Temperature
387 logical :: have_temp
388 real(kind=kind_real), pointer :: pt(:,:,:) ! Potential temperature
389 real(kind=kind_real), allocatable :: pkz(:,:,:) ! Pressure ^ kappa
390 real(kind=kind_real), allocatable :: t(:,:,:) ! Temperature
391 
392 ! Humidity
393 logical :: have_rhum
394 real(kind=kind_real), pointer :: q(:,:,:) ! Specific humidity
395 real(kind=kind_real), allocatable :: qsat(:,:,:) ! Saturation specific humidity
396 real(kind=kind_real), allocatable :: rh(:,:,:) ! Relative humidity
397 
398 ! Clouds
399 logical :: have_qiql, have_cfrc
400 real(kind=kind_real), allocatable :: qi(:,:,:) ! Cloud liquid ice
401 real(kind=kind_real), allocatable :: ql(:,:,:) ! Cloud liquid water
402 real(kind=kind_real), pointer :: qils(:,:,:) ! Cloud liquid ice large scale
403 real(kind=kind_real), pointer :: qicn(:,:,:) ! Cloud liquid ice convective
404 real(kind=kind_real), pointer :: qlls(:,:,:) ! Cloud liquid water large scale
405 real(kind=kind_real), pointer :: qlcn(:,:,:) ! Cloud liquid water convective
406 real(kind=kind_real), allocatable :: qilsf(:,:,:) ! Fraction ice large scale
407 real(kind=kind_real), allocatable :: qicnf(:,:,:) ! Fraction ice convective
408 
409 ! Virtual temperature
410 logical :: have_virt
411 real(kind=kind_real), allocatable :: tv(:,:,:) ! Virtual temperature
412 
413 ! Copy fields that are the same in both
414 ! -------------------------------------
415 call copy_subset(xana%fields, xctl%fields, fields_to_do)
416 
417 ! If variable change is the identity early exit
418 ! ---------------------------------------------
419 if (.not.allocated(fields_to_do)) return
420 
421 ! Wind variables
422 ! --------------
423 have_pcvd = .false.
424 if (xana%has_field('ud') .and. xana%has_field('vd')) then
425  call xana%get_field('ud', ud)
426  call xana%get_field('vd', vd)
427  allocate( psi(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
428  allocate( chi(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
429  allocate(vort(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
430  allocate(divg(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
431  call udvd_to_psichi(geom, self%grid, self%oprs, ud, vd, psi, chi, &
432  self%lprocs, self%lev_start, self%lev_final, vort, divg)
433  have_pcvd = .true.
434 endif
435 
436 ! Pressure
437 ! --------
438 have_pres = .false.
439 if (xana%has_field('delp')) then
440  call xana%get_field('delp', delp)
441  allocate(ps(geom%isc:geom%iec,geom%jsc:geom%jec,1))
442  ps(:,:,1) = sum(delp,3)
443  have_pres = .true.
444 elseif (xana%has_field('pe')) then
445  call xana%get_field('pe', pe )
446  allocate( ps(geom%isc:geom%iec,geom%jsc:geom%jec,1))
447  allocate(delp(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
448  call pe_to_delp(geom, pe, delp)
449  ps(:,:,1) = pe(:,:,geom%npz+1)
450  have_pres = .true.
451 endif
452 
453 ! Temperature
454 ! -----------
455 have_temp = .false.
456 if (xana%has_field('t')) then
457  call xana%get_field('t', t)
458  have_temp = .true.
459 elseif (xana%has_field('pt')) then
460  allocate(t(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
461  call xana%get_field('pt', pt)
462  if (xana%has_field('pkz')) then
463  call xana%get_field('pkz', pkz)
464  have_temp = .true.
465  elseif (have_pres) then
466  allocate( pkz(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
467  call ps_to_pkz(geom, ps, pkz)
468  have_temp = .true.
469  endif
470  if (have_temp) call pt_to_t(geom, pkz, pt, t)
471 endif
472 
473 ! Humidity
474 ! --------
475 have_rhum = .false.
476 if (xana%has_field('sphum') .and. have_temp .and. have_pres) then
477  allocate(qsat(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
478  allocate( rh(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
479  call xana%get_field('sphum', q)
480  call get_qsat(geom, delp, t, q, qsat)
481  call q_to_rh(geom, qsat, q, rh)
482  have_rhum = .true.
483 endif
484 
485 ! Clouds
486 ! ------
487 have_qiql = .false.
488 have_cfrc = .false.
489 if (xana%has_field('ice_wat') .and. xana%has_field('liq_wat')) then
490  call xana%get_field('ice_wat', qi)
491  call xana%get_field('liq_wat', ql)
492  have_qiql = .true.
493 elseif (xana%has_field('qils') .and. xana%has_field('qicn') .and. &
494  xana%has_field('qlls') .and. xana%has_field('qlcn')) then
495  call xana%get_field('qils', qils)
496  call xana%get_field('qicn', qicn)
497  call xana%get_field('qlls', qlls)
498  call xana%get_field('qlcn', qlcn)
499  allocate(qi(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
500  allocate(ql(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
501  allocate(qilsf(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
502  allocate(qicnf(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
503  call q4_to_q2(geom, qils, qicn, qlls, qlcn, qi, ql, qilsf, qicnf)
504  have_cfrc = .true.
505  have_qiql = .true.
506 endif
507 
508 ! Virtual temperature
509 ! -------------------
510 have_virt = .false.
511 if (have_temp .and. xana%has_field('sphum')) then
512  allocate(tv(geom%isc:geom%iec,geom%jsc:geom%jec,1:geom%npz))
513  call xana%get_field('sphum', q)
514  call t_to_tv(geom, t, q, tv)
515  have_virt = .true.
516 endif
517 
518 ! Loop over the fields not found in the input state and work through cases
519 ! ------------------------------------------------------------------------
520 do f = 1, size(fields_to_do)
521 
522  call xctl%get_field(trim(fields_to_do(f)), field_ptr)
523 
524  select case (trim(fields_to_do(f)))
525 
526  case ("psi")
527 
528  if (.not. have_pcvd) call field_fail(fields_to_do(f))
529  field_ptr = psi
530 
531  case ("chi")
532 
533  if (.not. have_pcvd) call field_fail(fields_to_do(f))
534  field_ptr = chi
535 
536  case ("vort")
537 
538  if (.not. have_pcvd) call field_fail(fields_to_do(f))
539  field_ptr = vort
540 
541  case ("divg")
542 
543  if (.not. have_pcvd) call field_fail(fields_to_do(f))
544  field_ptr = divg
545 
546  case ("ps")
547 
548  if (.not. have_pres) call field_fail(fields_to_do(f))
549  field_ptr = ps
550 
551  case ("delp")
552 
553  if (.not. have_pres) call field_fail(fields_to_do(f))
554  field_ptr = delp
555 
556  case ("t")
557 
558  if (.not. have_temp) call field_fail(fields_to_do(f))
559  field_ptr = t
560 
561  case ("tv")
562 
563  if (.not. have_virt) call field_fail(fields_to_do(f))
564  field_ptr = tv
565 
566  case ("rh")
567 
568  if (.not. have_rhum) call field_fail(fields_to_do(f))
569  field_ptr = rh
570 
571  case ("ice_wat")
572 
573  if (.not. have_qiql) call field_fail(fields_to_do(f))
574  field_ptr = qi
575 
576  case ("liq_wat")
577 
578  if (.not. have_qiql) call field_fail(fields_to_do(f))
579  field_ptr = ql
580 
581  case ("qilsf")
582 
583  if (.not. have_cfrc) call field_fail(fields_to_do(f))
584  field_ptr = qilsf
585 
586  case ("qicnf")
587 
588  if (.not. have_cfrc) call field_fail(fields_to_do(f))
589  field_ptr = qicnf
590 
591  end select
592 
593 enddo
594 
595 ! Copy calendar infomation
596 xctl%calendar_type = xana%calendar_type
597 xctl%date_init = xana%date_init
598 
599 end subroutine changevarinverse
600 
601 ! ------------------------------------------------------------------------------
602 
603 end module fv3jedi_varcha_c2a_mod
pressure_vt_mod::ps_to_delp
subroutine, public ps_to_delp(geom, ps, delp)
Definition: pressure_variables_mod.f90:205
fv3jedi_state_mod::fv3jedi_state
Fortran derived type to hold FV3JEDI state.
Definition: fv3jedi_state_mod.F90:30
fv3jedi_field_mod
Definition: fv3jedi_field_mod.f90:6
moisture_vt_mod
Definition: moisture_variables_mod.f90:6
fv3jedi_varcha_c2a_mod::changevarinverse
subroutine, public changevarinverse(self, geom, xana, xctl)
Definition: fv3jedi_varcha_c2a_mod.f90:359
fv3jedi_fieldfail_mod::field_fail
subroutine, public field_fail(field)
Definition: fv3jedi_fieldfail_mod.f90:14
pressure_vt_mod::ps_to_pkz
subroutine, public ps_to_pkz(geom, ps, pkz)
Definition: pressure_variables_mod.f90:188
temperature_vt_mod::t_to_tv
subroutine, public t_to_tv(geom, t, q, tv)
Definition: temperature_variables_mod.f90:27
fv3jedi_varcha_c2a_mod::changevar
subroutine, public changevar(self, geom, xctl, xana)
Definition: fv3jedi_varcha_c2a_mod.f90:167
moisture_vt_mod::q2_to_q4
subroutine, public q2_to_q4(geom, qi, ql, qilsf, qicnf, qils, qicn, qlls, qlcn)
Definition: moisture_variables_mod.f90:419
fv3jedi_state_mod
Definition: fv3jedi_state_mod.F90:6
fv3jedi_field_mod::copy_subset
subroutine, public copy_subset(field_in, field_ou, not_copied)
Definition: fv3jedi_field_mod.f90:236
fv3jedi_varcha_c2a_mod::delete
subroutine, public delete(self)
Definition: fv3jedi_varcha_c2a_mod.f90:153
fv3jedi_geom_mod
Fortran module handling geometry for the FV3 model.
Definition: fv3jedi_geom_mod.f90:8
fv3jedi_varcha_c2a_mod::fv3jedi_varcha_c2a
Definition: fv3jedi_varcha_c2a_mod.f90:45
wind_vt_mod::psichi_to_vortdivg
subroutine, public psichi_to_vortdivg(geom, grid, oprs, psi, chi, lsize, lev_start, lev_final, vor, div)
Definition: wind_variables_mod.f90:522
fv3jedi_geom_mod::fv3jedi_geom
Fortran derived type to hold geometry data for the FV3JEDI model.
Definition: fv3jedi_geom_mod.f90:46
fv3jedi_fieldfail_mod
Definition: fv3jedi_fieldfail_mod.f90:1
wind_vt_mod
Definition: wind_variables_mod.f90:6
wind_vt_mod::udvd_to_psichi
subroutine, public udvd_to_psichi(geom, grid, oprs, u_in, v_in, psi, chi, lsize, lev_start, lev_final, vor_out, div_out)
Definition: wind_variables_mod.f90:321
temperature_vt_mod
Definition: temperature_variables_mod.f90:6
moisture_vt_mod::q_to_rh
subroutine, public q_to_rh(geom, qsat, q, rh)
Definition: moisture_variables_mod.f90:444
fv3jedi_varcha_c2a_mod
Definition: fv3jedi_varcha_c2a_mod.f90:6
pressure_vt_mod
Definition: pressure_variables_mod.f90:6
fv3jedi_varcha_c2a_mod::create
subroutine, public create(self, geom, conf)
Definition: fv3jedi_varcha_c2a_mod.f90:59
temperature_vt_mod::t_to_pt
subroutine, public t_to_pt(geom, pkz, t, pt)
Definition: temperature_variables_mod.f90:179
moisture_vt_mod::get_qsat
subroutine, public get_qsat(geom, delp, t, q, qsat)
Definition: moisture_variables_mod.f90:486
pressure_vt_mod::pe_to_pkz
subroutine, public pe_to_pkz(geom, pe, pkz)
Definition: pressure_variables_mod.f90:68
fv3jedi_kinds_mod::kind_real
integer, parameter, public kind_real
Definition: fv3jedi_kinds_mod.f90:14
temperature_vt_mod::pt_to_t
subroutine, public pt_to_t(geom, pkz, pt, t)
Definition: temperature_variables_mod.f90:131
moisture_vt_mod::q4_to_q2
subroutine, public q4_to_q2(geom, qils, qicn, qlls, qlcn, qi, ql, qilsf, qicnf)
Definition: moisture_variables_mod.f90:387
pressure_vt_mod::delp_to_pe_p_logp
subroutine, public delp_to_pe_p_logp(geom, delp, pe, p, logp)
Definition: pressure_variables_mod.f90:33
fv3jedi_kinds_mod
Definition: fv3jedi_kinds_mod.f90:6
wind_vt_mod::psichi_to_udvd
subroutine, public psichi_to_udvd(geom, psi, chi, u, v)
Definition: wind_variables_mod.f90:153
pressure_vt_mod::pe_to_delp
subroutine, public pe_to_delp(geom, pe, delp)
Definition: pressure_variables_mod.f90:90
fv3jedi_field_mod::field_clen
integer, parameter, public field_clen
Definition: fv3jedi_field_mod.f90:31