FV3 Bundle
CRTM_LowFrequency_MWSSEM.f90
Go to the documentation of this file.
1 !
2 ! CRTM_LowFrequency_MWSSEM
3 !
4 ! Module containgin routines to compute microwave ocean emissivity components
5 ! (FWD, TL, and AD) for low frequencies.
6 !
7 !
8 ! CREATION HISTORY:
9 ! Written by: Masahiro Kazumori, JCSDA 31-Jul-2006
10 ! Masahiro.Kazumori@noaa.gov
11 ! Quanhua Liu, QSS Group Inc.
12 ! Quanhua.Liu@noaa.gov
13 !
14 ! Refactored by: Paul van Delst, CIMSS/SSEC, April 2007
15 ! paul.vandelst@ssec.wisc.edu
16 !
17 
19 
20  ! -----------------
21  ! Environment setup
22  ! -----------------
23  ! Module use
24  USE type_kinds , ONLY: fp
25  USE fresnel , ONLY: fvar_type => ivar_type , &
29  USE guillou , ONLY: gpvar_type => ivar_type, &
33  USE ellison , ONLY: epvar_type => ivar_type, &
37  USE crtm_parameters , ONLY: pi, degrees_to_radians, &
38  zero, one, two, point_5
39  USE crtm_interpolation, ONLY: npts , &
40  lpoly_type , &
41  clear_lpoly , &
42  find_index , &
43  lpoly , &
44  lpoly_tl , &
45  lpoly_ad , &
46  interp_2d , &
47  interp_2d_tl, &
49  ! Disable implicit typing
50  IMPLICIT NONE
51 
52 
53  ! ------------
54  ! Visibilities
55  ! ------------
56  ! Everything private by default
57  PRIVATE
58  ! Data types
59  PUBLIC :: ivar_type
60  ! Science routines
61  PUBLIC :: lowfrequency_mwssem
62  PUBLIC :: lowfrequency_mwssem_tl
63  PUBLIC :: lowfrequency_mwssem_ad
64 
65 
66  ! -----------------
67  ! Module parameters
68  ! -----------------
69  ! RCS Id for the module
70  CHARACTER(*), PARAMETER :: module_rcs_id = &
71  '$Id: CRTM_LowFrequency_MWSSEM.f90 60152 2015-08-13 19:19:13Z paul.vandelst@noaa.gov $'
72 
73  ! Various quantities
74  REAL(fp), PARAMETER :: low_f_threshold = 20.0_fp ! Frequency threshold for permittivity models(GHz)
75  REAL(fp), PARAMETER :: smallscale_f_threshold = 15.0_fp ! Frequency threshold for small scale limit(GHz)
76  REAL(fp), PARAMETER :: foam_threshold = 7.0_fp ! Wind speed threshold (m/s) for foam
77  REAL(fp), PARAMETER :: ghz_to_hz = 1.0e+09_fp ! Frequency units conversion
78 
79  ! LUT dimension arrays
80  REAL(fp), PARAMETER :: d_frequency = 2.0_fp
81  INTEGER , PARAMETER :: n_frequencies = 11
82  REAL(fp), PARAMETER :: frequency_sdd(n_frequencies) = (/ &
83  3.0_fp, 5.0_fp, 7.0_fp, 9.0_fp, 11.0_fp, &
84  13.0_fp, 15.0_fp, 17.0_fp, 19.0_fp, 21.0_fp, &
85  23.0_fp /)
86 
87  REAL(fp), PARAMETER :: d_wind_speed = 0.5_fp
88  INTEGER , PARAMETER :: n_wind_speeds = 40
89  REAL(fp), PARAMETER :: wind_speed_sdd(n_wind_speeds) = (/ &
90  0.5_fp, 1.0_fp, 1.5_fp, 2.0_fp, 2.5_fp, 3.0_fp, 3.5_fp, 4.0_fp, &
91  4.5_fp, 5.0_fp, 5.5_fp, 6.0_fp, 6.5_fp, 7.0_fp, 7.5_fp, 8.0_fp, &
92  8.5_fp, 9.0_fp, 9.5_fp,10.0_fp,10.5_fp,11.0_fp,11.5_fp,12.0_fp, &
93  12.5_fp,13.0_fp,13.5_fp,14.0_fp,14.5_fp,15.0_fp,15.5_fp,16.0_fp, &
94  16.5_fp,17.0_fp,17.5_fp,18.0_fp,18.5_fp,19.0_fp,19.5_fp,20.0_fp /)
95 
96  ! Foam coverage coefficients
97  ! See Eqn.(1) in
98  ! Liu, Q. et al. (1998) Monte Carlo simulations of the
99  ! microwave emissivity of the sea surface.
100  ! JGR, Volume 103, No.C11, Pages 24983-24989
101  REAL(fp), PARAMETER :: fc1 = 7.751e-6_fp ! FWD model
102  REAL(fp), PARAMETER :: fc2 = 3.231_fp ! FWD model
103  REAL(fp), PARAMETER :: fc3 = fc1*fc2 ! TL model
104  REAL(fp), PARAMETER :: fc4 = fc2-one ! TL model
105 
106  ! Reflectivity large scale correction coefficients
107  REAL(fp), PARAMETER :: lsc_coeff(8) = &
108  (/ 5.209e-04_fp, 1.582e-05_fp, -3.510e-07_fp, & ! Rv coefficients
109  5.209e-04_fp, -1.550e-05_fp, -1.356e-07_fp, & ! Rh coefficients
110  1.90896_fp, 0.120448_fp /) ! Frequency coefficients
111 
112 
113  ! --------------------------------------
114  ! Structure definition to hold internal
115  ! variables across FWD, TL, and AD calls
116  ! --------------------------------------
117  TYPE :: ivar_type
118  PRIVATE
119  ! Validity indicator
120  LOGICAL :: is_valid = .false.
121  ! Forward model input values
122  REAL(fp) :: frequency = zero
123  REAL(fp) :: zenith_angle = zero
124  REAL(fp) :: temperature = zero
125  REAL(fp) :: salinity = zero
126  REAL(fp) :: wind_speed = zero
127  ! The zenith angle terms
128  REAL(fp) :: cos_z = zero, cos2_z = zero
129  ! The interpolating polynomials
130  TYPE(lpoly_type) :: flp ! Frequency
131  TYPE(lpoly_type) :: wlp ! Wind speed
132  ! The LUT interpolation indices
133  INTEGER :: i1 = -1, i2 = -1 ! Frequency
134  INTEGER :: j1 = -1, j2 = -1 ! Wind speed
135  ! The interpolation input
136  REAL(fp) :: f_int = zero ! Frequency
137  REAL(fp) :: w_int = zero ! Wind speed
138  ! The interpolation output
139  REAL(fp) :: sdd_int = zero ! Ocean surface height variance
140  ! Ocean permittivity
141  COMPLEX(fp) :: permittivity = zero
142  ! Fresnel reflectivties
143  REAL(fp) :: rv_fresnel = zero
144  REAL(fp) :: rh_fresnel = zero
145  ! Foam reflectivties
146  REAL(fp) :: rv_foam = zero
147  REAL(fp) :: rh_foam = zero
148  ! Foam coverage
149  REAL(fp) :: foam_cover = zero
150  ! Small scale correction intermediate term
151  REAL(fp) :: r_sterm = zero
152  ! Small scale correction to reflectivities
153  REAL(fp) :: rv_small = zero
154  REAL(fp) :: rh_small = zero
155  ! Large scale correction intermediate terms
156  REAL(fp) :: rv_lterm = zero
157  REAL(fp) :: rh_lterm = zero
158  ! Final reflectivity
159  REAL(fp) :: rv = zero
160  REAL(fp) :: rh = zero
161  ! Internal variables for subcomponents
162  TYPE(fvar_type) :: fvar
163  TYPE(epvar_type) :: epvar
164  TYPE(gpvar_type) :: gpvar
165  END TYPE ivar_type
166 
167 
168  ! -----------------------------------------------------
169  ! LUT data. Ocean surface height variance. Units of m^2
170  ! -----------------------------------------------------
172  INTEGER :: n
173 
174  DATA (sdd(n,1),n=1,n_frequencies) / &
175  3.45700e-03_fp, 9.60280e-03_fp, 1.86610e-02_fp, 2.63970e-02_fp, 2.99600e-02_fp, &
176  3.11630e-02_fp, 3.11300e-02_fp, 3.02250e-02_fp, 2.87810e-02_fp, 2.73210e-02_fp, &
177  2.57200e-02_fp /
178 
179  DATA (sdd(n,2),n=1,n_frequencies) / &
180  4.04780e-02_fp, 4.26570e-02_fp, 4.30630e-02_fp, 4.22650e-02_fp, 4.07540e-02_fp, &
181  3.94740e-02_fp, 3.78900e-02_fp, 3.60050e-02_fp, 3.41810e-02_fp, 3.24220e-02_fp, &
182  3.04150e-02_fp /
183 
184  DATA (sdd(n,3),n=1,n_frequencies) / &
185  4.59590e-02_fp, 4.48090e-02_fp, 4.46420e-02_fp, 4.36330e-02_fp, 4.24890e-02_fp, &
186  4.15140e-02_fp, 4.01950e-02_fp, 3.88670e-02_fp, 3.74780e-02_fp, 3.60030e-02_fp, &
187  3.53260e-02_fp /
188 
189  DATA (sdd(n,4),n=1,n_frequencies) / &
190  4.59450e-02_fp, 4.58550e-02_fp, 4.47040e-02_fp, 4.46270e-02_fp, 4.41120e-02_fp, &
191  4.32700e-02_fp, 4.28010e-02_fp, 4.19510e-02_fp, 4.10440e-02_fp, 4.08170e-02_fp, &
192  4.11710e-02_fp /
193 
194  DATA (sdd(n,5),n=1,n_frequencies) / &
195  4.64550e-02_fp, 4.56360e-02_fp, 4.59130e-02_fp, 4.55520e-02_fp, 4.53160e-02_fp, &
196  4.51440e-02_fp, 4.53070e-02_fp, 4.51610e-02_fp, 4.52720e-02_fp, 4.59090e-02_fp, &
197  4.67220e-02_fp /
198 
199  DATA (sdd(n,6),n=1,n_frequencies) / &
200  4.66200e-02_fp, 4.59520e-02_fp, 4.58590e-02_fp, 4.60200e-02_fp, 4.59480e-02_fp, &
201  4.64690e-02_fp, 4.65140e-02_fp, 4.70470e-02_fp, 4.72740e-02_fp, 4.85060e-02_fp, &
202  4.90510e-02_fp /
203 
204  DATA (sdd(n,7),n=1,n_frequencies) / &
205  4.61430e-02_fp, 4.62460e-02_fp, 4.63810e-02_fp, 4.67430e-02_fp, 4.77040e-02_fp, &
206  4.85840e-02_fp, 4.95640e-02_fp, 5.12130e-02_fp, 5.25020e-02_fp, 5.29830e-02_fp, &
207  5.38460e-02_fp /
208 
209  DATA (sdd(n,8),n=1,n_frequencies) / &
210  4.66920e-02_fp, 4.60310e-02_fp, 4.66870e-02_fp, 4.77130e-02_fp, 4.86120e-02_fp, &
211  4.99180e-02_fp, 5.17960e-02_fp, 5.37230e-02_fp, 5.46800e-02_fp, 5.54400e-02_fp, &
212  5.60010e-02_fp /
213 
214  DATA (sdd(n,9),n=1,n_frequencies) / &
215  4.63630e-02_fp, 4.69850e-02_fp, 4.79640e-02_fp, 4.94210e-02_fp, 5.30620e-02_fp, &
216  5.50180e-02_fp, 5.59020e-02_fp, 5.75140e-02_fp, 5.88560e-02_fp, 5.99310e-02_fp, &
217  6.00960e-02_fp /
218 
219  DATA (sdd(n,10),n=1,n_frequencies) / &
220  4.65880e-02_fp, 4.71020e-02_fp, 4.87240e-02_fp, 5.08410e-02_fp, 5.46400e-02_fp, &
221  5.65520e-02_fp, 5.82880e-02_fp, 5.97860e-02_fp, 6.10270e-02_fp, 6.20060e-02_fp, &
222  6.20490e-02_fp /
223 
224  DATA (sdd(n,11),n=1,n_frequencies) / &
225  4.63770e-02_fp, 4.79790e-02_fp, 5.07890e-02_fp, 5.49620e-02_fp, 5.83330e-02_fp, &
226  6.02050e-02_fp, 6.19130e-02_fp, 6.33790e-02_fp, 6.45780e-02_fp, 6.55040e-02_fp, &
227  6.61630e-02_fp /
228 
229  DATA (sdd(n,12),n=1,n_frequencies) / &
230  4.66490e-02_fp, 4.87630e-02_fp, 5.48840e-02_fp, 5.66550e-02_fp, 6.01050e-02_fp, &
231  6.20090e-02_fp, 6.37430e-02_fp, 6.52270e-02_fp, 6.64330e-02_fp, 6.73580e-02_fp, &
232  6.80090e-02_fp /
233 
234  DATA (sdd(n,13),n=1,n_frequencies) / &
235  4.73370e-02_fp, 5.04070e-02_fp, 5.78990e-02_fp, 5.98650e-02_fp, 6.35580e-02_fp, &
236  6.55950e-02_fp, 6.74360e-02_fp, 6.90000e-02_fp, 7.02620e-02_fp, 7.12200e-02_fp, &
237  7.18850e-02_fp /
238 
239  DATA (sdd(n,14),n=1,n_frequencies) / &
240  4.79880e-02_fp, 5.43970e-02_fp, 6.08900e-02_fp, 6.31260e-02_fp, 6.56630e-02_fp, &
241  6.93250e-02_fp, 7.01740e-02_fp, 7.19480e-02_fp, 7.33710e-02_fp, 7.44530e-02_fp, &
242  7.52120e-02_fp /
243 
244  DATA (sdd(n,15),n=1,n_frequencies) / &
245  4.88910e-02_fp, 5.87050e-02_fp, 6.30810e-02_fp, 6.58100e-02_fp, 6.87010e-02_fp, &
246  7.13730e-02_fp, 7.37010e-02_fp, 7.56450e-02_fp, 7.72040e-02_fp, 7.83900e-02_fp, &
247  7.83620e-02_fp /
248 
249  DATA (sdd(n,16),n=1,n_frequencies) / &
250  5.08950e-02_fp, 6.31620e-02_fp, 6.58390e-02_fp, 7.13070e-02_fp, 7.45450e-02_fp, &
251  7.61300e-02_fp, 7.88420e-02_fp, 7.99710e-02_fp, 8.18550e-02_fp, 8.23400e-02_fp, &
252  8.34470e-02_fp /
253 
254  DATA (sdd(n,17),n=1,n_frequencies) / &
255  5.54760e-02_fp, 6.46340e-02_fp, 6.81130e-02_fp, 7.41470e-02_fp, 7.61420e-02_fp, &
256  7.96030e-02_fp, 8.12890e-02_fp, 8.38430e-02_fp, 8.48200e-02_fp, 8.54770e-02_fp, &
257  8.67520e-02_fp /
258 
259  DATA (sdd(n,18),n=1,n_frequencies) / &
260  5.87440e-02_fp, 6.91110e-02_fp, 7.31410e-02_fp, 7.77400e-02_fp, 8.03610e-02_fp, &
261  8.28910e-02_fp, 8.64080e-02_fp, 8.81030e-02_fp, 8.94090e-02_fp, 9.03320e-02_fp, &
262  9.08970e-02_fp /
263 
264  DATA (sdd(n,19),n=1,n_frequencies) / &
265  6.37450e-02_fp, 7.11440e-02_fp, 7.64340e-02_fp, 8.18560e-02_fp, 8.50240e-02_fp, &
266  8.79790e-02_fp, 9.05320e-02_fp, 9.26230e-02_fp, 9.31090e-02_fp, 9.43720e-02_fp, &
267  9.52240e-02_fp /
268 
269  DATA (sdd(n,20),n=1,n_frequencies) / &
270  6.75330e-02_fp, 7.14800e-02_fp, 7.78540e-02_fp, 8.39610e-02_fp, 8.75990e-02_fp, &
271  9.09100e-02_fp, 9.23600e-02_fp, 9.47960e-02_fp, 9.67050e-02_fp, 9.70420e-02_fp, &
272  9.80910e-02_fp /
273 
274  DATA (sdd(n,21),n=1,n_frequencies) / &
275  6.83640e-02_fp, 7.42600e-02_fp, 8.16470e-02_fp, 8.63360e-02_fp, 9.07880e-02_fp, &
276  9.47130e-02_fp, 9.66050e-02_fp, 9.94340e-02_fp, 1.00460e-01_fp, 1.01120e-01_fp, &
277  1.02470e-01_fp /
278 
279  DATA (sdd(n,22),n=1,n_frequencies) / &
280  7.34960e-02_fp, 7.95300e-02_fp, 8.45850e-02_fp, 9.00820e-02_fp, 9.51380e-02_fp, &
281  9.79040e-02_fp, 1.00300e-01_fp, 1.02240e-01_fp, 1.04930e-01_fp, 1.05860e-01_fp, &
282  1.06410e-01_fp /
283 
284  DATA (sdd(n,23),n=1,n_frequencies) / &
285  7.10970e-02_fp, 8.02870e-02_fp, 8.95720e-02_fp, 9.32400e-02_fp, 9.89780e-02_fp, &
286  1.02230e-01_fp, 1.03530e-01_fp, 1.05910e-01_fp, 1.07750e-01_fp, 1.09070e-01_fp, &
287  1.09960e-01_fp /
288 
289  DATA (sdd(n,24),n=1,n_frequencies) / &
290  7.36820e-02_fp, 8.40180e-02_fp, 9.12120e-02_fp, 9.81670e-02_fp, 1.02370e-01_fp, &
291  1.04470e-01_fp, 1.07850e-01_fp, 1.10580e-01_fp, 1.11430e-01_fp, 1.13090e-01_fp, &
292  1.14270e-01_fp /
293 
294  DATA (sdd(n,25),n=1,n_frequencies) / &
295  8.18210e-02_fp, 8.72020e-02_fp, 9.52010e-02_fp, 1.00400e-01_fp, 1.05370e-01_fp, &
296  1.08020e-01_fp, 1.11860e-01_fp, 1.13550e-01_fp, 1.16100e-01_fp, 1.16820e-01_fp, &
297  1.18310e-01_fp /
298 
299  DATA (sdd(n,26),n=1,n_frequencies) / &
300  8.28020e-02_fp, 8.99060e-02_fp, 9.88250e-02_fp, 1.02260e-01_fp, 1.08030e-01_fp, &
301  1.11260e-01_fp, 1.15580e-01_fp, 1.17650e-01_fp, 1.19200e-01_fp, 1.20280e-01_fp, &
302  1.22110e-01_fp /
303 
304  DATA (sdd(n,27),n=1,n_frequencies) / &
305  8.31730e-02_fp, 9.21910e-02_fp, 9.91400e-02_fp, 1.06140e-01_fp, 1.10390e-01_fp, &
306  1.16000e-01_fp, 1.17460e-01_fp, 1.20050e-01_fp, 1.22050e-01_fp, 1.23510e-01_fp, &
307  1.25680e-01_fp /
308 
309  DATA (sdd(n,28),n=1,n_frequencies) / &
310  8.30700e-02_fp, 9.41130e-02_fp, 1.02090e-01_fp, 1.09790e-01_fp, 1.14520e-01_fp, &
311  1.18740e-01_fp, 1.22290e-01_fp, 1.23670e-01_fp, 1.26010e-01_fp, 1.27770e-01_fp, &
312  1.29050e-01_fp /
313 
314  DATA (sdd(n,29),n=1,n_frequencies) / &
315  8.91170e-02_fp, 9.98440e-02_fp, 1.04780e-01_fp, 1.13210e-01_fp, 1.16370e-01_fp, &
316  1.21240e-01_fp, 1.25310e-01_fp, 1.27100e-01_fp, 1.29790e-01_fp, 1.30580e-01_fp, &
317  1.32220e-01_fp /
318 
319  DATA (sdd(n,30),n=1,n_frequencies) / &
320  8.81060e-02_fp, 1.01160e-01_fp, 1.07230e-01_fp, 1.13940e-01_fp, 1.20090e-01_fp, &
321  1.25390e-01_fp, 1.28140e-01_fp, 1.30340e-01_fp, 1.33410e-01_fp, 1.34500e-01_fp, &
322  1.35220e-01_fp /
323 
324  DATA (sdd(n,31),n=1,n_frequencies) / &
325  9.35810e-02_fp, 1.02220e-01_fp, 1.12600e-01_fp, 1.16940e-01_fp, 1.23650e-01_fp, &
326  1.27510e-01_fp, 1.30780e-01_fp, 1.33410e-01_fp, 1.35440e-01_fp, 1.38270e-01_fp, &
327  1.39280e-01_fp /
328 
329  DATA (sdd(n,32),n=1,n_frequencies) / &
330  9.35500e-02_fp, 1.04710e-01_fp, 1.13220e-01_fp, 1.21530e-01_fp, 1.26710e-01_fp, &
331  1.31310e-01_fp, 1.35180e-01_fp, 1.38300e-01_fp, 1.40760e-01_fp, 1.41310e-01_fp, &
332  1.42850e-01_fp /
333 
334  DATA (sdd(n,33),n=1,n_frequencies) / &
335  9.85130e-02_fp, 1.09590e-01_fp, 1.18240e-01_fp, 1.24170e-01_fp, 1.29910e-01_fp, &
336  1.34960e-01_fp, 1.37450e-01_fp, 1.41030e-01_fp, 1.43870e-01_fp, 1.44760e-01_fp, &
337  1.46610e-01_fp /
338 
339  DATA (sdd(n,34),n=1,n_frequencies) / &
340  9.61960e-02_fp, 1.09960e-01_fp, 1.19920e-01_fp, 1.26660e-01_fp, 1.32990e-01_fp, &
341  1.36530e-01_fp, 1.41320e-01_fp, 1.43620e-01_fp, 1.46860e-01_fp, 1.48090e-01_fp, &
342  1.50270e-01_fp /
343 
344  DATA (sdd(n,35),n=1,n_frequencies) / &
345  1.00720e-01_fp, 1.14610e-01_fp, 1.21450e-01_fp, 1.29020e-01_fp, 1.35940e-01_fp, &
346  1.39920e-01_fp, 1.43320e-01_fp, 1.47690e-01_fp, 1.49730e-01_fp, 1.51310e-01_fp, &
347  1.52550e-01_fp /
348 
349  DATA (sdd(n,36),n=1,n_frequencies) / &
350  1.06720e-01_fp, 1.16120e-01_fp, 1.24380e-01_fp, 1.32860e-01_fp, 1.38180e-01_fp, &
351  1.44960e-01_fp, 1.48770e-01_fp, 1.50260e-01_fp, 1.52880e-01_fp, 1.54980e-01_fp, &
352  1.56700e-01_fp /
353 
354  DATA (sdd(n,37),n=1,n_frequencies) / &
355  1.03610e-01_fp, 1.15990e-01_fp, 1.28950e-01_fp, 1.34940e-01_fp, 1.40860e-01_fp, &
356  1.46100e-01_fp, 1.50500e-01_fp, 1.54100e-01_fp, 1.57010e-01_fp, 1.57970e-01_fp, &
357  1.60030e-01_fp /
358 
359  DATA (sdd(n,38),n=1,n_frequencies) / &
360  1.09100e-01_fp, 1.21630e-01_fp, 1.31540e-01_fp, 1.38450e-01_fp, 1.45050e-01_fp, &
361  1.50820e-01_fp, 1.53830e-01_fp, 1.57980e-01_fp, 1.59840e-01_fp, 1.62740e-01_fp, &
362  1.63900e-01_fp /
363 
364  DATA (sdd(n,39),n=1,n_frequencies) / &
365  1.13210e-01_fp, 1.21190e-01_fp, 1.32510e-01_fp, 1.40300e-01_fp, 1.47510e-01_fp, &
366  1.53760e-01_fp, 1.57150e-01_fp, 1.59950e-01_fp, 1.63770e-01_fp, 1.65520e-01_fp, &
367  1.67040e-01_fp /
368 
369  DATA (sdd(n,40),n=1,n_frequencies) / &
370  1.06990e-01_fp, 1.22540e-01_fp, 1.33950e-01_fp, 1.41810e-01_fp, 1.49100e-01_fp, &
371  1.53340e-01_fp, 1.58850e-01_fp, 1.61700e-01_fp, 1.64030e-01_fp, 1.65960e-01_fp, &
372  1.68950e-01_fp /
373 
374 
375 CONTAINS
376 
377 
378 !################################################################################
379 !################################################################################
380 !## ##
381 !## ## PUBLIC MODULE ROUTINES ## ##
382 !## ##
383 !################################################################################
384 !################################################################################
385 
386 
387 !--------------------------------------------------------------------------------
388 !
389 ! NAME:
390 ! LowFrequency_MWSSEM
391 !
392 ! PURPOSE:
393 ! Subroutine to compute microwave sea surface emissivity for the
394 ! frequency range 5GHz < f < 20GHz
395 !
396 ! CALLING SEQUENCE:
397 ! CALL LowFrequency_MWSSEM( &
398 ! Frequency , & ! Input
399 ! Zenith_Angle, & ! Input
400 ! Temperature , & ! Input
401 ! Salinity , & ! Input
402 ! Wind_Speed , & ! Input
403 ! Emissivity , & ! Output
404 ! iVar ) ! Internal variable output
405 !
406 ! INPUT ARGUMENTS:
407 ! Frequency: Microwave frequency. Valid input range is
408 ! UNITS: GHz
409 ! TYPE: REAL(fp)
410 ! DIMENSION: Scalar
411 ! ATTRIBUTES: INTENT(IN)
412 !
413 ! Zenith_Angle: Satellite zenith angle at the
414 ! sea surface
415 ! UNITS: Degrees
416 ! TYPE: REAL(fp)
417 ! DIMENSION: Scalar
418 ! ATTRIBUTES: INTENT(IN)
419 !
420 ! Temperature: Sea surface temperature
421 ! UNITS: Kelvin, K
422 ! TYPE: REAL(fp)
423 ! DIMENSION: Scalar
424 ! ATTRIBUTES: INTENT(IN)
425 !
426 ! Salinity: Water salinity
427 ! UNITS: ppt (parts per thousand)
428 ! TYPE: REAL(fp)
429 ! DIMENSION: Scalar
430 ! ATTRIBUTES: INTENT(IN)
431 !
432 ! Wind_Speed: Sea surface wind speed
433 ! UNITS: m/s
434 ! TYPE: REAL(fp)
435 ! DIMENSION: Scalar
436 ! ATTRIBUTES: INTENT(IN)
437 !
438 ! OUTPUT ARGUMENTS:
439 ! Emissivity: The surface emissivity at vertical
440 ! and horizontal polarizations.
441 ! UNITS: N/A
442 ! TYPE: REAL(fp)
443 ! DIMENSION: Rank-1
444 ! ATTRIBUTES: INTENT(OUT)
445 !
446 ! iVar: Structure containing internal variables required for
447 ! subsequent tangent-linear or adjoint model calls.
448 ! The contents of this structure are NOT accessible
449 ! outside of this module.
450 ! UNITS: N/A
451 ! TYPE: iVar_type
452 ! DIMENSION: Scalar
453 ! ATTRIBUTES: INTENT(OUT)
454 !
455 ! COMMENTS:
456 ! For computational speed, no input checking is done.
457 ! - Valid frequency range is 5GHz < f < 20GHz
458 ! - The size of the output emissivity array is not checked on input.
459 ! It is assumed the caller has dimensioned it correctly.
460 !
461 !--------------------------------------------------------------------------------
462 
463  SUBROUTINE lowfrequency_mwssem( Frequency , & ! Input
464  Zenith_Angle, & ! Input
465  Temperature , & ! Input
466  Salinity , & ! Input
467  Wind_Speed , & ! Input
468  Emissivity , & ! Output
469  iVar ) ! Internal variable output
470  ! Arguments
471  REAL(fp), INTENT(IN) :: frequency
472  REAL(fp), INTENT(IN) :: zenith_angle
473  REAL(fp), INTENT(IN) :: temperature
474  REAL(fp), INTENT(IN) :: salinity
475  REAL(fp), INTENT(IN) :: wind_speed
476  REAL(fp), INTENT(OUT) :: emissivity(:)
477  TYPE(ivar_type), INTENT(IN OUT) :: ivar
478  ! Local variables
479  LOGICAL :: f_outbound, w_outbound
480  REAL(fp) :: rv_large, rh_large
481 
482  ! Setup
483  ! ...Save forward input variables for TL and AD calculations
484  ivar%Frequency = frequency
485  ivar%Zenith_Angle = zenith_angle
486  ivar%Temperature = temperature
487  ivar%Salinity = salinity
488  ivar%Wind_Speed = wind_speed
489  ! ...Save derived variables
490  ivar%cos_z = cos(zenith_angle*degrees_to_radians)
491  ivar%cos2_z = ivar%cos_z * ivar%cos_z
492 
493 
494  ! Find the wind speed and frequency indices for interpolation
495  ivar%f_int = max(min(frequency_sdd(n_frequencies),frequency),frequency_sdd(1))
496  CALL find_index(frequency_sdd, d_frequency, ivar%f_int, ivar%i1, ivar%i2, f_outbound)
497 
498  ivar%w_int = max(min(wind_speed_sdd(n_wind_speeds),wind_speed),wind_speed_sdd(1))
499  CALL find_index(wind_speed_sdd, d_wind_speed, ivar%w_int, ivar%j1, ivar%j2, w_outbound)
500 
501 
502  ! Calculate the interpolating polynomials
503  ! ...Frequency term
504  CALL lpoly( frequency_sdd(ivar%i1:ivar%i2), & ! Input
505  ivar%f_int , & ! Input
506  ivar%flp ) ! Output
507  ! ...Wind speed term
508  CALL lpoly( wind_speed_sdd(ivar%j1:ivar%j2), & ! Input
509  ivar%w_int , & ! Input
510  ivar%wlp ) ! Output
511 
512 
513  ! Perform the 2-D interpolation
514  CALL interp_2d( sdd(ivar%i1:ivar%i2, ivar%j1:ivar%j2), &
515  ivar%flp, ivar%wlp, ivar%sdd_int )
516 
517 
518  ! Permittivity Calculation
519  IF( frequency < low_f_threshold ) THEN
520  CALL guillou_ocean_permittivity( temperature, salinity, frequency, &
521  ivar%Permittivity, &
522  ivar%gpVar )
523  ELSE
524  CALL ellison_ocean_permittivity( temperature, frequency, &
525  ivar%Permittivity, &
526  ivar%epVar )
527  END IF
528 
529 
530  ! Fresnel reflectivity calculation
531  CALL fresnel_reflectivity( ivar%Permittivity,ivar%cos_z, &
532  ivar%Rv_Fresnel,ivar%Rh_Fresnel,&
533  ivar%fVar )
534 
535 
536  ! Foam reflectivity calculation
537  CALL foam_reflectivity( zenith_angle, ivar%Rv_Foam, ivar%Rh_Foam )
538 
539 
540  ! Foam Coverage calculation
541  CALL foam_coverage( wind_speed, ivar%Foam_Cover )
542 
543 
544  ! Small scale correction calculation
545  CALL small_scale_correction( ivar%sdd_int, ivar%cos2_z, frequency, &
546  ivar%Rv_Small, ivar%Rh_Small, &
547  ivar%R_STerm )
548  ivar%Rv = ivar%Rv_Fresnel * ivar%Rv_Small
549  ivar%Rh = ivar%Rh_Fresnel * ivar%Rh_Small
550 
551 
552  ! Large Scale Correction Calculation
553  CALL large_scale_correction( wind_speed, zenith_angle, frequency, &
554  rv_large, rh_large, &
555  ivar%Rv_Lterm, ivar%Rh_Lterm )
556  ivar%Rv = ivar%Rv + rv_large
557  ivar%Rh = ivar%Rh + rh_large
558 
559 
560  ! Emissivity Calculation
561  emissivity(1) = one - (one-ivar%Foam_Cover)*ivar%Rv - ivar%Foam_Cover*ivar%Rv_Foam
562  emissivity(2) = one - (one-ivar%Foam_Cover)*ivar%Rh - ivar%Foam_Cover*ivar%Rh_Foam
563 
564 
565  ! Flag the internal variable structure as valid
566  ivar%Is_Valid = .true.
567 
568  END SUBROUTINE lowfrequency_mwssem
569 
570 
571 !--------------------------------------------------------------------------------
572 !
573 ! NAME:
574 ! LowFrequency_MWSSEM_TL
575 !
576 ! PURPOSE:
577 ! Subroutine to compute the tangent-linear microwave sea surface
578 ! emissivity for the frequency range 5GHz < f < 20GHz
579 !
580 ! CALLING SEQUENCE:
581 ! CALL LowFrequency_MWSSEM_TL( &
582 ! Temperature_TL, & ! TL Input
583 ! Salinity_TL , & ! TL Input
584 ! Wind_Speed_TL , & ! TL Input
585 ! Emissivity_TL , & ! TL Output
586 ! iVar ) ! Internal variable input
587 !
588 ! INPUT ARGUMENTS:
589 ! Temperature_TL: The tangent-linear sea surface temperature
590 ! UNITS: Kelvin, K
591 ! TYPE: REAL(fp)
592 ! DIMENSION: Scalar
593 ! ATTRIBUTES: INTENT(IN)
594 !
595 ! Salinity_TL: The tangent-linear water salinity
596 ! UNITS: ppt (parts per thousand)
597 ! TYPE: REAL(fp)
598 ! DIMENSION: Scalar
599 ! ATTRIBUTES: INTENT(IN)
600 !
601 ! Wind_Speed_TL: The tangent-linear sea surface wind speed
602 ! UNITS: m/s
603 ! TYPE: REAL(fp)
604 ! DIMENSION: Scalar
605 ! ATTRIBUTES: INTENT(IN)
606 !
607 ! iVar: Structure containing internal forward variables
608 ! required for tangent-linear calculations. This
609 ! structure is output from the forward model.
610 ! The contents of this structure are NOT accessible
611 ! outside of this module.
612 ! UNITS: N/A
613 ! TYPE: iVar_type
614 ! DIMENSION: Scalar
615 ! ATTRIBUTES: INTENT(IN)
616 !
617 ! OUTPUT ARGUMENTS:
618 ! Emissivity_TL: The tangent-linear sea surface emissivity at
619 ! vertical and horizontal polarizations.
620 ! UNITS: N/A
621 ! TYPE: REAL(fp)
622 ! DIMENSION: Rank-1
623 ! ATTRIBUTES: INTENT(OUT)
624 !
625 ! COMMENTS:
626 ! For computational speed, no input checking is done.
627 ! - Valid frequency range is 5GHz < f < 20GHz
628 ! - The size of the output emissivity array is not checked on input.
629 ! It is assumed the caller has dimensioned it correctly.
630 !
631 !--------------------------------------------------------------------------------
632 
633  SUBROUTINE lowfrequency_mwssem_tl( &
634  Temperature_TL, & ! TL Input
635  Salinity_TL , & ! TL Input
636  Wind_Speed_TL , & ! TL Input
637  Emissivity_TL , & ! TL Output
638  iVar ) ! Internal variable input
639  ! Arguments
640  REAL(fp), INTENT(IN) :: temperature_tl
641  REAL(fp), INTENT(IN) :: salinity_tl
642  REAL(fp), INTENT(IN) :: wind_speed_tl
643  REAL(fp), INTENT(OUT) :: emissivity_tl(:)
644  TYPE(ivar_type), INTENT(IN) :: ivar
645  ! Local variables
646  REAL(fp) :: sdd_int_tl
647  COMPLEX(fp) :: permittivity_tl
648  REAL(fp) :: foam_cover_tl
649  REAL(fp) :: rv_fresnel_tl, rh_fresnel_tl
650  REAL(fp) :: rv_foam_tl , rh_foam_tl
651  REAL(fp) :: rv_small_tl , rh_small_tl
652  REAL(fp) :: rv_large_tl , rh_large_tl
653  REAL(fp) :: rv_tl , rh_tl
654  REAL(fp) :: w_int_tl
655  REAL(fp), DIMENSION(NPTS) :: w_tl
656  REAL(fp), DIMENSION(NPTS,NPTS) :: sdd_tl
657  TYPE(lpoly_type) :: flp_tl, wlp_tl
658 
659 
660  ! Setup
661  ! ...Check internal structure
662  IF ( .NOT. ivar%Is_Valid ) THEN
663  emissivity_tl = zero
664  RETURN
665  END IF
666  ! ...Initialise local TL variables
667  w_int_tl = wind_speed_tl
668  w_tl = zero
669  sdd_tl = zero
670 
671 
672  ! Calculate the TL interpolating polynomials
673  ! ...Frequency term.
674  ! The TL frequency term is always zero (so far at least)
675  ! so we just need to initialise the polynomial so it has
676  ! no impact on the TL interpolation below.
677  CALL clear_lpoly(flp_tl)
678  ! Wind speed term
679  CALL lpoly_tl( wind_speed_sdd(ivar%j1:ivar%j2), & ! FWD Input
680  ivar%w_int , & ! FWD Input
681  ivar%wlp , & ! FWD Input
682  w_tl , & ! TL Input
683  w_int_tl , & ! TL Input
684  wlp_tl ) ! TL Output
685 
686 
687  ! Perform the 2-D TL interpolation
688  CALL interp_2d_tl( sdd(ivar%i1:ivar%i2, ivar%j1:ivar%j2), & ! FWD Input
689  ivar%flp, ivar%wlp, & ! FWD Input
690  sdd_tl, & ! TL Input
691  flp_tl , wlp_tl, & ! TL Input
692  sdd_int_tl ) ! TL Output
693 
694 
695  ! Permittivity Calculation
696  IF( ivar%Frequency < low_f_threshold ) THEN
697  CALL guillou_ocean_permittivity_tl( temperature_tl, salinity_tl, ivar%Frequency, &
698  permittivity_tl, &
699  ivar%gpVar)
700  ELSE
701  CALL ellison_ocean_permittivity_tl( temperature_tl, &
702  permittivity_tl, &
703  ivar%epVar )
704  END IF
705 
706 
707  ! Fresnel reflectivity calculation
708  CALL fresnel_reflectivity_tl( permittivity_tl, ivar%cos_z, &
709  rv_fresnel_tl, rh_fresnel_tl, &
710  ivar%fVar )
711 
712 
713  ! Foam reflectivity "calculation"
714  rv_foam_tl = zero
715  rh_foam_tl = zero
716 
717 
718  ! Foam coverage calculation
719  CALL foam_coverage_tl( ivar%Wind_Speed, wind_speed_tl, foam_cover_tl )
720 
721 
722  ! Small scale correction calculation
723  CALL small_scale_correction_tl( sdd_int_tl, ivar%cos2_z, ivar%Frequency, &
724  rv_small_tl, rh_small_tl, &
725  ivar%R_STerm )
726  rv_tl = ivar%Rv_Fresnel*rv_small_tl + rv_fresnel_tl*ivar%Rv_Small
727  rh_tl = ivar%Rh_Fresnel*rh_small_tl + rh_fresnel_tl*ivar%Rh_Small
728 
729 
730  ! Large scale correction calculation
731  CALL large_scale_correction_tl( wind_speed_tl, rv_large_tl, rh_large_tl, &
732  ivar%Rv_Lterm, ivar%Rh_Lterm )
733  rv_tl = rv_tl + rv_large_tl
734  rh_tl = rh_tl + rh_large_tl
735 
736 
737  ! Emissivity calculation
738  emissivity_tl(1) = (ivar%Foam_Cover-one)*rv_tl + &
739  (ivar%Rv-ivar%Rv_Foam)*foam_cover_tl - &
740  ivar%Foam_Cover*rv_foam_tl
741  emissivity_tl(2) = (ivar%Foam_Cover-one)*rh_tl + &
742  (ivar%Rh-ivar%Rh_Foam)*foam_cover_tl - &
743  ivar%Foam_Cover*rh_foam_tl
744 
745  END SUBROUTINE lowfrequency_mwssem_tl
746 
747 
748 !--------------------------------------------------------------------------------
749 !
750 ! NAME:
751 ! LowFrequency_MWSSEM_AD
752 !
753 ! PURPOSE:
754 ! Subroutine to compute the adjoint microwave sea surface emissivity
755 ! for the frequency range 5GHz < f < 20GHz
756 !
757 ! CALLING SEQUENCE:
758 ! CALL LowFrequency_MWSSEM_AD( &
759 ! Emissivity_AD , & ! AD Input
760 ! Temperature_AD, & ! AD Output
761 ! Salinity_AD , & ! AD Output
762 ! Wind_Speed_AD , & ! AD Output
763 ! iVar ) ! Internal variable input
764 !
765 ! INPUT ARGUMENTS:
766 ! Emissivity_AD: The adjoint sea surface emissivity at
767 ! vertical and horizontal polarizations.
768 ! UNITS: N/A
769 ! TYPE: REAL(fp)
770 ! DIMENSION: Rank-1
771 ! ATTRIBUTES: INTENT(IN OUT)
772 !
773 ! iVar: Structure containing internal forward variables
774 ! required for adjoint calculations. This
775 ! structure is output from the forward model.
776 ! The contents of this structure are NOT accessible
777 ! outside of this module.
778 ! UNITS: N/A
779 ! TYPE: iVar_type
780 ! DIMENSION: Scalar
781 ! ATTRIBUTES: INTENT(IN)
782 !
783 ! OUTPUT ARGUMENTS:
784 ! Temperature_AD: The adjoint sea surface temperature
785 ! UNITS: (Kelvin, K)^-1
786 ! TYPE: REAL(fp)
787 ! DIMENSION: Scalar
788 ! ATTRIBUTES: INTENT(IN)
789 !
790 ! Salinity_AD: The adjoint water salinity
791 ! UNITS: (ppt)^-1
792 ! TYPE: REAL(fp)
793 ! DIMENSION: Scalar
794 ! ATTRIBUTES: INTENT(IN)
795 !
796 ! Wind_Speed_AD: The adjoint sea surface wind speed
797 ! UNITS: (m/s)^-1
798 ! TYPE: REAL(fp)
799 ! DIMENSION: Scalar
800 ! ATTRIBUTES: INTENT(IN)
801 !
802 ! SIDE EFFECTS:
803 ! The input adjoint variable, Emissivity_AD, is set to zero upon
804 ! exiting this routine.
805 !
806 ! COMMENTS:
807 ! For computational speed, no input checking is done.
808 ! - Valid frequency range is 5GHz < f < 20GHz
809 ! - The size of the input emissivity array is not checked on input.
810 ! It is assumed the caller has dimensioned it correctly.
811 !
812 !--------------------------------------------------------------------------------
813 
814  SUBROUTINE lowfrequency_mwssem_ad( &
815  Emissivity_AD , & ! AD Input
816  Temperature_AD, & ! AD Output
817  Salinity_AD , & ! AD Output
818  Wind_Speed_AD , & ! AD Output
819  iVar ) ! Internal variable input
820  ! Arguments
821  REAL(fp), INTENT(IN OUT) :: emissivity_ad(:)
822  REAL(fp), INTENT(IN OUT) :: temperature_ad
823  REAL(fp), INTENT(IN OUT) :: salinity_ad
824  REAL(fp), INTENT(IN OUT) :: wind_speed_ad
825  TYPE(ivar_type), INTENT(IN) :: ivar
826  ! Local variables
827  REAL(fp) :: sdd_int_ad
828  COMPLEX(fp) :: permittivity_ad
829  REAL(fp) :: foam_cover_ad
830  REAL(fp) :: rv_fresnel_ad, rh_fresnel_ad
831  REAL(fp) :: rv_foam_ad , rh_foam_ad
832  REAL(fp) :: rv_small_ad , rh_small_ad
833  REAL(fp) :: rv_large_ad , rh_large_ad
834  REAL(fp) :: rv_ad , rh_ad
835  REAL(fp) :: w_int_ad
836  REAL(fp), DIMENSION(NPTS) :: w_ad
837  REAL(fp), DIMENSION(NPTS,NPTS) :: sdd_ad
838  TYPE(lpoly_type) :: flp_ad, wlp_ad
839 
840  ! Setup
841  ! ...Check internal structure
842  IF ( .NOT. ivar%Is_Valid ) THEN
843  temperature_ad = zero
844  salinity_ad = zero
845  wind_speed_ad = zero
846  RETURN
847  END IF
848  ! ...Initialise local adjoint variables
849  w_int_ad = zero
850  w_ad = zero
851  sdd_ad = zero
852  sdd_int_ad = zero
853  permittivity_ad = cmplx(zero,zero,fp)
854  CALL clear_lpoly(flp_ad)
855  CALL clear_lpoly(wlp_ad)
856 
857 
858  ! Emissivity calculation
859  rh_foam_ad = -ivar%Foam_Cover *emissivity_ad(2)
860  foam_cover_ad = (ivar%Rh-ivar%Rh_Foam)*emissivity_ad(2)
861  rh_ad = (ivar%Foam_Cover-one) *emissivity_ad(2)
862 
863  rv_foam_ad = -ivar%Foam_Cover *emissivity_ad(1)
864  foam_cover_ad = (ivar%Rv-ivar%Rv_Foam)*emissivity_ad(1) + foam_cover_ad
865  rv_ad = (ivar%Foam_Cover-one) *emissivity_ad(1)
866 
867  emissivity_ad = zero
868 
869 
870  ! Large scale correction calculation
871  rh_large_ad = rh_ad
872  rv_large_ad = rv_ad
873  CALL large_scale_correction_ad( rv_large_ad, rh_large_ad, wind_speed_ad, &
874  ivar%Rv_Lterm, ivar%Rh_Lterm )
875 
876  ! Small scale correction calculation
877  rv_small_ad = ivar%Rv_Fresnel*rv_ad
878  rv_fresnel_ad = ivar%Rv_Small*rv_ad
879  rv_ad = zero
880  rh_small_ad = ivar%Rh_Fresnel*rh_ad
881  rh_fresnel_ad = ivar%Rh_Small*rh_ad
882  rh_ad = zero
883  CALL small_scale_correction_ad( rv_small_ad, rh_small_ad, ivar%cos2_z, ivar%Frequency, &
884  sdd_int_ad, &
885  ivar%R_STerm )
886 
887  ! Foam coverage calculation
888  CALL foam_coverage_ad( ivar%Wind_Speed, foam_cover_ad, wind_speed_ad)
889 
890  ! Fresnel reflectivity calculation
891  CALL fresnel_reflectivity_ad( rv_fresnel_ad, rh_fresnel_ad, ivar%cos_z, &
892  permittivity_ad, &
893  ivar%fVar)
894 
895  ! Permittivity Calculation
896  IF( ivar%Frequency < low_f_threshold ) THEN
897  CALL guillou_ocean_permittivity_ad( permittivity_ad, ivar%Frequency, &
898  temperature_ad, salinity_ad, &
899  ivar%gpVar)
900  ELSE
901  CALL ellison_ocean_permittivity_ad( permittivity_ad, &
902  temperature_ad, &
903  ivar%epVar)
904  END IF
905 
906 
907  ! Perform the adjoint interpolation
908  CALL interp_2d_ad( sdd(ivar%i1:ivar%i2, ivar%j1:ivar%j2), & ! FWD Input
909  ivar%flp, ivar%wlp, & ! FWD Input
910  sdd_int_ad, & ! AD Input
911  sdd_ad, & ! AD Output
912  flp_ad , wlp_ad ) ! AD Output
913 
914  ! Compute the adjoint of the interpolating polynomials
915  ! ...Wind speed term
916  CALL lpoly_ad( wind_speed_sdd(ivar%j1:ivar%j2), & ! FWD Input
917  ivar%w_int , & ! FWD Input
918  ivar%wlp , & ! FWD Input
919  wlp_ad , & ! AD Input
920  w_ad , & ! AD Input
921  w_int_ad ) ! AD Output
922  ! ...No frequency term LPoly_AD.
923 
924 
925  ! The AD outputs
926  wind_speed_ad = wind_speed_ad + w_int_ad
927 
928  END SUBROUTINE lowfrequency_mwssem_ad
929 
930 
931 !################################################################################
932 !################################################################################
933 !## ##
934 !## ## PRIVATE MODULE ROUTINES ## ##
935 !## ##
936 !################################################################################
937 !################################################################################
938 
939  ! =============================================================
940  ! Routine for forward model foam reflectivity
941  ! Function dependence is on zenith angle only
942  ! so no TL or AD routine.
943  ! See Eqns(18.44a) (18.44b) in
944  ! Ulaby, F.T. et al. (1986) Microwave Remote Sensing, Active
945  ! and Passive, vol.3, From Theory to Applications, pp1457.
946  ! =============================================================
947  SUBROUTINE foam_reflectivity(z, Rv, Rh)
948  ! Arguments
949  REAL(fp), INTENT(IN) :: z ! Zenith angle
950  REAL(fp), INTENT(OUT) :: Rv, Rh
951  ! Parameters
952  REAL(fp), PARAMETER :: FR_COEFF(9) = &
953  (/ -9.946e-4_fp, 3.218e-5_fp, -1.187e-6_fp, &
954  7.e-20_fp, 0.07_fp, -1.748e-3_fp, &
955  -7.336e-5_fp, 1.044e-7_fp, -0.93_fp /)
956  ! Local variables
957  REAL(fp) :: Fv, Fh
958  ! The vertical component
959  fv = one + z*(fr_coeff(1)+ z*(fr_coeff(2) + z*fr_coeff(3))) + fr_coeff(4)*z**10
960  rv = fr_coeff(5)
961  ! The horizontal component
962  fh = one + z*(fr_coeff(6) + z*(fr_coeff(7) + z*fr_coeff(8)))
963  rh = one + fr_coeff(9)*fh
964  END SUBROUTINE foam_reflectivity
965 
966 
967  ! ======================================================
968  ! Foam coverage routines
969  ! See Eqn.(1) in
970  ! Liu, Q. et al. (1998) Monte Carlo simulations of the
971  ! microwave emissivity of the sea surface.
972  ! JGR, Volume 103, No.C11, Pages 24983-24989
973  ! ======================================================
974  ! Forward model
975  SUBROUTINE foam_coverage(wind_speed, coverage)
976  REAL(fp), INTENT(IN) :: wind_speed
977  REAL(fp), INTENT(OUT) :: coverage
978  IF ( wind_speed < foam_threshold ) THEN
979  coverage = zero
980  ELSE
981  coverage = fc1 * (wind_speed**fc2)
982  END IF
983  END SUBROUTINE foam_coverage
984 
985  ! Tangent-linear model
986  SUBROUTINE foam_coverage_tl(wind_speed, wind_speed_TL, coverage_TL)
987  REAL(fp), INTENT(IN) :: wind_speed
988  REAL(fp), INTENT(IN) :: wind_speed_TL
989  REAL(fp), INTENT(OUT) :: coverage_TL
990  IF ( wind_speed < foam_threshold ) THEN
991  coverage_tl = zero
992  ELSE
993  coverage_tl = fc3 * (wind_speed**fc4) * wind_speed_tl
994  END IF
995  END SUBROUTINE foam_coverage_tl
996 
997  ! Adjoint model
998  SUBROUTINE foam_coverage_ad(wind_speed, coverage_AD, wind_speed_AD)
999  REAL(fp), INTENT(IN) :: wind_speed ! Input
1000  REAL(fp), INTENT(IN OUT) :: coverage_AD ! Input
1001  REAL(fp), INTENT(IN OUT) :: wind_speed_AD ! Output
1002  IF ( .NOT. (wind_speed < foam_threshold) ) THEN
1003  wind_speed_ad = wind_speed_ad + fc3*(wind_speed**fc4)*coverage_ad
1004  END IF
1005  coverage_ad = zero
1006  END SUBROUTINE foam_coverage_ad
1007 
1008 
1009  ! ======================================================
1010  ! Reflectivity small scale correction routines
1011  ! See Eqns.(17a) and (17b) in
1012  ! Liu, Q. et al. (1998) Monte Carlo simulations of the
1013  ! microwave emissivity of the sea surface.
1014  ! JGR, Volume 103, No.C11, Pages 24983-24989
1015  ! ======================================================
1016  ! Forward model
1017  SUBROUTINE small_scale_correction( sdd, cos2_z, f, Rv, Rh, &
1018  R_term )
1019  ! Arguments
1020  REAL(fp), INTENT(IN) :: sdd ! Ocean surface height variance
1021  REAL(fp), INTENT(IN) :: cos2_z ! cos^2 of zenith angle
1022  REAL(fp), INTENT(IN) :: f ! Frequency
1023  REAL(fp), INTENT(OUT) :: Rv, Rh ! Reflectivity correction
1024  REAL(fp), INTENT(OUT) :: R_term ! Intermediate variable output
1025 
1026  IF ( f > smallscale_f_threshold ) THEN
1027  r_term = exp(-sdd*cos2_z)
1028  ELSE
1029  r_term = one
1030  END IF
1031  rv = r_term
1032  rh = r_term
1033  END SUBROUTINE small_scale_correction
1034 
1035  ! Tangent-linear model
1036  SUBROUTINE small_scale_correction_tl( sdd_TL, cos2_z, f, Rv_TL, Rh_TL, &
1037  R_term )
1038  ! Arguments
1039  REAL(fp), INTENT(IN) :: sdd_TL
1040  REAL(fp), INTENT(IN) :: cos2_z
1041  REAL(fp), INTENT(IN) :: f
1042  REAL(fp), INTENT(OUT) :: Rv_TL, Rh_TL
1043  REAL(fp), INTENT(IN) :: R_term ! Intermediate variable input
1044  ! Local variables
1045  REAL(fp) :: R_term_TL
1046 
1047  IF ( f > smallscale_f_threshold ) THEN
1048  r_term_tl = -cos2_z * r_term * sdd_tl
1049  ELSE
1050  r_term_tl = zero
1051  END IF
1052  rv_tl = r_term_tl
1053  rh_tl = r_term_tl
1054  END SUBROUTINE small_scale_correction_tl
1055 
1056  ! Adjoint model
1057  SUBROUTINE small_scale_correction_ad( Rv_AD, Rh_AD, cos2_z, f, sdd_AD, &
1058  R_term )
1059  ! Arguments
1060  REAL(fp), INTENT(IN OUT) :: Rv_AD, Rh_AD ! Input
1061  REAL(fp), INTENT(IN) :: cos2_z ! Input
1062  REAL(fp), INTENT(IN) :: f ! Input
1063  REAL(fp), INTENT(IN OUT) :: sdd_AD ! Output
1064  REAL(fp), INTENT(IN) :: R_term ! Intermediate variable input
1065  ! Local variables
1066  REAL(fp) :: R_term_AD
1067 
1068  r_term_ad = rv_ad ; rv_ad = zero
1069  r_term_ad = r_term_ad + rh_ad; rh_ad = zero
1070  IF ( f > smallscale_f_threshold ) THEN
1071  sdd_ad = sdd_ad - cos2_z*r_term*r_term_ad
1072  ELSE
1073  sdd_ad = zero
1074  END IF
1075  END SUBROUTINE small_scale_correction_ad
1076 
1077 
1078  ! ============================================
1079  ! Reflectivity large scale correction routines
1080  ! ============================================
1081  ! Forward model
1082  SUBROUTINE large_scale_correction( v, z, f, Rv, Rh, &
1083  Rv_term, Rh_term )
1084  ! Arguments
1085  REAL(fp), INTENT(IN) :: v ! Wind speed
1086  REAL(fp), INTENT(IN) :: z ! Zenith angle
1087  REAL(fp), INTENT(IN) :: f ! Frequency
1088  REAL(fp), INTENT(OUT) :: Rv, Rh ! Reflectivity correction
1089  REAL(fp), INTENT(OUT) :: Rv_term, Rh_term ! Intermediate variable output
1090  ! Local variables
1091  REAL(fp) :: f_term
1092 
1093  f_term = f / (lsc_coeff(7) + f*lsc_coeff(8))
1094  rv_term = (lsc_coeff(1) + z*(lsc_coeff(2) + z*lsc_coeff(3))) * f_term
1095  rh_term = (lsc_coeff(4) + z*(lsc_coeff(5) + z*lsc_coeff(6))) * f_term
1096 
1097  rv = rv_term * v
1098  rh = rh_term * v
1099  END SUBROUTINE large_scale_correction
1100 
1101  ! Tangent-linear model
1102  SUBROUTINE large_scale_correction_tl( v_TL, Rv_TL, Rh_TL, &
1103  Rv_term, Rh_term )
1104  ! Arguments
1105  REAL(fp), INTENT(IN) :: v_TL
1106  REAL(fp), INTENT(OUT) :: Rv_TL, Rh_TL
1107  REAL(fp), INTENT(IN) :: Rv_term, Rh_term ! Intermediate variable input
1108 
1109  rv_tl = rv_term * v_tl
1110  rh_tl = rh_term * v_tl
1111  END SUBROUTINE large_scale_correction_tl
1112 
1113  ! Adjoint model
1114  SUBROUTINE large_scale_correction_ad( Rv_AD, Rh_AD, v_AD, &
1115  Rv_term, Rh_term )
1116  ! Arguments
1117  REAL(fp), INTENT(IN OUT) :: Rv_AD, Rh_AD ! Input
1118  REAL(fp), INTENT(IN OUT) :: v_AD ! Output
1119  REAL(fp), INTENT(IN) :: Rv_term, Rh_term ! Intermediate variable input
1120 
1121  v_ad = v_ad + rv_term*rv_ad
1122  v_ad = v_ad + rh_term*rh_ad
1123 
1124  rv_ad = zero
1125  rh_ad = zero
1126  END SUBROUTINE large_scale_correction_ad
1127 
1128 END MODULE crtm_lowfrequency_mwssem
real(fp), parameter smallscale_f_threshold
subroutine, public ellison_ocean_permittivity(Temperature, Frequency, Permittivity, iVar)
Definition: Ellison.f90:176
real(fp), parameter, public zero
integer, parameter, public fp
Definition: Type_Kinds.f90:124
subroutine large_scale_correction_tl(v_TL, Rv_TL, Rh_TL, Rv_term, Rh_term)
subroutine, public fresnel_reflectivity_ad(Rv_AD, Rh_AD, cos_i, permittivity_AD, iVar)
Definition: Fresnel.f90:324
subroutine, public guillou_ocean_permittivity_ad(Permittivity_AD, Frequency, Temperature_AD, Salinity_AD, iVar)
Definition: Guillou.f90:539
subroutine foam_reflectivity(z, Rv, Rh)
subroutine, public clear_lpoly(p)
subroutine, public lpoly_ad(x, x_int, p, p_AD, x_AD, x_int_AD)
subroutine, public lowfrequency_mwssem_ad(Emissivity_AD, Temperature_AD, Salinity_AD, Wind_Speed_AD, iVar)
subroutine small_scale_correction(sdd, cos2_z, f, Rv, Rh, R_term)
subroutine large_scale_correction_ad(Rv_AD, Rh_AD, v_AD, Rv_term, Rh_term)
real(fp), dimension(n_frequencies), parameter frequency_sdd
subroutine, public ellison_ocean_permittivity_tl(Temperature_TL, Permittivity_TL, iVar)
Definition: Ellison.f90:300
real(fp), parameter, public one
subroutine foam_coverage(wind_speed, coverage)
subroutine foam_coverage_tl(wind_speed, wind_speed_TL, coverage_TL)
real(fp), parameter, public two
real(fp), dimension(n_frequencies, n_wind_speeds) sdd
real(fp), parameter, public degrees_to_radians
subroutine, public interp_2d_ad(z, ulp, vlp, z_int_AD, z_AD, ulp_AD, vlp_AD)
subroutine, public lpoly(x, x_int, p)
subroutine, public lowfrequency_mwssem(Frequency, Zenith_Angle, Temperature, Salinity, Wind_Speed, Emissivity, iVar)
integer, parameter, public npts
subroutine, public lowfrequency_mwssem_tl(Temperature_TL, Salinity_TL, Wind_Speed_TL, Emissivity_TL, iVar)
subroutine large_scale_correction(v, z, f, Rv, Rh, Rv_term, Rh_term)
character(*), parameter module_rcs_id
real(fp), parameter, public point_5
subroutine foam_coverage_ad(wind_speed, coverage_AD, wind_speed_AD)
subroutine small_scale_correction_tl(sdd_TL, cos2_z, f, Rv_TL, Rh_TL, R_term)
#define max(a, b)
Definition: mosaic_util.h:33
real(fp), dimension(n_wind_speeds), parameter wind_speed_sdd
subroutine, public guillou_ocean_permittivity(Temperature, Salinity, Frequency, Permittivity, iVar)
Definition: Guillou.f90:209
subroutine, public ellison_ocean_permittivity_ad(Permittivity_AD, Temperature_AD, iVar)
Definition: Ellison.f90:424
real(fp), dimension(8), parameter lsc_coeff
#define min(a, b)
Definition: mosaic_util.h:32
subroutine small_scale_correction_ad(Rv_AD, Rh_AD, cos2_z, f, sdd_AD, R_term)
subroutine, public lpoly_tl(x, x_int, p, x_TL, x_int_TL, p_TL)
subroutine, public fresnel_reflectivity(permittivity, cos_i, Rv, Rh, iVar)
Definition: Fresnel.f90:127
subroutine, public interp_2d_tl(z, ulp, vlp, z_TL, ulp_TL, vlp_TL, z_int_TL)
subroutine, public fresnel_reflectivity_tl(permittivity_TL, cos_i, Rv_TL, Rh_TL, iVar)
Definition: Fresnel.f90:222
subroutine, public guillou_ocean_permittivity_tl(Temperature_TL, Salinity_TL, Frequency, Permittivity_TL, iVar)
Definition: Guillou.f90:371
real(fp), parameter, public pi
subroutine, public interp_2d(z, ulp, vlp, z_int)