kawin.precipitation.coupling

kawin.kawin.precipitation.coupling.Strength

def ignore_numpy_warnings(func):

The strength contributions may be undefined at rss=0 or ls=0, so this
is just an easy way to ignore any divide by 0 warnings

TODO: a better way would be to know the function limits for rss=0 and ls=0,
      then apply the function limits at those values instead of computing, but
      ignoring the warnings here is a lot easier for now

def GrainGrowthModel.wrapper(*args, **kwargs):

DislocationParameters

class DislocationParameters(self, G, b, nu = 1/3, ri = None, theta = 90, psi = 120):

Parameters for dislocation line tension, used for all precipitate strength models

Parameters
----------
G : float
    Shear modulus of matrix (Pa)
b : float
    Burgers vector (meters)
nu : float
    Poisson ratio
ri : float (optional)
    Dislocation core radius (meters)
    If None, ri will be set to Burgers vector
r0 : float (optional)
    Closest distance between parallel dislocations
    For shearable precipitates, r0 is average distance between particles on slip plane
    For non-shearable precipitates, r0 is average particle diameter on slip plane
    If None, r0 will be set such that ln(r0/ri) = 2*pi
theta : float (optional)
    Dislocation characteristic, 90 for edge, 0 for screw (default is 90)
psi : float (optional)
    Dislocation bending angle (default is 120)

def DislocationParameters.tension(self, r0, theta = None):

StrengthContributionBase(ABC)

class StrengthContributionBase(self, phase = ‘all’):

def StrengthContributionBase.r0Weak(self, Ls, dislocations: DislocationParameters):

def StrengthContributionBase.r0Strong(self, Ls, dislocations: DislocationParameters):

def StrengthContributionBase.computeCRSS(self, r, Ls, dislocations: DislocationParameters):

OrowanContribution(StrengthContributionBase)

Orowan strengthening contribution

def OrowanContribution.computeCRSS(self, r, Ls, dislocations: DislocationParameters):

CoherencyContribution(StrengthContributionBase)

class CoherencyContribution(self, eps, phase=‘all’):

Parameters for coherency effect

Parameters
----------
eps : float
    Lattice misfit strain

def CoherencyContribution.latticeMisfit(delta, dislocations: DislocationParameters):

Strain (eps) from lattice misfit (delta)

def CoherencyContribution.computeCRSS(self, r, Ls, dislocations: DislocationParameters):

ModulusContribution(StrengthContributionBase)

class ModulusContribution(self, Gp, w1=0.05, w2=0.85, phase=‘all’):

Parameters for modulus effect

Parameters
----------
Gp : float
    Shear modulus of precipitate
w1 : float (optional)
    First factor for Nembach model taking value between 0.0175 and 0.0722
    Default at 0.05
w2 : float (optional)
    Second factor for Nembach model taking value of 0.81 +/- 0.09
    Default at 0.85

def ModulusContribution.f(self, r, G, b):

Term for modulus effect

def ModulusContribution.computeCRSS(self, r, Ls, dislocations: DislocationParameters):

Modulus effect for mixed dislocation on weak and shearable particles

APBContribution(StrengthContributionBase)

class APBContribution(self, yAPB, s=2, beta=1, V=2.8, phase=‘all’):

Parameters for anti-phase boundary effect for ordered precipitates in a disordered matrix

Parameters
----------
yAPB : float
    Anti-phase boundary energy
s : int (optional)
    Number of leading + trailing dislocations to repair anti-phase boundary
    Default at 2
beta : float (optional)
    Factor representating dislocation shape (between 0 and 1)
    Default at 1
V : float (optional)
    Correction factor accounting for extra dislocations and uncertainties
    Default at 2.8

def APBContribution.computeCRSS(self, r, Ls, dislocations: DislocationParameters):

SFEContribution(StrengthContributionBase)

class SFEContribution(self, ySFM, ySFP, bp = None, phase=‘all’):

Parameters for stacking fault energy effect

Parameters
----------
ySFM : float
    Stacking fault energy of matrix
ySFP : float
    Stacking fault energy of precipitate
bp : float (optional)
    Burgers vector in precipitate
    If None, will be set to burgers vector in matrix

def SFEContribution.K(self, dislocations: DislocationParameters):

def SFEContribution.Weff(self, dislocations: DislocationParameters):

Effective stacking fault width for mixed dislocations

def SFEContribution.f(self, r, dislocations: DislocationParameters):

Stacking fault term

def SFEContribution.computeCRSS(self, r, Ls, dislocations: DislocationParameters):

InterfacialContribution(StrengthContributionBase)

class InterfacialContribution(self, gamma, phase=‘all’):

Parameters for interfacial effect

Parameters
----------
gamma : float
    Interfacial energy of matrix/precipitate surface

def InterfacialContribution.computeCRSS(self, r, Ls, dislocations: DislocationParameters):

SolidSolutionStrength

class SolidSolutionStrength(self, weights = {}, exp = {}):

Model for solid solution strengthening

Parameters
----------
weights : dictionary
    Dictionary mapping element (str) to weight (float). Weights are in units of (Pa)
exp : float
    Exponential factor for weights

def SolidSolutionStrength.compute(self, composition, elements):

def computeCRSS(rss, Ls, contributions: list[StrengthContributionBase], dislocations: DislocationParameters, phase):

Computes critical resolved shear stress from precipitate contributions over rss, Ls
Weak and strong values will be mapped to the specific contribution
Orowan is a single contribution that's added by default

def combineCRSS(weak, strong, owo, exp = 1.8, returnContributions = False):

Sums the different critical resolved shear stress contributions together

def SolidSolutionStrength.sumArray(contribution):

def SolidSolutionStrength.processArray(contribution):

StrengthModel

class StrengthModel(self, phases: list[PrecipitateParameters | str], contributions:StrengthContributionBase|list[StrengthContributionBase], dislocations:DislocationParameters, ssModel:SolidSolutionStrength=None, sigma0:float=0):

Strength model for coupling with precipitate model

Parameters
----------
phases: list[PrecipitateParameters|str]
    List of phases to model
    If coupling with a precipitate model, must be the same phases
contributions: StrengthContributionBase | list[StrengthContributionBase]
dislocations: DislocationParameters
ssModel: SolidSolutionStrength (optional)
    If None, then solid solution strength model will output 0
sigma0: float (optional)
    Base strength of all
    Default is 0

def StrengthModel.save(self, filename):

Saves strength model data

Note, this only saves solid solution strength, the rss and ls terms
Parameters should be free so user can load model and evaluate different parameters

def StrengthModel.load(self, filename):

def StrengthModel.updateCoupledModel(self, model: PrecipitateModel):

Computes rss and Ls from current state of the PrecipitateModel

rss - mean projected radius of particles
Ls - mean surface to surface distance between particles

r1 = first ordered moment of particle size distribution
r2 = second ordered moment of particle size distribution
rss = sqrt(2/3) * r2 / r1
ls = sqrt(ln(3)/(2*pi*r1) + (2*rss)^2) - 2*rss

Parameters
----------
model : PrecpitateModel

def StrengthModel.computePrecipitateStrength(self, model: PrecipitateModel):

Computes yield strength from precipitate contributions

Precipitate contributions compute CRSS, so we multiply by Taylor factor (M) to get YS

def StrengthModel.totalStrength(self, model: PrecipitateModel, returnContributions = False):

Compute total strength from precipitate, solid solution and base contributions
All strength are assumed to be yield strength

def _get_strength_units(units = ‘Pa’):

def _plotContributionOverX(x, r, Ls, contribution: StrengthContributionBase, dislocations: DislocationParameters, strengthUnits=‘MPa’, ax=None, *args, **kwargs):

Plots single strengthening contribution

Parameters
----------
x: np.array
    x coordinates
r: np.array
    Mean projected particle radius
Ls: np.array
    Surface-surface particle distance
contribution: StrengthContributionBase
dislocations: DislocationParameters
strengthUnits: str (optional)
    Pa, kPa, MPa or GPa
    Default is MPa
ax: matplotlib axis

def plotContribution(r, Ls, contribution: StrengthContributionBase, dislocations: DislocationParameters, strengthUnits=‘MPa’, ax=None, *args, **kwargs):

Plots single strengthening contribution
This plots the critical resolved shear stress vs. mean projected particle radius

Parameters
----------
r: np.array
    Mean projected particle radius
Ls: np.array
    Surface-surface particle distance
contribution: StrengthContributionBase
dislocations: DislocationParameters
strengthUnits: str (optional)
    Pa, kPa, MPa or GPa
    Default is MPa
ax: matplotlib axis

def plotContributionOverTime(model: PrecipitateModel, strengthModel: StrengthModel, contribution: StrengthContributionBase, timeUnits=’s', strengthUnits=‘MPa’, ax=None, *args, **kwargs):

Plots single strengthening contribution from precipitate evolution in a precipitation model
This plots the critical resolved shear stress vs. time

Parameters
----------
model: PrecipitateModel
strengthModel: StrengthModel
contribution: StrengthContributionBase
timeUnits: str (optional)
    s, min, or h
    Default is s
strengthUnits: str (optional)
    Pa, kPa, MPa or GPa
    Default is MPa
ax: matplotlib axis

def _plotPrecipitateStrengthOverX(x, r, Ls, strengthModel: StrengthModel, phase, plotContributions=False, strengthUnits=‘MPa’, ax=None, *args, **kwargs):

Plots total strengthening contribution of single precipitate

Parameters
----------
x: np.array
    x coordinates
r: np.array
    Mean projected particle radius
Ls: np.array
    Surface-surface particle distance
strengthModel: StrengthModel
phase: str (optional)
    If None, then first phase of strength model is used
plotContributions: bool (optional)
    If True, will plot the weak, strong, orowan and minimum contribution
    If False, will only plot the minimum contribution
strengthUnits: str (optional)
    Pa, kPa, MPa or GPa
    Default is MPa
ax: matplotlib axis

def plotPrecipitateStrength(r, Ls, strengthModel: StrengthModel, phase=None, plotContributions=False, strengthUnits=‘MPa’, ax=None, *args, **kwargs):

Plots total strengthening contribution of single precipitate
This plots the critical resolved shear stress vs. mean projected particle radius

Parameters
----------
r: np.array
    Mean projected particle radius
Ls: np.array
    Surface-surface particle distance
strengthModel: StrengthModel
phase: str (optional)
    If None, then first phase of strength model is used
plotContributions: bool (optional)
    If True, will plot the weak, strong, orowan and minimum contribution
    If False, will only plot the minimum contribution
strengthUnits: str (optional)
    Pa, kPa, MPa or GPa
    Default is MPa
ax: matplotlib axis

def plotPrecipitateStrengthOverTime(model: PrecipitateModel, strengthModel: StrengthModel, phase=None, plotContributions=False, timeUnits=’s', strengthUnits=‘MPa’, ax=None, *args, **kwargs):

Plots total strengthening contribution of single precipitate
This plots the critical resolved shear stress vs. time

Parameters
----------
model: PrecipitateModel
strengthModel: StrengthModel
phase: str (optional)
    If None, then first phase of strength model is used
plotContributions: bool (optional)
    If True, will plot the weak, strong, orowan and minimum contribution
    If False, will only plot the minimum contribution
timeUnits: str (optional)
    s, min, or h
    Default is s
strengthUnits: str (optional)
    Pa, kPa, MPa or GPa
    Default is MPa
ax: matplotlib axis

def plotAlloyStrength(model: PrecipitateModel, strengthModel: StrengthModel, plotContributions=False, timeUnits=’s', strengthUnits=‘MPa’, ax=None, *args, **kwargs):

Plots total predicted strength of an alloy during precipitate evolution
This plots the yield strength vs. time

Parameters
----------
model: PrecipitateModel
strengthModel: StrengthModel
plotContributions: bool (optional)
    If True, will plot the base, solid solution, precipitate and total contribution
    If False, will only plot the total contribution
timeUnits: str (optional)
    s, min, or h
    Default is s
strengthUnits: str (optional)
    Pa, kPa, MPa or GPa
    Default is MPa
ax: matplotlib axis

kawin.kawin.precipitation.coupling.GrainGrowth

GrainGrowthModel(GenericModel)

class GrainGrowthModel(self, gbe=0.5, M=1e-14, alpha=1, m={}, K={}, iterator=rk4Iterator, pbmKwargs={‘cMin’: 1e-10, ‘cMax’: 1e-8}):

Model for grain growth that can be coupled with the KWN model to account for Zener pinning

Following implentation described in
K. W, J. Jeppsson and P. Mason, J. Phase Equilib. Diffus. 43 (2022) 866-875

Parameters
----------
gbe: float (optional)
    Grain boundary energy
    Default = 0.5
M: float (optional)
    Grain boundary mobility
    Default = 1e-14
alpha: float (optional)
    Correction factor (for when fitting data to the model)
    Default = 1
m: dict[str, float] (optional)
    Factor related to spatial distribution of precipitates
    Maps a phase name to the factor
    Default = 1 for all phases
K: dict[str, float] (optional)
    Factor related to spatial distribution of precipitates
    Maps a phase name to the factor
    Default = 4/3 for all phases
iterator: iterator function (optional)
    Iterator function for model
    Default is rk4Iterator
pbmKwargs: dict[str, int | float] (optional)
    Arguments to define population balance model
    Default is {cMin = 1e-10, cMax = 1e-8}

def GrainGrowthModel.loadDistribution(self, data):

Creates a particle size distribution from a set of data

Parameters
----------
data : array of floats
    Array of data to be inserted into PSD

def GrainGrowthModel.loadDistributionFunction(self, function):

Creates a particle size distribution from a function

Parameters
----------
function : function
    Takes in R and returns density

def GrainGrowthModel.reset(self):

Resets model with initially loaded grain size distribution

def GrainGrowthModel.Rcr(self, x):

Critical radius, grains larger than Rcr will growth while smaller grains will shrink

Critical radius is defined so that the volume will be constant when applying the growth rate

Parameters
----------
x : np.array
    Grain size distribution corresponding to GrainGrowthModel.pbm.PSDbounds

def GrainGrowthModel.Rm(self, x):

Mean radius

Parameters
----------
x : np.array
    Grain size distribution corresponding to GrainGrowthModel.pbm.PSDbounds

def GrainGrowthModel.grainGrowth(self, x):

Grain growth model
dRi/dt = alpha * M * gbe * (1/Rcr - 1/Ri)

Parameters
----------
x : np.array
    Grain size distribution corresponding to GrainGrowthModel.pbm.PSDbounds

def GrainGrowthModel.normalize(self):

Normalize PSD to have a third moment of 1

Ideally, this isn't needed since the grain growth model accounts for constant volume
    But numerical errors will lead to small changes in volume over time

def GrainGrowthModel.constrainedGrowth(self, growthRate, z = 0):

Constrain growth rate due to zener pinning

The growth rate given the zener radius is defined by:
    dR/dt = alpha * M * gbe * ((1/Rcr - 1/Ri) +/- 1/Rz)
    Where 1/Rz is added if (1/Rcr - 1/Ri) + 1/Rz < 0 (inhibits grain dissolution)
    And   1/Rz is subtracted in (1/Rcr - 1/Ri) - 1/Rz) > 0 (inhibits grain growth)
    And   dR/dt is 0 for Ri between these two limits

Note: Rather than Rz (zener radius), we use z here which represents the drag force
    But these are related by z = 1/Rz

Parameters
----------
growthRate : array
    Growth rate for grain sizes
z : float (optional)
    Zener radius, default is 0, which will not change the growth rate

def GrainGrowthModel.getCurrentX(self):

Returns current time and grain size distribution

def GrainGrowthModel.getdXdt(self, t, x):

Returns dn_i/dt for the grain size distribution

Steps:
    1. Get grain growth rate and corrected it with zener drag force
    2. Get dn_i/dt from the PBM given the Eulerian implementation

def GrainGrowthModel.correctdXdt(self, dt, x, dXdt):

Corrects dn_i/dt with the new time step

def GrainGrowthModel.getDt(self, dXdt):

Calculated a suitable dt with the growth rate and new time step
We'll limit the max time step to the remaining time for solving

def GrainGrowthModel.postProcess(self, time, x):

Sets grain size distribution to x and record time and average grain size

Steps:
    1. Set grain size distribution
    2. Adjust PSD size classes
    3. Remove grains below the dissolution threshold
    4. Normalize grain size distribution to 1 (should be a tiny correction factor due to step 3)
    5. Record time and average grain size

def GrainGrowthModel.printHeader(self):

Header string before solving

def GrainGrowthModel.printStatus(self, iteration, modelTime, simTimeElapsed):

Status string that prints every n iteration

def GrainGrowthModel.computeZenerRadius(self, model: PrecipitateModel, x=None):

Gets zener radius/drag force from PrecipitateModel and PSD defined by x

Drag force is defined as z_j = f_j^m_j / (K_j * avgR_j)
    Where f_j is volume fraction for phase j
    And   avgR_j is average radius for phase j
The total drag force is the sum of z_j over all the phases

Parameters
----------
model : PrecpitateModel
x : list[np.array]
    List of particle size distributions in model

def GrainGrowthModel.updateCoupledModel(self, model):

Computes zener radius/drag force from the PrecipitateModel,
Then solves the grain growth model with the time step of the PrecipitateModel

Parameters
----------
model : PrecpitateModel

def _plot_grain_growth_generic(model: GrainGrowthModel, func, ax=None, *args, **kwargs):

def plotGrainPSD(model: GrainGrowthModel, ax=None, *args, **kwargs):

Plots grain size distribution

Parameters
----------
model: GrainGrowthModel
ax : matplotlib axes

def plotGrainPDF(model: GrainGrowthModel, ax=None, *args, **kwargs):

Plots grain size distribution density

Parameters
----------
model: GrainGrowthModel
ax : matplotlib axes

def plotGrainCDF(model: GrainGrowthModel, ax=None, *args, **kwargs):

Plots cumulative grain size distribution

Parameters
----------
model: GrainGrowthModel
ax : matplotlib axes

def plotRadiusvsTime(model: GrainGrowthModel, timeUnits=’s', ax=None, *args, **kwargs):

Plots average grain size vs time

Parameters
----------
model: GrainGrowthModel
timeUnits: str
    's', 'min' or 'hr'
ax : matplotlib axes