UFO
ufo_gnssro_bndropp1d_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 for gnssro bending angle ropp1d forward operator
7 !> following the ROPP (2018 Aug) implementation
8 
10 
11 use iso_c_binding
12 use kinds
13 use ufo_vars_mod
16 use ufo_basis_mod, only: ufo_basis
19 use obsspace_mod
20 use missing_values_mod
22 use fckit_log_module, only : fckit_log
23 
24 implicit none
25 public :: ufo_gnssro_bndropp1d
26 private
27 
28  !> Fortran derived type for gnssro trajectory
30  contains
31  procedure :: simobs => ufo_gnssro_bndropp1d_simobs
32 end type ufo_gnssro_bndropp1d
33 
34 contains
35 
36 ! ------------------------------------------------------------------------------
37 ! ------------------------------------------------------------------------------
38 subroutine ufo_gnssro_bndropp1d_simobs(self, geovals, hofx, obss)
39  use ropp_fm_types, only: state1dfm
40  use ropp_fm_types, only: obs1dbangle
41  use typesizes, only: wp => eightbytereal
42  use datetimetypes, only: dp
43 
44  implicit none
45  class(ufo_gnssro_bndropp1d), intent(in) :: self
46  type(ufo_geovals), intent(in) :: geovals
47  real(kind_real), intent(inout) :: hofx(:)
48  type(c_ptr), value, intent(in) :: obss
49  real(c_double) :: missing
50 
51  type(state1dfm) :: x
52  type(obs1dbangle) :: y
53 
54  character(len=*), parameter :: myname_="ufo_gnssro_bndropp1d_simobs"
55  real(kind=dp) :: ob_time
56  integer, parameter :: max_string = 800
57 
58  character(max_string) :: err_msg
59  integer :: nlev, nobs, iobs,nvprof
60  integer, allocatable, dimension(:) :: ichk
61  type(ufo_geoval), pointer :: t, q, prs, gph, gph_sfc
62  real(kind_real), allocatable :: obsLat(:), obsLon(:), obsImpP(:), obsLocR(:), obsGeoid(:)
63  integer :: iflip
64  write(err_msg,*) "TRACE: ufo_gnssro_bndropp1d_simobs: begin"
65  call fckit_log%info(err_msg)
66 
67 ! check if nlocs is consistent in geovals & hofx
68  if (geovals%nlocs /= size(hofx)) then
69  write(err_msg,*) myname_, ' error: nlocs inconsistent!'
70  call abor1_ftn(err_msg)
71  endif
72 
73 ! get variables from geovals
74  call ufo_geovals_get_var(geovals, var_ts, t) ! temperature
75  call ufo_geovals_get_var(geovals, var_q, q) ! specific humidity
76  call ufo_geovals_get_var(geovals, var_prs, prs) ! pressure
77  call ufo_geovals_get_var(geovals, var_z, gph) ! geopotential height
78  call ufo_geovals_get_var(geovals, var_sfc_geomz, gph_sfc) ! surface geopotential height
79 
80  missing = missing_value(missing)
81 
82  nlev = t%nval ! number of model levels
83  nobs = obsspace_get_nlocs(obss)
84 
85  if (nobs > 0) then
86  iflip = 0
87  if (prs%vals(1,1) .lt. prs%vals(prs%nval,1) ) then
88  iflip = 1
89  write(err_msg,'(a)') ' ufo_gnssro_bndropp1d_simobs:'//new_line('a')// &
90  ' Model vertical height profile is in descending order,'//new_line('a')// &
91  ' but ROPP requires it to be ascending order, need flip'
92  call fckit_log%info(err_msg)
93  end if
94 
95  ! set obs space struture
96  allocate(obslon(nobs))
97  allocate(obslat(nobs))
98  allocate(obsimpp(nobs))
99  allocate(obslocr(nobs))
100  allocate(obsgeoid(nobs))
101 
102  call obsspace_get_db(obss, "MetaData", "longitude", obslon)
103  call obsspace_get_db(obss, "MetaData", "latitude", obslat)
104  call obsspace_get_db(obss, "MetaData", "impact_parameter", obsimpp)
105  call obsspace_get_db(obss, "MetaData", "earth_radius_of_curvature", obslocr)
106  call obsspace_get_db(obss, "MetaData", "geoid_height_above_reference_ellipsoid", obsgeoid)
107 
108  nvprof = 1 ! number of vertical profiles (occultation points)
109  allocate(ichk(nvprof))
110  ichk(:) = 0 ! this will hold QC values for observation from QC flags
111 
112  write(err_msg,*) "TRACE: ufo_gnssro_bndropp1d_simobs: begin observation loop, nobs = ", nobs
113  call fckit_log%info(err_msg)
114 
115  obs_loop: do iobs = 1, nobs
116 
117  ob_time = 0.0
118  call init_ropp_1d_statevec(ob_time, &
119  obslon(iobs), &
120  obslat(iobs), &
121  t%vals(:,iobs), &
122  q%vals(:,iobs), &
123  prs%vals(:,iobs), &
124  gph%vals(:,iobs), &
125  nlev, &
126  gph_sfc%vals(1,iobs), &
127  x, iflip)
128 
129  call init_ropp_1d_obvec(nvprof, &
130  obsimpp(iobs), &
131  ichk, ob_time, &
132  obslat(iobs), &
133  obslon(iobs), &
134  obslocr(iobs), &
135  obsgeoid(iobs), &
136  y)
137 
138  call ropp_fm_bangle_1d(x,y)
139 
140  ! hack -- handling ropp missing value
141  if (y%bangle(nvprof) .lt. -900.0_wp ) then
142  hofx(iobs) = missing
143  y%bangle(nvprof) = missing
144  else
145  hofx(iobs) = y%bangle(nvprof) ! nvprof is just one point
146  end if
147  ! hack -- handling ropp missing value
148  call ropp_tidy_up_1d(x,y)
149 
150  end do obs_loop
151 
152  deallocate(ichk)
153  deallocate(obslat)
154  deallocate(obslon)
155  deallocate(obsimpp)
156  deallocate(obslocr)
157  deallocate(obsgeoid)
158  end if ! nobs > 0
159 
160  write(err_msg,*) "TRACE: ufo_gnssro_bndropp1d_simobs: completed"
161  call fckit_log%info(err_msg)
162 
163 end subroutine ufo_gnssro_bndropp1d_simobs
164 ! ------------------------------------------------------------------------------
165 
166 end module ufo_gnssro_bndropp1d_mod
Fortran module to prepare for Lagrange polynomial interpolation. based on GSI: lagmod....
Definition: lag_interp.F90:4
subroutine, public lag_interp_const(q, x, n)
Definition: lag_interp.F90:24
subroutine, public lag_interp_smthweights(x, xt, aq, bq, w, dw, n)
Definition: lag_interp.F90:174
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 for gnssro bending angle ropp1d forward operator following the ROPP (2018 Aug) impleme...
subroutine ufo_gnssro_bndropp1d_simobs(self, geovals, hofx, obss)
Fortran module to handle gnssro bending angle observations following the ROPP (2018 Aug) implementati...
subroutine, public init_ropp_1d_statevec(step_time, rlon, rlat, temp, shum, pres, phi, lm, phi_sfc, x, iflip)
subroutine, public init_ropp_1d_obvec(nvprof, obs_impact, ichk, ob_time, rlat, rlon, roc, undulat, y)
character(len=maxvarlen), parameter, public var_prs
character(len=maxvarlen), parameter, public var_q
character(len=maxvarlen), parameter, public var_sfc_geomz
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
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.