Code documentation

Module API

Module defining class Microstructure that contains the necessary methods and attributes to analyze experimental microstructures in form of EBSD maps to generate statistical descriptors for 3D microstructures, and to create synthetic RVE that fulfill the required statistical microstructure descriptors.

The methods of the class Microstructure for an API that can be used to generate Python workflows.

Authors: Alexander Hartmaier, Golsa Tolooei Eshlghi, Abhishek Biswas Institution: ICAMS, Ruhr University Bochum

class kanapy.api.Microstructure(descriptor=None, file=None, name='Microstructure')[source]

Bases: object

Define class for synthetic microstructures

Attributes:
namestr

Name of microstructure

nphasesint

Number of phases in microstructure

ngrainsndarray

Array of grain number in each phase

Ngrint

Total number of grains summed over all phases

nparticleslist

List with number of particles in each phase

descriptorlist

List of dictionaries describing the microstructure of each phase; Dict Keys: “Grains type”, “Equivalent diameter”, “Aspect ratio”, “Tilt Angle”, “RVE”, “Simulation”

precipitNone or float
Indicates microstructure with precipitates/pores/particles in continuous matrix. If type is float,

it gives volume the fraction of that precipitate phase

from_voxelsbool

Indicates whether microstructure object is imported from voxel file, not generated from particle simulation

particleslist

List of particle objects of class entities containing information object particle geometries

rveobject of class RVE_creator

Contains information about the RVE Attributes: dim, size, nparticles, periodic, units, packing_steps, particle_data, phase_names, phase_vf, ialloy

simboxObject of class Simulation_Box

Contains information about geometry of simulation box for particle simulation

meshobject of class mesh_creator
Attributes: dim, grain_dict, grain_ori_dict, grain_phase_dict, grains, ngrains_phase. nodes, nodes_smooth,

nvox, phases, prec_vf_voxels, vox_center_dict, voxel_dict

geometrydict

Dictionary of grain geometries; Dict keys: “Vertices”, “Points”, “Simplices”, “Facets”, “Grains”, “GBnodes”, GBarea” “GBfaces” “Grains” : dictionary with key grain_number Keys:”Vertices”, “Points”, “Center”, “Simplices”, “Volume”, “Area”, “Phase”, “eqDia”, “majDia”, “minDia”

rve_statslist

List of dictionaries containing statistical information on different RVE types: particles, voxels, polyhedral grains

rve_stats_labelslist

List of strings containing the labels for the RVE types, i.e. Partcls, Voxels, Grains

generate_grains()[source]

Writes out the particle- and grain diameter attributes for statistical comparison. Final RVE grain volumes and shared grain boundary surface areas info are written out as well.

generate_orientations(data, ang=None, omega=None, Nbase=5000, hist=None, shared_area=None, iphase=None)[source]

Calculates the orientations of grains to give a desired crystallographic texture.

Parameters

data ang omega Nbase hist shared_area

Returns

import_particles(file, path='./')[source]
init_RVE(descriptor=None, nsteps=1000)[source]

Creates particle distribution inside simulation box (RVE) based on the data provided in the data file.

Parameters

descriptor nsteps

Returns

init_stats(descriptor=None, gs_data=None, ar_data=None, porous=False, save_files=False)[source]

Legacy function for plot_stats_init.

output_abq(nodes=None, name=None, voxel_dict=None, grain_dict=None, faces=None, dual_phase=False, thermal=False, units=None)[source]

Legacy function for write_abq.

output_ang(ori=None, cut='xy', data=None, plot=True, cs=None, pos=None, fname=None, matname='XXXX', save_files=True, dual_phase=False, save_plot=False)[source]

Convert orientation information of microstructure into a .ang file, mimicking an EBSD map. If polygonalized microstructure is available, it will be used as data basis, otherwise or if data=’voxels’ the voxelized microstructure will be exported. If no orientations are provided, each grain will get a random Euler angle. Values in ANG file: phi1 Phi phi2 X Y imageQuality confidenseIndex phase semSignal Fit(/mad) Output of ang file can be deactivated if called for plotting of slice.

Parameters

ori(self.Ngr, )-array, optional

Euler angles of grains. The default is None.

cutstr, optional

Define cutting plane of slice as ‘xy’, ‘xz’ or ‘yz’. The default is ‘xy’.

datastr, optional

Define data basis for plotting as ‘voxels’ or ‘poly’. The default is None.

plotbool, optional

Indicate if slice is plotted. The default is True.

posstr or float

Position in which slice is taken, either as absolute value, or as one of ‘top’, ‘bottom’, ‘left’, ‘right’. The default is None.

csstr, Optional

Crystal symmetry. Default is None

fnamestr, optional

Filename of ang file. The default is None.

matnamestr, optional

Name of the material to be written in ang file. The default is ‘XXXX’

save_filesbool, optional

Indicate if ang file is saved, The default is True.

Returns

fnamestr

Name of ang file.

output_neper()[source]

Writes out particle position and weights files required for tessellation in Neper.

pack(particle_data=None, k_rep=0.0, k_att=0.0, fill_factor=None, poly=None, save_files=False, verbose=False)[source]

Packs the particles into a simulation box.

pckl(file=None, path='./')[source]

Write microstructure into pickle file. Usefull for to store complex structures.

Parameters

filestring (optional, default: None)

File name for pickled microstructure. The default is None, in which case the filename will be the microstructure name + ‘.pckl’.

pathstring

Path to location for pickles

Returns

None.

plot_ellipsoids(cmap='prism', dual_phase=False)[source]

Generates plot of particles

plot_grains(geometry=None, cmap='prism', alpha=0.4, ec=[0.5, 0.5, 0.5, 0.1], dual_phase=False)[source]

Plot polygonalized microstructure

plot_particles(cmap='prism', dual_phase=False, plot_hull=True)[source]

Generates plot of particles

plot_slice(cut='xy', data=None, pos=None, fname=None, dual_phase=False, save_files=False)[source]

Plot a slice through the microstructure.

If polygonalized microstructure is available, it will be used as data basis, otherwise or if data=’voxels’ the voxelized microstructure will be plotted.

This subroutine calls the output_ang function with plotting active and writing of ang file deactivated.

Parameters

cutstr, optional

Define cutting plane of slice as ‘xy’, ‘xz’ or ‘yz’. The default is ‘xy’.

datastr, optional

Define data basis for plotting as ‘voxels’ or ‘poly’. The default is None.

posstr or float

Position in which slice is taken, either as absolute value, or as one of ‘top’, ‘bottom’, ‘left’, ‘right’. The default is None.

fnamestr, optional

Filename of PDF file. The default is None.

save_filesbool, optional

Indicate if figure file is saved and PDF. The default is False.

Returns

None.

plot_stats(data=None, gs_data=None, gs_param=None, ar_data=None, ar_param=None, dual_phase=False, save_files=False, show_all=False, verbose=False, silent=False, enhanced_plot=False)[source]

Plots the particle- and grain diameter attributes for statistical comparison.

plot_stats_init(descriptor=None, gs_data=None, ar_data=None, porous=False, get_res=False, show_res=False, save_files=False, silent=False)[source]

Plots initial statistical microstructure descriptors .

plot_voxels(sliced=False, dual_phase=False, cmap='prism', ori=None, color_key=0, silent=False)[source]

Generates plot of voxel structure

Parameters

sliced dual_phase cmap ori color_key: int

selects the color key: 0: iphHSVKey, 1: BungeColorKey, 2: ipfHKLKey (optional, default: 0)

silent

Returns

smoothen(nodes_v=None, voxel_dict=None, grain_dict=None)[source]

Generates smoothed grain boundary from a voxelated mesh.

voxelize(particles=None, dim=None)[source]

Generates the RVE by assigning voxels to grains.

write_abq(nodes=None, file=None, path='./', voxel_dict=None, grain_dict=None, dual_phase=False, thermal=False, units=None, ialloy=None, nsdv=200)[source]

Writes out the Abaqus deck (.inp file) for the generated RVE. The parameter nodes should be a string indicating if voxel (“v”) or smoothened (“s”) mesh should be written. It can also provide an array of nodal positions. If dual_phase is true, the generated deck contains plain material definitions for each phase. Material parameters must be specified by the user. If ialloy is provided, the generated deck material definitions for each grain. For dual phase structures to be used with crystal plasticity, ialloy can be a list with all required material definitions. If the list ialloy is shorted than the number of phases in the RVE, plain material definitions for the remaining phases will be included in the input deck.

Parameters

nodes file path voxel_dict grain_dict dual_phase thermal units ialloy nsdv

Returns

write_abq_ori(ialloy=None, ori=None, file=None, path='./', nsdv=200)[source]
write_centers(file=None, path='./', grains=None)[source]

Write grain center positions into CSV file.

write_ori(angles=None, file=None, path='./')[source]
write_stl(data='grains', file=None, path='./')[source]

Write triangles of convex polyhedra forming grains in form of STL files in the format ‘ solid name

facet normal n1 n2 n3
outer loop

vertex p1x p1y p1z vertex p2x p2y p2z vertex p3x p3y p3z

endloop

endfacet

endsolid name ‘

Returns

None.

write_voxels(angles=None, script_name=None, file=None, path='./', mesh=True, source=None, system=False)[source]

Write voxel structure into JSON file.

Parameters

angles script_name file path mesh source system

Returns

Module CLI

kanapy.cli.chkVersion(matlab)[source]

Read the version of Matlab

kanapy.cli.setPaths()[source]

Starts matlab engine, after installation if required, and initializes MTEX.

kanapy.cli.start()[source]

Module initialization

Classes for RVE and mesh initialization

Authors: Mahesh Prassad, Abhishek Biswas, Alexander Hartmaier Ruhr University Bochum, Germany

March 2024

class kanapy.initializations.NodeSets(nodes)[source]

Bases: object

CreatePeriodicEdgeSets(Nodes, sortCoord, NodeSet)[source]
CreatePeriodicNodeSets(Nodes, sortCoord1, sortCoord2, NodeSet)[source]
RemoveItemInList(NodeList, ReplaceNodeList)[source]
class kanapy.initializations.RVE_creator(stats_dicts, nsteps=1000, from_voxels=False, poly=None)[source]

Bases: object

Creates an RVE based on user-defined statistics

Note

  1. Input parameters provided by the user in the input file are:

  • Standard deviation for ellipsoid equivalent diameter (Log-normal distribution)

  • Mean value of ellipsoid equivalent diameter (Log-normal distribution)

  • Minimum and Maximum cut-offs for ellipsoid equivalent diameters

  • Mean value for aspect ratio

  • Mean value for ellipsoid tilt angles (Normal distribution)

  • Standard deviation for ellipsoid tilt angles (Normal distribution)

  • Side dimension of the RVE

  • Discretization along the RVE sides

  1. Particle, RVE and simulation data are written as JSON files in a folder in the current working directory for later access.

  • Ellipsoid attributes such as Major, Minor, Equivalent diameters and its tilt angle.

  • RVE attributes such as RVE (Simulation domain) size, the number of voxels and the voxel resolution.

  • Simulation attributes such as periodicity and output unit scale (\(mm\) or \(\mu m\)) for ABAQUS .inp file.

Generates ellipsoid size distribution (Log-normal) based on user-defined statistics

Note

  1. Input parameters provided by the user in the input file are:

  • Standard deviation for ellipsoid equivalent diameter (Normal distribution)

  • Mean value of ellipsoid equivalent diameter (Normal distribution)

  • Minimum and Maximum cut-offs for ellipsoid equivalent diameters

  • Mean value for aspect ratio

  • Mean value for ellipsoid tilt angles (Normal distribution)

  • Standard deviation for ellipsoid tilt angles (Normal distribution)

  • Side dimension of the RVE

  • Discretization along the RVE sides

  1. Particle, RVE and simulation data are written as JSON files in a folder in the current working directory for later access.

  • Ellipsoid attributes such as Major, Minor, Equivalent diameters and its tilt angle.

  • RVE attributes such as RVE (Simulation domain) size, the number of voxels and the voxel resolution.

  • Simulation attributes such as periodicity and output unit scale (\(mm\) or \(\mu m\)) for ABAQUS .inp file.

class kanapy.initializations.mesh_creator(dim)[source]

Bases: object

create_voxels(sim_box)[source]

Generates voxels inside the defined RVE (Simulation box)

Parameters:

sim_box (entities.Cuboid) – Simulation box representing RVE dimensions

Returns:

  • Node list containing coordinates.

  • Element dictionary containing element IDs and nodal connectivities.

  • Voxel dictionary containing voxel ID and center coordinates.

Return type:

Tuple of Python dictionaries.

get_ind(addr)[source]

Return the index in an array from a generic address.

Parameters

addr

Returns

kanapy.initializations.set_stats(grains, ar=None, omega=None, deq_min=None, deq_max=None, asp_min=None, asp_max=None, omega_min=None, omega_max=None, size=100, voxels=60, gtype='Elongated', rveunit='um', periodicity=False, VF=None, phasename=None, phasenum=None, save_file=False)[source]

grains = [sigma, loc , scale of lognorm distrib. of equiv. diameter] ar = [sigma, loc, scale of logorm distrib of aspect ratio] omega = [kappa, loc of von Mises distrib. of tilt angle]

kanapy.initializations.stat_names(legacy=False)[source]

Module input_output

kanapy.input_output.export2abaqus(nodes, file, grain_dict, voxel_dict, units='mm', gb_area=None, dual_phase=False, thermal=False, ialloy=None, grain_phase_dict=None, periodic=False)[source]

Creates an ABAQUS input file with microstructure morphology information in the form of nodes, elements and element sets. If “dual_phase” is true, element sets with phase numbers will be defined and assigned to materials “PHASE_{phase_id}MAT” plain material definitions for phases will be included. Otherwise, it will be assumed that each grain refers to a material “GRAIN_{}grain_id}MAT. In this case, a “_mat.inp” file with the same name trunc will be included, in which the alloy number and Euler angles for each grain must be defined.

Note

The nodal coordinates are written out in units of 1 mm or 1 \(\mu\) m scale, as requested by the user in the input file.

kanapy.input_output.import_stats(file, path='./')[source]

Write statistical descriptors of microstructure to JSON file.

Parameters

filestring

File name of pickled microstructure to be read.

pathstring

Path under which pickle-files are stored (optional, default: ‘./’)

Returns

desclist or dict

(List of) dict with statistical microstructure descriptors

kanapy.input_output.import_voxels(file, path='./')[source]
kanapy.input_output.pickle2microstructure(file, path='./')[source]

Read pickled microstructure file.

Parameters

filestring

File name of pickled microstructure to be read.

pathstring

Path under which pickle-files are stored (optional, default: ‘./’)

Returns

pclMaterial object

unpickled microstructure

kanapy.input_output.read_dump(file)[source]

Reads the (.dump) file to extract information for voxelization (meshing) routine

Parameters:

file (document) – Contains information of ellipsoids generated in the packing routine.

Returns:

  • Cuboid object that represents the RVE.

  • List of ellipsoid objects that represent the grains.

Return type:

Tuple of python objects (Cuboid, Ellipsoid)

kanapy.input_output.writeAbaqusMat(ialloy, angles, file=None, path='./', grain_phase_dict=None, nsdv=200)[source]

Export Euler angles to Abaqus input deck that can be included in the _geom.inp file. If parameter “grain_phase_dict” is given, the phase number for each grain will be used to select the corresponding ialloy from a list. If the list ialloy is shorter than the number of phases in grain_phase_dict, no angles for phases with no corresponding ialloy will be written.

Parameters:

ialloyint or list

Identifier, alloy number in ICAMS CP-UMAT: mod_alloys.f

anglesdict or (N, 3)-ndarray
Dict with Euler angles for each grain or array with number of N rows (= number of grains) and

three columns phi1, Phi, phi2

filestr

Filename, optional (default: None)

pathstr

Path to save file, option (default: ‘./’)

grain_phase_dict: dict

Dict with phase for each grain, optional (default: None)

nsdvint

Number of state dependant variables, optional (default: 200)

kanapy.input_output.write_dump(Ellipsoids, sim_box)[source]

Writes the (.dump) files into a sub-directory “dump_files”, which can be read by visualization software OVITO or imported again into Kanapy to avoid the packing simulation.

Parameters:
  • Ellipsoids (list) – Contains information of ellipsoids such as its position, axes lengths and tilt angles

  • sim_box (Cuboid) – Contains information of the dimensions of the simulation box

Note

This function writes (.dump) files containing simulation domain and ellipsoid attribute information.

kanapy.input_output.write_stats(stats, file, path='./')[source]

Write statistical descriptors of microstructure to JSON file.

Parameters

filestring

File name of pickled microstructure to be read.

pathstring

Path under which pickle-files are stored (optional, default: ‘./’)

Returns

desclist or dict

(List of) dict with statistical microstructure descriptors

Module packing

kanapy.packing.calculateForce(Ellipsoids, sim_box, periodicity, k_rep=0.0, k_att=0.0)[source]

Calculate interaction force between ellipsoids.

Parameters:

Ellipsoids

:param sim_box : :type sim_box: :param periodicity: :type periodicity: :param k_rep: optional, default is 0.0 :type k_rep: float :param k_att: optional, default is 0.0 :type k_att: float

kanapy.packing.packingRoutine(particle_data, periodic, nsteps, sim_box, k_rep=0.0, k_att=0.0, fill_factor=None, poly=None, save_files=False, verbose=False)[source]
The main function that controls the particle packing routine using:

particle_grow() & particle_generator()

Note

Particle, RVE and simulation data are read from the JSON files generated by kanapy.input_output.particleStatGenerator(). They contain the following information:

  • Ellipsoid attributes such as Major, Minor, Equivalent diameters and its tilt angle.

  • RVE attributes such as RVE (Simulation domain) size, the number of voxels and the voxel resolution.

  • Simulation attributes such as total number of timesteps and periodicity.

kanapy.packing.particle_generator(particle_data, sim_box, poly)[source]

Initializes ellipsoids by assigning them random positions and speeds within the simulation box.

Parameters:
  • particle_data (Python dictionary) – Ellipsoid information such as Major, Minor, Equivalent diameters and its tilt angle.

  • sim_box (entities.Simulation_Box) – Simulation box representing RVE.

  • poly – Points defining primitive polygon inside ellipsoid (optional).

Returns:

Ellipsoids for the packing routine

Return type:

list

kanapy.packing.particle_grow(sim_box, Ellipsoids, periodicity, nsteps, k_rep=0.0, k_att=0.0, fill_factor=None, dump=False, verbose=False)[source]

Initializes the entities.Octree class and performs recursive subdivision with collision checks and response for the ellipsoids. At each time step of the simulation it increases the size of the ellipsoid by a factor, which depends on the user-defined value for total number of time steps.

Parameters:
  • sim_box (entities.Simulation_Box) – Simulation box representing RVE.

  • Ellipsoids (list) – Ellipsoids for the packing routine.

  • periodicity (boolean) – Status of periodicity.

  • nsteps (int) – Total simulation steps to fill box volume with particle volume.

  • k_rep (float) – Repulsion factor for particles

  • k_att (float) – Attraction factor for particles

  • fill_factor (float) – Target volume fraction for particle filling

  • dump (boolean) – Indicate if dump files for particles are written.

  • verbose (bool) – Indicate if detailed output in iteration steps occurs

Note

kanapy.input_output.write_dump() function is called at each time step of the simulation to write output (.dump) files. By default, periodic images are written to the output file, but this option can be disabled within the function.

Module entities

class kanapy.entities.Cuboid(left, top, right, bottom, front, back)[source]

Bases: object

Creates Cuboid objects for ellipsoids and Octree sub-branches.

Parameters:
  • left (float) – Bounding box minimum along x

  • top (float) – Bounding box minimum along y

  • right (float) – Bounding box maximum along x

  • bottom (float) – Bounding box maximum along y

  • front (float) – Bounding box minimum along z

  • back (float) – Bounding box maximum along z

intersect(other)[source]

Evaluates whether the Cuboid object of the ellipsoid intersects with the Cuboid object of the Octree sub-branch.

Parameters:

other (object of the class Cuboid) – Sub-branch cuboid object of the octree

Returns:

if intersection - True, else False

Return type:

boolean

class kanapy.entities.Ellipsoid(iden, x, y, z, a, b, c, quat, phasenum=0, dup=None, points=None)[source]

Bases: object

Creates Ellipsoid objects for each ellipsoid generated from input statistics.

Parameters:
  • iden (integer) – ID of the ellipsoid

  • center (floats) – Position \((x, y, z)\) of the ellipsoid center in the simulation domain

  • coefficient (floats) – Semi-major and semin-minor axes lengths \((a, b, c)\) of the ellipsoid

  • quat (numpy array) – Quaternion representing ellipsoid’s axis and tilt angle with respect to the positive x-axis

Note

  1. The orientations of ellipsoid \(i\) in the global coordinate space is defined by its tilt angle and axis vector and expressed in quaternion notation as,

    _images/quaternion_ell.png
  2. Ellipsoids are initilaized without a value for its velocity, and is later assigned a random value by kanapy.packing.particle_generator.

  3. An empty list for storing voxels belonging to the ellipsoid is initialized.

Bbox()[source]

Generates the bounding box limits along x, y, z directions using the surface points from surfacePointsGen()

Return type:

numpy array

create_poly(points)[source]

Creates a polygon inside the ellipsoid

get_coeffs()[source]

Returns the coefficients array of the ellipsoid

Return type:

numpy array

get_cub()[source]

Returns the cuboid object of the ellipsoid

Return type:

object of the class Cuboid

get_pos()[source]

Returns the position array of the ellipsoid

Return type:

numpy array

get_volume()[source]

Returns the volume of the ellipsoid

Return type:

float

gravity_effect(value)[source]

Moves the ellipsoid downwards to mimick the effect of gravity acting on it

Parameters:

value (float) – User defined value for downward movement

Note

The Cuboid object of the ellipsoid has to be updated everytime it moves

growth(fac)[source]

Increases the size of the ellipsoid along its axes governed by the volume increment

Parameters:

dv (float) – Volume increment per time step

move(dt)[source]

Moves the ellipsoid by updating its position vector according to the Verlet integration method

Note

The Cuboid object of the ellipsoid has to be updated everytime it moves

rotationMatrixGen()[source]

Evaluates the rotation matrix for the ellipsoid using the quaternion

Return type:

numpy array

set_cub()[source]

Initializes an object of the class Cuboid using the bounding box limits from Bbox()

surfacePointsGen(nang=20)[source]

Generates points on the outer surface of the ellipsoid using the rotation matrix from rotationMatrixGen()

Return type:

numpy array

sync_poly(scale=None)[source]

Moves the center of the polygon to the center of the ellipsoid and scales the hull to fit inside the ellipsoid

wallCollision(sim_box, periodicity)[source]

Evaluates whether the ellipsoid collides with the boundaries of the simulation box.

  • If periodicity is enabled -> Creates duplicates of the ellipsoid on opposite faces of the box

  • If periodicity is disabled -> Mimicks the bouncing back effect.

Parameters:
  • sim_box (object of the class Simulation_Box) – Simulation box

  • periodicity (boolean) – Status of periodicity

Returns:

if periodic - ellipsoid duplicates, else None

Return type:

list

Note

The Cuboid object of the ellipsoid has to be updated everytime it moves

class kanapy.entities.Octree(level, cub, particles=[])[source]

Bases: object

Creates Octree objects for tree trunk and its sub-branches.

Parameters:
  • level (int) – Current level of the Octree

  • cub (object of the class Cuboid) – Cuboid object of the tree trunk / sub-branches

  • particles (list) – Particles within the tree trunk / sub-branches

Note

  1. level is set to zero for the trunk of the Octree.

  2. cub should be entire simulation box for the tree trunk.

  3. particles list contains all the ellipsoids in the simulation domain for the tree trunk.

add_particle(particle)[source]

Update the particle list of the Octree sub-branch

collisionsTest()[source]

Tests for collision between all ellipsoids in the particle list of octree

get_cub()[source]

Returns the cuboid object of the octree sub-branch

Return type:

object Cuboid

make_neighborlist()[source]

Finds the neighborlist for each particle

subdivide()[source]

Divides the given Octree sub-branch into eight further sub-branches and initializes each newly created sub-branch as an Octree object

subdivide_particles()[source]

Evaluates whether ellipsoids belong to a particular Octree sub-branch by calling intersect() on each ellipsoid.

update()[source]

Updates the Octree and begins the recursive process of subdividing or collision testing

class kanapy.entities.Simulation_Box(size)[source]

Bases: object

Creates Simulation_Box objects for the defined simulation domain.

kanapy.entities.cub_oct_split(cub)[source]

Splits cuboid object of the class Cuboid into eight smaller cuboid objects

Parameters:

cub (object of the class Cuboid) – Branch cuboid object containing ellipsoids

Returns:

Eight new sub-branch cuboid objects in a list

Return type:

List

Module collisions

kanapy.collisions.collide_detect(coef_i, coef_j, r_i, r_j, A_i, A_j)[source]

Implementation of Algebraic separation condition developed by W. Wang et al. 2001 for overlap detection between two static ellipsoids. Kudos to ChatGPT for support with translation from C++ code.

Parameters:
  • coef_i (numpy array) – Coefficients of ellipsoid \(i\)

  • coef_j (numpy array) – Coefficients of ellipsoid \(j\)

  • r_i (numpy array) – Position of ellipsoid \(i\)

  • r_j (numpy array) – Position of ellipsoid \(j\)

  • A_i (numpy array) – Rotation matrix of ellipsoid \(i\)

  • A_j (numpy array) – Rotation matrix of ellipsoid \(j\)

Returns:

True if ellipsoids \(i, j\) overlap, else False

Return type:

boolean

kanapy.collisions.collision_react(ell1, ell2)[source]

Evaluates and modifies the magnitude and direction of the ellipsoid’s velocity after collision.

Parameters:
  • ell1 (object Ellipsoid) – Ellipsoid \(i\)

  • ell2 (object Ellipsoid) – Ellipsoid \(j\)

Note

Consider two ellipsoids \(i, j\) at collision. Let them occupy certain positions in space defined by the position vectors \(\mathbf{r}^{i}, \mathbf{r}^{j}\) and have certain velocities represented by \(\mathbf{v}^{i}, \mathbf{v}^{j}\) respectively. The objective is to find the velocity vectors after collision. The elevation angle \(\phi\) between the ellipsoids is determined by,

_images/elevation_ell.png
where \(dx, dy, dz\) are defined as the distance between the two ellipsoid centers along

\(x, y, z\) directions given by,

_images/dist_ell.png

Depending on the magnitudes of \(dx, dz\) as projected on the \(x-z\) plane, the angle \(\Theta\) is computed. The angles \(\Theta\) and \(\phi\) determine the in-plane and out-of-plane directions along which the ellipsoid \(i\) would bounce back after collision. Thus, the updated velocity vector components along the \(x, y, z\) directions are determined by,

_images/velocities_ell.png
kanapy.collisions.collision_routine(E1, E2)[source]

Calls the method :meth:’collide_detect’ to determine whether the given two ellipsoid objects overlap using the Algebraic separation condition developed by W. Wang et al. A detailed description is provided therein.

Also calls the collision_react() to evaluate the response after collision.

Parameters:
  • E1 (object Ellipsoid) – Ellipsoid \(i\)

  • E2 (object Ellipsoid) – Ellipsoid \(j\)

Note

  1. If both the particles to be tested for overlap are spheres, then the bounding sphere hierarchy is sufficient to determine whether they overlap.

  2. Else, if either of them is an ellipsoid, then their coefficients, positions & rotation matrices are used to determine whether they overlap.

Module voxelization

kanapy.voxelization.assign_voxels_to_ellipsoid(cooDict, Ellipsoids, voxel_dict, vf_target=None)[source]

Determines voxels belonging to each ellipsoid

Parameters:
  • cooDict (Python dictionary) – Voxel dictionary containing voxel IDs and center coordinates.

  • Ellipsoids (list) – Ellipsoids from the packing routine.

  • voxel_dict (Python dictionary) – Element dictionary containing element IDs and nodal connectivities.

  • vf_target – target value for volume fraction of particles (optional, default: None)

Type:

float

kanapy.voxelization.points_in_convexHull(Points, hull)[source]

Determines whether the given array of points lie inside the convex hull or outside.

Parameters:
  • Points (numpy array) – Array of points to be tested whether they lie inside the hull or not.

  • hull (Scipy’s ConvexHull object) – Ellipsoid represented by a convex hull created from its outer surface points.

Returns:

Boolean values representing the status. If inside: True, else False

Return type:

numpy array

kanapy.voxelization.reassign_shared_voxels(cooDict, Ellipsoids, voxel_dict)[source]

Assigns shared voxels between ellipsoids to the ellipsoid with the closest center.

Parameters:
  • cooDict (Python dictionary) – Voxel dictionary containing voxel IDs and center coordinates.

  • Ellipsoids (list) – Ellipsoids from the packing routine.

  • voxel_dict (dict) – Dictionary of element definitions

kanapy.voxelization.voxelizationRoutine(Ellipsoids, mesh, nphases, prec_vf=None)[source]

The main function that controls the voxelization routine using: kanapy.input_output.read_dump(), create_voxels(), assign_voxels_to_ellipsoid(), reassign_shared_voxels()

Note

  1. The RVE attributes such as RVE (Simulation domain) size, the number of voxels and the voxel resolution is read by loading the JSON file that is generated by kanapy.input_output.read_dump().

  2. The following dictionaries are written as json files into a folder in the current working directory.

  • Node list containing coordinates.

  • Element dictionary containing element ID and nodal connectivities.

  • Element set dictionary containing element set ID and group of elements each representing a grain of the RVE.

Module grains

Subroutines for analysis of grains in microstructure

@author: Alexander Hartmaier ICAMS, Ruhr University Bochum, Germany March 2024

kanapy.grains.calc_polygons(rve, mesh, tol=0.001)[source]

Evaluates the grain volume and the grain boundary shared surface area between neighbouring grains in voxelized microstrcuture. Generate vertices as contact points between 3 or more grains. Generate polyhedral convex hull for vertices.

Parameters

rvekanapy object

Contains information about RVE geometry

meshkanapy object

Contains information about the voxel mesh and grain definitions for each voxel

tolfloat, optional

Tolerance value. The default is 1.e-3.

Returns

geometry : dict

ISSUES: for periodic structures, large grains with segments in both halves of the box and touching one boundary are split wrongly

Module texture

Tools for analysis of EBSD maps in form of .ang files

@author: Alexander Hartmaier, Abhishek Biswas, ICAMS, Ruhr-Universität Bochum

class kanapy.textures.EBSDmap(fname, matname=None, gs_min=3, vf_min=0.03, plot=True, hist=None)[source]

Bases: object

Class to store attributes and methods to import EBSD maps and filter out their statistical data needed to generate synthetic RVEs

calcORI(Ng, iphase=0, shared_area=None, nbins=12)[source]

Estimate optimum kernel half-width and produce reduced set of orientations for given number of grains.

Parameters

Ngint

Numbr of grains for which orientation is requested.

iphaseint, optional

Phase for which data is requested. The default is 0.

shared_areaarray, optional

Grain boundary data. The default is None.

nbinsint, optional

number of bins for GB texture. The default is 12.

Returns

ori(Ng, 3)-array

Array with Ng Euler angles.

showIPF()[source]

Plot IPF key.

Returns

None.

kanapy.textures.createOriset(num, ang, omega, hist=None, shared_area=None, cs='m-3m', Nbase=10000)[source]

Create a set of num Euler angles according to the ODF defined by the set of Euler angles ang and the kernel half-width omega. Example: Goss texture: ang = [0, 45, 0], omega = 5

Parameters

numint

Number of Euler angles in set to be created.

ang(3, ) or (M, 3) array

Set of Euler angles (in degrees) defining the ODF.

omegafloat

Half-width of kernel in degrees.

histarray, optional

Histogram of MDF. The default is None.

shared_area: array, optional

The shared area between pairs of grains. The default in None.

csstr, optional

Crystal symmetry group. The default is ‘m3m’.

Nbaseint, optional

Base number of orientations for random texture. The default is 10000

Returns

ori(num, 3) array

Set of Euler angles

kanapy.textures.createOrisetRandom(num, omega=7.5, hist=None, shared_area=None, cs='m-3m', Nbase=5000)[source]

Create a set of num Euler angles for Random texture. Other than knpy.createOriset() this method does not create an artificial EBSD which is reduced in a second step to num discrete orientations but directly samples num randomly distributed orientations.s

Parameters

numint

Number of Euler angles in set to be created.

omegafloat

Halfwidth of kernel in degrees (optional, default: 7.5)

histarray, optional

Histogram of MDF. The default is None.

shared_area: array, optional

The shared area between pairs of grains. The default in None.

csstr, optional

Crystal symmetry group. The default is ‘m3m’.

Nbaseint, optional

Base number of orientations for random texture. The default is 5000

Returns

ori(num, 3) array

Set of Euler angles

kanapy.textures.get_ipf_colors(ori_list, color_key=0)[source]

Get colors of list of orientations (in radians). Assumes cubic crystal symmetry and cubic specimen symmetry.

Parameters

ori_list: (N, 3) ndarray

List of N Euler angles in radians

Returns

colors: (N, 3) ndarray

List of RGB values

Module plotting

Subroutines for plotting of microstructures

Created on Mon Oct 4 07:52:55 2021

@author: Alexander Hartmaier, Golsa Tolooei Eshlaghi, Ronak Shoghi

kanapy.plotting.plot_ellipsoids_3D(particles, cmap='prism', dual_phase=False, silent=False, asp_arr=None)[source]

Display ellipsoids after packing procedure Parameters ———- particles : Class particles

Ellipsoids in microstructure before voxelization.

cmapcolor map, optional

Color map for ellipsoids. The default is ‘prism’.

dual_phasebool, optional

Whether to display the ellipsoids in red/green contrast or in colors

kanapy.plotting.plot_init_stats(stats_dict, gs_data=None, ar_data=None, gs_param=None, ar_param=None, save_files=False, silent=False)[source]

Plot initial microstructure descriptors, including cut-offs, based on user defined statistics

kanapy.plotting.plot_output_stats(data_list, labels, iphase=None, gs_data=None, gs_param=None, ar_data=None, ar_param=None, save_files=False, silent=False, enhanced_plot=False)[source]

Evaluates particle- and output RVE grain statistics with respect to Major, Minor & Equivalent diameters and plots the distributions.

kanapy.plotting.plot_particles_3D(particles, cmap='prism', dual_phase=False, plot_hull=True, silent=False, asp_arr=None)[source]

Display inner polyhedra of ellipsoids after packing procedure

Parameters

particlesClass particles

Ellipsoids in microstructure before voxelization.

cmapcolor map, optional

Color map for ellipsoids. The default is ‘prism’.

dual_phasebool, optional

Whether to display the ellipsoids in red/green contrast or in colors

plot_hullbool, optional

Whether to display the ellipsoids as hulls. Defaults to True

Returns

None.

kanapy.plotting.plot_polygons_3D(geometry, cmap='prism', alpha=0.4, ec=[0.5, 0.5, 0.5, 0.1], dual_phase=False, silent=False, asp_arr=None)[source]

Plot triangularized convex hulls of grains, based on vertices, i.e. connection points of 4 up to 8 grains or the end points of triple or quadruple lines between grains.

Parameters

silent geometry : dict

Dictionary with ‘vertices’ (node numbers) and ‘simplices’ (triangles)

cmapcolor map, optional

Color map for triangles. The default is ‘prism’.

alphafloat, optional

Transparency of plotted objects. The default is 0.4.

eccolor, optional

Color of edges. The default is [0.5,0.5,0.5,0.1].

dual_phasebool, optional

Whether to plot red/green contrast for dual phase microstructure or colored grains

Returns

None.

kanapy.plotting.plot_stats_dict(sdict, title=None, save_files=False)[source]

Plot statistical data on semi axes of effective ellipsoids in RVE as histogram.

Parameters

save_files sdict title

Returns

kanapy.plotting.plot_voxels_3D(data, Ngr=1, sliced=False, dual_phase=False, mask=None, cmap='prism', alpha=1.0, silent=False, clist=None, asp_arr=None)[source]

Plot voxels in microstructure, each grain with a different color. Sliced indicates whether one eighth of the box should be removed to see internal structure. With alpha, the transparency of the grains can be adjusted.

Parameters

dataint array

Grain number or phase number associated to each voxel

Ngrint, optional

Number of grains. The default is 1.

slicedBoolean, optional

Indicate if one eighth of box is invisible. The default is False.

dual_phaseBoolean, optional

Indicate if microstructure is dual phase. The default is False.

maskbool array, optional

Mask for voxels to be displayed. The default is None, in which case all voxels will be plotted (except sliced).

cmapcolor map, optional

Color map for voxels. The default is ‘prism’.

alphafloat, optional

Adjust transparency of voxels in alpha channel of color map. The default is 1.0.

showbool

Indicate whether to show the plot or to return the axis for further use

clist(Ngr, 3)-ndarray

List of RGB colors for each grain

Returns

axmatplotlib.axes

Axes handle for 3D subplot

Module smoothing

class kanapy.smoothingGB.Node(iden, px, py, pz)[source]

Bases: object

compute_acc(fx, fy, fz, mass)[source]
get_Oripos()[source]
get_pos()[source]
get_vel()[source]
update_pos(dt)[source]
update_vel(dt)[source]
kanapy.smoothingGB.initalizeSystem(nodes_v, grain_facesDict)[source]
kanapy.smoothingGB.readGrainFaces(nodes_v, elmtDict, elmtSetDict)[source]
kanapy.smoothingGB.relaxSystem(allNodes, anchDict, dt, N, k, c, RVE_xmin, RVE_xmax, RVE_ymin, RVE_ymax, RVE_zmin, RVE_zmax)[source]
kanapy.smoothingGB.smoothingRoutine(nodes_v, elmtDict, elmtSetDict)[source]