Models

The physics models used in the analyses

Isentrope

class F_UNCLE.Models.Isentrope.Isentrope(name=u'Isentrope', *args, **kwargs)

Abstract class for an isentrope

The equation of state for an isentropic expansion of high explosive is modeled by this class. The experiments for which this model is used occur at such short timescales that the process can be considered adiabatic

Units

Isentropes are assumed to be in CGS units

Diagram

_images/eos.png

The assumed shape of the equation of state isentrope

Options

Name Type Def Min Max Units Description
spline_N (int) 50 7 None ‘’ Number of knots in the EOS spline
spline_min (float) 0.1 0.0 None ‘cm**3/g’ Minimum value of volume modeled by EOS
spline_max (float) 1.0 0.0 None ‘cm**3/g’ Maximum value of volume modeled by EOS
spline_end (float) 4 0 None ‘’ Number of zero nodes at end of spline
__init__(name=u'Isentrope', *args, **kwargs)
Parameters:
  • *args – Variable length argument list.
  • **kwargs – Arbitrary keyword arguments.
Keyword Arguments:
 

name (str) – Name if the isentrope Def ‘Isentrope’

Returns:

None

_get_cj_point(vol_0, pres_0=0.0, debug=False)

Find CJ conditions using two nested line searches.

The CJ point is the location on the EOS where a Rayleigh line originating at the pre-detonation volume and pressure is tangent to the equation of state isentrope.

This method uses two nested line searches implemented by the scipy.optimize.brentq() algorithm to locate the velocity corresponding to this tangent Rayleigh line

Parameters:

vol_0 (float) – The specific volume of the equation of state before the shock arrives

Keyword Arguments:
 
  • pres_0 (float) – The pressure of the reactants
  • debug (bool) – Flag to print debug information
Returns:

Length 3 elements are:

  1. (float): The detonation velocity
  2. (float): The specific volume at the CJ point
  3. (float): The pressure at the CJ point
  4. (function): A function defining the Rayleigh line which passes through the CJ point

Return type:

(tuple)

get_constraints(scale=False)

Returns the G and h matrices corresponding to the model

Parameters:model (GaussianModel) – The physics model subject to physical constraints
Returns:():
Return type:()

Method

Calculate constraint matrix G and vector h. The constraint enforced by cvxopt.solvers.qp is

\[G*x \leq h\]

Equivalent to \(max(G*x - h) \leq 0\)

Since

\[c_{f_{new}} = c_f+x,\]
\[G(c_f+x) \leq 0\]

is the same as

\[G*x \leq -G*c_f,\]

and \(h = -G*c_f\)

Here are the constraints for \(p(v)\):

p’’ positive for all v p’ negative for v_max p positive for v_max

For cubic splines between knots, f’’ is constant and f’ is affine. Consequently, \(f''rho + 2f'\) is affine between knots and it is sufficient to check eq:star at the knots.

plot(axes=None, figure=None, linestyles=[u'-k'], labels=[u'Isentrope'], vrange=None, log=False, draw_rayl=False, *args, **kwargs)

Plots the EOS

Overloads the F_UNCLE.Utils.Struc.Struc.plot() method to plot the EOS over the range of volumes.

Parameters:
  • axes (plt.Axes) – The axes on which to plot the figure, if None, creates a new figure object on which to plot.
  • figure (plt.Figure) – The figure on which to plot ignored
  • linstyles (list) – Strings for the linestyles 0. ‘-k’, The Isentrope
  • labels (list) – Strings for the plot labels 0. ‘Isentrope’
  • vrange (tuple) – Specific volume range to plot
  • log (bool) – Flag to plot on semilogy
  • draw_rayl (bool) – Flag to draw rayleigh line
Returns:

A reference to the figure containing the plot

Return type:

(plt.Figure)

shape()

Overloaded class to get isentrope DOF’s

Overloads F_UNCLE.Utils.PhysModel.PhysModel.shape()

Returns:(n,1) where n is the number of dofs
Return type:(tuple)

BumpEOS

class F_UNCLE.Models.Isentrope.EOSBump(name=u'Bump EOS', *args, **kwargs)

Model of an ideal isentrope with Gaussian bumps

This is treated as the true EOS

Options

Name Type Def Min Max Units Description
const_C (float) 2.56e9 0.0 None ‘Pa’ ‘Constant p = C/v**3’
bumps (list)

[(0.4 0.1 0.25)

(0.5 0.1 -0.3)]

None None ‘’ ‘Gaussian bumps to the EOS’
__call__(vol)

Solve the EOS

Calculates the pressure for a given volume, is called the same way as the EOS model but uses underlying equation rather than the spline

Parameters:vol (np.ndarray) – Specific volume
Returns:Pressure
Return type:pr(np.ndarray)
__init__(name=u'Bump EOS', *args, **kwargs)

Instantiate the bump EOS

Parameters:
  • *args – Variable length argument list.
  • **kwargs – Arbitrary keyword arguments.
Keyword Arguments:
 

name (str) – Name if the isentrope Def ‘Bump EOS’

derivative(order=1)

Returns the nth order derivative

Keyword Arguments:
 order (int) – The order of the derivative. Def 1
Returns:Function object yielded first derivative of pressure w.r.t volume
Return type:d1_fun(function)

EOSModel

class F_UNCLE.Models.Isentrope.EOSModel(p_fun, name=u'Equation of State Spline', *args, **kwargs)

Spline based EOS model

Multiply inherited structure from both Isentrope and Spline

Options

Name Type Def Min Max Units Description
spline_sigma float 5e-3 0.0 None ‘??’ ‘Multiplicative uncertainty of the prior (1/2%)
precondition bool False None None ‘’ Precondition flag
__init__(p_fun, name=u'Equation of State Spline', *args, **kwargs)

Instantiates the object

Parameters:
  • p_fun (function) – A function defining the prior EOS
  • *args – Variable length argument list.
  • **kwargs – Arbitrary keyword arguments.
Keyword Arguments:
 

name (str) – Name of the isentrope Def ‘Equation of State Spline’

_get_knot_spacing()

Returns a list of knot locations based on the spline parameters

If the option spacing is ‘lin’, uses linear spacing
‘log’, uses log spacing

Places ‘spline_N’ knots between ‘spline_min’ and ‘spline_max’

_on_str()

Added information on the EOS

_on_update_dof(model)

An extra method to perform special post-processing tasks when the DOF has been updated

Parameters:model (EOSModel) – The new physics model
Returns:The post-processed model
Return type:(EOSModel)
get_dof(*args, **kwargs)

Returns the spline coefficients as the model degrees of freedom

Returns:The degrees of freedom of the model
Return type:(np.ndarray)
get_gram()

Creates the gramm matrix

Parameters:None
Returns:(np.ndarray) A nxn gram matrix where n is the numer of knots
get_kernel(alpha, beta)
get_scaling()

Returns a scaling matrix to make the dofs of the same scale

The scaling matrix is a diagonal matrix with off diagonal terms zero the terms along the diagonal are the prior DOFs times the variance in the DOF values.

Returns:A nxn matrix where n is the number of model DOFs.
Return type:(np.ndarray)
get_sigma()

Returns the co-variance matrix of the spline

Returns:Co-variance matrix for the EOS shape is (nxn) where n is the dof of the model
Return type:(np.ndarray)
plot_basis(axes=None, fig=None, labels=[], linstyles=[])

Plots the basis function and their first and second derivatives

Parameters:
  • fig (plt.Figure) – A valid figure object on which to plot
  • axes (plt.Axes) – A valid axes, Ignored
  • labels (list) – The labels, Ignored
  • linestyles (list) – The linestyles, Ignored
Returns:

The figure

Return type:

(plt.Figure)

plot_diff(isentropes, axes=None, figure=None, linestyles=[u'-k', u'--k', u'-.k'], labels=None, vrange=None)

Plots the difference between the current EOS and other Isentropes

Plots the difference vs the prior

Parameters:
  • isentropes (list) – A list of isentrope objects to compare
  • axes (plt.Axes) – The Axes on which to plot
  • figure (plt.Figure) – The figure object Ignored
  • linestyles (list) – A list of styles to plot
  • labels (list) – A list of labels for the Isentropes
  • vrange (tuple) – the specific volume range to plot
Returns:

(plt.Axes)

static plot_sens_matrix(sens_matrix, simid, models, mkey, fig=None)

Prints the sensitivity matrix

Parameters:
  • sens_matrix (dict) – Dictionary of sensitivity matrices
  • simid (str) – Key for simulation
  • models (OrderedDict) – Ordered dictionary of models
  • mkey (str) – Key in models corresponding to the EOSModel
Keyword Args
fig(plt.Figure): A valid matplotlib figure on which to plot.
If None, creates a new figure
Returns:The figure
Return type:(plt.Figure)
update_dof(c_in, *args, **kwargs)

Sets the spline coefficients

Parameters:c_in (Iterable) – The knot positions of the spline
Returns:A copy of self with the new dofs
Return type:(EOSModel)

PTW

class F_UNCLE.Models.Ptw.Ptw(name=u'Ptw flow stress', *args, **kwargs)

PTW Flow stress model

Usage

  1. Instantiate a Ptw object with desired options
  2. optional set the options as desired
  3. Call the Ptw object with a valid temperature, strain rate and material specification
  4. Call the object again, material must be specified each time
__call__(temp, strain_rate, material, **overrides)

Solves for the yield stress and flow stress at the given condition

Args:
temp(float): Temperature, in degrees Kelvin

strain_rate(float): Strain rate, in sec**-1 material(str): Key for material type

Keyword Args:
**overrides(dict): Passed as a chain of keyword arguments. These
arguments override any material property
Return:
flow_stress(float): Flow stress in ??Pa?? yield_stress(float): Yield stress in ??Pa??
__init__(name=u'Ptw flow stress', *args, **kwargs)

Instantiates the structure

Parameters:
  • name (str) – Name of the structure
  • *args – Variable length argument list.
  • **kwargs – Arbitrary keyword arguments.
Returns:

None

__pre__(temp, strain_rate, material, **overrides)

Performs data conditioning

Performs the following tasks

  1. Ensure the temperature input is in a valid range
  2. Ensures the strain rate input is in a valid range
  3. Ensure a valid material has been specified
  4. Populates the materials database
  5. Applies the overrides
  6. Returns the normalizing factors
Parameters:
  • temp (float) – Temperature, in degrees Kelvin
  • strain_rate (float) – Strain rate, in sec**-1
  • material (str) – Key for material type
Keyword Arguments:
 

**overrides (dict) – Passed as a chain of keyword arguments. These arguments override any material property

Returns
(float): The shear modulus (float): The characteristic temperature. Temp normalized by melt temp (float): The characteristic timescale. Time for a vibration to pass through one atom
get_mat_data(material)

Gets the material data corresponding to the given material

Parameters:material (str) – Key for material type
get_melt_temperature(temp, strain_rate)

Gets the melt temperature for the operating conditions

Parameters:
  • temp (float) – Temperature, in degrees Kelvin
  • strain_rate (float) – Strain rate, in sec**-1
Returns:

The melt temperature in Kelvin

Return type:

(float)

get_shear_modulus(temp, strain_rate, t_melt)

Gets the shear modulus for the operating conditions

Parameters:
  • temp (float) – Temperature, in degrees Kelvin
  • strain_rate (float) – Strain rate, in sec**-1
  • t_melt (float) – Melt temperature of the material in Kelvin
get_stress(strain, strain_rate, temp, material, **overrides)

Returns the stress in the material material

Parameters:
  • strain (float or np.array) – The strain in the material
  • strain_rate (float) – The train rate in the material
  • temp (float) – The temperature of the material
  • mat (str) – A valid material specifier
Keyword Arguments:
 
  • materials properties can be passed as kwargs to (valid) –
  • default values (override) –
Returns:

The stress in the material

Return type:

(float or np.array)

get_stress_strain(temp, strain_rate, material, min_strain=0.05, max_strain=0.7, **overrides)

Returns the stress strain relationship for the material

Random Varialbe

class F_UNCLE.Models.RandomVariable.RandomVariable(prior_mean, name='Random Variable', *args, **kwargs)

A model of a normally distributed random variable of known variance

__call__(*args, **kwargs)

Returns the mean of the distribution

Parameters:None
Returns:The mean
Return type:(float)
__init__(prior_mean, name='Random Variable', *args, **kwargs)
_on_str()
get_dof()
get_sigma(*args, **kwargs)
plot(axes=None, labels=['Coeff'], linestyles=['ok'], *args, **kwargs)
shape()
update_dof(dof_in)

Random Process

class F_UNCLE.Models.RandomProcess.RandomProcess(prior_mean, prior_var, name='Random Process', *args, **kwargs)

A model of a normally distributed random process of known variance

__call__(*args, **kwargs)

Returns the mean of the distribution

Parameters:None
Returns:The mean
Return type:(float)
__init__(prior_mean, prior_var, name='Random Process', *args, **kwargs)
_on_str()
get_dof()
get_sigma(*args, **kwargs)
plot(axes=None, labels=['Coeff'], linestyles=['ok'], *args, **kwargs)
shape()
update_dof(dof_in, var_in=None)