ppcpy.cloudmask#
ppcpy.cloudmask.cloudscreen#
- ppcpy.cloudmask.cloudscreen.smooth_signal(signal: ndarray, window_len: int) ndarray[source]#
Uniformly smooth the input signal
- Parameters:
- singalndarray
Signal to be smooth
- window_lenint
Width of the applied uniform filter
- Returns:
- ndarray
Smoothed signal
Notes
History
2026-02-04: Changed from scipy.ndimage.uniform_filter1d to ppcpy.misc.helper.uniform_filter
2026-03-17: Changed back to scipy.ndimage.uniform_filter1d due to NaN-related issues in Zhao’s cloud screening algorithm.
- ppcpy.cloudmask.cloudscreen.cloudscreen(data_cube, wv: int = 532, collect_debug: bool = False) ndarray[source]#
Preform cloud screening.
- Parameters:
- data_cubeobject
Main PicassoProc object.
- wvint, optional
Wavelength of the channel to perfom cloud screening [nm]. Default is 532.
- collect_debugbool, optional
If true, collects debug information. Default is False.
- Returns:
- falgCloudFreendarray
1 dimensional temporal boolean array. 0 = cloudy, 1 = cloud free.
Notes
The cloud screening is done for both the far range and near range total channels of wavelength wv if the channels exist.
History
xxxx-xx-xx: First edition by …
2026-03-18: Updated to include necessary configurations for cloudScreen_Zhao.
2026-04-09: Added ‘maxCloudSearchHeight’ config parameters.
- ppcpy.cloudmask.cloudscreen.cloudScreen_MSG(height: ndarray, signal: ndarray, slope_thres: float, search_region: list, smooth_win: int = 9) float[source]#
Cloud screen with maximum signal gradient.
- Parameters:
- heightndarray
Height in meters.
- signalndarray (time, height)
Range corrected photon count rate signal [MHz].
- slope_thresfloat
Threshold of the slope to determine whether there is strong backscatter signal. [MHz*m]
- search_regionlist or array (2 elements)
[baseHeight, topHeight] in meters.
- smooth_winint, optional
Width of the applied smoothing filter [bins]. Default is 9.
- Returns:
- flagCloudFreendarray
A boolean array indicates whether the profile is cloud free.
- layerStatus: ndarray (time x height)
Layer status for each bin (0: unknown, 1: cloud, 2: aerosol).
Notes
This function does not yield any layerStatus information.
History
2021-05-18: First edition by Zhenping.
2025-03-20: Translated into python.
- ppcpy.cloudmask.cloudscreen.cloudScreen_Zhao(height: ndarray, signal: ndarray, bg: ndarray, search_region: list = [0, 10000], minDepth: float = 100, heightFullOverlap: float = 600, smoothWin: int = 7, minSNR: float = 1, showDetails: bool = False) tuple[source]#
Cloud layer detection based on Zhao’s algorithm.
- Parameters:
- heightndarray
Height above ground [m].
- signalndarray (time, height)
Background corrected photon count rate signal [MHz].
- bgndarray
Background of the photon count rate signal [MHz].
- search_regionlist
Bottom and top height for cloud detection [m].
- minDepthfloat
Minimum layer depth [m]. Default is 100.
- heightFullOverlapfloat
Minimum height with full overlap [m]. Default is 600.
- smoothWinint
Smoothing window width [bins]. Default is 7.
- minSNRfloat
Minimum layer mean signal-noise-ratio. Default is 1.
- Returns:
- flagCloudFreendarray
A boolean array indicates whether the profile is cloud free.
- layerStatus: ndarray (time x height)
Layer status for each bin (0: unknown, 1: cloud, 2: aerosol).
References
Zhao, C., Y. Wang, Q. Wang, Z. Li, Z. Wang, and D. Liu (2014), A new cloud and aerosol layer detection method based on micropulse lidar measurements, Journal of Geophysical Research: Atmospheres, 119(11), 6788-6802.
History
2021-05-18: First edition by Zhengping.
2026-03-11: Translated into python.
- ppcpy.cloudmask.cloudscreen.VDE_cld(signal: ndarray, height: ndarray, BG: float, minLayerDepth: float = 0.2, minHeight: float = 0.4, smoothWin: int = 3, savgolSmoothWin: int = 9, minSNR: int = 5, showDetails: bool = False) tuple[source]#
Cloud layer detection with VDE method. This method only required elstic signal.
- Parameters:
- signalndarray (height)
raw signal without background [photon count].
- heightndarray
height above the ground [km].
- BGfloat
background signal [photon count].
- minLayerDepthfloat
minimum layer geometrical depth [km]. Default is 0.2.
- minHeightfloat
minimum height to start the searching with [km]. Default is 0.4.
- smoothWinint
smoothing window bins. Default is 3.
- savgolSmoothWinint
smoothing window bins of the Savitzky-Golay filter. Default is 9.
- minSNRint
minimum layer SNR to filter the fake features. Default is 5.
- Returns:
- layerInfolist of dicts
- idint
identity of the layer.
- baseHeightfloat
the layer base height [km].
- topHeightfloat
the layer top height [km].
- peakHeightfloat
the layer height with maximum backscatter signal [km].
- layerDepthint
geometrical depth of the layer [km].
- flagCloud: bool
cloud flag.
- PDndarray
SDP signal.
- PNndarray
VDE signal
Notes
Some of the for-loops could be vectorised for enhanced performance. However, not in a pretty way.
Also could consider using numba and its @njit decorator. That would also enhance the performance but the code needs to be rewritten to fit numba’s requirements.
This function does not work with NaN values in the signal.
History
2021-06-13: First edition by Zhenping.
2026-03-11: Translated into pyhton.
2026-04-09: Enhanced robustness against NaN values.
References
Zhao, C., Y. Wang, Q. Wang, Z. Li, Z. Wang, and D. Liu (2014), A new cloud and aerosol layer detection method based on micropulse lidar measurements, Journal of Geophysical Research: Atmospheres, 119(11), 6788-6802.
ppcpy.cloudmask.profilesegment#
- ppcpy.cloudmask.profilesegment.segment(data_cube) ndarray[source]#
Perform cloud free profile segmentation
- Parameters:
- data_cubeobject
Main PicassoProc object
- Returns:
- clFreeGrpsndarray
Time index of the detected cloud free regions [start, end].
- History
- xxxx-xx-xx: First edition by …
- 2026-03-18: Updated flagValPrf to include ‘shutterOnMask’ and ‘fogMask’.
- ppcpy.cloudmask.profilesegment.clFreeSeg(prfFlag, nIntPrf, minNIntPrf)[source]#
splits continuous cloud-free profiles into small sections.
- Parameters:
- prfFlag: array-like (boolean)
Cloud-free flags for each profile.
- nIntPrf: int
Number of integral profiles.
- minNIntPrf: int
Minimum number of integral profiles.
- Returns:
- clFreSegs: 2D numpy array
Start and stop indexes for each cloud-free section. [[start1, stop1], [start2, stop2], …]
Notes
History
2021-05-22: First edition by Zhenping
2025-03-20: Translated to python