FV3 Bundle
ufo_conventional_profile_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 
7 
8  use iso_c_binding
9  use kinds
10  use ufo_vars_mod
11  use ufo_geovals_mod
13  use vert_interp_mod
14  use ufo_basis_mod, only: ufo_basis
15  use obsspace_mod
16 
17  integer, parameter :: max_string=800
18 
20  private
21  integer :: nval, nobs
22  real(kind_real), allocatable :: wf(:)
23  integer, allocatable :: wi(:)
24  contains
25  procedure :: simobs => conventional_profile_simobs_
27 contains
28 
29 ! ------------------------------------------------------------------------------
30 
31  subroutine conventional_profile_simobs_(self, geovals, hofx, obss)
32 
33  implicit none
34  class(ufo_conventional_profile), intent(in) :: self
35  type(ufo_geovals), intent(in) :: geovals
36  real(c_double), intent(inout) :: hofx(:)
37  type(c_ptr), value, intent(in) :: obss
38 
39  character(len=*), parameter :: myname_="ufo_conventional_profile_simobs"
40  character(max_string) :: err_msg
41 
42  integer :: iobs
43  real(kind_real) :: wf
44  integer :: wi, ierr, nobs
45  real(kind_real), allocatable :: pressure(:)
46  type(ufo_geoval), pointer :: prsl, tv
47 
48  ! check if nobs is consistent in geovals & hofx
49  if (geovals%nobs /= size(hofx)) then
50  write(err_msg,*) myname_, ' error: nobs inconsistent!'
51  call abor1_ftn(err_msg)
52  endif
53 
54  ! check if prsl variable is in geovals and get it
55  call ufo_geovals_get_var(geovals, var_prsl, prsl,status=ierr)
56  if (ierr/=0) then
57  write(err_msg,*) myname_, trim(var_prsl), ' doesnt exist'
58  call abor1_ftn(err_msg)
59  endif
60 
61  ! check if tv variable is in geovals and get it
62  call ufo_geovals_get_var(geovals, var_tv, tv,status=ierr)
63  if (ierr/=0) then
64  write(err_msg,*) myname_, trim(var_tv), ' doesnt exist'
65  call abor1_ftn(err_msg)
66  endif
67 
68  ! observation of pressure (for vertical interpolation)
69  nobs = obsspace_get_nobs(obss)
70  allocate(pressure(nobs))
71  call obsspace_get_db(obss, "ObsValue", "air_pressure", pressure)
72 
73  ! obs operator
74  do iobs = 1, geovals%nobs
75  call vert_interp_weights(prsl%nval,log(pressure(iobs)/10.),prsl%vals(:,iobs),wi,wf)
76  call vert_interp_apply(tv%nval, tv%vals(:,iobs), hofx(iobs), wi, wf)
77  enddo
78 
79  ! cleanup
80  deallocate(pressure)
81  end subroutine conventional_profile_simobs_
82 
83 ! ------------------------------------------------------------------------------
84 
subroutine, public ufo_geovals_get_var(self, varname, geoval, status)
subroutine vert_interp_apply(nlev, fvec, f, wi, wf)
Definition: vert_interp.F90:67
integer, parameter max_string
integer function, public obsspace_get_nobs(c_dom)
Return the number of observations.
Fortran module to perform linear interpolation.
Definition: vert_interp.F90:8
character(len=maxvarlen), public var_prsl
subroutine conventional_profile_simobs_(self, geovals, hofx, obss)
type(registry_t), public ufo_geovals_registry
Linked list interface - defines registry_t type.
Fortran interface to ObsSpace.
Definition: obsspace_mod.F90:9
character(len=maxvarlen), public var_tv
subroutine vert_interp_weights(nlev, obl, vec, wi, wf)
Definition: vert_interp.F90:20