This page was generated from wta/elasticity3D.ipynb.

3D Solid Mechanics

[1]:
from netgen.occ import *
from netgen.webgui import Draw as DrawGeo
import ngsolve
[2]:
box = Box((0,0,0), (3,0.6,1))
box.faces.name="outer"
cyl = sum( [Cylinder((0.5+i,0,0.5), Y, 0.25,0.8) for i in range(3)] )
cyl.faces.name="cyl"
geo = box-cyl

DrawGeo(geo);

find edges between box and cylinder, and build chamfers (requires OCC 7.4 or newer):

[3]:
cylboxedges = geo.faces["outer"].edges * geo.faces["cyl"].edges
cylboxedges.name = "cylbox"
geo = geo.MakeChamfer(cylboxedges, 0.03)

name faces for boundary conditions:

[4]:
geo.faces.Min(X).name = "fix"
geo.faces.Max(X).name = "force"

DrawGeo(geo);
[5]:
from ngsolve import *
from ngsolve.webgui import Draw
mesh = Mesh(OCCGeometry(geo).GenerateMesh(maxh=0.1)).Curve(3)
Draw (mesh);

Linear elasticity

Displacement: \(u : \Omega \rightarrow R^3\)

Linear strain:

\[\varepsilon(u) := \tfrac{1}{2} ( \nabla u + (\nabla u)^T )\]

Stress by Hooke’s law:

\[\sigma = 2 \mu \varepsilon + \lambda \operatorname{tr} \varepsilon I\]

Equilibrium of forces:

\[\operatorname{div} \sigma = f\]

Displacement boundary conditions:

\[u = u_D \qquad \text{on} \, \Gamma_D\]

Traction boundary conditions:

\[\sigma n = g \qquad \text{on} \, \Gamma_N\]

Variational formulation:

Find: \(u \in H^1(\Omega)^3\) such that \(u = u_D\) on \(\Gamma_D\)

\[\int_\Omega \sigma(\varepsilon(u)) : \varepsilon(v) \, dx = \int_\Omega f v dx + \int_{\Gamma_N} g v ds\]

holds for all \(v = 0\) on \(\Gamma_D\).

[6]:
E, nu = 210, 0.2
mu  = E / 2 / (1+nu)
lam = E * nu / ((1+nu)*(1-2*nu))

def Stress(strain):
    return 2*mu*strain + lam*Trace(strain)*Id(3)
[7]:
fes = VectorH1(mesh, order=3, dirichlet="fix")
u,v = fes.TnT()
gfu = GridFunction(fes)

with TaskManager():
    a = BilinearForm(InnerProduct(Stress(Sym(Grad(u))), Sym(Grad(v))).Compile()*dx)
    pre = Preconditioner(a, "bddc")
    a.Assemble()
[8]:
force = CF( (1e-3,0,0) )
f = LinearForm(force*v*ds("force")).Assemble()
[9]:
from ngsolve.krylovspace import CGSolver
inv = CGSolver(a.mat, pre, printrates=True, tol=1e-8)
gfu.vec.data = inv * f.vec
CG iteration 1, residual = 0.00017989045056298667
CG iteration 2, residual = 7.59308964480146e-05
CG iteration 3, residual = 9.553190911073734e-05
CG iteration 4, residual = 6.915155168636555e-05
CG iteration 5, residual = 6.408034527301735e-05
CG iteration 6, residual = 4.695723193578787e-05
CG iteration 7, residual = 4.074247131805333e-05
CG iteration 8, residual = 2.964351045495066e-05
CG iteration 9, residual = 2.3189596705094926e-05
CG iteration 10, residual = 1.7242540984822378e-05
CG iteration 11, residual = 1.384025284756599e-05
CG iteration 12, residual = 1.1371290859109167e-05
CG iteration 13, residual = 7.557236098508645e-06
CG iteration 14, residual = 5.7724640041776386e-06
CG iteration 15, residual = 4.233341381713356e-06
CG iteration 16, residual = 3.465784460724113e-06
CG iteration 17, residual = 2.3625074117511388e-06
CG iteration 18, residual = 1.9603349624429085e-06
CG iteration 19, residual = 1.4687989703180474e-06
CG iteration 20, residual = 1.1178053345310562e-06
CG iteration 21, residual = 7.786087109406076e-07
CG iteration 22, residual = 5.912354552243856e-07
CG iteration 23, residual = 4.871534255435881e-07
CG iteration 24, residual = 3.1939698996758273e-07
CG iteration 25, residual = 2.6448994915477074e-07
CG iteration 26, residual = 1.7131785458603245e-07
CG iteration 27, residual = 1.5039339583039758e-07
CG iteration 28, residual = 1.0180139619024654e-07
CG iteration 29, residual = 7.508935450463062e-08
CG iteration 30, residual = 5.4855918275149467e-08
CG iteration 31, residual = 4.239828755821046e-08
CG iteration 32, residual = 2.7875622339808062e-08
CG iteration 33, residual = 2.192466867521554e-08
CG iteration 34, residual = 1.612943181026612e-08
CG iteration 35, residual = 1.1134506193370303e-08
CG iteration 36, residual = 8.848793390786646e-09
CG iteration 37, residual = 5.853704493673767e-09
CG iteration 38, residual = 4.7199351069567645e-09
CG iteration 39, residual = 3.0692428259226866e-09
CG iteration 40, residual = 2.282139549259245e-09
CG iteration 41, residual = 1.5489433668891869e-09
CG iteration 42, residual = 1.40237026809983e-09
CG iteration 43, residual = 1.1043216587386159e-09
CG iteration 44, residual = 8.426004121436543e-10
CG iteration 45, residual = 5.038484098351175e-10
CG iteration 46, residual = 3.9807508211627856e-10
CG iteration 47, residual = 3.024443341396366e-10
CG iteration 48, residual = 2.0539088405776478e-10
CG iteration 49, residual = 1.6563424655935663e-10
CG iteration 50, residual = 1.5250682851819616e-10
CG iteration 51, residual = 9.588624673016277e-11
CG iteration 52, residual = 7.065009588616109e-11
CG iteration 53, residual = 5.81416427325876e-11
CG iteration 54, residual = 4.690803905100855e-11
CG iteration 55, residual = 2.75487272834282e-11
CG iteration 56, residual = 1.9844703640900036e-11
CG iteration 57, residual = 1.4157578901846039e-11
CG iteration 58, residual = 1.0499534465754418e-11
CG iteration 59, residual = 7.383207006908612e-12
CG iteration 60, residual = 4.934955596283276e-12
CG iteration 61, residual = 4.169350317398152e-12
CG iteration 62, residual = 2.7542810775044715e-12
CG iteration 63, residual = 1.928754978451529e-12
CG iteration 64, residual = 2.0000899682503005e-12
CG iteration 65, residual = 1.054079573433895e-12
[10]:
with TaskManager():
    fesstress = MatrixValued(H1(mesh,order=3), symmetric=True)
    gfstress = GridFunction(fesstress)
    gfstress.Interpolate (Stress(Sym(Grad(gfu))))
[11]:
Draw (5e4*gfu, mesh);
[12]:
Draw (Norm(gfstress), mesh, deformation=1e4*gfu, draw_vol=False, order=3);
[ ]: