UFO
ufo_gnssro_ref_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 !> Fortran module to handle gnssro refractivity observations
7 
9  use fckit_configuration_module, only: fckit_configuration
10  use kinds
11  use ufo_vars_mod
12  use ufo_geovals_mod
14  use vert_interp_mod
15  use ufo_basis_mod, only: ufo_basis
16  use obsspace_mod
17  use gnssro_mod_conf
19 
20  implicit none
21  public :: ufo_gnssro_ref
22  private
23 
24  !> Fortran derived type for gnssro trajectory
25  type, extends(ufo_basis) :: ufo_gnssro_ref
26  private
27  type(gnssro_conf) :: roconf
28  contains
29  procedure :: setup => ufo_gnssro_ref_setup
30  procedure :: simobs => ufo_gnssro_ref_simobs
31  end type ufo_gnssro_ref
32 
33 contains
34 ! ------------------------------------------------------------------------------
35  subroutine ufo_gnssro_ref_setup(self, f_conf)
36  implicit none
37  class(ufo_gnssro_ref), intent(inout) :: self
38  type(fckit_configuration), intent(in) :: f_conf
39 
40  call gnssro_conf_setup(self%roconf,f_conf)
41 
42  end subroutine ufo_gnssro_ref_setup
43 
44  subroutine ufo_gnssro_ref_simobs(self, geovals, hofx, obss)
46  implicit none
47  class(ufo_gnssro_ref), intent(in) :: self
48  type(ufo_geovals), intent(in) :: geovals
49  real(kind_real), intent(inout) :: hofx(:)
50  type(c_ptr), value, intent(in) :: obss
51 
52  character(len=*), parameter :: myname_="ufo_gnssro_ref_simobs"
53  character(max_string) :: err_msg
54  integer :: iobs,nlocs
55  real(kind_real) :: wf
56  integer :: wi
57  type(ufo_geoval), pointer :: t,q,prs,gph
58  real(kind_real) :: refr1, refr2,refr3
59  real(kind_real), allocatable :: obsZ(:), obsLat(:)
60  real(kind_real) :: obsH, gesT,gesQ, gesTv, gesTv0,gesP
61 
62  ! check if nlocs is consistent in geovals & hofx
63  if (geovals%nlocs /= size(hofx)) then
64  write(err_msg,*) myname_, ' error: nlocs inconsistent!'
65  call abor1_ftn(err_msg)
66  endif
67  ! get variables from geovals
68  call ufo_geovals_get_var(geovals, var_prs, prs)
69  call ufo_geovals_get_var(geovals, var_ts,t)
70  call ufo_geovals_get_var(geovals, var_q, q)
71  call ufo_geovals_get_var(geovals, var_z, gph)
72 
73  nlocs = obsspace_get_nlocs(obss)
74 
75  allocate(obsz(nlocs))
76  allocate(obslat(nlocs))
77 
78  call obsspace_get_db(obss, "MetaData", "altitude", obsz)
79  call obsspace_get_db(obss, "MetaData", "latitude", obslat)
80 
81  call gnssro_ref_constants(self%roconf%use_compress)
82 
83 
84  ! obs operator
85  do iobs = 1, geovals%nlocs
86 
87  ! Convert geometric height at observation to geopotential height
88  call geometric2geop(obslat(iobs), obsz(iobs), obsh)
89  call vert_interp_weights(gph%nval,obsh, gph%vals(:,iobs),wi,wf) ! calculate weights
90  call vert_interp_apply(t%nval, t%vals(:,iobs), gest, wi, wf)
91  call vert_interp_apply(q%nval, q%vals(:,iobs), gesq, wi, wf)
92 
93  ! use hypsometric equation to calculate pressure
94  gestv = 0.0
95  gestv0 = 0.0
96  gestv = gest*(one + (rv_over_rd-one)* (gesq/(1-gesq) ) )
97  gestv0 = t%vals(wi,iobs)*(one + (rv_over_rd-one) * (q%vals(wi,iobs)/(1-q%vals(wi,iobs)) ))
98  gesp = prs%vals(wi,iobs)/exp(two*grav*(obsh-gph%vals(wi,iobs))/(rd*(gestv+gestv0)))
99  refr1 = n_a*gesp/gest
100  refr2 = n_b*gesp*gesq/ ( gest**2 * (rd_over_rv+(1-rd_over_rv)*gesq) )
101  refr3 = n_c*gesp*gesq/ ( gest * (rd_over_rv+(1-rd_over_rv)*gesq) )
102  hofx(iobs) = refr1 + refr2 + refr3
103  enddo
104 
105 
106  ! cleanup
107  deallocate(obsz)
108  deallocate(obslat)
109  end subroutine ufo_gnssro_ref_simobs
110 
111 ! ------------------------------------------------------------------------------
112 end module ufo_gnssro_ref_mod
subroutine, public gnssro_conf_setup(roconf, f_conf)
subroutine, public gnssro_ref_constants(use_compress)
real(kind_real), public n_a
real(kind_real), public n_b
real(kind_real), public n_c
subroutine geometric2geop(Latitude, geometricZ, geopotentialH)
type(registry_t), public ufo_geovals_registry
Linked list interface - defines registry_t type.
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Fortran module to handle gnssro refractivity observations.
subroutine ufo_gnssro_ref_setup(self, f_conf)
subroutine ufo_gnssro_ref_simobs(self, geovals, hofx, obss)
character(len=maxvarlen), parameter, public var_prs
character(len=maxvarlen), parameter, public var_q
character(len=maxvarlen), parameter, public var_z
character(len=maxvarlen), parameter, public var_ts
Fortran module to perform linear interpolation.
Definition: vert_interp.F90:8
subroutine vert_interp_weights(nlev, obl, vec, wi, wf)
Definition: vert_interp.F90:22
subroutine vert_interp_apply(nlev, fvec, f, wi, wf)
Definition: vert_interp.F90:73
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators
Fortran derived type for gnssro trajectory.