UFO
ufo_atmvertinterp_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 ! ------------------------------------------------------------------------------
11 
12  type, public :: ufo_atmvertinterp
13  type(oops_variables), public :: geovars
14  type(oops_variables), public :: obsvars
15  character(len=MAXVARLEN), public :: v_coord ! GeoVaL to use to interpolate in vertical
16  logical, public :: use_ln ! if T, use ln(v_coord) not v_coord
17  contains
18  procedure :: setup => atmvertinterp_setup_
19  procedure :: simobs => atmvertinterp_simobs_
20  end type ufo_atmvertinterp
21 
22 ! ------------------------------------------------------------------------------
23 contains
24 ! ------------------------------------------------------------------------------
25 
26 subroutine atmvertinterp_setup_(self, grid_conf)
27  use iso_c_binding
28  use fckit_configuration_module, only: fckit_configuration
29  implicit none
30  class(ufo_atmvertinterp), intent(inout) :: self
31  type(fckit_configuration), intent(in) :: grid_conf
32 
33  character(kind=c_char,len=:), allocatable :: coord_name
34  integer :: ivar, nvars
35 
36  !> Size of variables
37  nvars = self%obsvars%nvars()
38  !> Fill in geovars: variables we need from the model
39  ! need additional slot to hold vertical coord.
40  do ivar = 1, nvars
41  call self%geovars%push_back(self%obsvars%variable(ivar))
42  enddo
43  !> grab what vertical coordinate/variable to use from the config
44  self%use_ln = .false.
45  if( grid_conf%has("vertical coordinate") ) then
46  call grid_conf%get_or_die("vertical coordinate",coord_name)
47  self%v_coord = coord_name
48  if( trim(self%v_coord) .eq. var_prs ) self%use_ln = .true.
49  else ! default
50  self%v_coord = var_prs
51  self%use_ln = .true.
52  endif
53  call self%geovars%push_back(self%v_coord)
54 
55 end subroutine atmvertinterp_setup_
56 
57 ! ------------------------------------------------------------------------------
58 
59 subroutine atmvertinterp_simobs_(self, geovals, obss, nvars, nlocs, hofx)
60  use kinds
61  use obsspace_mod
62  use vert_interp_mod
63  use ufo_geovals_mod
64  implicit none
65  class(ufo_atmvertinterp), intent(in) :: self
66  integer, intent(in) :: nvars, nlocs
67  type(ufo_geovals), intent(in) :: geovals
68  real(c_double), intent(inout) :: hofx(nvars, nlocs)
69  type(c_ptr), value, intent(in) :: obss
70 
71  integer :: iobs, ivar
72  real(kind_real), dimension(:), allocatable :: obsvcoord
73  type(ufo_geoval), pointer :: vcoordprofile, profile
74  real(kind_real), allocatable :: wf(:)
75  integer, allocatable :: wi(:)
76  character(len=MAXVARLEN) :: geovar
77 
78  real(kind_real), allocatable :: tmp(:)
79  real(kind_real) :: tmp2
80 
81  ! Get pressure profiles from geovals
82  call ufo_geovals_get_var(geovals, self%v_coord, vcoordprofile)
83 
84  ! Get the observation vertical coordinates
85  allocate(obsvcoord(nlocs))
86  call obsspace_get_db(obss, "MetaData", self%v_coord, obsvcoord)
87 
88  ! Allocate arrays for interpolation weights
89  allocate(wi(nlocs))
90  allocate(wf(nlocs))
91 
92  ! Calculate the interpolation weights
93  allocate(tmp(vcoordprofile%nval))
94  do iobs = 1, nlocs
95  if (self%use_ln) then
96  tmp = log(vcoordprofile%vals(:,iobs))
97  tmp2 = log(obsvcoord(iobs))
98  else
99  tmp = vcoordprofile%vals(:,iobs)
100  tmp2 = obsvcoord(iobs)
101  end if
102  call vert_interp_weights(vcoordprofile%nval, tmp2, tmp, wi(iobs), wf(iobs))
103  enddo
104 
105  do ivar = 1, nvars
106  ! Get the name of input variable in geovals
107  geovar = self%geovars%variable(ivar)
108 
109  ! Get profile for this variable from geovals
110  call ufo_geovals_get_var(geovals, geovar, profile)
111 
112  ! Interpolate from geovals to observational location into hofx
113  do iobs = 1, nlocs
114  call vert_interp_apply(profile%nval, profile%vals(:,iobs), &
115  & hofx(ivar,iobs), wi(iobs), wf(iobs))
116  enddo
117  enddo
118  ! Cleanup memory
119  deallocate(obsvcoord)
120  deallocate(wi)
121  deallocate(wf)
122 
123  deallocate(tmp)
124 
125 end subroutine atmvertinterp_simobs_
126 
127 ! ------------------------------------------------------------------------------
128 
129 end module ufo_atmvertinterp_mod
ufo_atmvertinterp_mod::atmvertinterp_setup_
subroutine atmvertinterp_setup_(self, grid_conf)
Definition: ufo_atmvertinterp_mod.F90:27
ufo_atmvertinterp_mod
Definition: ufo_atmvertinterp_mod.F90:6
ufo_atmvertinterp_mod::ufo_atmvertinterp
Definition: ufo_atmvertinterp_mod.F90:12
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_vars_mod
Definition: ufo_variables_mod.F90:8
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_atmvertinterp_mod::atmvertinterp_simobs_
subroutine atmvertinterp_simobs_(self, geovals, obss, nvars, nlocs, hofx)
Definition: ufo_atmvertinterp_mod.F90:60
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
ufo_vars_mod::var_prs
character(len=maxvarlen), parameter, public var_prs
Definition: ufo_variables_mod.F90:25