defdap.hrdic module¶
- class defdap.hrdic.Map(*args, **kwargs)[source]¶
Bases:
Map
Class to encapsulate DIC map data and useful analysis and plotting methods.
- corrVal¶
Correlation value.
- Type:
- ebsd_map¶
EBSD map linked to DIC map.
- Type:
- crop_dists¶
Crop distances (default all zeros).
- Type:
- data¶
- Must contain after loading data (maps):
- coordinatenumpy.ndarray
X and Y coordinates
- displacementnumpy.ndarray
X and Y displacements
- Generated data:
- fnumpy.ndarray
Components of the deformation gradient (0=x, 1=y).
- enumpy.ndarray
Components of the green strain (0=x, 1=y).
- max_shearnumpy.ndarray
Max shear component np.sqrt(((e11 - e22) / 2.)**2 + e12**2).
- Derived data:
Grain list data to map data from all grains
- Type:
- MAPNAME = 'hrdic'¶
- property original_shape¶
- property crystal_sym¶
- load_data(file_name, data_type=None)[source]¶
Load DIC data from file.
- Parameters:
file_name (pathlib.Path) – Name of file including extension.
data_type (str, {'Davis', 'OpenPIV'}) – Type of data file.
- load_corr_val_data(file_name, data_type=None)[source]¶
Load correlation value for DIC data
- Parameters:
file_name (pathlib.Path or str) – Path to file.
data_type (str, {'DavisImage'}) – Type of data file.
- set_scale(scale)[source]¶
Sets the scale of the map.
- Parameters:
scale (float) – Length of pixel in original BSE image in micrometres.
- property scale¶
Returns the number of micrometers per pixel in the DIC map.
- print_stats_table(percentiles, components)[source]¶
Print out a statistics table for a DIC map
- Parameters:
percentiles (list of float) – list of percentiles to print i.e. 0, 50, 99.
components (list of str) – list of map components to print i.e. e, f, max_shear.
- set_crop(*, left=None, right=None, top=None, bottom=None, update_homog_points=False)[source]¶
Set a crop for the DIC map.
- Parameters:
left (int) – Distance to crop from left in pixels (formally xMin)
right (int) – Distance to crop from right in pixels (formally xMax)
top (int) – Distance to crop from top in pixels (formally yMin)
bottom (int) – Distance to crop from bottom in pixels (formally yMax)
update_homog_points (bool, optional) – If true, change homologous points to reflect crop.
- crop(map_data, binning=None)[source]¶
Crop given data using crop parameters stored in map i.e. cropped_data = DicMap.crop(DicMap.data_to_crop).
- Parameters:
map_data (numpy.ndarray) – Bap data to crop.
binning (int) – True if mapData is binned i.e. binned BSE pattern.
- link_ebsd_map(ebsd_map, transform_type='affine', **kwargs)[source]¶
Calculates the transformation required to align EBSD dataset to DIC.
- Parameters:
ebsd_map (defdap.ebsd.Map) – EBSD map object to link.
transform_type (str, optional) – affine, piecewiseAffine or polynomial.
kwargs – All arguments are passed to estimate method of the transform.
- warp_to_dic_frame(map_data, **kwargs)[source]¶
Warps a map to the DIC frame.
- Parameters:
map_data (numpy.ndarray) – Data to warp.
kwargs – All other arguments passed to
defdap.experiment.Experiment.warp_map()
.
- Returns:
Map (i.e. EBSD map data) warped to the DIC frame.
- Return type:
- generate_threshold_mask(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 dic_map where max_shear is above 0.8, use:
>>> mask = dic_map.data.max_shear > 0.8
To remove data points in dic_map where e11 is above 1 or less than -1, use:
>>> mask = (dic_map.data.e[0, 0] > 1) | (dic_map.data.e[0, 0] < -1)
To remove data points in dic_map where corrVal is less than 0.4, use:
>>> mask = dic_map.corr_val < 0.4
Note: correlation value data needs to be loaded seperately from the DIC map, see
defdap.hrdic.load_corr_val_data()
- plot_grain_av_max_shear(**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.plot_grain_data_map()
.
- find_grains(algorithm=None, min_grain_size=10)[source]¶
Finds grains in the DIC map.
- Parameters:
algorithm (str {'warp', 'floodfill'}) – Use floodfill or warp algorithm.
min_grain_size (int) – Minimum grain area in pixels for floodfill algorithm.
- grain_inspector(vmax=0.1, correction_angle=0, rdr_line_length=3)[source]¶
Run the grain inspector interactive tool.
- Parameters:
vmax (float) – Maximum value of the colour map.
correction_angle (float) – Correction angle in degrees to subtract from measured angles to account for small rotation between DIC and EBSD frames. Approximately the rotation component of affine transform.
rdr_line_length (int) – Length of lines perpendicular to slip trace used to calculate RDR.
- class defdap.hrdic.Grain(grain_id, dicMap, group_id)[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:
- ebsd_grain¶
EBSD grain ID that this DIC grain corresponds to.
- Type:
- ebsd_map¶
EBSD map that this DIC grain belongs to.
- Type:
- points_list¶
Start and end points for lines drawn using defdap.inspector.GrainInspector.
- Type:
- groups_list¶
Groups, angles and slip systems detected for lines drawn using defdap.inspector.GrainInspector.
- data¶
- Must contain after creating:
- pointlist of tuples
(x, y) in cropped map
Generated data:
- Derived data:
Map data to list data from the map the grain is part of
- Type:
- plot_max_shear(**kwargs)[source]¶
Plot a maximum shear map for a grain.
- Parameters:
kwargs – All arguments are passed to
defdap.base.plot_grain_data()
.- Return type:
- property ref_ori¶
Returns average grain orientation.
- Return type:
- property slip_traces¶
Returns list of slip trace angles based on EBSD grain orientation.
- Return type:
- calc_slip_traces(slip_systems=None)[source]¶
Calculates list of slip trace angles based on EBSD grain orientation.
- Parameters:
slip_systems (defdap.crystal.SlipSystem, optional) –
- calc_slip_bands(grain_map_data, thres=None, min_dist=None)[source]¶
Use Radon transform to detect slip band angles.
- Parameters:
grain_map_data (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: