import kanapy as knpy
import numpy as np
knpy.__version__
'6.1.6'
2. Analyse EBSD map with MTEX module and determine statistical microstructure descriptors¶
Microstructure geometry parameters are calculated based on the given EBSD map.
fname = '../EBSD_analysis/ebsd_316L_500x500.ang' # name of ang file to be imported
if knpy.MTEX_AVAIL:
# read EBSD map and evaluate statistics of microstructural features
ebsd = knpy.EBSDmap(fname)
ms_data = ebsd.ms_data[0] # analyse only data for majority phase with order parameter "0"
matname = ms_data['name'].lower().replace(' ', '_') # material name = MTEX Phase name (lower case and w/o blanks)
gs_param = ms_data['gs_param'] # parameters for lognorm distr of grain size: [std dev., location, scale]
ar_param = ms_data['ar_param'] # parameters for lognorm distr. of aspect ratios [std dev., location, scale]
om_param = ms_data['om_param'] # parameters for normal distribution of ellipsoid tilt angles [std dev., mean]
print('\n*** Statistical information on microstructure ***')
print(f'Found {len(ebsd.ms_data)} phases in EBSD map.')
for i, ms in enumerate(ebsd.ms_data):
print('\n----------------------------')
print(f'=== Phase {i+1}: ', ms['name'], '===')
print('---------------------------')
print(f'Volume fraction: {ms["vf"]}')
print('==== Grain size (equivalent grain diameter) ====')
print('-- Moments --')
print(f'median: {ms["gs_moments"][0].round(3)} micron, ' +
f'std. deviation: {ms["gs_moments"][1].round(3)} micron')
print('-- Fit parameters --')
print(f'scale: {ms["gs_param"][2].round(3)}, ' +
f'location: {ms["gs_param"][1].round(3)}, ' +
f'std. deviation: {ms["gs_param"][0].round(3)}')
print('==== Aspect ratio ====')
print('-- Moments --')
print(f'median: {ms["ar_moments"][0].round(3)}, ' +
f'std. deviation: {ms["ar_moments"][1].round(3)}')
print('-- Fit parameters --')
print(f'scale: {ms["ar_param"][2].round(3)}, ' +
f'location: {ms["ar_param"][1].round(3)}, ' +
f'std. deviation: {ms["ar_param"][0].round(3)}')
print('==== Tilt angle ====')
print('-- Moments --')
print(f'median: {(ms["om_moments"][0]*180/np.pi).round(3)}°, ' +
f'std. deviation: {(ms["om_moments"][1]*180/np.pi).round(3)}°')
print('-- Fit parameters --')
print(f'mean value: {ms["om_param"][1].round(3)}, ' +
f'std. deviation: {ms["om_param"][0].round(3)}')
else:
# MTEX is not available fall back to predefined values
print('***Warning: Anaysis of EBSD maps is only possible with an ' +
'existing MTEX installation in Matlab.')
print('\nWill continue with predefined settings for microstructure.')
matname = 'iron_fcc' # material name
gs_param = np.array([1.06083584, 6.24824603, 13.9309554])
ar_param = np.array([0.67525027, 0.76994992, 1.69901906])
om_param = np.array([0.53941709, 1.44160447])
print('*** Statistical information on microstructure ***')
print(f'=== Phase:', matname, '===')
print('==== Grain size (equivalent grain diameter) ====')
print(f'scale: {gs_param[2].round(3)}, ' +
f'location: {gs_param[1].round(3)}, ' +
f'std. deviation: {gs_param[0].round(3)}')
print('==== Aspect ratio ====')
print(f'scale: {ar_param[2].round(3)}, ' +
f'location: {ar_param[1].round(3)}, ' +
f'std. deviation: {ar_param[0].round(3)}')
print('==== Tilt angle ====')
print(f'mean value: {(om_param[1] * 180 / np.pi).round(3)}, ' +
f'std. deviation: {(om_param[0] * 180 / np.pi).round(3)}')
I found another version of MTEX and remove it from the current search path! initialize MTEX 5.5.2 ....Warning: "/Volumes/Extern/Codes/kanapy/libs/mtex/extern/nfft_openMP" not found in path. > In rmpath>doRemoveFromPath (line 102) In rmpath (line 59) In check_installation (line 7) In startup_mtex (line 79) In startup (line 9) ******************************************************************************** ******************************************************************************** done! For compatibility reasons MTEX is not using OpenMP. You may want to switch on OpenMP in the file <a href="matlab: edit mtex_settings">mtex_settings.m</a> I'm going to colorize the orientation data with the standard MTEX ipf colorkey. To view the colorkey do: ipfKey = ipfColorKey(ori_variable_name) plot(ipfKey)
Analyzed microstructure with 183 grains. Median values: equiv. diameter: 20.189 micron, aspect ratio: 2.469, tilt angle: 82.82° Std. dev: equivalent diameter: 35.189 micron, aspect ratio: 1.622, tilt angle: 30.966° *** Statistical information on microstructure *** Found 1 phases in EBSD map. ---------------------------- === Phase 1: Iron fcc === --------------------------- Volume fraction: 1.0 ==== Grain size (equivalent grain diameter) ==== -- Moments -- median: 20.189 micron, std. deviation: 35.189 micron -- Fit parameters -- scale: 13.962, location: 6.227, std. deviation: 1.059 ==== Aspect ratio ==== -- Moments -- median: 2.469, std. deviation: 1.622 -- Fit parameters -- scale: 1.697, location: 0.771, std. deviation: 0.676 ==== Tilt angle ==== -- Moments -- median: 82.82°, std. deviation: 30.966° -- Fit parameters -- mean value: 1.445, std. deviation: 4.066
3. Generate representative volume element (RVE)¶
In this part, a representative volume element (RVE) is generated based on the experimental input data.
Step 1: To generate the statistics dictionary required by Kanapy, the set_stats()
method is called with arguments defineing the statistical descriptors for grain size, aspect ratio and omega parameters as determined from the EBSD map. Other parameters need to be set manually. They include lower and upper cutoff values for distributions of equivalent diameters (deq_min
, deq_max
), aspect ratio (asp_min
, asp_max
) and tilt angle (omega_min
, omega_max
). Furthermore, the number of voxels (voxels
), side length of RVE (size
), its periodicity (periodicity
), and volume fraction, phase name and number (VF
, phasename
, and phasenum
) need to be specified.
# create dictionary with statistical information obtained from EBSD map
Nv = 30 # number of voxels per side
size = 50 # side length of generated RVE in micron
periodic = True # create periodic RVE
ms_stats = knpy.set_stats(gs_param, ar_param, om_param,
deq_min=8.0, deq_max=19.0, asp_min=0.95, asp_max=3.5,
omega_min=0.0, omega_max=2*np.pi, voxels=Nv, size=size,
periodicity=periodic,
VF=1.0, phasename=matname, phasenum=0)
print('*** Statistical microstructure descriptors ***')
ms_stats
*** Statistical microstructure descriptors ***
{'Grain type': 'Elongated', 'Equivalent diameter': {'sig': 1.0590239046147008, 'loc': 6.227370026015077, 'scale': 13.961517101996103, 'cutoff_min': 8.0, 'cutoff_max': 19.0}, 'RVE': {'sideX': 50, 'sideY': 50, 'sideZ': 50, 'Nx': 30, 'Ny': 30, 'Nz': 30}, 'Simulation': {'periodicity': True, 'output_units': 'um'}, 'Phase': {'Name': 'iron_fcc', 'Number': 0, 'Volume fraction': 1.0}, 'Aspect ratio': {'sig': 0.675672502695043, 'loc': 0.7714908342947019, 'scale': 1.6972691341554906, 'cutoff_min': 0.95, 'cutoff_max': 3.5}, 'Tilt angle': {'kappa': 4.06644119350967, 'loc': 1.4454743191447381, 'cutoff_min': 0.0, 'cutoff_max': 6.283185307179586}}
Step 2: Using the Microstructure
class of Kanapy, a microstructure object is generated, by passing the microstructure descriptors 'ms_stats' obtained from the EBSD analysis as input data to this class. The Python dictionary 'ms_stats' also contains the geometrical definition of the RVE box. After this, the method plot_stats_init()
plots the statistical distributions that will be used for the particle simulation. The experimental histogram of grain size and aspect ratio distributions are passed as argument 'gs_data' and 'ar_data'. In this case, the size distribution of the grains in the microstructure (orange bars) is fitted by a log-normal distribution (blue line). Finally, the method init_RVE()
creates the particle set that fulfills the required size distribution within the given lower (dashed vertical line) and upper (solid vertical line) bound.
# Create box with ellipsoids that represent future grains
ms = knpy.Microstructure(descriptor=ms_stats, name=matname+'_RVE')
if knpy.MTEX_AVAIL:
ms.plot_stats_init(gs_data=ms_data['gs_data'], ar_data=ms_data['ar_data'])
else:
ms.plot_stats_init()
ms.init_RVE()
Creating an RVE based on user defined statistics Analyzed statistical data for phase iron_fcc (0) Total number of particles = 55 RVE characteristics: RVE side lengths (X, Y, Z) = (50, 50, 50) (um) Number of voxels (X, Y, Z) = (30, 30, 30) Voxel resolution (X, Y, Z) = [1.6667 1.6667 1.6667](um) Total number of voxels = 27000 Considered phases (volume fraction): 0: iron_fcc (100.0%)
Step 3: Perform the particle simulation by invoking the method pack()
. To accomplish this, the particles generated in the last step are initally downscaled in size by a factor of 100 such that they can move freely through the simulation box but still fulfill the same charactistics of the size distribution with larger and smaller particles in the desired quantities. To avoid overlap of particles, collsions are considered. In this way, the particles form a realistic distribution within the space of the simulation box, as can be verified in the plot by plot_ellipsoids()
.
Note: As the particle simulations can be quite time consuming, particle configurations can be saved by passing the optional argument "save_files=True" to the pack()
method. Then, the last particle file can be imported with the import_particles()
method. The particles are saved a "dump" files that can also be imported and visualized with the Software package Ovito, see also the section "Visualize the packing simulation" in the Kanapy documentation.
ms.pack()
# **Alternative:** read particles from dump file, instead of executing ms.pack()
# Requires prior exection with particle save option activated: "ms.pack(save_files=True)"
# Import with: "ms.import_particles('dump_files/particle.699.dump')"
ms.plot_ellipsoids()
Starting particle simulation Creating particles from distribution statistics Particle packing by growth simulation
100%|█████████████████████████████████████████████████████████████████████████████████████████████████████████| 700/700 [00:05<00:00, 121.78it/s]
Completed particle packing 138 overlapping particles detected after packing Kinetic energy of particles after packing: 13.603902180400363 Initial kinetic energy: 13.603902180400363
Step 4: After the generation of the particle distribution of the desired volume fraction and fulfilling the experimentally found statistical distribution of particles, the method voxelize()
generates a voxelated structure in form of a regular 3D mesh of the simulation box. Each voxel lying inside a particle will be assigned to this phase. The remaining voxels will be assigned to the grains according to a growth and distribution algorithm.
In the figure created by the plot_voxels()
method, the individual grains are represented in arbitrary colors. With the optional argument "sliced" either the full RVE can be plotted (default) or one octand can be removed from the plot (sliced=True) to get an view into th einternal of the RVE.
# Create voxel structure and assign each voxel to a grain according to ellipsoids
# at this point, the mesh can be re-defined with the optional argument "dim", e.g. "ms.voxelize(dim=(20,20,20))"
ms.voxelize()
ms.plot_voxels(sliced=True)
Generating voxels inside RVE Starting RVE voxelization Assigning voxels to grains
100%|████████████████████████████████████████████████████████████████████████████████████████████████| 27000/27000 [00:00<00:00, 44340723.57it/s]
Completed RVE voxelization
Step 5: The grain geometry in the voxelized 3D microstructure (representative volume element, RVE) will be analysed statistically with respect to the grain size distribution and the aspect ratio distribution. This is accomplished by running the generate_grains()
method of the microstructure object, which performs a Delaunay tesseleation of the voxel structure by which a ployhedral hull around each grain is constructed. The resulting polyhedral representation of the microstructure can be plotted with the plot_grains()
method.
Furthermore, during the generation of grains, statistical data on the 3D microstructure is gathered and analyzed with respect to grain size and aspect ratio. The results are plotted with the plot_stats()
method, which compares the resulting grain microstructure statistics with that of the particles.
In this example, also the microstructural data contained gathered from the ebsd object is passed to the plot_stats()
method together with the fit parameters describing the statistics of that data, such that both information items can be plotted together with the grain statistics of the RVE.
# Create polygonized structure for voxelized grains and analyse statistics of synthetic microstructure
ms.generate_grains()
ms.plot_grains()
if knpy.MTEX_AVAIL:
ms.plot_stats(gs_param=ms_data['gs_param'], gs_data=ms_data['gs_data'],
ar_param=ms_data['ar_param'], ar_data=ms_data['ar_data'])
else:
ms.plot_stats()
Generated Delaunay tesselation of grain vertices. Assigning 1849 tetrahedra to grains ...
1849it [00:29, 63.14it/s] WARNING:root:No tet assigned to grain 14. WARNING:root:No tet assigned to grain 16. WARNING:root:Inconsistent volume of voxelized grains: 124523.14814814816, Reference volume: 125000. Grians missing in polyhedral structure. WARNING:root:5 grains are not represented in polyhedral geometry.
Finished generating polyhedral hulls for grains. Total volume of RVE: 125000 um^3 Total volume of polyhedral grains: 124999.99999999999 um^3 Mean relative error of polyhedral vs. voxel volume of individual grains: 0.194 for mean volume of grains: 2272.727 um^3. Computing the L1-error between input and output diameter distributions. L1 error phase 0 between particle and grain geometries: 0.69297
Plotting input & output statistics for phase 0
4. Generate grain orientations and export Abaqus input files¶
In this step, the texture module of Kanapy is invoked to assign a crystallographic orientation to each grain. The grain orientations are generated in form of Euler angles with the method ms.generate_orientations()
, which is based on the Matlab package MTEX.
The grain orientations are generated to follow the orientation distribution function (ODF) of the EBSD map, contained in the ebsd object given as argument. After the orientations for each grain have been generated, the grain structure can be plotted in colors of the IPF key, mimicking an EBSD map.
The geometrical information of the RVE is exported as Abaqus input file with the ending "_geom.inp" by the method write_abq()
. This input deck defines a finite element for each voxel. Furthermore, for each grain a set of elements is defined and assigned to a user defined material with the name "GRAIN{grain_id}_MAT". By the same method, the information on the grain orientations is exported as Abaqus input file with the ending "_mat.inp". This file contains the information on the material number to be used in the ICAMS Crystal Plasticity UMAT and the Euler angles for each grain, identified by the name "GRAIN{grain_id}_MAT". The _geom.inp file contains an *INCLUDE command to read the material definitions in the _mat.inp file.
ms.generate_orientations(ebsd) # downsample ODF of EBSD map to actual number of grains
ms.plot_voxels(ori=True) # plot voxelized grains in color code of IPF key
# write Abaqus input deck {name}_geom.inp for voxelated structure
# file {name}_mat.inp contains material definitions for ICAMS CP-UMAT together with Euler angles
ptag = 'pbc' if periodic else 'non-pbc'
ftrunc = f'from_ebsd_vox{Nv}_gr{ms.Ngr}_{ptag}'
ms.write_abq(nodes='v', file=ftrunc+'_geom.inp', ialloy=4) # parameter "ialloy" refers to material number of austenite in ICAMS CP-UMAT
Optimum kernel estimated from mean orientations of grains estimating optimal kernel halfwidth: 0% estimate optimal kernel halfwidth: 10% estimate optimal kernel halfwidth: 20% estimate optimal kernel halfwidth: 30% estimate optimal kernel halfwidth: 40% estimate optimal kernel halfwidth: 50% estimate optimal kernel halfwidth: 60% estimate optimal kernel halfwidth: 70% estimate optimal kernel halfwidth: 80% estimate optimal kernel halfwidth: 90% estimate optimal kernel halfwidth: 100% time = 10.7174 I found another version of MTEX and remove it from the current search path! initialize MTEX 5.5.2 ....Warning: "/Volumes/Extern/Codes/kanapy/libs/mtex/extern/nfft_openMP" not found in path. > In rmpath>doRemoveFromPath (line 102) In rmpath (line 59) In check_installation (line 7) In startup_mtex (line 79) In startup (line 9) ******************************************************************************** ******************************************************************************** done! For compatibility reasons MTEX is not using OpenMP. You may want to switch on OpenMP in the file <a href="matlab: edit mtex_settings">mtex_settings.m</a>
Writing RVE as ABAQUS file "./from_ebsd_vox30_gr55_pbc_geom.inp" Using brick element type C3D8. ---->DONE!
'./from_ebsd_vox30_gr55_pbc_geom.inp'
Note 1: To generate a valid Abaqus model the {}_geom.inp file can be imported directly into Abaqus CAE or included into a generic Abaqus deck (.inp file) via "*INCLUDE" commands. After that, the load step and boundary conditions must be specified.
Note 2: An Abaqus model with periodic boundary conditions can be generated by executing the Python script "Create_PeriodicBC_EDGE_3D.py" in the form
python Create_PeriodicBC_EDGE_3D.py {file_trunc}_geom.inp
This generates a subfolder-structure with different files. The resulting periodic geometry can by included into an Abaqus deck via the command "*Include, input=PeriodicData/geometry_Periodic.inp". An example for such an input deck for a periodic structure is avalaible in the file "Fatigue_BC_Amp_05.inp".
5. Safe voxelized structure¶
With the method write_voxels()
the information about geometry of the voxel structure and the grain orientations are written into a JSON file, which can be re-imported into Kanapy to safe the packing and voxelization steps, as demonstrated below.
# Save Kanapy microstructure as voxel file
ms.write_voxels(file=f'{matname}_voxels.json', script_name='ebsd2rve.ipynb', mesh=False, system=False)
Writing voxel information of microstructure to ./iron_fcc_voxels.json.
5.1 Import voxel file¶
JSON files with voxel structures can be imported instead of RVE generation and processed further. Note that also grain orientations are stored in the JSON file and can be recaptured.
ms2 = knpy.import_voxels(f'{matname}_voxels.json')
ms2.plot_voxels(ori=True)
ms2.generate_grains()
ms2.plot_stats()
Creating an RVE from voxel input RVE characteristics: RVE side lengths (X, Y, Z) = (50, 50, 50) (um) Number of voxels (X, Y, Z) = (30, 30, 30) Voxel resolution (X, Y, Z) = [1.6667 1.6667 1.6667](um) Total number of voxels = 27000 Considered phases (volume fraction): 0: iron_fcc (100.0%) Generating voxels inside RVE Voxel structure imported. I found another version of MTEX and remove it from the current search path! initialize MTEX 5.5.2 ....Warning: "/Volumes/Extern/Codes/kanapy/libs/mtex/extern/nfft_openMP" not found in path. > In rmpath>doRemoveFromPath (line 102) In rmpath (line 59) In check_installation (line 7) In startup_mtex (line 79) In startup (line 9) ******************************************************************************** ******************************************************************************** done! For compatibility reasons MTEX is not using OpenMP. You may want to switch on OpenMP in the file <a href="matlab: edit mtex_settings">mtex_settings.m</a>
Generated Delaunay tesselation of grain vertices. Assigning 1849 tetrahedra to grains ...
1849it [00:28, 64.63it/s] WARNING:root:No tet assigned to grain 14. WARNING:root:No tet assigned to grain 16. WARNING:root:Inconsistent volume of voxelized grains: 124523.14814814816, Reference volume: 125000. Grians missing in polyhedral structure. WARNING:root:5 grains are not represented in polyhedral geometry.
Finished generating polyhedral hulls for grains. Total volume of RVE: 125000 um^3 Total volume of polyhedral grains: 124999.99999999999 um^3 Mean relative error of polyhedral vs. voxel volume of individual grains: 0.194 for mean volume of grains: 2272.727 um^3. Plotting input & output statistics for phase 0