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.0001798273150819834
CG iteration 2, residual = 7.527505823604908e-05
CG iteration 3, residual = 8.667360647586901e-05
CG iteration 4, residual = 6.611763953432723e-05
CG iteration 5, residual = 6.0095304950833445e-05
CG iteration 6, residual = 4.24015118164131e-05
CG iteration 7, residual = 3.1746596641092696e-05
CG iteration 8, residual = 2.7007176859764126e-05
CG iteration 9, residual = 2.0093499070754567e-05
CG iteration 10, residual = 1.6202784652567277e-05
CG iteration 11, residual = 1.1164241035834792e-05
CG iteration 12, residual = 8.707361209404838e-06
CG iteration 13, residual = 6.2633002189683115e-06
CG iteration 14, residual = 4.827662227893142e-06
CG iteration 15, residual = 3.5048301823911805e-06
CG iteration 16, residual = 2.674487250730767e-06
CG iteration 17, residual = 1.9959626253864133e-06
CG iteration 18, residual = 1.417393961500008e-06
CG iteration 19, residual = 1.0735929420646392e-06
CG iteration 20, residual = 7.342807879188788e-07
CG iteration 21, residual = 6.364131834900322e-07
CG iteration 22, residual = 4.13030782908918e-07
CG iteration 23, residual = 3.6675837210598077e-07
CG iteration 24, residual = 2.5942494636472225e-07
CG iteration 25, residual = 1.8056292944575385e-07
CG iteration 26, residual = 1.3642541682194336e-07
CG iteration 27, residual = 9.96538334015704e-08
CG iteration 28, residual = 7.638450216779979e-08
CG iteration 29, residual = 5.2824180300395875e-08
CG iteration 30, residual = 4.0166965995025714e-08
CG iteration 31, residual = 2.853050965153184e-08
CG iteration 32, residual = 2.1713860987773947e-08
CG iteration 33, residual = 1.695593122041132e-08
CG iteration 34, residual = 1.1309617391624032e-08
CG iteration 35, residual = 8.114589692957366e-09
CG iteration 36, residual = 5.851886240215636e-09
CG iteration 37, residual = 4.026669357050303e-09
CG iteration 38, residual = 3.0617221485914917e-09
CG iteration 39, residual = 2.5243980277360883e-09
CG iteration 40, residual = 1.6419934217423393e-09
CG iteration 41, residual = 1.1584813633399818e-09
CG iteration 42, residual = 9.166054512379719e-10
CG iteration 43, residual = 6.364022799733947e-10
CG iteration 44, residual = 4.198907294689068e-10
CG iteration 45, residual = 4.0766964806976395e-10
CG iteration 46, residual = 2.7020087958916477e-10
CG iteration 47, residual = 1.917679894951738e-10
CG iteration 48, residual = 1.478350469255369e-10
CG iteration 49, residual = 1.2578997411319497e-10
CG iteration 50, residual = 9.219849373017894e-11
CG iteration 51, residual = 6.5451592838057e-11
CG iteration 52, residual = 4.127598947610903e-11
CG iteration 53, residual = 3.0961514955158123e-11
CG iteration 54, residual = 2.0256104285748734e-11
CG iteration 55, residual = 1.5070684923169194e-11
CG iteration 56, residual = 1.1829538676432467e-11
CG iteration 57, residual = 7.406496945825333e-12
CG iteration 58, residual = 5.29412372172693e-12
CG iteration 59, residual = 3.45466822179098e-12
CG iteration 60, residual = 2.4050730427341846e-12
CG iteration 61, residual = 1.983243780240082e-12
CG iteration 62, residual = 1.6164852296036066e-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);
[ ]: