A square patch made of two triangular plate elements

Basic implementation test with applied loads.

Testing the tangent stiffness computation for a Triangle6() (using quadratic shape functions).

Using

free   free
 ^     ^
 |     |
 3--6--2 -> free
 |\  b | >
 | \   | >
 7  8  5 > (w = 1.0)
 |   \ | >
 | a  \| >
 0--4--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: [ 10./6, 0.0]
    node 2: [ 10./6, 0.0]
    node 3: [ 0.0, 0.0]
    node 4: [ 0.0, 0.0]
    node 5: [ 10.*2/3, 0.0]
    node 6: [ 0.0, 0.0]
    node 7: [ 0.0, 0.0]
    node 8: [ 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  ]
    node 4: [ 0.477,  0.000  ]
    node 5: [ 0.954, -0.1525 ]
    node 6: [ 0.477, -0.305  ]
    node 7: [ 0.954, -0.1525 ]
    node 8: [ 0.477, -0.1525 ]

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 Triangle6
from femedu.materials import PlaneStress


class ExamplePlate02b(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)
        nd4 = Node( a/2, 0.0)
        nd5 = Node(   a, b/2)
        nd6 = Node( a/2,   b)
        nd7 = Node( 0.0, b/2)
        nd8 = Node( a/2, b/2)

        nd0.fixDOF('ux', 'uy')
        nd1.fixDOF('uy')
        nd3.fixDOF('ux')

        model.addNode(nd0, nd1, nd2, nd3, nd4, nd5, nd6, nd7, nd8)

        elemA = Triangle6(nd0, nd1, nd3, nd4, nd8, nd7, PlaneStress(params))
        elemB = Triangle6(nd2, nd3, nd1, nd6, nd8, nd5, PlaneStress(params))

        model.addElement(elemA, elemB)

        elemA.setSurfaceLoad(face=2, pn=1.0)
        elemB.setSurfaceLoad(face=2, pn=1.0)

        model.setLoadFactor(1.0)

        nd0.setDisp([0.0,  0.0])
        nd1.setDisp([5.0,  0.0])
        nd2.setDisp([5.0, -5.0])
        nd3.setDisp([0.0, -5.0])
        nd4.setDisp([2.5,  0.0])
        nd5.setDisp([5.0, -2.5])
        nd6.setDisp([2.5, -5.0])
        nd7.setDisp([0.0, -2.5])
        nd8.setDisp([2.5, -2.5])

        elemA.updateState()
        elemB.updateState()

        model.report()

        model.plot(factor=1.0, filename="plate02b_deformed.png")




        model.setLoadFactor(0.0)

        # model.solver.assemble()
        # model.solver.showKt()
        #
        model.solve()

        model.report()  # activate this line for lots of debug info
        model.plot(factor=0.0, title="Undeformed system", filename="plate02b_undeformed.png", show_bc=1)





        model.setLoadFactor(1.0)
        model.solve(verbose=1)
        model.plot(factor=1.0, filename="plate02b_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 = ExamplePlate02b()
    ex.run()
  • Deformed System (magnification=1.00)
  • Undeformed system
  • Deformed System (magnification=1.00)
System Analysis Report
=======================

Nodes:
---------------------
  Node_0:
      x:    [0. 0.]
      fix:  ['ux', 'uy']
      u:    [0. 0.]
  Node_1:
      x:    [10.  0.]
      fix:  ['uy']
      u:    [5. 0.]
  Node_2:
      x:    [10. 10.]
      u:    [ 5. -5.]
  Node_3:
      x:    [ 0. 10.]
      fix:  ['ux']
      u:    [ 0. -5.]
  Node_4:
      x:    [5. 0.]
      u:    [2.5 0. ]
  Node_5:
      x:    [10.  5.]
      u:    [ 5.  -2.5]
  Node_6:
      x:    [ 5. 10.]
      u:    [ 2.5 -5. ]
  Node_7:
      x:    [0. 5.]
      u:    [ 0.  -2.5]
  Node_8:
      x:    [5. 5.]
      u:    [ 2.5 -2.5]

Elements:
---------------------
  Triangle6_0: nodes ( Node_0 Node_1 Node_3 Node_4 Node_8 Node_7 )
      material: PlaneStress
      strain 0: xx=6.250e-01 yy=-3.750e-01 xy=0.000e+00 zz=-7.500e-02
      stress 0: xx=5.632e+00 yy=-2.060e+00 xy=0.000e+00 zz=0.000e+00
      strain 1: xx=6.250e-01 yy=-3.750e-01 xy=0.000e+00 zz=-7.500e-02
      stress 1: xx=5.632e+00 yy=-2.060e+00 xy=0.000e+00 zz=0.000e+00
      strain 2: xx=6.250e-01 yy=-3.750e-01 xy=0.000e+00 zz=-7.500e-02
      stress 2: xx=5.632e+00 yy=-2.060e+00 xy=0.000e+00 zz=0.000e+00
  Triangle6_1: nodes ( Node_2 Node_3 Node_1 Node_6 Node_8 Node_5 )
      material: PlaneStress
      strain 0: xx=6.250e-01 yy=-3.750e-01 xy=-1.865e-15 zz=-7.500e-02
      stress 0: xx=5.632e+00 yy=-2.060e+00 xy=-7.174e-15 zz=0.000e+00
      strain 1: xx=6.250e-01 yy=-3.750e-01 xy=-2.665e-16 zz=-7.500e-02
      stress 1: xx=5.632e+00 yy=-2.060e+00 xy=-1.025e-15 zz=0.000e+00
      strain 2: xx=6.250e-01 yy=-3.750e-01 xy=-7.994e-16 zz=-7.500e-02
      stress 2: xx=5.632e+00 yy=-2.060e+00 xy=-3.074e-15 zz=0.000e+00

+

System Analysis Report
=======================

Nodes:
---------------------
  Node_0:
      x:    [0. 0.]
      fix:  ['ux', 'uy']
      u:    [0. 0.]
  Node_1:
      x:    [10.  0.]
      fix:  ['uy']
      u:    [-2.2004262e-12  0.0000000e+00]
  Node_2:
      x:    [10. 10.]
      u:    [-2.19815779e-12  1.08940151e-11]
  Node_3:
      x:    [ 0. 10.]
      fix:  ['ux']
      u:    [0.00000000e+00 1.08945491e-11]
  Node_4:
      x:    [5. 0.]
      u:    [-1.10086125e-12  3.39908785e-17]
  Node_5:
      x:    [10.  5.]
      u:    [-2.20151180e-12  5.44645326e-12]
  Node_6:
      x:    [ 5. 10.]
      u:    [-1.09955281e-12  1.08931920e-11]
  Node_7:
      x:    [0. 5.]
      u:    [3.17666985e-16 5.44624201e-12]
  Node_8:
      x:    [5. 5.]
      u:    [-1.10032167e-12  5.44614435e-12]

Elements:
---------------------
  Triangle6_0: nodes ( Node_0 Node_1 Node_3 Node_4 Node_8 Node_7 )
      material: PlaneStress
      strain 0: xx=-2.198e-13 yy=1.089e-12 xy=6.353e-17 zz=-2.609e-13
      stress 0: xx=1.176e-12 yy=1.125e-11 xy=2.444e-16 zz=0.000e+00
      strain 1: xx=-2.198e-13 yy=1.089e-12 xy=7.054e-33 zz=-2.609e-13
      stress 1: xx=1.176e-12 yy=1.125e-11 xy=2.713e-32 zz=0.000e+00
      strain 2: xx=-2.200e-13 yy=1.090e-12 xy=-6.353e-17 zz=-2.609e-13
      stress 2: xx=1.174e-12 yy=1.125e-11 xy=-2.444e-16 zz=0.000e+00
  Triangle6_1: nodes ( Node_2 Node_3 Node_1 Node_6 Node_8 Node_5 )
      material: PlaneStress
      strain 0: xx=-2.199e-13 yy=1.089e-12 xy=3.553e-16 zz=-2.608e-13
      stress 0: xx=1.174e-12 yy=1.124e-11 xy=1.366e-15 zz=0.000e+00
      strain 1: xx=-2.198e-13 yy=1.089e-12 xy=-3.553e-16 zz=-2.609e-13
      stress 1: xx=1.176e-12 yy=1.125e-11 xy=-1.366e-15 zz=0.000e+00
      strain 2: xx=-2.203e-13 yy=1.089e-12 xy=-3.553e-16 zz=-2.607e-13
      stress 2: xx=1.171e-12 yy=1.124e-11 xy=-1.366e-15 zz=0.000e+00

/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)
norm of the out-of-balance force:   1.2019e+01
norm of the out-of-balance force:   2.4414e+00
norm of the out-of-balance force:   5.2114e-02
norm of the out-of-balance force:   2.8292e-05
norm of the out-of-balance force:   9.9149e-12
+

System Analysis Report
=======================

Nodes:
---------------------
  Node_0:
      x:    [0. 0.]
      fix:  ['ux', 'uy']
      u:    [0. 0.]
  Node_1:
      x:    [10.  0.]
      fix:  ['uy']
      u:    [0.93634551 0.        ]
  Node_2:
      x:    [10. 10.]
      u:    [ 0.94185675 -0.30050703]
  Node_3:
      x:    [ 0. 10.]
      fix:  ['ux']
      u:    [ 0.         -0.31608555]
  Node_4:
      x:    [5. 0.]
      u:    [0.4643119  0.01196108]
  Node_5:
      x:    [10.  5.]
      u:    [ 0.99272035 -0.16590176]
  Node_6:
      x:    [ 5. 10.]
      u:    [ 0.47573439 -0.31903817]
  Node_7:
      x:    [0. 5.]
      u:    [-0.05024981 -0.14188329]
  Node_8:
      x:    [5. 5.]
      u:    [ 0.47172318 -0.15330258]

Elements:
---------------------
  Triangle6_0: nodes ( Node_0 Node_1 Node_3 Node_4 Node_8 Node_7 )
      material: PlaneStress
      strain 0: xx=1.011e-01 yy=-2.840e-02 xy=-8.894e-03 zz=-2.181e-02
      stress 0: xx=1.017e+00 yy=2.126e-02 xy=-3.421e-02 zz=0.000e+00
      strain 1: xx=1.028e-01 yy=-3.297e-02 xy=-8.777e-04 zz=-2.095e-02
      stress 1: xx=1.021e+00 yy=-2.341e-02 xy=-3.376e-03 zz=0.000e+00
      strain 2: xx=1.138e-01 yy=-3.464e-02 xy=8.745e-03 zz=-2.375e-02
      stress 2: xx=1.137e+00 yy=-5.469e-03 xy=3.363e-02 zz=0.000e+00
  Triangle6_1: nodes ( Node_2 Node_3 Node_1 Node_6 Node_8 Node_5 )
      material: PlaneStress
      strain 0: xx=1.012e-01 yy=-2.751e-02 xy=-8.782e-03 zz=-2.211e-02
      stress 0: xx=1.022e+00 yy=3.139e-02 xy=-3.378e-02 zz=0.000e+00
      strain 1: xx=1.033e-01 yy=-3.359e-02 xy=-9.099e-04 zz=-2.092e-02
      stress 1: xx=1.025e+00 yy=-2.850e-02 xy=-3.500e-03 zz=0.000e+00
      strain 2: xx=1.133e-01 yy=-3.356e-02 xy=8.834e-03 zz=-2.393e-02
      stress 2: xx=1.135e+00 yy=4.784e-03 xy=3.398e-02 zz=0.000e+00

Total running time of the script: (0 minutes 1.385 seconds)

Gallery generated by Sphinx-Gallery