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
ngsolve.__version__
[1]:
'6.2.2105-64-g8568b5f22'
[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
find edges between box and cylinder, and build chamfers:
[3]:
cylboxedges = geo.faces["outer"].edges * geo.faces["cyl"].edges
cylboxedges.name = "cylbox"
geo = geo.MakeChamfer(cylboxedges, 0.03)
geo.faces.Min(X).name = "fix"
geo.faces.Max(X).name = "force"
[4]:
DrawGeo(geo);
[5]:
from ngsolve import *
from ngsolve.webgui import Draw
mesh = Mesh(OCCGeometry(geo).GenerateMesh(maxh=0.1)).Curve(3)
Draw (mesh);
[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]:
f = LinearForm(1e-3*v[0]*ds("force")).Assemble()
[9]:
from ngsolve.krylovspace import CGSolver
inv = CGSolver(a.mat, pre, printrates='\r')
gfu.vec.data = inv * f.vec
CG converged in 93 iterations to residual 1.690213829965995e-16
[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=2*gfu, draw_vol=False, order=3);
[ ]: