defdap.hrdic module¶
- class defdap.hrdic.Map(path, fname, dataType=None)[source]¶
Bases:
Map
Class to encapsulate DIC map data and useful analysis and plotting methods.
- xc¶
X coordinates.
- Type:
- yc¶
Y coordinates.
- Type:
- xd¶
X displacement.
- Type:
- yd¶
Y displacement.
- Type:
- corrVal¶
Correlation value.
- Type:
- ebsdMap¶
EBSD map linked to DIC map.
- Type:
- ebsdTransform¶
Transform from EBSD to DIC coordinates.
- Type:
various
- ebsdTransformInv¶
Transform from DIC to EBSD coordinates.
- Type:
various
- plotHomog¶
Map to use for defining homologous points (defaults to plotMaxShear).
- patScale¶
Size of pixel in loaded pattern relative to pixel size of dic data i.e 1 means they are the same size and 2 means the pixels in the pattern are half the size of the dic data.
- Type:
- self.x_map¶
Map of u displacement component along x.
- Type:
- self.y_map¶
Map of v displacement component along x.
- Type:
- f11, f22, f12, f21 ; numpy.ndarray
Components of the deformation gradient, where 1=x and 2=y.
- e11, e22, e12
Components of the green strain , where 1=x and 2=y.
- Type:
- eMaxShear¶
Max shear component np.sqrt(((e11 - e22) / 2.)**2 + e12**2).
- Type:
- cropDists¶
Crop distances (default all zeros).
- Type:
- property crystalSym¶
- setScale(micrometrePerPixel)[source]¶
Sets the scale of the map.
- Parameters:
micrometrePerPixel (float) – Length of pixel in original BSE image in micrometres.
- property scale¶
Returns the number of micrometers per pixel in the DIC map.
- setCrop(xMin=None, xMax=None, yMin=None, yMax=None, updateHomogPoints=False)[source]¶
Set a crop for the DIC map.
- Parameters:
- crop(mapData, binned=True)[source]¶
Crop given data using crop parameters stored in map i.e. cropped_data = DicMap.crop(DicMap.data_to_crop).
- Parameters:
mapData (numpy.ndarray) – Bap data to crop.
binned (bool) – True if mapData is binned i.e. binned BSE pattern.
- setHomogPoint(points=None, display=None, **kwargs)[source]¶
Set homologous points. Uses interactive GUI if points is None.
- Parameters:
points (list, optional) – homologous points to set.
display (string, optional) – Use max shear map if set to ‘maxshear’ or pattern if set to ‘pattern’.
- linkEbsdMap(ebsdMap, transformType='affine', **kwargs)[source]¶
Calculates the transformation required to align EBSD dataset to DIC.
- Parameters:
ebsdMap (defdap.ebsd.Map) – EBSD map object to link.
transformType (str, optional) – affine, piecewiseAffine or polynomial.
kwargs – All arguments are passed to estimate method of the transform.
- warpToDicFrame(mapData, cropImage=True, order=1, preserve_range=False)[source]¶
Warps a map to the DIC frame.
- Parameters:
mapData (numpy.ndarray) – Data to warp.
cropImage (bool, optional) – Crop to size of DIC map if true.
order (int, optional) – Order of interpolation (0: Nearest-neighbor, 1: Bi-linear…).
preserve_range (bool, optional) – Keep the original range of values.
- Returns:
Map (i.e. EBSD map data) warped to the DIC frame.
- Return type:
- warp_lines_to_dic_frame(lines)[source]¶
Warp a set of lines to the DIC reference frame.
- Parameters:
lines (list of tuples) – Lines to warp. Each line is represented as a tuple of start and end coordinates (x, y).
- Returns:
List of warped lines with same representation as input.
- Return type:
list of tuples
- property boundaries¶
Returns EBSD map grain boundaries warped to DIC frame.
- property boundaryLines¶
- property phaseBoundaryLines¶
- generateThresholdMask(mask, dilation=0, preview=True)[source]¶
Generate a dilated mask, based on a boolean array and previews the appication of this mask to the max shear map.
- Parameters:
mask (numpy.array(bool)) – A boolean array where points to be removed are True
dilation (int, optional) – Number of pixels to dilate the mask by. Useful to remove anomalous points around masked values. No dilation applied if not specified.
preview (bool) – If true, show the mask and preview the masked effective shear strain map.
Examples
To remove data points in dicMap where eMaxShear is above 0.8, use:
>>> mask = dicMap.eMaxShear > 0.8
To remove data points in dicMap where e11 is above 1 or less than -1, use:
>>> mask = (dicMap.e11 > 1) | (dicMap.e11 < -1)
To remove data points in dicMap where corrVal is less than 0.4, use:
>>> mask = dicMap.corrVal < 0.4
Note: correlation value data needs to be loaded seperately from the DIC map, see
defdap.hrdic.loadCorrValData()
- plotPattern(**kwargs)[source]¶
Plot BSE image of Map. For use with setting homog points.
- Parameters:
kwargs – All arguments are passed to
defdap.plotting.MapPlot.create()
.- Return type:
- plotMap(component, **kwargs)[source]¶
Plot a map from the DIC data.
- Parameters:
component – Map component to plot i.e. e11, f11, eMaxShear.
kwargs – All arguments are passed to
defdap.plotting.MapPlot.create()
.
- Returns:
Plot containing map.
- Return type:
- plotMaxShear(**kwargs)[source]¶
Plot a map of maximum shear strain.
- Parameters:
kwargs – All arguments are passed to
defdap.hrdic.plotMap()
.- Returns:
Plot containing map.
- Return type:
- plotGrainAvMaxShear(**kwargs)[source]¶
Plot grain map with grains filled with average value of max shear. This uses the max shear values stored in grain objects, to plot other data use
plotGrainAv()
.- Parameters:
kwargs – All arguments are passed to
defdap.base.Map.plotGrainDataMap()
.
- findGrains(algorithm=None, minGrainSize=10)[source]¶
Finds grains in the DIC map.
- Parameters:
algorithm (str {'warp', 'floodfill'}) – Use floodfill or warp algorithm.
minGrainSize (int) – Minimum grain area in pixels for floodfill algorithm.
- floodFill(x, y, grainIndex, points_left)[source]¶
Flood fill algorithm that uses the combined x and y boundary array to fill a connected area around the seed point. The points are inserted into a grain object and the grain map array is updated.
- Parameters:
x (int) – Seed point x for flood fill
y (int) – Seed point y for flood fill
grainIndex (int) – Value to fill in grain map
points_left (numpy.ndarray) – Boolean map of the points that have not been assigned a grain yet
- Returns:
currentGrain – New grain object with points added
- Return type:
- class defdap.hrdic.Grain(grainID, dicMap)[source]¶
Bases:
Grain
Class to encapsulate DIC grain data and useful analysis and plotting methods.
- dicMap¶
DIC map this grain is a member of
- Type:
- ownerMap¶
DIC map this grain is a member of
- Type:
- ebsdGrain¶
EBSD grain ID that this DIC grain corresponds to.
- Type:
- ebsdMap¶
EBSD map that this DIC grain belongs to.
- Type:
- pointsList¶
Start and end points for lines drawn using defdap.inspector.GrainInspector.
- Type:
- groupsList¶
Groups, angles and slip systems detected for lines drawn using defdap.inspector.GrainInspector.
- property plotDefault¶
- plotMaxShear(**kwargs)[source]¶
Plot a maximum shear map for a grain.
- Parameters:
kwargs – All arguments are passed to
defdap.base.plotGrainData()
.- Return type:
- property refOri¶
Returns average grain orientation.
- Return type:
- property slipTraces¶
Returns list of slip trace angles based on EBSD grain orientation.
- Return type:
- calcSlipTraces(slipSystems=None)[source]¶
Calculates list of slip trace angles based on EBSD grain orientation.
- Parameters:
slipSystems (defdap.crystal.SlipSystem, optional) –
- calcSlipBands(grainMapData, thres=None, min_dist=None)[source]¶
Use Radon transform to detect slip band angles.
- Parameters:
grainMapData (numpy.ndarray) – Data to find bands in.
thres (float, optional) – Normalised threshold for peaks.
min_dist (int, optional) – Minimum angle between bands.
- Returns:
Detected slip band angles
- Return type: