UFO
ufo_identity_tlad_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 identity tl/ad observation operator
7 
9 
10  use oops_variables_mod
11  use ufo_vars_mod
12  implicit none
13 
14  ! ------------------------------------------------------------------------------
15 
16  !> Fortran derived type for the tl/ad observation operator
17  type, public :: ufo_identity_tlad
18  private
19  type(oops_variables), public :: obsvars
20  type(oops_variables), public :: geovars
21  contains
22  procedure :: setup => identity_tlad_setup_
23  procedure :: settraj => identity_tlad_settraj_
24  procedure :: simobs_tl => identity_simobs_tl_
25  procedure :: simobs_ad => identity_simobs_ad_
26  end type ufo_identity_tlad
27 
28 contains
29 
30 ! ------------------------------------------------------------------------------
31 subroutine identity_tlad_setup_(self)
32  implicit none
33  class(ufo_identity_tlad), intent(inout) :: self
34 
35  integer :: ivar, nvars
36 
37  !> copy observed variables to variables requested from the model
38  nvars = self%obsvars%nvars()
39  do ivar = 1, nvars
40  call self%geovars%push_back(self%obsvars%variable(ivar))
41  enddo
42 
43 end subroutine identity_tlad_setup_
44 
45 !------------------------------------------------------------------------------
46 subroutine identity_tlad_settraj_(self, geovals, obss)
47  use iso_c_binding
48  use ufo_geovals_mod, only: ufo_geovals
49  implicit none
50  class(ufo_identity_tlad), intent(inout) :: self
51  type(ufo_geovals), intent(in) :: geovals
52  type(c_ptr), value, intent(in) :: obss
53 
54 ! since observation operator is linear, don't care about trajectory itself
55 
56 end subroutine identity_tlad_settraj_
57 
58 ! ------------------------------------------------------------------------------
59 subroutine identity_simobs_tl_(self, geovals, hofx, obss)
60  use iso_c_binding
61  use kinds
62  use ufo_geovals_mod, only: &
63  ufo_geovals, &
64  ufo_geoval, &
66  use obsspace_mod
67  implicit none
68  class(ufo_identity_tlad), intent(in) :: self
69  type(ufo_geovals), intent(in) :: geovals
70  real(c_double), intent(inout) :: hofx(:)
71  type(c_ptr), value, intent(in) :: obss
72 
73  integer :: iobs, ivar, nlocs, nvars
74  type(ufo_geoval), pointer :: point
75  character(len=MAXVARLEN) :: geovar
76 
77  !> Get the observation vertical coordinates
78  nlocs = obsspace_get_nlocs(obss)
79  nvars = self%obsvars%nvars()
80  do ivar = 1, nvars
81  !> Get the name of input variable in geovals
82  geovar = self%geovars%variable(ivar)
83 
84  !> Get profile for this variable from geovals
85  call ufo_geovals_get_var(geovals, geovar, point)
86 
87  !> Here we just apply a identity hofx
88  do iobs = 1, nlocs
89  hofx(ivar + (iobs-1)*nvars) = point%vals(1,iobs)
90  enddo
91  enddo
92 
93 end subroutine identity_simobs_tl_
94 
95 ! ------------------------------------------------------------------------------
96 subroutine identity_simobs_ad_(self, geovals, hofx, obss)
97  use iso_c_binding
98  use kinds
99  use ufo_geovals_mod, only: &
100  ufo_geovals, &
101  ufo_geoval, &
103  use obsspace_mod
104  use missing_values_mod
105  implicit none
106  class(ufo_identity_tlad), intent(in) :: self
107  type(ufo_geovals), intent(inout) :: geovals
108  real(c_double), intent(in) :: hofx(:)
109  type(c_ptr), value, intent(in) :: obss
110 
111  integer :: iobs, ivar, nlocs, nvars
112  type(ufo_geoval), pointer :: point
113  character(len=MAXVARLEN) :: geovar
114  real(c_double) :: missing
115 
116  !> Set missing value
117  missing = missing_value(missing)
118 
119  if (.not. geovals%linit ) geovals%linit=.true.
120 
121  nlocs = obsspace_get_nlocs(obss)
122  nvars = self%obsvars%nvars()
123  do ivar = 1, nvars
124  !> Get the name of input variable in geovals
125  geovar = self%geovars%variable(ivar)
126 
127  !> Get profile for this variable from geovals
128  call ufo_geovals_get_var(geovals, geovar, point)
129 
130  if (.not.(allocated(point%vals))) then
131  point%nval=1
132  allocate(point%vals(1,size(hofx,1)))
133  point%vals = 0.0
134  end if
135 
136  ! backward obs operator
137  do iobs = 1, nlocs
138  if (hofx(ivar + (iobs-1)*nvars) /= missing) then
139  point%vals(1,iobs) = hofx(ivar + (iobs-1)*nvars)
140  endif
141  enddo
142  enddo
143 
144 end subroutine identity_simobs_ad_
145 
146 end module ufo_identity_tlad_mod
ufo_identity_tlad_mod::identity_simobs_ad_
subroutine identity_simobs_ad_(self, geovals, hofx, obss)
Definition: ufo_identity_tlad_mod.F90:97
ufo_geovals_mod
Definition: ufo_geovals_mod.F90:7
ufo_identity_tlad_mod::identity_tlad_settraj_
subroutine identity_tlad_settraj_(self, geovals, obss)
Definition: ufo_identity_tlad_mod.F90:47
ufo_identity_tlad_mod::ufo_identity_tlad
Fortran derived type for the tl/ad observation operator.
Definition: ufo_identity_tlad_mod.F90:17
ufo_vars_mod
Definition: ufo_variables_mod.F90:8
ufo_identity_tlad_mod::identity_simobs_tl_
subroutine identity_simobs_tl_(self, geovals, hofx, obss)
Definition: ufo_identity_tlad_mod.F90:60
ufo_geovals_mod::ufo_geovals_get_var
subroutine, public ufo_geovals_get_var(self, varname, geoval)
Definition: ufo_geovals_mod.F90:128
ufo_identity_tlad_mod::identity_tlad_setup_
subroutine identity_tlad_setup_(self)
Definition: ufo_identity_tlad_mod.F90:32
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_identity_tlad_mod
Fortran module for identity tl/ad observation operator.
Definition: ufo_identity_tlad_mod.F90:8