UFO
ufo_rttovonedvarcheck_profindex_mod.f90
Go to the documentation of this file.
1 ! (C) British Crown Copyright 2017-2018 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 containing profile index
7 
9 
10 use kinds
11 use fckit_log_module, only : fckit_log
14 
15 implicit none
16 private
17 
19 
20  ! general
21  integer :: nprofelements !< number of profile elements being used
22  integer :: nlevels !< number of model levels
23 
24  ! atmosphere (locate start and end points of relevant fields)
25  integer :: t(2) !< temperature profile
26  integer :: q(2) !< water vapour profile (specific humidity)
27  integer :: ql(2) !< liquid water profile
28  integer :: qt(2) !< total water profile
29  integer :: qi(2) !< frozen ice profile
30  integer :: cf(2) !< cloud fraction profile
31  integer :: o3total !< total column ozone
32  integer :: o3profile(2) !< ozone profile
33  integer :: lwp !< liquid water path
34 
35  ! surface
36  integer :: t2 !< screen temperature
37  integer :: q2 !< screen specific humidity
38  integer :: rh2 !< screen relative humidity
39  integer :: tstar !< skin temperature
40  integer :: pstar !< surface pressure
41  integer :: mwemiss(2) !< microwave emissivity
42  integer :: emisspc(2) !< emissivity principal components
43  integer :: windspeed !< surface windspeed
44 
45  ! cloud (single-level grey cloud only)
46  integer :: cloudtopp !< single-level cloud top pressure
47  integer :: cloudfrac !< effective cloud fraction
48 
49  ! other
50  integer :: t70hpa !< temperature at 70hpa - used for ozone profile, not currently implemented
51  integer :: t700hpa !< temperature at 700hpa - used for ozone profile, not currently implemented
52  integer :: t950hpa !< temperature at 950hpa - used for ozone profile, not currently implemented
53  integer :: t1000hpa !< temperature at 1000hpa - used for ozone profile, not currently implemented
54  integer :: qsurf !< surface specific humidity
55 
56 contains
59 
60 end type
61 
62 contains
63 
64 !-------------------------------------------------------------------------------
65 !> Profile index setup
66 !!
67 !! \details Heritage: Ops_SatRad_MapProfileToB
68 !!
69 !! Setup the profile index which requires the bmatrix object.
70 !!
71 !! \author Met Office
72 !!
73 !! \date 09/06/2020: Created
74 !!
75 subroutine ufo_rttovonedvarcheck_profindex_setup(self, bmatrix, nlevels)
76 
77 implicit none
78 
79 ! subroutine arguments:
80 class(ufo_rttovonedvarcheck_profindex), intent(inout) :: self !< profindex structure
81 type(ufo_metoffice_bmatrixstatic), intent(in) :: bmatrix !< state error covariances
82 integer, intent(in) :: nlevels !< number of model levels
83 
84 ! local constants:
85 character(len=*), parameter :: routinename = "ufo_rttovonedvarcheck_profindex_setup"
86 
87 ! local variables:
88 integer :: i,j
89 integer :: firstelement
90 integer :: lastelement
91 integer :: nelements
92 
93 ! Initialise to zeros
94 call self % delete()
95 nelements = 0
96 self % nlevels = nlevels
97 
98 !loop through fields in bmatrix. the number of elements that the field is composed of
99 !is held in bmatrix % fields(i,2).
100 do j = 1, bmatrix % nfields
101  firstelement = nelements + 1
102  lastelement = nelements + bmatrix % fields(j,2)
103 
104  !assign start and end points. if the field wasn't found then assign a value of
105  !zero, which indicates absence.
106 
107  select case( bmatrix % fields(j,1) )
108 
109  !----------
110  !atmosphere (set start and end points for multi-level fields)
111  !----------
112 
114  self % t(1) = firstelement
115  self % t(2) = lastelement
116 
118  self % q(1) = firstelement
119  self % q(2) = lastelement
120 
122  call abor1_ftn("ufo_metoffice_fieldtype_ql: Not currently implemented aborting")
123 
125  call abor1_ftn("ufo_metoffice_fieldtype_qi: Not currently implemented aborting")
126 
128  call abor1_ftn("ufo_metoffice_fieldtype_cf: Not currently implemented aborting")
129 
131  self % qt(1) = firstelement
132  self % qt(2) = lastelement
133 
135  call abor1_ftn("ufo_metoffice_fieldtype_o3profile: Not currently implemented aborting")
136 
138  call abor1_ftn("ufo_metoffice_fieldtype_o3total: Not currently implemented aborting")
139 
140  !-------
141  !surface
142  !-------
143 
145  self % t2 = firstelement
146 
148  self % q2 = firstelement
149 
151  self % tstar = firstelement
152 
154  self % pstar = firstelement
155 
157  self % windspeed = firstelement
158 
160  call abor1_ftn("ufo_metoffice_fieldtype_mwemiss: Not currently implemented aborting")
161 
163  call abor1_ftn("ufo_metoffice_fieldtype_emisspc: Not currently implemented aborting")
164 
165  !------------------------------------
166  !cloud (single-level grey cloud only)
167  !------------------------------------
168 
170  call abor1_ftn("ufo_metoffice_fieldtype_cloudtopp: Not currently implemented aborting")
171 
173  call abor1_ftn("ufo_metoffice_fieldtype_cloudfrac: Not currently implemented aborting")
174 
175  case( ufo_metoffice_fieldtype_not_used ) ! currently unused
176  continue
177 
178  case default
179  write(*,*) 'invalid field type in b matrix file: ',j
180  cycle
181 
182  end select
183 
184  if ( firstelement /= 0 ) nelements = nelements + bmatrix % fields(j,2)
185 
186 end do
187 
188 self % nprofelements = nelements
189 
191 
192 !-------------------------------------------------------------------------------
193 !> Delete profile index
194 !!
195 !! \details Heritage: Ops_SatRad_InitProfInfo.f90
196 !!
197 !! Reset profile index
198 !!
199 !! \author Met Office
200 !!
201 !! \date 09/06/2020: Created
202 !!
204 
205 implicit none
206 class(ufo_rttovonedvarcheck_profindex), intent(inout) :: self !< profile index structure
207 
208 ! Zero all values
209 self % nprofelements = 0
210 self % t = 0
211 self % q = 0
212 self % ql = 0
213 self % qt = 0
214 self % qi = 0
215 self % cf = 0
216 self % o3total = 0
217 self % o3profile = 0
218 self % t2 = 0
219 self % q2 = 0
220 self % rh2 = 0
221 self % tstar = 0
222 self % pstar = 0
223 self % windspeed = 0
224 self % t70hpa = 0
225 self % t700hpa = 0
226 self % t950hpa = 0
227 self % t1000hpa = 0
228 self % qsurf = 0
229 self % lwp = 0
230 self % mwemiss = 0
231 self % cloudtopp = 0
232 self % cloudfrac = 0
233 self % emisspc = 0
234 
236 
237 ! ------------------------------------------------------------------------------
238 
Fortran module containing the full b-matrix data type and methods for the 1D-Var.
integer, parameter, public ufo_metoffice_fieldtype_cloudtopp
single-level cloud top pressure
integer, parameter, public ufo_metoffice_fieldtype_q
specific humidity profile
integer, parameter, public ufo_metoffice_fieldtype_q2
surface spec humidity
integer, parameter, public ufo_metoffice_fieldtype_qi
ice profile - not currently setup
integer, parameter, public ufo_metoffice_fieldtype_t
temperature
integer, parameter, public ufo_metoffice_fieldtype_t2
surface air temperature
integer, parameter, public ufo_metoffice_fieldtype_emisspc
emissivity prinipal components - not currently setup
integer, parameter, public ufo_metoffice_fieldtype_qt
total water profile
integer, parameter, public ufo_metoffice_fieldtype_ql
liquid water profile - not currently setup
integer, parameter, public ufo_metoffice_fieldtype_not_used
not currently in use - not currently setup
integer, parameter, public ufo_metoffice_fieldtype_windspeed
surface wind speed
integer, parameter, public ufo_metoffice_fieldtype_cf
cloud fraction profile - not currently setup
integer, parameter, public ufo_metoffice_fieldtype_o3profile
ozone - not currently setup
integer, parameter, public ufo_metoffice_fieldtype_mwemiss
microwave emissivity - not currently setup
integer, parameter, public ufo_metoffice_fieldtype_o3total
total column ozone - not currently setup
integer, parameter, public ufo_metoffice_fieldtype_cloudfrac
effective cloud fraction
integer, parameter, public ufo_metoffice_fieldtype_pstar
surface pressure
integer, parameter, public ufo_metoffice_fieldtype_tstar
surface skin temperature
Fortran module constants used throughout the rttovonedvarcheck filter.
Fortran module containing profile index.
subroutine ufo_rttovonedvarcheck_profindex_delete(self)
Delete profile index.
subroutine ufo_rttovonedvarcheck_profindex_setup(self, bmatrix, nlevels)
Profile index setup.