Note
Go to the end to download the full example code
A square patch made of two triangular plate elements
- Basic implementation test with applied loads.
Testing the tangent stiffness computation for a
Triangle()
(using linear shape functions).
Using
elements.linear.Triangle
(see Triangle class)materials.PlaneStress
(see PlaneStress material class)
free free
^ ^
| |
3-----2 -> free
|\ b | >
| \ | >
| \ | > (w = 1.0)
| \ | >
| a \| >
0-----1 -> free
width: 10.
height: 10.
Material parameters: St. Venant-Kirchhoff, plane stress
E = 10.0
nu = 0.30
t = 1.0
Element loads:
node 0: [ 0.0, 0.0]
node 1: [ 5.0, 0.0]
node 2: [ 5.0, 0.0]
node 3: [ 0.0, 0.0]
2nd Piola-Kirchhoff stress:
S_XX = w = 1.000
S_YY = S_XY = S_YX = S_ZZ = 0.000
Green Lagrange strain:
eps_XX = (1/E) ((1.000) - (0.30)(0.000)) = 0.100
eps_YY = (1/E) ((0.000) - (0.30)(1.000)) = -0.030
eps_XY = eps_YX = 0.000
eps_ZZ = -nu * (eps_XX + eps_YY) = -0.021
Stretches:
lam_X = sqrt(1 + 2 eps_XX) = 1.095
lam_Y = sqrt(1 + 2 eps_YY) = 0.9695
Displacements:
ux = (lam_X - 1) * x, uy = (lam_Y - 1) * y
node 0: [ 0.000, 0.000 ]
node 1: [ 0.954, 0.000 ]
node 2: [ 0.954, -0.305 ]
node 3: [ 0.000, -0.305 ]
Author: Peter Mackenzie-Helnwein
from femedu.examples import Example
from femedu.domain import System, Node
from femedu.solver import NewtonRaphsonSolver
from femedu.elements.linear import Triangle
from femedu.materials import PlaneStress
class ExamplePlate02(Example):
def problem(self):
params = dict(
E = 10., # Young's modulus
nu = 0.3, # Poisson's ratio
t = 1.0, # thickness of the plate
fy = 1.e30 # yield stress
)
a = 10. # length of the plate in the x-direction
b = 10. # length of the plate in the y-direction
model = System()
model.setSolver(NewtonRaphsonSolver())
nd0 = Node( 0.0, 0.0)
nd1 = Node( a, 0.0)
nd2 = Node( a, b)
nd3 = Node( 0.0, b)
nd0.fixDOF('ux', 'uy')
nd1.fixDOF('uy')
nd3.fixDOF('ux')
model.addNode(nd0, nd1, nd2, nd3)
elemA = Triangle(nd0, nd1, nd3, PlaneStress(params))
elemB = Triangle(nd2, nd3, nd1, PlaneStress(params))
model.addElement(elemA, elemB)
elemB.setSurfaceLoad(face=2, pn=1.0)
model.setLoadFactor(0.0)
model.solve()
#model.report() # activate this line for lots of debug info
model.plot(factor=0.0, title="Undeformed system", filename="plate02_undeformed.png", show_bc=1)
model.setLoadFactor(1.0)
model.solve()
model.plot(factor=1.0, filename="plate02_deformed.png")
model.report()
Run the example by creating an instance of the problem and executing it by calling Example.run()
if __name__ == "__main__":
ex = ExamplePlate02()
ex.run()
+
/usr/local/lib/python3.11/site-packages/matplotlib/quiver.py:632: RuntimeWarning: Mean of empty slice.
amean = a.mean()
/usr/local/lib/python3.11/site-packages/numpy/core/_methods.py:129: RuntimeWarning: invalid value encountered in scalar divide
ret = ret.dtype.type(ret / rcount)
+
System Analysis Report
=======================
Nodes:
---------------------
Node_633:
x: [0. 0.]
fix: ['ux', 'uy']
u: [0. 0.]
Node_634:
x: [10. 0.]
fix: ['uy']
u: [1. 0.]
Node_635:
x: [10. 10.]
u: [ 1. -0.3]
Node_636:
x: [ 0. 10.]
fix: ['ux']
u: [ 0. -0.3]
Elements:
---------------------
Triangle_917: nodes ( Node_633 Node_634 Node_636 )
material: PlaneStress
strain: xx=1.000e-01 yy=-3.000e-02 xy=0.000e+00 zz=-2.100e-02
stress: xx=1.000e+00 yy=0.000e+00 xy=0.000e+00 zz=0.000e+00
Triangle_918: nodes ( Node_635 Node_636 Node_634 )
material: PlaneStress
strain: xx=1.000e-01 yy=-3.000e-02 xy=1.776e-16 zz=-2.100e-02
stress: xx=1.000e+00 yy=2.442e-15 xy=6.832e-16 zz=0.000e+00
Total running time of the script: (0 minutes 0.252 seconds)