UFO
ufo_radarreflectivity_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 radarreflectivity observation operator
7 
9 
10  use ufo_vars_mod
11  use oops_variables_mod
12 
13  implicit none
14  private
15 
16 !> Fortran derived type for the observation type
17 
18  type, public :: ufo_radarreflectivity
19  private
20  type(oops_variables), public :: obsvars
21  type(oops_variables), public :: geovars
22  character(len=MAXVARLEN), public :: v_coord ! GeoVaL to use to interpolate in vertical
23  contains
24  procedure :: setup => ufo_radarreflectivity_setup
25  procedure :: simobs => ufo_radarreflectivity_simobs
26  end type ufo_radarreflectivity
27 
28  character(len=maxvarlen), dimension(1), parameter :: geovars_default = (/var_refl/)
29 
30 contains
31 
32 ! ------------------------------------------------------------------------------
33 subroutine ufo_radarreflectivity_setup(self, yaml_conf)
34  use fckit_configuration_module, only: fckit_configuration
35  use iso_c_binding
36  implicit none
37  class(ufo_radarreflectivity), intent(inout) :: self
38  type(fckit_configuration), intent(in) :: yaml_conf
39 
40  character(kind=c_char,len=:), allocatable :: coord_name
41 
42  call self%geovars%push_back(geovars_default)
43 
44  if( yaml_conf%has("VertCoord") ) then
45  call yaml_conf%get_or_die("VertCoord",coord_name)
46  self%v_coord = coord_name
47  if( trim(self%v_coord) .ne. var_z ) then
48  call abor1_ftn("ufo_radarreflectivity: incorrect vertical coordinate specified")
49  endif
50  else ! default
51  self%v_coord = var_z
52  endif
53 
54  call self%geovars%push_back(self%v_coord)
55 
56 end subroutine ufo_radarreflectivity_setup
57 
58 ! ------------------------------------------------------------------------------
59 ! Code in this routine is for radarreflectivity only
60 subroutine ufo_radarreflectivity_simobs(self, geovals, obss, nvars, nlocs, hofx)
61  use kinds
62  use vert_interp_mod
64  use obsspace_mod
65  implicit none
66  class(ufo_radarreflectivity), intent(in) :: self
67  integer, intent(in) :: nvars, nlocs
68  type(ufo_geovals), intent(in) :: geovals
69  real(c_double), intent(inout) :: hofx(nvars, nlocs)
70  type(c_ptr), value, intent(in) :: obss
71 
72  ! Local variables
73  integer :: iobs, ivar
74  real(kind_real), dimension(:), allocatable :: obsvcoord
75  type(ufo_geoval), pointer :: vcoordprofile, profile
76  real(kind_real), allocatable :: wf(:)
77  integer, allocatable :: wi(:)
78 
79  character(len=MAXVARLEN) :: geovar
80 
81  real(kind_real), allocatable :: tmp(:)
82  real(kind_real) :: tmp2
83 
84  !@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@
85 
86  ! Get height profiles from geovals
87  call ufo_geovals_get_var(geovals, self%v_coord, vcoordprofile)
88 
89 
90  ! Get the observation vertical coordinates
91  allocate(obsvcoord(nlocs))
92  call obsspace_get_db(obss, "MetaData", "height", obsvcoord)
93 
94  ! Allocate arrays for interpolation weights
95 
96  allocate(wi(nlocs))
97  allocate(wf(nlocs))
98 
99  ! Calculate the interpolation weights
100  allocate(tmp(vcoordprofile%nval))
101  do iobs = 1, nlocs
102  tmp = vcoordprofile%vals(:,iobs)
103  tmp2 = obsvcoord(iobs)
104  call vert_interp_weights(vcoordprofile%nval, tmp2, tmp, wi(iobs), wf(iobs))
105  enddo
106 
107  do ivar = 1, nvars
108  ! Get the name of input variable in geovals
109  geovar = self%geovars%variable(ivar)
110 
111  ! Get profile for this variable from geovals
112  call ufo_geovals_get_var(geovals, geovar, profile)
113 
114  ! Interpolate from geovals to observational location into hofx
115  do iobs = 1, nlocs
116  call vert_interp_apply(profile%nval, profile%vals(:,iobs), &
117  & hofx(ivar,iobs), wi(iobs), wf(iobs))
118  enddo
119  enddo
120 
121  ! Cleanup memory
122  deallocate(obsvcoord)
123  deallocate(wi)
124  deallocate(wf)
125 
126 end subroutine ufo_radarreflectivity_simobs
127 
128 
129 ! ------------------------------------------------------------------------------
130 
131 end module ufo_radarreflectivity_mod
ufo_radarreflectivity_mod::ufo_radarreflectivity
Fortran derived type for the observation type.
Definition: ufo_radarreflectivity_mod.F90:18
ufo_radarradialvelocity_mod::geovars_default
character(len=maxvarlen), dimension(3), parameter geovars_default
Definition: ufo_radarradialvelocity_mod.F90:27
ufo_vars_mod::var_refl
character(len=maxvarlen), parameter, public var_refl
Definition: ufo_variables_mod.F90:86
vert_interp_mod
Fortran module to perform linear interpolation.
Definition: vert_interp.F90:8
vert_interp_mod::vert_interp_weights
subroutine vert_interp_weights(nlev, obl, vec, wi, wf)
Definition: vert_interp.F90:22
ufo_geovals_mod
Definition: ufo_geovals_mod.F90:7
ufo_radarreflectivity_mod::ufo_radarreflectivity_setup
subroutine ufo_radarreflectivity_setup(self, yaml_conf)
Definition: ufo_radarreflectivity_mod.F90:34
ufo_vars_mod
Definition: ufo_variables_mod.F90:8
ufo_radarreflectivity_mod::ufo_radarreflectivity_simobs
subroutine ufo_radarreflectivity_simobs(self, geovals, obss, nvars, nlocs, hofx)
Definition: ufo_radarreflectivity_mod.F90:61
ufo_vars_mod::var_z
character(len=maxvarlen), parameter, public var_z
Definition: ufo_variables_mod.F90:29
vert_interp_mod::vert_interp_apply
subroutine vert_interp_apply(nlev, fvec, f, wi, wf)
Definition: vert_interp.F90:73
ufo_geovals_mod::ufo_geovals_get_var
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Definition: ufo_geovals_mod.F90:128
ufo_radarreflectivity_mod
Fortran module for radarreflectivity observation operator.
Definition: ufo_radarreflectivity_mod.F90:8
ufo_geovals_mod::ufo_geovals
type to hold interpolated fields required by the obs operators
Definition: ufo_geovals_mod.F90:47
ufo_geovals_mod::ufo_geoval
type to hold interpolated field for one variable, one observation
Definition: ufo_geovals_mod.F90:40