constitutive : The Constitutive Library


In [1]:
%matplotlib inline

import numpy as np
import matplotlib.pyplot as plt
from simmit import smartplus as sim
import os

L_iso

Provides the elastic stiffness tensor for an isotropic material. The two first arguments are a couple of elastic properties. The third argument specifies which couple has been provided and the nature and order of coefficients. Exhaustive list of possible third argument : ‘Enu’,’nuE,’Kmu’,’muK’, ‘KG’, ‘GK’, ‘lambdamu’, ‘mulambda’, ‘lambdaG’, ‘Glambda’. Return a numpy ndarray. Example :


In [2]:
E = 70000.0
nu = 0.3
L = sim.L_iso(E,nu,"Enu")
print np.array_str(L, precision=4, suppress_small=True)

d = sim.check_symetries(L)
print(d['umat_type'])
print(d['props'])

x = sim.L_iso_props(L)
print(x)


[[ 94230.7692  40384.6154  40384.6154      0.          0.          0.    ]
 [ 40384.6154  94230.7692  40384.6154      0.          0.          0.    ]
 [ 40384.6154  40384.6154  94230.7692      0.          0.          0.    ]
 [     0.          0.          0.      26923.0769      0.          0.    ]
 [     0.          0.          0.          0.      26923.0769      0.    ]
 [     0.          0.          0.          0.          0.      26923.0769]]
ELISO
[  7.00000000e+04   3.00000000e-01]
[  7.00000000e+04   3.00000000e-01]

M_iso

Provides the elastic compliance tensor for an isotropic material. The two first arguments are a couple of elastic properties. The third argument specify which couple has been provided and the nature and order of coefficients. Exhaustive list of possible third argument : ‘Enu’,’nuE,’Kmu’,’muK’, ‘KG’, ‘GK’, ‘lambdamu’, ‘mulambda’, ‘lambdaG’, ‘Glambda’.


In [3]:
E = 70000.0
nu = 0.3
M = sim.M_iso(E,nu,"Enu")
print np.array_str(M, precision=2)

L_inv = np.linalg.inv(M)
d = sim.check_symetries(L_inv)
print(d['umat_type'])
print(d['props'])

x = sim.M_iso_props(M)
print(x)


[[  1.43e-05  -4.29e-06  -4.29e-06   0.00e+00   0.00e+00   0.00e+00]
 [ -4.29e-06   1.43e-05  -4.29e-06   0.00e+00   0.00e+00   0.00e+00]
 [ -4.29e-06  -4.29e-06   1.43e-05   0.00e+00   0.00e+00   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   3.71e-05   0.00e+00   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   0.00e+00   3.71e-05   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   0.00e+00   0.00e+00   3.71e-05]]
ELISO
[  7.00000000e+04   3.00000000e-01]
[  7.00000000e+04   3.00000000e-01]

L_cubic

Provides the elastic stiffness tensor for a cubic material. Arguments are the stiffness coefficients C11, C12 and C44, or the elastic constants E, nu, G Exhaustive list of possible third argument : ‘Cii’,’EnuG, the by-default argument is 'Cii'


In [4]:
E = 70000.0
nu = 0.3
G = 23000.0
L = sim.L_cubic(E,nu,G,"EnuG")
print np.array_str(L, precision=2)

d = sim.check_symetries(L)
print(d['umat_type'])
print(d['props'])

x = sim.L_cubic_props(L)
print(x)


[[ 94230.77  40384.62  40384.62      0.        0.        0.  ]
 [ 40384.62  94230.77  40384.62      0.        0.        0.  ]
 [ 40384.62  40384.62  94230.77      0.        0.        0.  ]
 [     0.        0.        0.    23000.        0.        0.  ]
 [     0.        0.        0.        0.    23000.        0.  ]
 [     0.        0.        0.        0.        0.    23000.  ]]
ELCUB
[  7.00000000e+04   3.00000000e-01   2.30000000e+04]
[  7.00000000e+04   3.00000000e-01   2.30000000e+04]

M_cubic

Provides the elastic compliance tensor for a cubic material. Arguments are the stiffness coefficients C11, C12 and C44, or the elastic constants E, nu, G Exhaustive list of possible third argument : ‘Cii’,’EnuG, the by-default argument is 'Cii'


In [5]:
E = 70000.0
nu = 0.3
G = 23000.0
M = sim.M_cubic(E,nu,G,"EnuG")
print np.array_str(M, precision=2)

L = np.linalg.inv(M)
d = sim.check_symetries(L)
print(d['umat_type'])
print(d['props'])

x = sim.L_cubic_props(L)
print(x)


[[  1.43e-05  -4.29e-06  -4.29e-06   0.00e+00   0.00e+00   0.00e+00]
 [ -4.29e-06   1.43e-05  -4.29e-06   0.00e+00   0.00e+00   0.00e+00]
 [ -4.29e-06  -4.29e-06   1.43e-05   0.00e+00   0.00e+00   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   4.35e-05   0.00e+00   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   0.00e+00   4.35e-05   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   0.00e+00   0.00e+00   4.35e-05]]
ELCUB
[  7.00000000e+04   3.00000000e-01   2.30000000e+04]
[  7.00000000e+04   3.00000000e-01   2.30000000e+04]

L_isotrans

Provides the elastic stiffness tensor for an isotropic transverse material. Arguments are longitudinal Young modulus EL, transverse young modulus, Poisson’s ratio for loading along the longitudinal axis nuTL, Poisson’s ratio for loading along the transverse axis nuTT, shear modulus GLT and the axis of symmetry.


In [6]:
EL = 70000.0
ET = 20000.0
nuTL = 0.08
nuTT = 0.3
GLT = 12000.0
axis = 3
L = sim.L_isotrans(EL,ET,nuTL,nuTT,GLT,axis)
print np.array_str(L, precision=2)

d = sim.check_symetries(L)
print(d['umat_type'])
print(d['axis'])
print np.array_str(d['props'], precision=2)

x = sim.L_isotrans_props(L,axis)
print np.array_str(x, precision=2)


[[ 22954.82   7570.21   8547.01      0.        0.        0.  ]
 [  7570.21  22954.82   8547.01      0.        0.        0.  ]
 [  8547.01   8547.01  74786.32      0.        0.        0.  ]
 [     0.        0.        0.     7692.31      0.        0.  ]
 [     0.        0.        0.        0.    12000.        0.  ]
 [     0.        0.        0.        0.        0.    12000.  ]]
ELIST
3
[  7.00e+04   2.00e+04   8.00e-02   3.00e-01   1.20e+04]
[  7.00e+04   2.00e+04   8.00e-02   3.00e-01   1.20e+04]
bp::def("L_iso", L_iso);
bp::def("M_iso", M_iso);
bp::def("L_cubic", L_cubic);
bp::def("M_cubic", M_cubic);
bp::def("L_ortho", L_ortho);
bp::def("M_ortho", M_ortho);
bp::def("L_isotrans", L_isotrans);
bp::def("M_isotrans", M_isotrans);

bp::def("check_symetries", check_symetries);
bp::def("L_iso_props", L_iso_props);
bp::def("M_iso_props", M_iso_props);
bp::def("L_isotrans_props", L_isotrans_props);
bp::def("M_isotrans_props", M_isotrans_props);
bp::def("L_cubic_props", L_cubic_props);
bp::def("M_cubic_props", M_cubic_props);
bp::def("L_ortho_props", L_ortho_props);
bp::def("M_ortho_props", M_ortho_props);
bp::def("M_aniso_props", M_aniso_props);

M_isotrans

Provides the elastic compliance tensor for an isotropic transverse material. Arguments are longitudinal Young modulus EL, transverse young modulus, Poisson’s ratio for loading along the longitudinal axis nuTL, Poisson’s ratio for loading along the transverse axis nuTT, shear modulus GLT and the axis of symmetry.


In [7]:
EL = 70000.0
ET = 20000.0
nuTL = 0.08
nuTT = 0.3
GLT = 12000.0
axis = 3
M = sim.M_isotrans(EL,ET,nuTL,nuTT,GLT,axis)
print np.array_str(M, precision=2)

x = sim.M_isotrans_props(M,axis)
print np.array_str(x, precision=2)


[[  5.00e-05  -1.50e-05  -4.00e-06   0.00e+00   0.00e+00   0.00e+00]
 [ -1.50e-05   5.00e-05  -4.00e-06   0.00e+00   0.00e+00   0.00e+00]
 [ -4.00e-06  -4.00e-06   1.43e-05   0.00e+00   0.00e+00   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   1.30e-04   0.00e+00   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   0.00e+00   8.33e-05   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   0.00e+00   0.00e+00   8.33e-05]]
[  7.00e+04   2.00e+04   8.00e-02   3.00e-01   1.20e+04]

L_ortho

Provides the elastic stiffness tensor for an orthotropic material. Arguments are either (convention 'EnuG'):

  1. The Young modulus of axis 1 $E_1$,
  2. The Young modulus of axis 2 $E_2$,
  3. The Young modulus of axis 3 $E_3$,
  4. The Poisson ratio of direction 1 with respect to 2 $\nu_{12}$,
  5. The Poisson ratio of direction 1 with respect to 3 $\nu_{13}$,
  6. The Poisson ratio of direction 2 with respect to 3 $\nu_{13}$,
  7. The shear modulus of direction 1 with respect to 2 $G_{12}$,
  8. The shear modulus of direction 1 with respect to 3 $G_{13}$,
  9. The shear modulus of direction 2 with respect to 3 $G_{23}$,
  10. The axial coefficient of thermal expansion in direction 1 $\alpha_1$,
  11. The axial coefficient of thermal expansion in direction 1 $\alpha_2$,
  12. The axial coefficient of thermal expansion in direction 1 $\alpha_3$,

or the list of Cii (C11, C12, C13, C22, C23, C33, C44, C55, C66) (convention 'Cii')


In [8]:
E_1 = 4500.0
E_2 = 2300.0
E_3 = 2700.0
nu_12 = 0.06
nu_13 = 0.08
nu_23 = 0.3
G_12 = 2200.0
G_13 = 2100.0
G_23 = 2400.0

L = sim.L_ortho(E_1,E_2,E_3,nu_12,nu_13,nu_23,G_12,G_13,G_23,'EnuG')
print np.array_str(L, precision=2)

d = sim.check_symetries(L)
print(d['umat_type'])
print(d['axis'])
print np.array_str(d['props'], precision=2)

x = sim.L_ortho_props(L)
print np.array_str(x, precision=2)


[[ 4537.59   228.65   298.33     0.       0.       0.  ]
 [  228.65  2583.23   920.72     0.       0.       0.  ]
 [  298.33   920.72  3038.57     0.       0.       0.  ]
 [    0.       0.       0.    2200.       0.       0.  ]
 [    0.       0.       0.       0.    2100.       0.  ]
 [    0.       0.       0.       0.       0.    2400.  ]]
ELORT
0
[  4.50e+03   2.30e+03   2.70e+03   6.00e-02   8.00e-02   3.00e-01
   2.20e+03   2.10e+03   2.40e+03]
[  4.50e+03   2.30e+03   2.70e+03   6.00e-02   8.00e-02   3.00e-01
   2.20e+03   2.10e+03   2.40e+03]

M_ortho

Provides the elastic compliance tensor for an orthotropic material. Arguments are either (convention 'EnuG'):

  1. The Young modulus of axis 1 $E_1$,
  2. The Young modulus of axis 2 $E_2$,
  3. The Young modulus of axis 3 $E_3$,
  4. The Poisson ratio of direction 1 with respect to 2 $\nu_{12}$,
  5. The Poisson ratio of direction 1 with respect to 3 $\nu_{13}$,
  6. The Poisson ratio of direction 2 with respect to 3 $\nu_{13}$,
  7. The shear modulus of direction 1 with respect to 2 $G_{12}$,
  8. The shear modulus of direction 1 with respect to 3 $G_{13}$,
  9. The shear modulus of direction 2 with respect to 3 $G_{23}$,
  10. The axial coefficient of thermal expansion in direction 1 $\alpha_1$,
  11. The axial coefficient of thermal expansion in direction 1 $\alpha_2$,
  12. The axial coefficient of thermal expansion in direction 1 $\alpha_3$,

or the list of Cii (C11, C12, C13, C22, C23, C33, C44, C55, C66) (convention 'Cii')


In [9]:
E_1 = 4500.0
E_2 = 2300.0
E_3 = 2700.0
nu_12 = 0.06
nu_13 = 0.08
nu_23 = 0.3
G_12 = 2200.0
G_13 = 2100.0
G_23 = 2400.0

M = sim.M_ortho(E_1,E_2,E_3,nu_12,nu_13,nu_23,G_12,G_13,G_23,'EnuG')
print np.array_str(M, precision=2)

x = sim.M_ortho_props(M)
print np.array_str(x, precision=2)


[[  2.22e-04  -1.33e-05  -1.78e-05   0.00e+00   0.00e+00   0.00e+00]
 [ -1.33e-05   4.35e-04  -1.30e-04   0.00e+00   0.00e+00   0.00e+00]
 [ -1.78e-05  -1.30e-04   3.70e-04   0.00e+00   0.00e+00   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   4.55e-04   0.00e+00   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   0.00e+00   4.76e-04   0.00e+00]
 [  0.00e+00   0.00e+00   0.00e+00   0.00e+00   0.00e+00   4.17e-04]]
[  4.50e+03   2.30e+03   2.70e+03   6.00e-02   8.00e-02   3.00e-01
   2.20e+03   2.10e+03   2.40e+03]

In [ ]: