fmask

Implement the cloud and shadow algorithms known collectively as Fmask, as published in

Zhu, Z. and Woodcock, C.E. (2012). Object-based cloud and cloud shadow detection in Landsat imagery Remote Sensing of Environment 118 (2012) 83-94.

and

Zhu, Z., Wang, S. and Woodcock, C.E. (2015). Improvement and expansion of the Fmask algorithm: cloud, cloud shadow, and snow detection for Landsats 4-7, 8, and Sentinel 2 images Remote Sensing of Environment 159 (2015) 269-277.

Taken from Neil Flood’s implementation by permission.

The notation and variable names are largely taken from the paper. Equation numbers are also from the paper.

Input is a top of atmosphere (TOA) reflectance file.

The output file is a single thematic raster layer with codes representing null, clear, cloud, shadow, snow and water. These are the values 0-5 respectively, but there are constants defined for the different codes, as fmask.fmask.OUTCODE_*

fmask.fmask.accumHist(counts, vals)[source]

Accumulate the given values into the given (partial) counts

fmask.fmask.calcBTthresholds(otherargs)[source]

Calculate some global thresholds based on the results of the first pass

fmask.fmask.calcCDI(ref, fmaskConfig, refBands)[source]

Calculate the Cloud Displacement Index, as per Frantz et al (2018).

fmask.fmask.cloudFinalPass(info, inputs, outputs, otherargs)[source]

Called from RIOS

Final pass of cloud mask layer

fmask.fmask.cloudShapeFunc(info, inputs, outputs, otherargs)[source]

Called from RIOS.

Calculate the 3d cloud shape image. Requires that RIOS be run on the whole image at once, as substantial spatial structure is required. Needed to use RIOS because of the possibility of the thermal input being on a different resolution to the cloud input

Returns the result as whole arrays in otherargs, rather than writing them to output files, as they would just be read in again as whole arrays immediately.

fmask.fmask.clumpClouds(cloudmaskfile)[source]

Clump cloud pixels to make a layer of cloud objects. Currently assumes that the cloud mask contains only zeros and ones.

fmask.fmask.deleteRaster(filename)[source]

Use GDAL’s Driver.Delete() method to fully delete the given raster file

fmask.fmask.doCloudLayerFinalPass(fmaskFilenames, fmaskConfig, pass1file, pass2file, landThreshold, Tlow, missingThermal)[source]

Final pass

fmask.fmask.doFmask(fmaskFilenames, fmaskConfig)[source]

Main routine for whole Fmask algorithm. Calls all other routines in sequence. Parameters:

If fmask.config.FmaskConfig.setKeepIntermediates() has been called with True, then a dictionary of intermediate files will be returned. Otherwise None is returned.

fmask.fmask.doPotentialCloudFirstPass(fmaskFilenames, fmaskConfig, missingThermal)[source]

Run the first pass of the potential cloud layer. Also finds the temperature thresholds which will be needed in the second pass, because it has the relevant data handy.

fmask.fmask.doPotentialCloudSecondPass(fmaskFilenames, fmaskConfig, pass1file, Twater, Tlow, Thigh, missingThermal, nonNullCount)[source]

Second pass for potential cloud layer

fmask.fmask.doPotentialShadows(fmaskFilenames, fmaskConfig, NIR_17)[source]

Make potential shadow layer, as per section 3.1.3 of Zhu&Woodcock.

fmask.fmask.finalizeAll(fmaskFilenames, fmaskConfig, interimCloudmask, interimShadowmask, pass1file)[source]

Use the cloud and shadow masks to mask the snow layer (as per Zhu & Woodcock). Apply the optional extra buffer to the cloud mask, and write to final file.

fmask.fmask.focalVariance(img, winSize)[source]

Calculate the focal variance of the given 2-d image, over a moving window of size winSize pixels.

fmask.fmask.getIntersectionCoords(filelist)[source]

Use the RIOS utilities to get the correct area of intersection for a set of files, although we are not going to read the files using RIOS itself, but with GDAL directly.

fmask.fmask.make3Dclouds(fmaskFilenames, fmaskConfig, clumps, numClumps, missingThermal)[source]

Create 3-dimensional cloud objects from the cloud mask, and the thermal information. Assumes a constant lapse rate to convert temperature into height. Resulting cloud heights are relative to cloud base.

Returns an image of relative cloud height (relative to cloud base for each cloud object), and a dictionary of cloud base temperature, for each cloud object, and valueindexes.ValueIndexes object for use in extracting the location of every pixel for a given cloud object.

fmask.fmask.makeBufferKernel(buffsize)[source]

Make a 2-d array for buffering. It represents a circle of radius buffsize pixels, with 1 inside the circle, and zero outside.

fmask.fmask.makeCloudShadowShapes(fmaskFilenames, fmaskConfig, cloudShape, cloudClumpNdx)[source]

Project the 3d cloud shapes onto horizontal surface, along the sun vector, to make the 2d shape of the shadow.

fmask.fmask.maskAndBuffer(info, inputs, outputs, otherargs)[source]

Called from RIOS

Apply cloud and shadow masks to snow layer, and buffer cloud layer

The main aims of all this re-masking are:
  1. A pixel should be either cloud, shadow, snow or not, but never more than one

  2. Areas which are null in the input imagery should be null in the mask, even after buffering, etc.

fmask.fmask.matchOneShadow(cloudmask, shadowEntry, potentialShadow, Tcloudbase, Tlow, Thigh, xRes, yRes, cloudID, nullmask)[source]

Given the temperatures and sun angles for a single cloud object, and a shadow shape, search along the sun vector for a matching shadow object.

fmask.fmask.matchShadows(fmaskConfig, interimCloudmask, potentialShadowsFile, shadowShapesDict, cloudBaseTemp, Tlow, Thigh, pass1file)[source]

Match the cloud shadow shapes to the potential cloud shadows. Write an output file of the resulting shadow layer. Includes a 3-pixel buffer on the final shadows.

fmask.fmask.potentialCloudFirstPass(info, inputs, outputs, otherargs)[source]

Called from RIOS.

Calculate the first pass potential cloud layer (equation 6)

fmask.fmask.potentialCloudSecondPass(info, inputs, outputs, otherargs)[source]

Called from RIOS

Second pass of potential cloud layer

fmask.fmask.refDNtoUnits(refDN, fmaskConfig)[source]

Convert the given reflectance pixel value array to physical units, using parameters given in fmaskConfig.

Scaling is ref = (dn+offset)/scaleVal

fmask.fmask.scoreatpcnt(counts, pcnt)[source]

Given histogram counts (binned on the range 0-255), find the value which corresponds to the given percentile value (0-100).

fmask.fmask.singleRefDNtoUnits(refDN, scaleVal, offset)[source]

Apply the given scale and offset to transform a single band of reflectance from digital number (DN) to reflectance units.

Calculation is

ref = (refDN + offset) / scaleVal

fmask.fmask.B4_SCALE = 500.0

Gain to scale b4 reflectances to 0-255 for histograms

fmask.fmask.BT_OFFSET = 176

An offset so we can scale brightness temperature (BT, in deg C) to the range 0-255, for use in histograms.

fmask.fmask.OUTCODE_CLEAR = 1

Output pixel value for clear land

fmask.fmask.OUTCODE_CLOUD = 2

Output pixel value for cloud

fmask.fmask.OUTCODE_NULL = 0

Output pixel value for null

fmask.fmask.OUTCODE_SHADOW = 3

Output pixel value for cloud shadow

fmask.fmask.OUTCODE_SNOW = 4

Output pixel value for snow

fmask.fmask.OUTCODE_WATER = 5

Output pixel value for water

fmask.fmask.PROB_SCALE = 100.0

For scaling probability values so I can store them in 8 bits

fmask.fmask.RIOS_WINDOW_SIZE = 512

Global RIOS window size