Skip to content

Commit 0514a26

Browse files
update error calculation for extinction to smart mode (part1)
1 parent 38953ec commit 0514a26

3 files changed

Lines changed: 113 additions & 4 deletions

File tree

lib/interface/picassoProcV3.m

Lines changed: 3 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -1567,10 +1567,11 @@
15671567
sig387 = squeeze(sum(data.signal(flag387FR, :, flagClFre), 3));
15681568
bg387 = squeeze(sum(data.bg(flag387FR, :, flagClFre), 3));
15691569

1570-
thisAerExt355_raman = pollyRamanExt_smart(data.distance0, sig387, 355, 387, mExt355(iGrp,:), mExt387(iGrp,:), number_density(iGrp, :), PollyConfig.angstrexp, PollyConfig.smoothWin_raman_355, 'moving');
1570+
[thisAerExt355_raman, thisAerExtStd355_raman] = pollyRamanExt_smart_MC(data.distance0, sig387, 355, 387, mExt355(iGrp,:), mExt387(iGrp,:), number_density(iGrp, :), PollyConfig.angstrexp, PollyConfig.smoothWin_raman_355, 'moving',15,bg387);
1571+
% thisAerExt355_raman = pollyRamanExt_smart(data.distance0, sig387, 355, 387, mExt355(iGrp,:), mExt387(iGrp,:), number_density(iGrp, :), PollyConfig.angstrexp, PollyConfig.smoothWin_raman_355, 'moving');
15711572
% thisAerExt355_raman = pollyRamanExt(data.distance0, sig387, 355, 387, PollyConfig.angstrexp, data.pressure(iGrp, :), data.temperature(iGrp, :) + 273.17, PollyConfig.smoothWin_raman_355, 380, 70, 'moving');
15721573

1573-
thisAerExtStd355_raman = pollyRamanExtStd(data.distance0, sig387, bg387, 355, 387, PollyConfig.angstrexp, data.pressure(iGrp, :), data.temperature(iGrp, :) + 273.17, PollyConfig.smoothWin_raman_355, 380, 70, 15);
1574+
%thisAerExtStd355_raman = pollyRamanExtStd(data.distance0, sig387, bg387, 355, 387, PollyConfig.angstrexp, data.pressure(iGrp, :), data.temperature(iGrp, :) + 273.17, PollyConfig.smoothWin_raman_355, 380, 70, 15);
15741575
data.aerExt355_raman(iGrp, :) = thisAerExt355_raman;
15751576
data.aerExtStd355_raman(iGrp, :) = thisAerExtStd355_raman;
15761577

lib/retrievals/pollyRamanExt_smart.m

Lines changed: 1 addition & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -31,7 +31,7 @@
3131
% molecular scattering coefficient at Raman wavelength in m^-1 sr^-1
3232
% number_density: array
3333
% molecular number density
34-
%
34+
% Mc_countNumber of MC iterations
3535
% OUTPUTS:
3636
% ext_aer: array
3737
% aerosol extinction coefficient [m^{-1}]
@@ -74,5 +74,4 @@
7474

7575
ext_aer = (deriv_ratio - alpha_molecular_elastic - alpha_molecular_Raman) ./ ...
7676
(1 + (lambda_emit ./ lambda_Raman) .^ angstrom);
77-
7877
end
Lines changed: 109 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,109 @@
1+
function [ ext_aer, ext_error] = pollyRamanExt_smart_MC(height, sig, lambda_emit, ...
2+
lambda_Raman,alpha_molecular_elastic, alpha_molecular_Raman, number_density, angstrom, window_size, method, MC_count,bg)
3+
% POLLYRAMANEXT etrieve the aerosol extinction coefficient with Raman method
4+
%
5+
% USAGE:
6+
% [ ext_aer ] = pollyRamanExt(height, sig, lambda_emit, lambda_Raman, angstrom, pressure, temperature, window_size, C, rh, method, measure_error)
7+
%
8+
% INPUTS:
9+
% height: array
10+
% height[m]
11+
% sig: array
12+
% measured raman signal. Unit: Photon Count
13+
% lambda_emit: float
14+
% the wavelength of the emitted laser beam.[nm]
15+
% lambda_Raman: float
16+
% the wavelength of raman sigal.[nm]
17+
% angstrom: float
18+
% the angstrom exponent for aerosol extinction coefficient
19+
% window_size: integer
20+
% window_size for smoothing the signal with sgolay filter.
21+
% order: integer
22+
% order of the implemented sgolay filter.
23+
% method: char
24+
% specify the method to calculate the slope of the signal. You can
25+
% choose from 'moving', 'smoothing' and 'chi2'.
26+
% measure_error: array
27+
% measurement error for each bin.
28+
% alpha_molecular_elastic: array
29+
% molecular scattering coefficient at emitted wavelength in m^-1 sr^-1
30+
% alpha_molecular_Raman: array
31+
% molecular scattering coefficient at Raman wavelength in m^-1 sr^-1
32+
% number_density: array
33+
% molecular number density
34+
% Mc_countNumber of MC iterations
35+
% OUTPUTS:
36+
% ext_aer: array
37+
% aerosol extinction coefficient [m^{-1}]
38+
%
39+
% REFERENCES:
40+
% https://bitbucket.org/iannis_b/lidar_processing
41+
%
42+
% Ansmann, A. et al. Independent measurement of extinction and backscatter profiles in cirrus clouds by using a combined Raman elastic-backscatter lidar. Applied Optics Vol. 31, Issue 33, pp. 7113-7131 (1992)
43+
%
44+
% HISTORY:
45+
% - 2021-05-31: first edition by Zhenping
46+
%
47+
% .. Authors: - zhenping@tropos.de
48+
49+
% default method is movingslope
50+
if ~ exist('method', 'var')
51+
method = 'movingslope';
52+
end
53+
54+
if ~ exist('measure_error', 'var')
55+
measure_error = zeros(size(sig));
56+
end
57+
58+
temp = number_density ./ (sig .* height.^2);
59+
temp(temp <= 0) = NaN;
60+
ratio = log(temp);
61+
62+
if strcmpi(method, 'moving') || strcmpi(method, 'movingslope')
63+
deriv_ratio = movingslope_variedWin(ratio, window_size) ./ ...
64+
[height(2) - height(1), diff(height)];
65+
elseif strcmpi(method, 'smoothing') || strcmpi(method, 'smooth')
66+
deriv_ratio = movingsmooth_variedWin(ratio, window_size) ./ [height(2) - ...
67+
height(1), diff(height)];
68+
elseif strcmpi(method, 'chi2')
69+
deriv_ratio = movingLinfit_variedWin(height, ratio, measure_error, ...
70+
window_size);
71+
else
72+
error('Please set a valid method for calculate the extinction coefficient.');
73+
end
74+
75+
ext_aer = (deriv_ratio - alpha_molecular_elastic - alpha_molecular_Raman) ./ ...
76+
(1 + (lambda_emit ./ lambda_Raman) .^ angstrom);
77+
78+
79+
%calculation of error
80+
if MC_count>1
81+
ext_aer_MC = NaN(MC_count, length(sig));
82+
83+
signalGen = sigGenWithNoise(sig, sqrt(sig + bg), MC_count, 'norm');
84+
85+
for iProfile = 1:MC_count
86+
%sig_MC = signalGen(:, iProfile)';
87+
temp = number_density ./ ((signalGen(:, iProfile)') .* height.^2);
88+
temp(temp <= 0) = NaN;
89+
ratio = log(temp);
90+
if strcmpi(method, 'moving') || strcmpi(method, 'movingslope')
91+
deriv_ratio = movingslope_variedWin(ratio, window_size) ./ ...
92+
[height(2) - height(1), diff(height)];
93+
elseif strcmpi(method, 'smoothing') || strcmpi(method, 'smooth')
94+
deriv_ratio = movingsmooth_variedWin(ratio, window_size) ./ [height(2) - ...
95+
height(1), diff(height)];
96+
elseif strcmpi(method, 'chi2')
97+
deriv_ratio = movingLinfit_variedWin(height, ratio, measure_error, ...
98+
window_size);
99+
else
100+
error('Please set a valid method for calculate the extinction coefficient.');
101+
end
102+
ext_aer_MC(iProfile, :) = (deriv_ratio - alpha_molecular_elastic - alpha_molecular_Raman) ./ ...
103+
(1 + (lambda_emit ./ lambda_Raman) .^ angstrom);
104+
end
105+
ext_error = nanstd(ext_aer_MC, 1, 0);
106+
else
107+
ext_error = NaN(length(ext_aer));
108+
end
109+
end

0 commit comments

Comments
 (0)