UFO
ufo_atmvertinterp_tlad_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2017-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  use oops_variables_mod
9  use ufo_vars_mod
10  use ufo_geovals_mod
11  use vert_interp_mod
12  use missing_values_mod
13 
14 
15 ! ------------------------------------------------------------------------------
16 
17  type, public :: ufo_atmvertinterp_tlad
18  private
19  type(oops_variables), public :: obsvars ! Variables to be simulated
20  integer, allocatable, public :: obsvarindices(:) ! Indices of obsvars in the list of all
21  ! simulated variables in the ObsSpace
22  type(oops_variables), public :: geovars
23  integer :: nval, nlocs
24  real(kind_real), allocatable :: wf(:)
25  integer, allocatable :: wi(:)
26  character(len=MAXVARLEN), public :: v_coord ! GeoVaL to use to interpolate in vertical
27  character(len=MAXVARLEN), public :: o_v_coord ! Observation vertical coordinate
28  logical, public :: use_ln ! if T, use ln(v_coord) not v_coord
29  contains
30  procedure :: setup => atmvertinterp_tlad_setup_
31  procedure :: cleanup => atmvertinterp_tlad_cleanup_
32  procedure :: settraj => atmvertinterp_tlad_settraj_
33  procedure :: simobs_tl => atmvertinterp_simobs_tl_
34  procedure :: simobs_ad => atmvertinterp_simobs_ad_
35  final :: destructor
36  end type ufo_atmvertinterp_tlad
37 
38 ! ------------------------------------------------------------------------------
39 contains
40 ! ------------------------------------------------------------------------------
41 
42 subroutine atmvertinterp_tlad_setup_(self, grid_conf)
43  use fckit_configuration_module, only: fckit_configuration
44  implicit none
45  class(ufo_atmvertinterp_tlad), intent(inout) :: self
46  type(fckit_configuration), intent(in) :: grid_conf
47 
48  character(kind=c_char,len=:), allocatable :: coord_name
49  integer :: ivar, nvars
50 
51  !> Fill in variables requested from the model
52  nvars = self%obsvars%nvars()
53  do ivar = 1, nvars
54  call self%geovars%push_back(self%obsvars%variable(ivar))
55  enddo
56  !> grab what vertical coordinate/variable to use from the config
57  self%use_ln = .false.
58 
59  if( grid_conf%has("vertical coordinate") ) then
60  call grid_conf%get_or_die("vertical coordinate",coord_name)
61  self%v_coord = coord_name
62  if( (self%v_coord .eq. var_prs) .or. (self%v_coord .eq. var_prsi) ) then
63  self%use_ln = .true.
64  endif
65  else ! default
66  self%v_coord = var_prs
67  self%use_ln = .true.
68  endif
69 
70  !> Determine observation vertical coordinate.
71  ! Use the model vertical coordinate unless the option
72  ! 'observation vertical coordinate' is specified.
73  if ( grid_conf%has("observation vertical coordinate") ) then
74  call grid_conf%get_or_die("observation vertical coordinate",coord_name)
75  self%o_v_coord = coord_name
76  else
77  self%o_v_coord = self%v_coord
78  endif
79 
80 end subroutine atmvertinterp_tlad_setup_
81 
82 ! ------------------------------------------------------------------------------
83 
84 subroutine atmvertinterp_tlad_settraj_(self, geovals, obss)
85  use obsspace_mod
86  implicit none
87  class(ufo_atmvertinterp_tlad), intent(inout) :: self
88  type(ufo_geovals), intent(in) :: geovals
89  type(c_ptr), value, intent(in) :: obss
90 
91  real(kind_real), allocatable :: obsvcoord(:)
92  type(ufo_geoval), pointer :: vcoordprofile
93  integer :: iobs
94  real(kind_real), allocatable :: tmp(:)
95  real(kind_real) :: tmp2
96 
97  ! Make sure nothing already allocated
98  call self%cleanup()
99 
100  ! Get pressure profiles from geovals
101  call ufo_geovals_get_var(geovals, self%v_coord, vcoordprofile)
102  self%nval = vcoordprofile%nval
103 
104  ! Get the observation vertical coordinates
105  self%nlocs = obsspace_get_nlocs(obss)
106  allocate(obsvcoord(self%nlocs))
107  call obsspace_get_db(obss, "MetaData", self%o_v_coord, obsvcoord)
108 
109  ! Allocate arrays for interpolation weights
110  allocate(self%wi(self%nlocs))
111  allocate(self%wf(self%nlocs))
112 
113  ! Calculate the interpolation weights
114  allocate(tmp(vcoordprofile%nval))
115  do iobs = 1, self%nlocs
116  if (self%use_ln) then
117  tmp = log(vcoordprofile%vals(:,iobs))
118  tmp2 = log(obsvcoord(iobs))
119  else
120  tmp = vcoordprofile%vals(:,iobs)
121  tmp2 = obsvcoord(iobs)
122  end if
123  call vert_interp_weights(vcoordprofile%nval, tmp2, tmp, self%wi(iobs), self%wf(iobs))
124  enddo
125 
126  ! Cleanup memory
127  deallocate(obsvcoord)
128  deallocate(tmp)
129 
130 end subroutine atmvertinterp_tlad_settraj_
131 
132 ! ------------------------------------------------------------------------------
133 
134 subroutine atmvertinterp_simobs_tl_(self, geovals, obss, nvars, nlocs, hofx)
135  implicit none
136  class(ufo_atmvertinterp_tlad), intent(in) :: self
137  type(ufo_geovals), intent(in) :: geovals
138  integer, intent(in) :: nvars, nlocs
139  real(c_double), intent(inout) :: hofx(nvars, nlocs)
140  type(c_ptr), value, intent(in) :: obss
141 
142  integer :: iobs, iobsvar, ivar
143  type(ufo_geoval), pointer :: profile
144  character(len=MAXVARLEN) :: geovar
145 
146  do iobsvar = 1, size(self%obsvarindices)
147  ! Get the index of the row of hofx to fill
148  ivar = self%obsvarindices(iobsvar)
149 
150  ! Get the name of input variable in geovals
151  geovar = self%geovars%variable(iobsvar)
152 
153  ! Get profile for this variable from geovals
154  call ufo_geovals_get_var(geovals, geovar, profile)
155 
156  ! Interpolate from geovals to observational location into hofx
157  do iobs = 1, nlocs
158  call vert_interp_apply_tl(profile%nval, profile%vals(:,iobs), &
159  & hofx(ivar,iobs), self%wi(iobs), self%wf(iobs))
160  enddo
161  enddo
162 end subroutine atmvertinterp_simobs_tl_
163 
164 ! ------------------------------------------------------------------------------
165 
166 subroutine atmvertinterp_simobs_ad_(self, geovals, obss, nvars, nlocs, hofx)
167  implicit none
168  class(ufo_atmvertinterp_tlad), intent(in) :: self
169  type(ufo_geovals), intent(inout) :: geovals
170  integer, intent(in) :: nvars, nlocs
171  real(c_double), intent(in) :: hofx(nvars, nlocs)
172  type(c_ptr), value, intent(in) :: obss
173 
174  integer :: iobs, iobsvar, ivar
175  type(ufo_geoval), pointer :: profile
176  character(len=MAXVARLEN) :: geovar
177  real(c_double) :: missing
178 
179  missing = missing_value(missing)
180 
181  do iobsvar = 1, size(self%obsvarindices)
182  ! Get the index of the row of hofx to fill
183  ivar = self%obsvarindices(iobsvar)
184 
185  ! Get the name of input variable in geovals
186  geovar = self%geovars%variable(iobsvar)
187 
188  ! Get pointer to profile for this variable in geovals
189  call ufo_geovals_get_var(geovals, geovar, profile)
190 
191  ! Adjoint of interpolate, from hofx into geovals
192  do iobs = 1, self%nlocs
193  if (hofx(ivar,iobs) /= missing) then
194  call vert_interp_apply_ad(profile%nval, profile%vals(:,iobs), &
195  & hofx(ivar,iobs), self%wi(iobs), self%wf(iobs))
196  endif
197  enddo
198  enddo
199 end subroutine atmvertinterp_simobs_ad_
200 
201 ! ------------------------------------------------------------------------------
202 
204  implicit none
205  class(ufo_atmvertinterp_tlad), intent(inout) :: self
206  self%nval = 0
207  self%nlocs = 0
208  if (allocated(self%wi)) deallocate(self%wi)
209  if (allocated(self%wf)) deallocate(self%wf)
210 end subroutine atmvertinterp_tlad_cleanup_
211 
212 ! ------------------------------------------------------------------------------
213 
214 subroutine destructor(self)
215  type(ufo_atmvertinterp_tlad), intent(inout) :: self
216 
217  call self%cleanup()
218 
219 end subroutine destructor
220 
221 ! ------------------------------------------------------------------------------
222 
subroutine atmvertinterp_simobs_tl_(self, geovals, obss, nvars, nlocs, hofx)
subroutine atmvertinterp_tlad_setup_(self, grid_conf)
subroutine atmvertinterp_simobs_ad_(self, geovals, obss, nvars, nlocs, hofx)
subroutine atmvertinterp_tlad_settraj_(self, geovals, obss)
subroutine, public ufo_geovals_get_var(self, varname, geoval)
character(len=maxvarlen), parameter, public var_prsi
character(len=maxvarlen), parameter, public var_prs
Fortran module to perform linear interpolation.
Definition: vert_interp.F90:8
subroutine vert_interp_apply_ad(nlev, fvec_ad, f_ad, wi, wf)
subroutine vert_interp_weights(nlev, obl, vec, wi, wf)
Definition: vert_interp.F90:22
subroutine vert_interp_apply_tl(nlev, fvec_tl, f_tl, wi, wf)
Definition: vert_interp.F90:92
type to hold interpolated field for one variable, one observation
type to hold interpolated fields required by the obs operators