UFO
ufo_marinevertinterp_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 
6 !> Fortran marinevertinterp module for observation operator
7 
9 
10  use ufo_vars_mod
11  use oops_variables_mod
12 
13  implicit none
14  private
15 
16  integer, parameter :: max_string=800
17 
18 !> Fortran derived type for the observation type
19  type, public :: ufo_marinevertinterp
20  type(oops_variables), public :: geovars
21  type(oops_variables), public :: obsvars
22  contains
23  procedure :: setup => ufo_marinevertinterp_setup
24  procedure :: simobs => ufo_marinevertinterp_simobs
25  end type ufo_marinevertinterp
26 
27 contains
28 
29 ! ------------------------------------------------------------------------------
31 implicit none
32 class(ufo_marinevertinterp), intent(inout) :: self
33 character(max_string) :: err_msg
34 
35 integer :: ivar, nvars
36 
37 nvars = self%obsvars%nvars()
38 if (nvars /= 1) then
39  write(err_msg,*) 'ufo_marinevertinterp_setup error: only variables size 1 supported!'
40  call abor1_ftn(err_msg)
41 endif
42 
43 ! Set variables requested from the model
44 do ivar = 1, nvars
45  call self%geovars%push_back(self%obsvars%variable(ivar))
46 enddo
47 call self%geovars%push_back("sea_water_cell_thickness")
48 
49 end subroutine ufo_marinevertinterp_setup
50 
51 ! ------------------------------------------------------------------------------
52 subroutine ufo_marinevertinterp_simobs(self, geovals, hofx, obss)
56 use iso_c_binding
57 use kinds
59 use obsspace_mod
60 implicit none
61 class(ufo_marinevertinterp), intent(in) :: self
62 type(ufo_geovals), intent(in) :: geovals
63 real(c_double), intent(inout) :: hofx(:)
64 type(c_ptr), value, intent(in) :: obss
65 
66  character(len=*), parameter :: myname_="ufo_marinevertinterp_simobs"
67  character(max_string) :: err_msg
68 
69  integer :: iobs, ilev, nlev, nlocs
70  type(ufo_geoval), pointer :: var, h
71  real (kind_real), allocatable :: depth(:,:)
72  real(kind_real) :: deptho
73  real(kind_real), allocatable :: obs_depth(:)
74  integer :: obss_nlocs
75  real(kind_real) :: wf, sp, prs
76  integer :: wi
77 
78  ! check if nlocs is consistent in geovals & hofx
79  if (geovals%nlocs /= size(hofx,1)) then
80  write(err_msg,*) myname_, ' error: nlocs inconsistent!'
81  call abor1_ftn(err_msg)
82  endif
83 
84  ! Associate geoval pointers
85  call ufo_geovals_get_var(geovals, self%geovars%variable(1), var)
86  call ufo_geovals_get_var(geovals, var_ocn_lay_thick, h)
87 
88  ! Read in obs data
89  obss_nlocs = obsspace_get_nlocs(obss)
90  allocate(obs_depth(obss_nlocs))
91  call obsspace_get_db(obss, "MetaData", "depth", obs_depth)
92 
93  nlev = var%nval
94  nlocs = var%nlocs
95  allocate(depth(nlev,nlocs))
96  do iobs = 1,size(hofx,1)
97  !< Depth from layer thickness
98  depth(1,iobs)=0.5*h%vals(1,iobs)
99  do ilev = 2, nlev
100  depth(ilev,iobs)=sum(h%vals(1:ilev-1,iobs))+0.5*h%vals(ilev,iobs)
101  end do
102  end do
103 
104  hofx = 0.0
105  ! Vertical interpolation
106  do iobs = 1,size(hofx,1)
107 
108  deptho = obs_depth(iobs)
109 
110  !< Interpolation weight
111  call vert_interp_weights(nlev, deptho, depth(:,iobs), wi, wf)
112 
113  !Apply vertical interpolation
114  call vert_interp_apply(nlev, var%vals(:,iobs), hofx(iobs), wi, wf)
115 
116  enddo
117 
118  deallocate(depth)
119  deallocate(obs_depth)
120 
121  end subroutine ufo_marinevertinterp_simobs
122 
123 end module ufo_marinevertinterp_mod
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Fortran marinevertinterp module for observation operator.
subroutine ufo_marinevertinterp_simobs(self, geovals, hofx, obss)
character(len=maxvarlen), public var_ocn_lay_thick
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.