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.00017989045056298722
CG iteration 2, residual = 7.59308964480149e-05
CG iteration 3, residual = 9.553190911073876e-05
CG iteration 4, residual = 6.915155168636542e-05
CG iteration 5, residual = 6.408034527301734e-05
CG iteration 6, residual = 4.6957231935788064e-05
CG iteration 7, residual = 4.0742471318053424e-05
CG iteration 8, residual = 2.964351045495067e-05
CG iteration 9, residual = 2.3189596705095004e-05
CG iteration 10, residual = 1.7242540984822564e-05
CG iteration 11, residual = 1.3840252847566092e-05
CG iteration 12, residual = 1.1371290859109187e-05
CG iteration 13, residual = 7.557236098508677e-06
CG iteration 14, residual = 5.772464004177656e-06
CG iteration 15, residual = 4.233341381713331e-06
CG iteration 16, residual = 3.465784460724063e-06
CG iteration 17, residual = 2.3625074117511024e-06
CG iteration 18, residual = 1.960334962442812e-06
CG iteration 19, residual = 1.468798970318046e-06
CG iteration 20, residual = 1.1178053345310793e-06
CG iteration 21, residual = 7.786087109406196e-07
CG iteration 22, residual = 5.912354552243819e-07
CG iteration 23, residual = 4.871534255435908e-07
CG iteration 24, residual = 3.1939698996758453e-07
CG iteration 25, residual = 2.644899491547712e-07
CG iteration 26, residual = 1.71317854586033e-07
CG iteration 27, residual = 1.5039339583039652e-07
CG iteration 28, residual = 1.0180139619024776e-07
CG iteration 29, residual = 7.5089354504632e-08
CG iteration 30, residual = 5.485591827515106e-08
CG iteration 31, residual = 4.239828755821844e-08
CG iteration 32, residual = 2.7875622339846714e-08
CG iteration 33, residual = 2.192466867596213e-08
CG iteration 34, residual = 1.612943181509729e-08
CG iteration 35, residual = 1.1134506230504945e-08
CG iteration 36, residual = 8.848794024774418e-09
CG iteration 37, residual = 5.853707439816875e-09
CG iteration 38, residual = 4.719993644160393e-09
CG iteration 39, residual = 3.069585047196085e-09
CG iteration 40, residual = 2.2878553041236014e-09
CG iteration 41, residual = 1.6665810402956316e-09
CG iteration 42, residual = 1.6829294119738026e-09
CG iteration 43, residual = 1.0264091347491698e-09
CG iteration 44, residual = 8.222899514625089e-10
CG iteration 45, residual = 5.030388799378177e-10
CG iteration 46, residual = 3.9802271320265304e-10
CG iteration 47, residual = 3.022798438839758e-10
CG iteration 48, residual = 2.041878504699321e-10
CG iteration 49, residual = 1.5113239479211898e-10
CG iteration 50, residual = 1.1537262074546107e-10
CG iteration 51, residual = 9.966447212371691e-11
CG iteration 52, residual = 7.22036346443098e-11
CG iteration 53, residual = 4.518676860938005e-11
CG iteration 54, residual = 3.974562314988915e-11
CG iteration 55, residual = 3.178909911058241e-11
CG iteration 56, residual = 2.0096556396900874e-11
CG iteration 57, residual = 1.4162845004822439e-11
CG iteration 58, residual = 1.0499841755272722e-11
CG iteration 59, residual = 7.383181982235099e-12
CG iteration 60, residual = 4.93456527014089e-12
CG iteration 61, residual = 4.166765777420578e-12
CG iteration 62, residual = 2.743044534172076e-12
CG iteration 63, residual = 1.7989282291377607e-12
CG iteration 64, residual = 1.4304344059296405e-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);
[ ]: