.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/plates/plot_plate05_quad_test2.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code. .. rst-class:: sphx-glr-example-title .. _sphx_glr_auto_examples_plates_plot_plate05_quad_test2.py: =============================================================== A square patch made of one quadrilateral plate elements =============================================================== Basic implementation test with applied loads. Testing the tangent stiffness computation. .. code:: free free ^ ^ | | 3-----2 -> free | | > | a | > (w = 1.0) | | > 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] Author: Peter Mackenzie-Helnwein .. GENERATED FROM PYTHON SOURCE LINES 36-131 .. code-block:: Python import numpy as np from femedu.examples import Example from femedu.domain import System, Node from femedu.solver import NewtonRaphsonSolver from femedu.elements.linear import Quad from femedu.materials import PlaneStress class ExamplePlate08(Example): # sphinx_gallery_thumbnail_number = 2 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 c = 1.5 model = System() model.setSolver(NewtonRaphsonSolver()) nd0 = Node( 0.0, 0.0) nd1 = Node(a/2+c,0.0) nd2 = Node( a, 0.0) nd3 = Node( 0.0, b) nd4 = Node( a/2-c, b) nd5 = Node( a, b) # nd0.fixDOF('ux', 'uy') # nd1.fixDOF('uy') # nd2.fixDOF('uy') # nd3.fixDOF('ux') model.addNode(nd0, nd1, nd2, nd3, nd4, nd5) elemA = Quad(nd0, nd1, nd4, nd3, PlaneStress(params)) elemB = Quad(nd1, nd2, nd5, nd4, PlaneStress(params)) model.addElement(elemA, elemB) #elemA.setSurfaceLoad(face=3, pn=1.0) elemB.setSurfaceLoad(face=1, pn=1.0) model.plot(factor=0.0, title="Undeformed system", filename="plate08_undeformed.png", show_bc=1) # %% # We can have a quick look at the stiffness mode shapes using the # buckling-mode plotter. These are simply eigenvalues and eigenvectors of Kt # at the current load level (0.0) # model.setLoadFactor(0.0) model.solve(tol=1000.) for k in range(8): name = f"plate08_mode{k:2d}.png" model.plotBucklingMode(mode=k,filename=name,factor=1.0) # %% # Note the three rigid body modes (lam=0.0). It can be shown that all three # are limear combinations of translations in x and y-directions and a # rigid body rotation. # # %% # Now it is time to add boundary conditions, apply loads # and check the convergence behavior. # nd0.fixDOF('ux', 'uy') nd1.fixDOF('uy') nd2.fixDOF('uy') nd3.fixDOF('ux') model.setLoadFactor(1.0) model.solve() # %% # The output shows that we do have a quadratic rate of convergence. # %% # Let's finish off with a nice plot of the deformed system. model.plot(factor=1.0, filename="plate08_deformed.png") #model.report() .. GENERATED FROM PYTHON SOURCE LINES 162-164 Run the example by creating an instance of the problem and executing it by calling :py:meth:`Example.run()` .. GENERATED FROM PYTHON SOURCE LINES 164-170 .. code-block:: Python if __name__ == "__main__": ex = ExamplePlate08() ex.run() .. rst-class:: sphx-glr-horizontal * .. image-sg:: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_001.png :alt: Undeformed system :srcset: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_001.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_002.png :alt: Mode Shape for $ \lambda = 0.00 $ :srcset: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_002.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_003.png :alt: Mode Shape for $ \lambda = 0.00 $ :srcset: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_003.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_004.png :alt: Mode Shape for $ \lambda = 0.00 $ :srcset: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_004.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_005.png :alt: Mode Shape for $ \lambda = 3.56 $ :srcset: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_005.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_006.png :alt: Mode Shape for $ \lambda = 4.57 $ :srcset: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_006.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_007.png :alt: Mode Shape for $ \lambda = 5.06 $ :srcset: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_007.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_008.png :alt: Mode Shape for $ \lambda = 5.83 $ :srcset: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_008.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_009.png :alt: Mode Shape for $ \lambda = 8.55 $ :srcset: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_009.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_010.png :alt: Deformed System (magnification=1.00) :srcset: /auto_examples/plates/images/sphx_glr_plot_plate05_quad_test2_010.png :class: sphx-glr-multi-img .. rst-class:: sphx-glr-script-out .. code-block:: none + + .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 0.237 seconds) .. _sphx_glr_download_auto_examples_plates_plot_plate05_quad_test2.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: plot_plate05_quad_test2.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_plate05_quad_test2.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: plot_plate05_quad_test2.zip ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_