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.0001798904505629855
CG iteration 2, residual = 7.593089644801498e-05
CG iteration 3, residual = 9.553190911073947e-05
CG iteration 4, residual = 6.915155168636519e-05
CG iteration 5, residual = 6.408034527301731e-05
CG iteration 6, residual = 4.695723193578819e-05
CG iteration 7, residual = 4.074247131805321e-05
CG iteration 8, residual = 2.964351045495062e-05
CG iteration 9, residual = 2.3189596705094828e-05
CG iteration 10, residual = 1.724254098482253e-05
CG iteration 11, residual = 1.3840252847566067e-05
CG iteration 12, residual = 1.1371290859109194e-05
CG iteration 13, residual = 7.55723609850864e-06
CG iteration 14, residual = 5.772464004177674e-06
CG iteration 15, residual = 4.233341381713366e-06
CG iteration 16, residual = 3.465784460724156e-06
CG iteration 17, residual = 2.3625074117511862e-06
CG iteration 18, residual = 1.96033496244296e-06
CG iteration 19, residual = 1.4687989703180572e-06
CG iteration 20, residual = 1.117805334531053e-06
CG iteration 21, residual = 7.786087109406105e-07
CG iteration 22, residual = 5.912354552243833e-07
CG iteration 23, residual = 4.871534255435869e-07
CG iteration 24, residual = 3.193969899675837e-07
CG iteration 25, residual = 2.6448994915477207e-07
CG iteration 26, residual = 1.7131785458603234e-07
CG iteration 27, residual = 1.5039339583039724e-07
CG iteration 28, residual = 1.0180139619024727e-07
CG iteration 29, residual = 7.508935450463137e-08
CG iteration 30, residual = 5.485591827514987e-08
CG iteration 31, residual = 4.2398287558210285e-08
CG iteration 32, residual = 2.7875622339803115e-08
CG iteration 33, residual = 2.1924668675121065e-08
CG iteration 34, residual = 1.6129431809651494e-08
CG iteration 35, residual = 1.1134506188643242e-08
CG iteration 36, residual = 8.848793310071993e-09
CG iteration 37, residual = 5.853704118524823e-09
CG iteration 38, residual = 4.7199276521340486e-09
CG iteration 39, residual = 3.069199227542277e-09
CG iteration 40, residual = 2.2814084553332596e-09
CG iteration 41, residual = 1.5323165169287548e-09
CG iteration 42, residual = 1.217548757393583e-09
CG iteration 43, residual = 9.284977148519436e-10
CG iteration 44, residual = 8.927200969933592e-10
CG iteration 45, residual = 5.267412598273425e-10
CG iteration 46, residual = 3.99390888088927e-10
CG iteration 47, residual = 3.061290837350734e-10
CG iteration 48, residual = 2.2723263709788606e-10
CG iteration 49, residual = 2.1436108560640088e-10
CG iteration 50, residual = 1.3552250248529552e-10
CG iteration 51, residual = 9.326552328670466e-11
CG iteration 52, residual = 7.09057242736711e-11
CG iteration 53, residual = 5.984101563816124e-11
CG iteration 54, residual = 4.608242963856395e-11
CG iteration 55, residual = 2.7518452721281784e-11
CG iteration 56, residual = 1.984352449382903e-11
CG iteration 57, residual = 1.4157550362007406e-11
CG iteration 58, residual = 1.0499528393240154e-11
CG iteration 59, residual = 7.383169262609538e-12
CG iteration 60, residual = 4.93453854109137e-12
CG iteration 61, residual = 4.166590316396411e-12
CG iteration 62, residual = 2.742274009782412e-12
CG iteration 63, residual = 1.7889978972486464e-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);
[ ]: