UFO
ufo_radarradialvelocity_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 radarradialvelocity observation operator
7 
9 
10  use oops_variables_mod
11  use ufo_vars_mod
12 
13  implicit none
14  private
15 
16 !> Fortran derived type for the observation type
17  type, public :: ufo_radarradialvelocity
18  private
19  type(oops_variables), public :: obsvars
20  type(oops_variables), public :: geovars
21  character(len=MAXVARLEN), public :: v_coord ! GeoVaL to use to interpolate in vertical
22  contains
23  procedure :: setup => ufo_radarradialvelocity_setup
24  procedure :: simobs => ufo_radarradialvelocity_simobs
26 
27  character(len=maxvarlen), dimension(3), parameter :: geovars_default = (/var_u, &
28  var_v, &
29  var_w /)
30 
31 
32 contains
33 
34 ! ------------------------------------------------------------------------------
35 subroutine ufo_radarradialvelocity_setup(self, yaml_conf)
36 use fckit_configuration_module, only: fckit_configuration
37 use iso_c_binding
38 implicit none
39 class(ufo_radarradialvelocity), intent(inout) :: self
40 type(fckit_configuration), intent(in) :: yaml_conf
41 
42 character(kind=c_char,len=:), allocatable :: coord_name
43 
44  call self%geovars%push_back(geovars_default)
45 
46  if( yaml_conf%has("VertCoord") ) then
47  call yaml_conf%get_or_die("VertCoord",coord_name)
48  self%v_coord = coord_name
49  if( trim(self%v_coord) .ne. var_z ) then
50  call abor1_ftn("ufo_radarradialvelocity: incorrect vertical coordinate specified")
51  endif
52  else ! default
53  self%v_coord = var_zm
54  endif
55 
56  call self%geovars%push_back(self%v_coord)
57 
58 end subroutine ufo_radarradialvelocity_setup
59 
60 ! ------------------------------------------------------------------------------
61 ! Code in this routine is for radar radialvelocity only
62 subroutine ufo_radarradialvelocity_simobs(self, geovals, obss, nvars, nlocs, hofx)
63  use kinds
64  use vert_interp_mod
66  use obsspace_mod
67  implicit none
68  class(ufo_radarradialvelocity), intent(in) :: self
69  integer, intent(in) :: nvars, nlocs
70  type(ufo_geovals), intent(in) :: geovals
71  real(c_double), intent(inout) :: hofx(nvars, nlocs)
72  type(c_ptr), value, intent(in) :: obss
73 
74  ! Local variables
75  integer :: iobs, ivar, nvars_geovars
76  real(kind_real), dimension(:), allocatable :: obsvcoord
77  real(kind_real), dimension(:), allocatable :: cosazm_costilt, sinazm_costilt, sintilt, vterminal
78  type(ufo_geoval), pointer :: vcoordprofile, profile
79  real(kind_real), allocatable :: wf(:)
80  integer, allocatable :: wi(:)
81 
82  character(len=MAXVARLEN) :: geovar
83 
84  real(kind_real), allocatable :: tmp(:)
85  real(kind_real) :: tmp2
86 
87  real(kind_real), allocatable :: vfields(:,:) ! background fields (u,v,w) interplated vertically to the observation height
88 
89 !@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
90 
91 ! Get heigh profiles from geovals
92  call ufo_geovals_get_var(geovals, self%v_coord, vcoordprofile)
93 
94 ! Get the observation vertical coordinates
95  allocate(obsvcoord(nlocs))
96  allocate(cosazm_costilt(nlocs))
97  allocate(sinazm_costilt(nlocs))
98  allocate(sintilt(nlocs))
99  allocate(vterminal(nlocs))
100 
101  call obsspace_get_db(obss, "MetaData", "geometric_height", obsvcoord)
102  call obsspace_get_db(obss, "MetaData", "cosazm_costilt", cosazm_costilt)
103  call obsspace_get_db(obss, "MetaData", "sinazm_costilt", sinazm_costilt)
104  call obsspace_get_db(obss, "MetaData", "sintilt", sintilt)
105 
106 ! put observation operator code here
107 ! Allocate arrays for interpolation weights
108 
109  allocate(wi(nlocs))
110  allocate(wf(nlocs))
111 
112 ! Calculate the interpolation weights
113 
114  allocate(tmp(vcoordprofile%nval))
115  do iobs = 1, nlocs
116  tmp = vcoordprofile%vals(:,iobs)
117  tmp2 = obsvcoord(iobs)
118  call vert_interp_weights(vcoordprofile%nval, tmp2, tmp, wi(iobs), wf(iobs))
119  enddo
120 
121 ! Number of variables in geovars (without the vertical coordinate)
122  nvars_geovars = self%geovars%nvars() - 1
123  allocate(vfields(nvars_geovars,nlocs))
124 
125  do ivar = 1, nvars_geovars
126 ! Get the name of input variable in geovals
127  geovar = self%geovars%variable(ivar)
128 
129 ! Get profile for this variable from geovals
130  call ufo_geovals_get_var(geovals, geovar, profile)
131 
132 ! Interpolate from geovals to observational location into hofx
133  do iobs = 1, nlocs
134  call vert_interp_apply(profile%nval, profile%vals(:,iobs), &
135  & vfields(ivar,iobs), wi(iobs), wf(iobs))
136  enddo
137  enddo
138 
139  vterminal=0.0
140  do ivar = 1, nvars
141  do iobs=1,nlocs
142  hofx(ivar,iobs) = vfields(1,iobs)*cosazm_costilt(iobs) &
143  + vfields(2,iobs)*sinazm_costilt(iobs) &
144  + (vfields(3,iobs)-vterminal(iobs))*sintilt(iobs)
145  enddo
146  end do
147 
148 ! Cleanup memory
149  deallocate(obsvcoord)
150  deallocate(cosazm_costilt)
151  deallocate(sinazm_costilt)
152  deallocate(sintilt )
153  deallocate(vterminal)
154  deallocate(wi)
155  deallocate(wf)
156 
157  deallocate(vfields)
158 
159 end subroutine ufo_radarradialvelocity_simobs
160 
161 
162 ! ------------------------------------------------------------------------------
163 
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Fortran module for radarradialvelocity observation operator.
character(len=maxvarlen), dimension(3), parameter geovars_default
subroutine ufo_radarradialvelocity_setup(self, yaml_conf)
subroutine ufo_radarradialvelocity_simobs(self, geovals, obss, nvars, nlocs, hofx)
character(len=maxvarlen), parameter, public var_zm
character(len=maxvarlen), parameter, public var_v
character(len=maxvarlen), parameter, public var_z
character(len=maxvarlen), parameter, public var_u
character(len=maxvarlen), parameter, public var_w
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 the observation type.