kawin.precipitation.parameters
On this page
- kawin.kawin.precipitation.parameters.Volume
- kawin.kawin.precipitation.parameters.ElasticFactors
- kawin.kawin.precipitation.parameters.LebedevNodes
- kawin.kawin.precipitation.parameters.ShapeFactors
- kawin.kawin.precipitation.NucleationRate
- kawin.kawin.precipitation.parameters.Nucleation
- NucleationDescriptionBase (ABC)
- BulkDescription(NucleationDescriptionBase)
- DislocationDescription(BulkDescription)
- GrainBoundaryDescription(NucleationDescriptionBase)
- GrainEdgeDescription(NucleationDescriptionBase)
- GrainCornerDescription(NucleationDescriptionBase)
- NucleationBarrierParameters
- NucleationSiteParameters
- kawin.kawin.precipitation.parameters.EffectiveDiffusion
kawin.kawin.precipitation.parameters.Volume
VolumeParameter
class VolumeParameter(self, value=None, volumeType=None, atomsPerCell=None):
def VolumeParameter.setVolume(self, value, volumeType, atomsPerCell):
Function to set lattice parameter, atomic volume and molar volume
Parameters
----------
value : float
Value for volume parameters (lattice parameter, atomic (unit cell) volume or molar volume)
valueType : VolumeParameter or str
States what volume term that value is
atomsPerCell : int
Number of atoms in the unit cell
kawin.kawin.precipitation.parameters.ElasticFactors
def convert2To4rankTensor(c2):
Converts 2nd rank elastic tensor to 4th rank
Parameters
----------
c : ndarray
2nd rank elastic tensor
Returns
-------
c4 : ndarray
4th rank elastic tensor
def convert4To2rankTensor(c4):
Converts 4th rank elastic tensor to 4nd rank
Parameters
----------
c4 : ndarray
4th rank elastic tensor
Returns
-------
c : ndarray
2nd rank elastic tensor
def invert4rankTensor(c4):
Inverts 4th rank tensor to give stiffness tensor
This is done by converting to 2nd rank, inverting, then converting back to 4th rank
def convertVecTo2rankTensor(v):
Converts strain/stress vector to 2nd rank tensor
Parameters
----------
v : 1d array
Strain/stress vector
Returns
-------
e : ndarray
2nd rank elastic tensor
def convert2rankToVec(c):
Converts 2nd rank tensor to vector
Parameter
---------
c : ndarray
3x3 tensor
Returns
-------
v : 1darray
Strain/stress vector
def rotateRank2Tensor(rot, tensor):
Rotates a 2nd rank tensor
T_ij = r_il * r_jk * T_lk
Parameters
----------
tensor : ndarray
2nd rank tensor to rotate (3x3 array)
def rotateRank4Tensor(rot, tensor):
Rotates a 4th rank tensor
T_ijkl = r_im * r_jn * r_ok * r_lp * T_mnop
Parameters
----------
tensor : ndarray
4th rank tensor to rotate (3x3x3x3 array)
def elasticConstantToC(c11, c12, c44):
Creates elastic tensor from c11, c12 and c44 constants assuming isotropic system
Parameters
----------
c11 : float
c12 : float
c44 : float
def moduliToC(E = None, nu = None, G = None, lam = None, K = None, M = None):
Set elastic tensor by defining at least two of the moduli
Everything will be converted to E, nu and G
If more than two parameters are defined, then the following priority will be used:
E - Youngs modulus
nu - Poisson's ratio
G - shear modulus
lam - Lame's first parameter
K - bulk modulus
M - P-wave modulus
NOTE: There's gotta be a better way to implement the conversions
EffectiveDiffusionFunctions
StrainEnergyDescriptionBase (ABC)
class StrainEnergyDescriptionBase (self, params: StrainEnergyParameters = None):
def StrainEnergyDescriptionBase .computeStrainEnergy(self, radius):
ConstantEnergyDescription(StrainEnergyDescriptionBase)
def ConstantEnergyDescription.computeStrainEnergy(self, radius):
SphericalEnergyDescription(StrainEnergyDescriptionBase)
def SphericalEnergyDescription._Khachaturyan(self, I1, I2, radius):
Khachaturyan's approximation for strain energy of spherical and cuboidal precipitates
def SphericalEnergyDescription.computeStrainEnergy(self, radius):
Strain energy of perfect sphere (cubic factor = 1)
CuboidalEnergyDescription(SphericalEnergyDescription)
def CuboidalEnergyDescription.computeStrainEnergy(self, radius):
Strain energy of perfect cube (cubic factor = sqrt(2))
EllipsoidalEnergyDescription(StrainEnergyDescriptionBase)
class EllipsoidalEnergyDescription(self):
def EllipsoidalEnergyDescription.setOhmInverseFunction(self, method = ‘quick’):
Sets method to invert the ohm term in calculating eshelby's tensor
Parameters
----------
method : str
'numpy' - uses np.linalg.inv, which can be slower for batch, but runs through
multiple checks for whether values are real/complex or if inverse exists
'quick' - quick inverse using Cramer's rule assuming that values are real and inverse exists
- this is the recommended method since strain energy will be computed many times
def EllipsoidalEnergyDescription.setIntegrationIntervals(self, phiInt, thetaInt, assumeSymmetric=True):
Number of intervals to split domain along phi and theta for integration
Parameters
----------
phiInt : int
Number of intervals to divide along phi
thetaInt : int
Number of intervals to divide along theta
assumeSymmetric : bool (optional)
If True (default), will only integrate Eshelby's tensor on a single quadrant and
multiply the results by 8
def EllipsoidalEnergyDescription.setLebedevIntegration(self, order = ‘high’):
Creates Lebedev quadrature points and nodes for integrating Eshebly's tensor
This is preferred over discretizing phi and theta
Parameters
----------
order : str
'low' - uses quadrature order or 53 (974 points)
'mid' - uses quadrature order or 83 (2354 points)
'high' (default) - uses quadrature order or 131 (5810 points)
def EllipsoidalEnergyDescription._n(self, phi, theta):
Unit normal vector of sphere
Parameters
----------
phi - azimuthal angle
theta - polar angle
def EllipsoidalEnergyDescription._beta(self, a, b, c, phi, theta):
Distance from center to surface of ellipsoid
Parameters
----------
a, b, c - radii of ellipsoid along x,y,z axes
phi, theta - azimuthal, polar angle
def EllipsoidalEnergyDescription._OhmGeneral(self, n, c4):
Ohm term for general system
Ohm_ij = inverse(C_iklj * n_k * n_l)
def EllipsoidalEnergyDescription._ohm_quickInverse(self, m):
Hard coded inverse of m which is of shape (3,3,n)
numpy inv is more optimized for larger matrices, but can be slower for small
matrices such as a 2x2 or 3x3. We can take advantage of 3x3 matrices having a computable
inverse to make it faster
NOTE: this only works since we know that m has a shape of (3,3,n) and is only composed of real numbers
This function can probably be a bit more efficient, but quick
profiling on sphInt gives around a 35x speedup compared to doing
np.transpose(np.linalg.inv(np.transpose(m, (2,0,1))), (1,2,0))
where the slowdown was in np.linalg.inv
For matrix:
| a b c |
| d e f |
| g h i |
Inverse is defined as:
| ei-fh fg-di dh-eg | | A B C |
| ch-bi ai-cg bg-ah | / det -> | D E F | / det
| bf-ce cd-af ae-bd | | G H I |
Where det = aA + bB + cC
def EllipsoidalEnergyDescription._ohm_npinv(self, m):
Inverts ohm term using np.linalg.inv
numpy inverse function takes in an array of shape (m,n,n) and inverts each nxn matrix
So we have to transpose m from (3,3,n) -> (n,3,3), then invert, then transpose (n,3,3) ->(3,3,n)
def EllipsoidalEnergyDescription.sphInt(self, radius, c4):
Faster version of calculating the Dijkl, which avoids
using a for loop to iterate across phi and theta
def EllipsoidalEnergyDescription.Dijkl(self, radius, c4):
Dijkl term for Eshelby's theory
def EllipsoidalEnergyDescription.Sijmn(self, D):
S_ijmn = -0.5 * C_lkmn * (D_iklj + D_jkli)
def EllipsoidalEnergyDescription._multiply(self, a, b):
Multiplies 2 tensors
4th x 2nd -> c_ij = a_ijkl * b_kl
4th x 4th -> c_ijkl = a_ijmn * b_mnkl
def EllipsoidalEnergyDescription._strainEnergy(self, stress, strain, V):
u = -0.5 * V * sigma_ij * strain_ij
def EllipsoidalEnergyDescription.strainEnergyEllipsoidWithStress(self, radius):
Strain energy of ellipsoidal particle with applied stress
def EllipsoidalEnergyDescription.strainEnergyEllipsoid(self, radius):
Strain energy of ellipsoidal particle
def EllipsoidalEnergyDescription.strainEnergyEllipsoid2ndRank(self, radius):
Alternative method of strain energy on ellipsoidal particle using 2nd rank tensors
def EllipsoidalEnergyDescription.strainEnergyBohm(self, radius):
Strain energy of particle for when matrix and precipitate phases have different elastic tensors
def EllipsoidalEnergyDescription.strainEnergyBohm2ndRank(self, radius):
Strain energy of particle for when matrix and precipitate phases have different elastic tensors using 2nd rank tensors
def EllipsoidalEnergyDescription.computeStrainEnergy(self, radius):
def EllipsoidalEnergyDescription._OhmCubic(self, n, c4):
Ohm term for cubic crystal symmetry
def EllipsoidalEnergyDescription._D(self, n):
Needed for the Ohm term with cubic crystal symmetry
def EllipsoidalEnergyDescription._xi(self):
Needed for the Ohm term with cubic crystal symmetry
StrainEnergy
class StrainEnergy(self, shape=‘constant’):
Defines class for calculating strain energy of a precipitate
Ellipsoidal precipitates will use the Eshelby's tensor
Spherical and Cuboidal precipitates will use the Khachaturyan's approximation
def StrainEnergy.description(self):
def StrainEnergy.description(self, value):
def StrainEnergy.unrotated_cMatrix_4th(self):
def StrainEnergy.unrotated_cMatrix_4th(self, value):
def StrainEnergy.unrotated_cPrec_4th(self):
def StrainEnergy.unrotated_cPrec_4th(self, value):
def StrainEnergy.setShape(self, shape):
def StrainEnergy.setSpherical(self):
Assumes spherical geometry for strain energy calculation
Uses Khachaturyan's approximation
def StrainEnergy.setCuboidal(self):
Assumes cuboidal geometry for strain energy calculation
Uses Khachaturyan's approximation
def StrainEnergy.setEllipsoidal(self):
Assumes ellipsoidal geometry for strain energy calculation
Uses Eshelby's tensor
def StrainEnergy.setConstantElasticEnergy(self, energy):
If elastic strain energy is known to be a constant, this can be use to greatly
simplify calculations
Parameters
----------
energy - strain energy in J/m3
def StrainEnergy.setElasticTensor(self, tensor):
Sets elastic tensor of matrix using 2nd rank tensor
Parameters
----------
tensor : 6x6 array
2nd rank elastic tensor
def StrainEnergy.setElasticConstants(self, c11, c12, c44):
Sets elastic tensor of matrix by elastic constants, assuming cubic symmetry
Parameters
----------
c11 : float
Modulus for compression
c11 = E(1-nu) / (1+nu)(1-2nu)
c12 : float
Modulus for dilation (accounts for compression and Poisson's ratio)
c12 = E nu / (1+nu)(1-2nu)
c44 : float
Modulus for shear
c44 = (c11-c12)/2
def StrainEnergy.setModuli(self, E = None, nu = None, G = None, lam = None, K = None, M = None):
Sets elastic tensor of matrix by 2 moduli
Parameters (only 2 has to be defined)
----------
E : float
Elastic modulus
nu : float
Poisson's ratio
G : float
Shear modulus
lam : float
Lame's first parameter
K : float
Bulk modulus
M : float
P-wave modulus
def StrainEnergy.setElasticTensorPrecipitate(self, tensor):
Sets elastic tensor of precipitate using 2nd rank tensor
Parameters
----------
tensor : 6x6 array
2nd rank elastic tensor
def StrainEnergy.setElasticConsantsPrecipitate(self, c11, c12, c44):
Sets elastic tensor of precipitate by elastic constants, assuming cubic symmetry
Parameters
----------
c11 : float
Modulus for compression
c11 = E(1-nu) / (1+nu)(1-2nu)
c12 : float
Modulus for dilation (accounts for compression and Poisson's ratio)
c12 = E nu / (1+nu)(1-2nu)
c44 : float
Modulus for shear
c44 = (c11-c12)/2
def StrainEnergy.setModuliPrecipitate(self, E = None, nu = None, G = None, lam = None, K = None, M = None):
Sets elastic tensor of precipitate by 2 moduli
Parameters (only 2 has to be defined)
----------
E : float
Elastic modulus
nu : float
Poisson's ratio
G : float
Shear modulus
lam : float
Lame's first parameter
K : float
Bulk modulus
M : float
P-wave modulus
def StrainEnergy.setRotationMatrix(self, rot):
Sets rotation matrix to be applied to the matrix
This is for cases where the axes of the precipitate does not align with the axes of the matrix
(e.g., the long/short axes of the precipitate is not parallel to the <100> directions of the matrix)
Parameters
----------
rot : matrix
3x3 rotation matrix
def StrainEnergy.setRotationPrecipitate(self, rot):
def StrainEnergy.setEigenstrain(self, strain):
Sets eigenstrain of precipitate
Parameters
----------
strain : float, array or matrix
float - assume strain is the same along all 3 axis
array - each index corresponds to strain in a given axis
matrix - full 2nd rank strain tensor
NOTE: when using in conjunction with ShapeFactors, the axis are order specific
For needle-like precipitates, x1, x2 and x3 correspond to (short axis, short axis, long axis)
For plate-like precipitates, x1, x2 and x3 correspond to (long axis, long axis, short axis)
def StrainEnergy.setAppliedStress(self, stress):
Sets applied stress tensor
Axes of stress tensor should be the same as the matrix
Parameters
----------
stress : float, array or matrix
float - assume stress is the same along all 3 axis
array - each index corresponds to stress in a given axis
matrix - full 2nd rank stress tensor
NOTE: The applied stress is in reference to the coordinate axes of the precipitates
Thus, this is only valid if the following conditions are met:
The matrix phase has a single orientation (either as a single crystal or highly textured)
Precipitates will form only in a single orientation with respect to the matrix
NOTE: this is only available in EllipsoidalEnergyDescription.strainEnergyEllipsoidWithStress and is currently not used in the KWNModel
TODO: It will be nice to expand on this.
For polycrystalline matrix and randomly oriented precipitates, it should be possible
to average the strain energy contributions over all matrix/precipitate orientations
def StrainEnergy._computeAppliedStrain(self, cM2, stress):
Calculates applied strain tensor from applied stress tensor and elastic tensor
def StrainEnergy.update(self):
def StrainEnergy.compute(self, r):
def StrainEnergy.setAspectRatioResolution(self, resolution = 0.01, cachedRange = 5):
Sets resolution to which equilibrium aspect ratios are calculated
Equilibrium aspect ratios are found by calculated strain energy for a range of aspect ratios
and finding the aspect ratio giving the minimum energy (strain + interfacial energy)
If aspect ratio does not vary much in a given system, then the default parameters may lead to poor
prediction of the aspect ratios
Parameters
----------
resolution : float (optional)
Minimum distance between aspect ratios when calculating cache (default at 0.01)
cachedRange : float (optional)
Range of aspect ratio to calculate strain energy when updated cache (default at 5)
def StrainEnergy.setInterfacialEnergyMethod(self, method = ‘thermo’):
Sets method for calculating interfacial energy as a function of aspect ratio
Parameters
----------
method : str
'eqradius' - interfacial energy is determined using the equivalent spherical radius
'thermo' - interfacial energy is determined using dG/dSA (default)
def StrainEnergy.updateCache(self, normR):
Update cached calculations
def StrainEnergy.clearCache(self):
Clear cached calculations
def StrainEnergy.eqAR_byGR(self, Rsph, gamma, shpFactor : ShapeFactor, a=1.001, b=100):
Equilibrium aspect ratio using golden ratio search
Parameters
----------
Rsph : float or array
Equivalent spherical radius
gamma : float
Interfacial energy
shpFactor : ShapeFactor object
a, b : float
Min and max bounds
Default is 1.001 and 100
def StrainEnergy._GRsearch(self, Rsph, gamma, interfacial, normR, a, b):
Golden ratio search for a single radius
def StrainEnergy.eqAR_bySearch(self, Rsph, gamma, shpFactor : ShapeFactor):
Equilibrium aspect ratio by cached search
Parameters
----------
Rsph : float or array
Equivalent spherical radius
gamma : float
Interfacial energy
shpFactor : ShapeFactor object
def StrainEnergy._cachedSearch(self, Rsph, gamma, interfacial, normR):
Cached search for a single radius
kawin.kawin.precipitation.parameters.LebedevNodes
def loadPoints(order):
Load lebedov nodes of specific order
Parameters
----------
order : int
Can be 53, 83 or 131
Determines the order of polynomial at which an exact integral solution can be found
Returns
-------
phi coordinates, theta coordinates, weights
kawin.kawin.precipitation.parameters.ShapeFactors
ShapeDescriptionBase (ABC)
class ShapeDescriptionBase (self):
Defines functions to describe a precipitate shape
Must implement (as a function of aspect ratio)
_eqRadius - equivalent spherical radius (radius of sphere that gives the same volume)
_normalRadii - radius along normals to give a spherical volume of 1
_thermoFactor - factor that modifies the Gibbs-Thomson contribution
_kineticFactor - factor that modifies the growth rate
def ShapeDescriptionBase ._processAspectRatio(self, ar):
def ShapeDescriptionBase .eccentricity(self, ar):
Eccentricity given the aspect ratio (for needle and plate shaped precipitates)
def ShapeDescriptionBase .normalRadii(self, ar):
Radius along the 3 axis to give a spherical volume of 1
Parameters
----------
R : float or array
Equivalent spherical radius
Returns
-------
3 length array for each axis if R is scalar
n x 3 array if R is an array
def ShapeDescriptionBase .eqRadiusFactor(self, ar):
Equivalent spherical radius factor for generic precipitates
Parameters
----------
R : float or array
Equivalent spherical radius
Returns
-------
Eq. radius factor with same shape as R
def ShapeDescriptionBase .kineticFactor(self, ar):
Kinetic factor for generic precipitates
Parameters
----------
R : float or array
Equivalent spherical radius
Returns
-------
Kinetic factor with same shape as R
def ShapeDescriptionBase .thermoFactor(self, ar):
Thermodynamic factor for generic precipitates
Parameters
----------
R : float or array
Equivalent spherical radius
Returns
-------
Thermodynamic factor with same shape as R
def ShapeDescriptionBase ._eqRadius(self, ar):
def ShapeDescriptionBase ._normalRadii(self, ar):
def ShapeDescriptionBase ._kineticFactor(self, ar):
def ShapeDescriptionBase ._thermoFactor(self, ar):
SphereDescription(ShapeDescriptionBase)
def SphereDescription._eqRadius(self, ar):
Equivalent radius for a sphere (returns 1)
def SphereDescription._normalRadii(self, ar):
Returns radius along the 3-axis for a volume of 1
def SphereDescription._kineticFactor(self, ar):
Kinetic factor for a sphere (returns 1)
def SphereDescription._thermoFactor(self, ar):
Thermodynamic factor for a sphere (returns 1)
NeedleDescription(ShapeDescriptionBase)
def NeedleDescription._eqRadius(self, ar):
Equivalent radius for needle shaped precipitate
def NeedleDescription._normalRadii(self, ar):
Returns radius along the 3-axis for a volume of 1
def NeedleDescription._kineticFactor(self, ar):
Kinetic factor for needle shaped precipitate
def NeedleDescription._thermoFactor(self, ar):
Thermodynamic factor for needle shaped precipitate
PlateDescription(ShapeDescriptionBase)
def PlateDescription._eqRadius(self, ar):
Equivalent radius for plate shaped precipitate
def PlateDescription._normalRadii(self, ar):
Returns radius along the 3-axis for a volume of 1
def PlateDescription._kineticFactor(self, ar):
Kinetic factor for plate shaped precipitate
def PlateDescription._thermoFactor(self, ar):
Thermodynamic factor for plate shaped precipitate
CuboidalDescription(ShapeDescriptionBase)
class CuboidalDescription(self):
def CuboidalDescription._eqRadius(self, ar):
Equivalent radius for cuboidal shaped precipitate
def CuboidalDescription._normalRadii(self, ar):
Returns radius along the 3-axis for a volume of 1
def CuboidalDescription._kineticFactor(self, ar):
Kinetic factor for cuboidal shaped precipitate
def CuboidalDescription._thermoFactor(self, ar):
Thermodynamic factor for cuboidal shaped precipitate
ShapeFactor
class ShapeFactor(self, precipitateShape=‘sphere’, ar=1):
Defines functions for shape factors of needle, plate and cuboidal precipitates
Shape factors involve the following:
Aspect Ratio - ratio of long axis to short axis
Eccentricity - characteristic value of the precipitate shape that depends on the aspect ratio
Equivalent radius factor - correction factor to the radius to give the radius of a sphere with equivalent volume
Kinetic factor - correction factor for the growth rate of the particle
Thermodynamic factor / Area factor - correction factor for the critical driving force/radius for nucleation
NOTE:
normalRadaii, eqRadiusFactor, kineticFactor and thermoFactor are functions of radius
ShapeFactor.function(radius) -> factor
The equivalent functions in the description are functions of aspect ratio
ShapeFactor.description(aspect ratio) -> factor
def ShapeFactor.description(self):
def ShapeFactor.description(self, value):
def ShapeFactor.setPrecipitateShape(self, precipitateShape, ar = 1):
General shape setting function
Defaults to spherical
Parameters
----------
precipitateShape: str
'sphere', 'needle', 'plate' or 'cubic'
ar: float or callable
Aspect ratio. If callable, then it must be a function of equivalent spherical radius
def ShapeFactor.setSpherical(self, ar = 1):
Sets factors for spherical precipitates
def ShapeFactor.setNeedleShape(self, ar = 1):
Factors for needle shaped precipitates
def ShapeFactor.setPlateShape(self, ar = 1):
Factors for plate shaped precipitates
def ShapeFactor.setCuboidalShape(self, ar = 1):
Factors for cuboidal shaped precipitates
def ShapeFactor.setAspectRatio(self, ar):
Sets aspect ratio as a function of R (equivalent spherical radius)
Parameters
----------
ar : float or function
Aspect ratio
If function, must take in radius and return aspect ratio
def ShapeFactor._scalarAspectRatioEquation(self, R):
Aspect ratio as a function of radius if the aspect ratio is set as a scalar
def ShapeFactor.normalRadii(self, R):
Radius along the 3 axis to give a spherical volume of 1
Parameters
----------
R : float or array
Equivalent spherical radius
Returns
-------
3 length array for each axis if R is scalar
n x 3 array if R is an array
def ShapeFactor.eqRadiusFactor(self, R):
Equivalent spherical radius factor for generic precipitates
Parameters
----------
R : float or array
Equivalent spherical radius
Returns
-------
Eq. radius factor with same shape as R
def ShapeFactor.kineticFactor(self, R):
Kinetic factor for generic precipitates
Parameters
----------
R : float or array
Equivalent spherical radius
Returns
-------
Kinetic factor with same shape as R
def ShapeFactor.thermoFactor(self, R):
Thermodynamic factor for generic precipitates
Parameters
----------
R : float or array
Equivalent spherical radius
Returns
-------
Thermodynamic factor with same shape as R
def ShapeFactor._findRcritScalar(self, RcritSphere, Rmax):
Critical radius given a scalar aspect ratio
def ShapeFactor._findRcrit(self, RcritSphere, Rmax):
Critical radius given aspect ratio as a function of radius
Found by bisection method
kawin.kawin.precipitation.NucleationRate
def volumetricDrivingForce(therm: GeneralThermodynamics, x, T, precipitate: PrecipitateParameters, aspectRatio = 1, removeCache = False):
Computes volumetric driving force (chemical DG / VM - strain energy)
Strain energy will always reduce the driving force since the precipitate will add strain to the matrix
In the case where the matrix is prestrained and the precipitate will relax the matrix, then the strain
energy is negative
def nucleationBarrier(volumeDrivingForce, precipitate : PrecipitateParameters, aspectRatio = 1):
Critical Gibbs free energy and radius at the nucleation barrier
For bulk and dislocation nucleation
Rcrit = 2*f*gamma / dG (where f is the thermodynamic correction factor from the precipitate shape)
Gcrit = (4*pi/3)*gamma*Rcrit^2
For grain boundary, edge or corner nucleation, critical G and R is computed to in NucleationBarrierParameters to account for grain boundary energy
def zeldovich(T, Rcrit, precipitate : PrecipitateParameters):
Zeldovich factor
Z = sqrt(3*fv/4*pi) * Vm * sqrt(gamma/kB*T) / (2*pi*Nv*Rcrit^2)
def betaBinary1(therm : BinaryThermodynamics, x, T, Rcrit, matrix : MatrixParameters, precipitate : PrecipitateParameters, removeCache = False):
Impingement rate for binary systems
beta = fa * Rcrit^2 * x * D / a**4
def betaBinary2(therm : BinaryThermodynamics, x, T, Rcrit, matrix : MatrixParameters, precipitate : PrecipitateParameters, xEqAlpha = None, xEqBeta = None, removeCache = False):
Impingement rate for binary systems similar to how multicomponent systems are computed
beta = fa * Rcrit^2 * Dterm / a**4
D = [(xB - xA)^2 / (xA*D) + (xB - xA)^2 / ((1-xA)*D)]^-1
def betaMulti(therm : MulticomponentThermodynamics, x, T, Rcrit, matrix : MatrixParameters, precipitate : PrecipitateParameters, removeCache = False, searchDir = None):
Impingement rate for multicomponent systems
beta = fa * Rcrit^2 * Dterm / a**4
Dterm = [sum_i((xB_i - xA_i)^2 / (xA_i*D_i))]^-1
def incubationTime(beta, Z, matrix : MatrixParameters):
Returns incubation time (tau) and a time offset (this is to be compatible with the nonisothermal calculation)
tau = 1 / (theta * beta * Z^2)
def incubationTimeNonIsothermal(Z, currBeta, currTime, currTemp, betas, times, temperatures, matrix : MatrixParameters):
Note: beta, times and temperature is a subslide from startIndex:n+1
Start index is the incubationOffset term
Solve tau for int_0^tau (beta(t)dt) = 1 / (theta*Z(tau)^2)
def nucleationRate(Z, beta, Gcrit, T, tau, time = np.inf):
Nucleation rate
d#/dt = Z * beta * exp(-Gcrit / kB*T) * exp(-tau / t)
Units are 1/t
def nucleationRadius(T, Rcrit, precipitateParameters: PrecipitateParameters):
Adds 1/2 * sqrt(kb T / pi gamma) to critical radius to ensure they grow when growth rates are calculated
def computeSteadyStateNucleation(therm : GeneralThermodynamics, x, T, precipitate: PrecipitateParameters, matrix : MatrixParameters, betaFunc = None, aspectRatio = 1, removeCache = False):
kawin.kawin.precipitation.parameters.Nucleation
NucleationDescriptionBase (ABC)
def NucleationDescriptionBase ._createArrays(self, gbk):
Creates arrays for grain boundary factors, this is done here to allow limits
to be placed on the grain boundary ratio
def NucleationDescriptionBase ._formatArray(self, values, indices, setInvalidToNan = True):
By default, setInvalidToNan will be true, will will set all invalid indices to nan
This is for the user API so that it's easy to tell when a gb energy / interfacial energy
ratio is invalid
When setting factors internally, this will be set to false so it's easier to compare
when a value is invalid (since np.nan == np.nan will return False)
def NucleationDescriptionBase .gbRatio(self, gbEnergy, gamma):
Grain boundary to interfacial energy ratio
def NucleationDescriptionBase .gbRemoval(self, gbk = None, setInvalidToNan = True):
def NucleationDescriptionBase ._gbRemoval(self, gbk):
def NucleationDescriptionBase .areaFactor(self, gbk = None, setInvalidToNan = True):
def NucleationDescriptionBase ._areaFactor(self, gbk):
def NucleationDescriptionBase .volumeFactor(self, gbk = None, setInvalidToNan = True):
def NucleationDescriptionBase ._volumeFactor(self, gbk):
def NucleationDescriptionBase .areaRemoval(self, gbk = None, setInvalidToNan = True):
def NucleationDescriptionBase ._areaRemoval(self, gbk):
def NucleationDescriptionBase .isGrainBoundaryNucleation(self):
BulkDescription(NucleationDescriptionBase)
def BulkDescription._gbRemoval(self, gbk):
def BulkDescription._areaFactor(self, gbk):
def BulkDescription._volumeFactor(self, gbk):
def BulkDescription._areaRemoval(self, gbk):
def BulkDescription.isGrainBoundaryNucleation(self):
DislocationDescription(BulkDescription)
GrainBoundaryDescription(NucleationDescriptionBase)
def GrainBoundaryDescription._gbRemoval(self, gbk):
def GrainBoundaryDescription._areaFactor(self, gbk):
def GrainBoundaryDescription._volumeFactor(self, gbk):
GrainEdgeDescription(NucleationDescriptionBase)
def GrainEdgeDescription.alpha(self, gbk):
def GrainEdgeDescription.beta(self, gbk):
def GrainEdgeDescription._gbRemoval(self, gbk):
def GrainEdgeDescription._areaFactor(self, gbk):
def GrainEdgeDescription._volumeFactor(self, gbk):
GrainCornerDescription(NucleationDescriptionBase)
def GrainCornerDescription.K(self, gbk):
def GrainCornerDescription.phi(self, gbk):
def GrainCornerDescription.delta(self, gbk):
def GrainCornerDescription._gbRemoval(self, gbk):
def GrainCornerDescription._areaFactor(self, gbk):
def GrainCornerDescription._volumeFactor(self, gbk):
NucleationBarrierParameters
class NucleationBarrierParameters(self, site=‘dislocations’, gamma = None, gbEnergy=0.3):
Defines nucleation factors for bulk, dislocation,
surface, edge and corner nucleation
This includes: surface area, volume, GB area removal, Rcrit and Gcrit
Attributes
----------
gbRemoval - factor to multiply by r**2 to get area of grain boundary removed
areaFactor - factor to multiply by r**2 to get surface area of precipitate
volumeFactor - factor to multiply by r**3 to get volume of precipitate
areaRemoval - factor to multiply by r to get radius of eliminated grain boundary area
def NucleationBarrierParameters.description(self):
def NucleationBarrierParameters.description(self, value):
def NucleationBarrierParameters.gamma(self):
def NucleationBarrierParameters.gamma(self, value):
def NucleationBarrierParameters.gbEnergy(self):
def NucleationBarrierParameters.gbEnergy(self, value):
def NucleationBarrierParameters._resetFactors(self):
def NucleationBarrierParameters.setNucleationType(self, site):
Sets nucleation site type
Parameters
----------
site - int
Index for site type based off list on top of this file
def NucleationBarrierParameters._validateInputs(self):
def NucleationBarrierParameters._validateGBk(self):
def NucleationBarrierParameters.GBk(self):
def NucleationBarrierParameters.areaFactor(self):
def NucleationBarrierParameters.volumeFactor(self):
def NucleationBarrierParameters.gbRemoval(self):
def NucleationBarrierParameters.areaRemoval(self):
def NucleationBarrierParameters.isGrainBoundaryNucleation(self):
def NucleationBarrierParameters.Rcrit(self, dG):
Critical radius for nucleation
This is only done for GB nucleation
Bulk and dislocation nucleation needs to account for shape factors
R* = 2 * (b * \gamma_{\alpha \alpha} - a * \gamma_{\alpha \beta}) / (3 * c * dG)
Parameters
----------
dG - float
Volumetric driving force
def NucleationBarrierParameters.Gcrit(self, dG, Rcrit):
Critical driving force for nucleation
G* = 4/27 * (b * \gamma_{\alpha \alpha} - a * \gamma_{\alpha \beta})^3 / (c * dG)^2
or in terms of R*
G* = -R*^2 * (b * \gamma_{\alpha \alpha} - a * \gamma_{\alpha \beta} + c * r* * dG)
This is calculated in terms of R* since there is a check in the nucleation calculation that will
increase R* to minimum radius if R* is too small (very low nucleation barriers)
Parameters
----------
dG - float
Volumetric driving force
Rcrit - float
Critical radius for nucleation
NucleationSiteParameters
class NucleationSiteParameters(self, grainSize = 100, aspectRatio = 1, dislocationDensity = 5e12):
def NucleationSiteParameters.grainSize(self):
def NucleationSiteParameters.grainSize(self, value):
def NucleationSiteParameters.grainAspectRatio(self):
def NucleationSiteParameters.grainAspectRatio(self, value):
def NucleationSiteParameters.dislocationDensity(self):
def NucleationSiteParameters.dislocationDensity(self, value):
def NucleationSiteParameters.setNucleationDensity(self, grainSize = 100, aspectRatio = 1, dislocationDensity = 5e12, bulkN0 = None):
Sets grain size and dislocation density which determines the available nucleation sites
Parameters
----------
grainSize : float (optional)
Average grain size in microns (default at 100um if this function is not called)
aspectRatio : float (optional)
Aspect ratio of grains (default at 1)
dislocationDensity : float (optional)
Dislocation density (m/m3) (default at 5e12)
bulkN0 : float (optional)
This allows for the use to override the nucleation site density for bulk precipitation
By default (None), this is calculated by the number of lattice sites containing a solute atom
However, for calibration purposes, it may be better to set the nucleation site density manually
def NucleationSiteParameters.setGrainSize(self, grainSize = 100, aspectRatio = 1):
def NucleationSiteParameters.setDislocationDensity(self, dislocationDensity):
def NucleationSiteParameters.setBulkDensity(self, bulkN0):
def NucleationSiteParameters.setBulkDensityFromComposition(self, x0):
def NucleationSiteParameters.bulkN0(self):
def NucleationSiteParameters.bulkN0(self, value):
def NucleationSiteParameters.dislocationN0(self):
def NucleationSiteParameters.GBareaN0(self):
def NucleationSiteParameters.GBedgeN0(self):
def NucleationSiteParameters.GBcornerN0(self):
def NucleationSiteParameters.dislocationSites(self, VmAlpha):
def NucleationSiteParameters.grainBoundaryDensity(self, grainSize, grainAspectRatio):
def NucleationSiteParameters.grainBoundarySites(self, grainSize, grainAspectRatio, VmAlpha):
def NucleationSiteParameters.grainEdgeDensity(self, grainSize, grainAspectRatio):
def NucleationSiteParameters.grainEdgeSites(self, grainSize, grainAspectRatio, VmAlpha):
def NucleationSiteParameters.grainCornerDensity(self, grainSize, grainAspectRatio):
def NucleationSiteParameters.grainCornerSites(self, grainSize, grainAspectRatio, VmAlpha):
def NucleationSiteParameters._validateVolume(self, term):
kawin.kawin.precipitation.parameters.EffectiveDiffusion
EffectiveDiffusionFunctions
class EffectiveDiffusionFunctions(self, isEnabled=True):
Stores variables and functions to calculate effective diffusion distance
def EffectiveDiffusionFunctions.setupInterpolation(self, n = 250, lmin = -5, lmax = 1):
Sets up interpolation points for effective diffusion distance
Parameters
----------
n : int
Number of points to interpolate between, two endpoints will be added
at supersaturation of 0 and 1
lmin : float
Minimum lambda value in log scale
lmax : float
Maximum lambda value in log scale
def EffectiveDiffusionFunctions.effectiveDiffusionDistance(self, supersaturation):
Effective diffusion distance given supersaturation
This gives a constant to be multiplied by the particle radius
The effective diffusion distance is given by
eps = Q/2*lambda
Where Q is the super saturation and lambda is given by
Q = 2*lambda^2 - 2*lambda^3 * sqrt(pi) * exp(lambda^2) * (1 - erf(lambda))
Parameters
----------
supersaturation : float or array of floats
(x - x_P) / (x_M - x_P) where x is matrix composition, x_M is parent composition and x_P is precipitate composition
Returns
-------
Effective diffusion distance (eps) to be used as (eps*R) in the growth rate equation
def EffectiveDiffusionFunctions.lambdaLow(self, supersaturation):
Lambda when Q approaches 0
This is done to prevent precision errors when multiplying exp*(1-erf)
def EffectiveDiffusionFunctions.lambdaHigh(self, supersaturation):
Lambda when Q approaches 1
This is done to prevent precision errors when multiplying exp*(1-erf)
def EffectiveDiffusionFunctions.effectiveDiffusionDistanceApprox(self, supersaturation):
Effective diffusion distance given supersaturation
This gives a constant to be multiplied by the particle radius
When Q approaches 1, this equation oscillates due to floating point errors
Therefore, interpolation is performed between the two limits (Q -> 0 and Q -> 1)
This may be removed in later versions
Parameters
----------
supersaturation : float or array of floats
(x - x_P) / (x_M - x_P) where x is matrix composition, x_M is parent composition and x_P is precipitate composition
Returns
-------
Effective diffusion distance (eps) to be used as (eps*R) in the growth rate equation
def EffectiveDiffusionFunctions.noDiffusionDistance(self, supersaturation):
If neglecting the effective diffusion distance, return 1
Parameters
----------
supersaturation : float or array of floats
(x - x_P) / (x_M - x_P) where x is matrix composition, x_M is parent composition and x_P is precipitate composition
Returns
-------
1 or array of 1s
def EffectiveDiffusionFunctions.call(self, supersaturation):