UFO
ufo_rttovonedvarcheck_ob_mod.f90
Go to the documentation of this file.
1 ! (C) copyright 2020 Met Office
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 which contains the observation metadata for a single observation
7 
9 
10 use kinds
11 use missing_values_mod
12 use ufo_constants_mod, only: zero
17 
18 implicit none
19 private
20 
21 ! Ob info type definition
22 type, public :: ufo_rttovonedvarcheck_ob
23 
24  character(len=max_string) :: forward_mod_name !< forward model name (RTTOV only one at the moment)
25  integer :: nlocs !< number of locations = 1
26  integer :: surface_type !< surface type of observation
27  integer :: niter
28  integer, allocatable :: channels_used(:) !< channels used for this observation
29  integer, allocatable :: channels_all(:) !< all channels used for output
30  real(kind_real) :: latitude !< latitude of observation
31  real(kind_real) :: longitude !< longitude of observation
32  real(kind_real) :: elevation !< elevation above sea level of observation
33  real(kind_real) :: sensor_zenith_angle !< sensor zenith of observation
34  real(kind_real) :: sensor_azimuth_angle !< sensor azimuth of observation
35  real(kind_real) :: solar_zenith_angle !< solar zenith of observation
36  real(kind_real) :: solar_azimuth_angle !< solar azimuth of observation
37  real(kind_real) :: cloudtopp !< cloud top pressure (used in if cloudy retrieval used)
38  real(kind_real) :: cloudfrac !< cloud fraction (used in if cloudy retrieval used)
39  real(kind_real) :: final_cost !< final cost at solution
40  real(kind_real) :: lwp !< retrieved liquid water path. This is output for future filters
41  real(kind_real), allocatable :: yobs(:) !< satellite BTs
42  real(kind_real), allocatable :: final_bt_diff(:) !< final bt difference if converged
43  real(kind_real), allocatable :: emiss(:) !< surface emissivity
44  real(kind_real), allocatable :: background_t(:) !< background temperature used by qsplit
45  real(kind_real), allocatable :: output_profile(:) !< retrieved state at converge as profile vector
46  real(kind_real), allocatable :: output_bt(:) !< Brightness temperatures using retrieved state
47  real(kind_real), allocatable :: background_bt(:) !< Brightness temperatures from 1st itreration
48  logical :: retrievecloud !< flag to turn on retrieve cloud
49  logical :: mwscatt !< flag to use rttov-scatt model through the interface
50  logical :: mwscatt_totalice !< flag to use total ice (rather then ciw) for rttov-scatt simulations
51  logical, allocatable :: calc_emiss(:) !< flag to decide if RTTOV calculates emissivity
52  type(ufo_rttovonedvarcheck_pcemis), pointer :: pcemis !< pointer to the IR pc emissivity object
53 
54 contains
55  procedure :: setup => ufo_rttovonedvarcheck_initob
56  procedure :: delete => ufo_rttovonedvarcheck_deleteob
57  procedure :: info => ufo_rttovonedvarcheck_printob
58 
60 
61 contains
62 
63 !-------------------------------------------------------------------------------
64 !> Initialize observation object
65 !!
66 !! \author Met Office
67 !!
68 !! \date 09/06/2020: Created
69 !!
70 subroutine ufo_rttovonedvarcheck_initob(self, & ! out
71  nchans, & ! in
72  nlevels, & ! in
73  nprofelements, & ! in
74  nchans_all ) ! in
75 
76 implicit none
77 
78 ! subroutine arguments:
79 class(ufo_rttovonedvarcheck_ob), intent(out) :: self !< observation metadata type
80 integer, intent(in) :: nchans !< number of channels used for this particular observation
81 integer, intent(in) :: nlevels !< number of levels in the profile
82 integer, intent(in) :: nprofelements !< number of profile elements used
83 integer :: nchans_all !< Size of all channels in ObsSpace
84 
85 character(len=*), parameter :: routinename = "ufo_rttovonedvarcheck_InitOb"
86 real(kind_real) :: missing
87 
88 missing = missing_value(missing)
89 
90 call self % delete()
91 
92 allocate(self % yobs(nchans))
93 allocate(self % final_bt_diff(nchans))
94 allocate(self % channels_used(nchans))
95 allocate(self % channels_all(nchans_all))
96 allocate(self % emiss(nchans_all))
97 allocate(self % background_T(nlevels))
98 allocate(self % output_profile(nprofelements))
99 allocate(self % output_BT(nchans_all))
100 allocate(self % background_BT(nchans_all))
101 allocate(self % calc_emiss(nchans_all))
102 
103 self % yobs(:) = missing
104 self % final_bt_diff(:) = missing
105 self % emiss(:) = zero
106 self % background_T(:) = missing
107 self % output_profile(:) = missing
108 self % output_BT(:) = missing
109 self % background_BT(:) = missing
110 self % calc_emiss(:) = .true.
111 
112 end subroutine ufo_rttovonedvarcheck_initob
113 
114 !------------------------------------------------------------------------------
115 !> Delete the single observation object
116 !!
117 !! \author Met Office
118 !!
119 !! \date 09/06/2020: Created
120 !!
121 subroutine ufo_rttovonedvarcheck_deleteob(self) ! inout
122 
123 implicit none
124 
125 ! subroutine arguments:
126 class(ufo_rttovonedvarcheck_ob), intent(inout) :: self !< observation metadata type
127 
128 character(len=*), parameter :: routinename = "ufo_rttovonedvarcheck_DeleteOb"
129 real(kind_real) :: missing
130 
131 missing = missing_value(missing)
132 
133 self % nlocs = 1
134 self % latitude = missing
135 self % longitude = missing
136 self % elevation = missing
137 self % surface_type = 0
138 self % niter = 0
139 self % sensor_zenith_angle = missing
140 self % sensor_azimuth_angle = missing
141 self % solar_zenith_angle = missing
142 self % solar_azimuth_angle = missing
143 self % cloudtopp = 500.0_kind_real
144 self % cloudfrac = zero
145 self % final_cost = missing
146 self % LWP = missing
147 self % retrievecloud = .false.
148 self % mwscatt = .false.
149 self % mwscatt_totalice = .false.
150 
151 if (allocated(self % yobs)) deallocate(self % yobs)
152 if (allocated(self % final_bt_diff)) deallocate(self % final_bt_diff)
153 if (allocated(self % channels_used)) deallocate(self % channels_used)
154 if (allocated(self % channels_all)) deallocate(self % channels_all)
155 if (allocated(self % emiss)) deallocate(self % emiss)
156 if (allocated(self % background_T)) deallocate(self % background_T)
157 if (allocated(self % output_profile)) deallocate(self % output_profile)
158 if (allocated(self % output_BT)) deallocate(self % output_BT)
159 if (allocated(self % background_BT)) deallocate(self % background_BT)
160 if (allocated(self % calc_emiss)) deallocate(self % calc_emiss)
161 
162 self % pcemis => null()
163 
164 end subroutine ufo_rttovonedvarcheck_deleteob
165 
166 !------------------------------------------------------------------------------
167 !> Print information about the single observation object
168 !!
169 !! \author Met Office
170 !!
171 !! \date 09/06/2020: Created
172 !!
174 
175 implicit none
176 
177 ! subroutine arguments:
178 class(ufo_rttovonedvarcheck_ob), intent(inout) :: self !< observation metadata type
179 
180 character(len=20) :: surface_type
181 
182 if (self % surface_type == rtland) then
183  surface_type = "land"
184 else if (self % surface_type == rtsea) then
185  surface_type = "sea"
186 else if (self % surface_type == rtice) then
187  surface_type = "ice"
188 else
189  surface_type = "unknown"
190 end if
191 
192 write(*,"(A,2F8.2)") "Lat,Long:",self % latitude, self % longitude
193 write(*,*) "Surface type for RTTOV: ",surface_type
194 write(*,"(A,F8.2)") "Surface height:",self % elevation
195 write(*,"(A,F8.2)") "Satellite zenith angle: ",self % sensor_zenith_angle
196 write(*,"(A,F8.2)") "Solar zenith angle: ",self % solar_zenith_angle
197 write(*,"(A)") "Background T profile: "
198 write(*,"(10F8.2)") self % background_T
199 
200 end subroutine
201 
202 !-------------------------------------------------------------------------------
203 
real(kind_real), parameter, public zero
Fortran module constants used throughout the rttovonedvarcheck filter.
integer, parameter, public rtsea
integer for sea surface type
integer, parameter, public rtland
integer for land surface type
integer, parameter, public rtice
integer for seaice surface type
Fortran module which contains the observation metadata for a single observation.
subroutine ufo_rttovonedvarcheck_deleteob(self)
Delete the single observation object.
subroutine ufo_rttovonedvarcheck_initob(self, nchans, nlevels, nprofelements, nchans_all)
Initialize observation object.
subroutine ufo_rttovonedvarcheck_printob(self)
Print information about the single observation object.
Fortran module which contains the methods for Infrared Principal Component Emissivity The science can...
Fortran module containing main type, setup and utilities for the main rttovonedvarcheck object.