2D_elasticity.py.html 25.9 KB
Newer Older
Jeremy BLEYER's avatar
Jeremy BLEYER committed
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240

<!DOCTYPE html PUBLIC "-//W3C//DTD XHTML 1.0 Transitional//EN"
  "http://www.w3.org/TR/xhtml1/DTD/xhtml1-transitional.dtd">

<html xmlns="http://www.w3.org/1999/xhtml">
  <head>
    <meta http-equiv="Content-Type" content="text/html; charset=utf-8" />
    <title>2D linear elasticity &#8212; Numerical tours of continuum mechanics using FEniCS master documentation</title>
    <link rel="stylesheet" href="../../_static/agogo.css" type="text/css" />
    <link rel="stylesheet" href="../../_static/pygments.css" type="text/css" />
    <script type="text/javascript">
      var DOCUMENTATION_OPTIONS = {
        URL_ROOT:    '../../',
        VERSION:     'master',
        COLLAPSE_INDEX: false,
        FILE_SUFFIX: '.html',
        HAS_SOURCE:  true,
        SOURCELINK_SUFFIX: '.txt'
      };
    </script>
    <script type="text/javascript" src="../../_static/jquery.js"></script>
    <script type="text/javascript" src="../../_static/underscore.js"></script>
    <script type="text/javascript" src="../../_static/doctools.js"></script>
    <script type="text/javascript" src="https://cdnjs.cloudflare.com/ajax/libs/mathjax/2.7.1/MathJax.js?config=TeX-AMS-MML_HTMLorMML"></script>
    <link rel="index" title="Index" href="../../genindex.html" />
    <link rel="search" title="Search" href="../../search.html" />
    <link rel="next" title="Orthotropic linear elasticity" href="orthotropic_elasticity.py.html" />
    <link rel="prev" title="Linear problems in solid mechanics" href="../../linear_problems.html" /> 
  </head>
  <body>
    <div class="header-wrapper" role="banner">
      <div class="header">
        <div class="headertitle"><a
          href="../../index.html">Numerical tours of continuum mechanics using FEniCS</a></div>
        <div class="rel" role="navigation" aria-label="related navigation">
          <a href="../../linear_problems.html" title="Linear problems in solid mechanics"
             accesskey="P">previous</a> |
          <a href="orthotropic_elasticity.py.html" title="Orthotropic linear elasticity"
             accesskey="N">next</a> |
          <a href="../../genindex.html" title="General Index"
             accesskey="I">index</a>
        </div>
       </div>
    </div>

    <div class="content-wrapper">
      <div class="content">
        <div class="document">
            
      <div class="documentwrapper">
        <div class="bodywrapper">
          <div class="body" role="main">
            
  <div class="section" id="d-linear-elasticity">
<span id="linearelasticity2d"></span><h1>2D linear elasticity<a class="headerlink" href="#d-linear-elasticity" title="Permalink to this headline"></a></h1>
<div class="section" id="introduction">
<h2>Introduction<a class="headerlink" href="#introduction" title="Permalink to this headline"></a></h2>
<p>In this first numerical tour, we will show how to compute a small strain solution for
a 2D isotropic linear elastic medium, either in plane stress or in plane strain,
in a tradtional displacement-based finite element formulation. The corresponding
file can be obtained from <a class="reference download internal" href="../../_downloads/2D_elasticity.py" download=""><code class="xref download docutils literal"><span class="pre">2D_elasticity.py</span></code></a>.</p>
<div class="admonition seealso">
<p class="first admonition-title">See also</p>
<p class="last">Extension to 3D is straightforward and an example can be found in the <a class="reference internal" href="../modal_analysis_dynamics/cantilever_modal.py.html#modalanalysis"><span class="std std-ref">Modal analysis of an elastic structure</span></a> example.</p>
</div>
<p>We consider here the case of a cantilever beam modeled as a 2D medium of dimensions
<span class="math">\(L\times  H\)</span>. Geometrical parameters and mesh density are first defined
and the rectangular domain is  generated using the <code class="docutils literal"><span class="pre">RectangleMesh</span></code> function.
We also choose a criss-crossed structured mesh:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="kn">from</span> <span class="nn">__future__</span> <span class="k">import</span> <span class="n">print_function</span>
<span class="kn">from</span> <span class="nn">fenics</span> <span class="k">import</span> <span class="o">*</span>

<span class="n">L</span> <span class="o">=</span> <span class="mf">25.</span>
<span class="n">H</span> <span class="o">=</span> <span class="mf">1.</span>
<span class="n">Nx</span> <span class="o">=</span> <span class="mi">250</span>
<span class="n">Ny</span> <span class="o">=</span> <span class="mi">10</span>
<span class="n">mesh</span> <span class="o">=</span> <span class="n">RectangleMesh</span><span class="p">(</span><span class="n">Point</span><span class="p">(</span><span class="mf">0.</span><span class="p">,</span> <span class="mf">0.</span><span class="p">),</span> <span class="n">Point</span><span class="p">(</span><span class="n">L</span><span class="p">,</span> <span class="n">H</span><span class="p">),</span> <span class="n">Nx</span><span class="p">,</span> <span class="n">Ny</span><span class="p">,</span> <span class="s2">&quot;crossed&quot;</span><span class="p">)</span>
</pre></div>
</div>
</div>
<div class="section" id="constitutive-relation">
<h2>Constitutive relation<a class="headerlink" href="#constitutive-relation" title="Permalink to this headline"></a></h2>
<p>We now define the material parameters which are here given in terms of a Young’s
modulus <span class="math">\(E\)</span> and a Poisson coefficient <span class="math">\(\nu\)</span>. In the following, we will
need to define the constitutive relation between the stress tensor <span class="math">\(\boldsymbol{\sigma}\)</span>
and the strain tensor <span class="math">\(\boldsymbol{\varepsilon}\)</span>. Let us recall
that the general expression of the linear elastic isotropic constitutive relation
for a 3D medium is given by:</p>
<div class="math" id="equation-constitutive_3D">
<span class="eqno">(1)<a class="headerlink" href="#equation-constitutive_3D" title="Permalink to this equation"></a></span>\[\boldsymbol{\sigma} = \lambda \text{tr}(\boldsymbol{\varepsilon})\mathbf{1} + 2\mu\boldsymbol{\varepsilon}\]</div>
<p>for a natural (no prestress) initial state where the Lamé coefficients are given by:</p>
<div class="math" id="equation-Lame_coeff">
<span class="eqno">(2)<a class="headerlink" href="#equation-Lame_coeff" title="Permalink to this equation"></a></span>\[\lambda = \dfrac{E\nu}{(1+\nu)(1-2\nu)}, \quad \mu = \dfrac{E}{2(1+\nu)}\]</div>
<p>In this demo, we consider a 2D model either in plane strain or in plane stress conditions.
Irrespective of this choice, we will work only with a 2D displacement vector <span class="math">\(\boldsymbol{u}=(u_x,u_y)\)</span>
and will subsequently define the strain operator <code class="docutils literal"><span class="pre">eps</span></code> as follows:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="k">def</span> <span class="nf">eps</span><span class="p">(</span><span class="n">v</span><span class="p">):</span>
    <span class="k">return</span> <span class="n">sym</span><span class="p">(</span><span class="n">grad</span><span class="p">(</span><span class="n">v</span><span class="p">))</span>
</pre></div>
</div>
<p>which computes the 2x2 plane components of the symmetrized gradient tensor of
any 2D vectorial field. In the plane strain case, the full 3D strain tensor is defined as follows:</p>
<div class="math">
\[\begin{split}\boldsymbol{\varepsilon} = \begin{bmatrix} \varepsilon_{xx} &amp; \varepsilon_{xy} &amp; 0\\
\varepsilon_{xy} &amp; \varepsilon_{yy} &amp; 0 \\ 0 &amp; 0 &amp; 0\end{bmatrix}\end{split}\]</div>
<p>so that the 2x2 plane part of the stress tensor is defined in the same way as for the 3D case
(the out-of-plane stress component being given by <span class="math">\(\sigma_{zz}=\lambda(\varepsilon_{xx}+\varepsilon_{yy})\)</span>.</p>
<p>In the plane stress case, an out-of-plane strain component <span class="math">\(\varepsilon_{zz}\)</span>
must be considered so that <span class="math">\(\sigma_{zz}=0\)</span>. Using this condition in the
3D constitutive relation, one has <span class="math">\(\varepsilon_{zz}=-\dfrac{\lambda}{\lambda+2\mu}(\varepsilon_{xx}+\varepsilon_{yy})\)</span>.
Injecting into <a class="reference internal" href="#equation-constitutive_3D">(1)</a>, we have for the 2D plane stress relation:</p>
<div class="math">
\[\boldsymbol{\sigma} = \lambda^* \text{tr}(\boldsymbol{\varepsilon})\mathbf{1} + 2\mu\boldsymbol{\varepsilon}\]</div>
<p>where <span class="math">\(\boldsymbol{\sigma}, \boldsymbol{\varepsilon}, \mathbf{1}\)</span> are 2D tensors and with
<span class="math">\(\lambda^* = \dfrac{2\lambda\mu}{\lambda+2\mu}\)</span>. Hence, the 2D constitutive relation
is identical to the plane strain case by changing only the value of the Lamé coefficient <span class="math">\(\lambda\)</span>.
We can then have:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="n">E</span> <span class="o">=</span> <span class="n">Constant</span><span class="p">(</span><span class="mf">1e5</span><span class="p">)</span>
<span class="n">nu</span> <span class="o">=</span> <span class="n">Constant</span><span class="p">(</span><span class="mf">0.3</span><span class="p">)</span>
<span class="n">model</span> <span class="o">=</span> <span class="s2">&quot;plane_stress&quot;</span>

<span class="n">mu</span> <span class="o">=</span> <span class="n">E</span><span class="o">/</span><span class="mi">2</span><span class="o">/</span><span class="p">(</span><span class="mi">1</span><span class="o">+</span><span class="n">nu</span><span class="p">)</span>
<span class="n">lmbda</span> <span class="o">=</span> <span class="n">E</span><span class="o">*</span><span class="n">nu</span><span class="o">/</span><span class="p">(</span><span class="mi">1</span><span class="o">+</span><span class="n">nu</span><span class="p">)</span><span class="o">/</span><span class="p">(</span><span class="mi">1</span><span class="o">-</span><span class="mi">2</span><span class="o">*</span><span class="n">nu</span><span class="p">)</span>
<span class="k">if</span> <span class="n">model</span> <span class="o">==</span> <span class="s2">&quot;plane_stress&quot;</span><span class="p">:</span>
    <span class="n">lmbda</span> <span class="o">=</span> <span class="mi">2</span><span class="o">*</span><span class="n">mu</span><span class="o">*</span><span class="n">lmbda</span><span class="o">/</span><span class="p">(</span><span class="n">lmbda</span><span class="o">+</span><span class="mi">2</span><span class="o">*</span><span class="n">mu</span><span class="p">)</span>

<span class="k">def</span> <span class="nf">sigma</span><span class="p">(</span><span class="n">v</span><span class="p">):</span>
    <span class="k">return</span> <span class="n">lmbda</span><span class="o">*</span><span class="n">tr</span><span class="p">(</span><span class="n">eps</span><span class="p">(</span><span class="n">v</span><span class="p">))</span><span class="o">*</span><span class="n">Identity</span><span class="p">(</span><span class="mi">2</span><span class="p">)</span> <span class="o">+</span> <span class="mf">2.0</span><span class="o">*</span><span class="n">mu</span><span class="o">*</span><span class="n">eps</span><span class="p">(</span><span class="n">v</span><span class="p">)</span>
</pre></div>
</div>
<div class="admonition note">
<p class="first admonition-title">Note</p>
<p>Note that we used the variable name <code class="docutils literal"><span class="pre">lmbda</span></code> to avoid any confusion with the
lambda functions of Python</p>
<p class="last">We also used an intrinsic formulation of the constitutive relation. Example of
constitutive relation implemented with a matrix/vector engineering notation
will be provided in the <a class="reference internal" href="orthotropic_elasticity.py.html#orthotropicelasticity"><span class="std std-ref">Orthotropic linear elasticity</span></a> example.</p>
</div>
</div>
<div class="section" id="variational-formulation">
<h2>Variational formulation<a class="headerlink" href="#variational-formulation" title="Permalink to this headline"></a></h2>
<p>For this example, we consider a continuous polynomial interpolation of degree 2
and a uniformly distributed loading <span class="math">\(\boldsymbol{f}=(0,-f)\)</span> corresponding
to the beam self-weight. The continuum mechanics variational formulation (obtained
from the virtual work principle) is given by:</p>
<div class="math">
\[\text{Find } \boldsymbol{u}\in V \text{ s.t. } \int_{\Omega}
\boldsymbol{\sigma}(\boldsymbol{u}):\boldsymbol{\varepsilon}(\boldsymbol{v}) d\Omega
= \int_{\Omega} \boldsymbol{f}\cdot\boldsymbol{v}  d\Omega \quad \forall\boldsymbol{v} \in V\]</div>
<p>which translates into the following FEniCS code:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="n">rho_g</span> <span class="o">=</span> <span class="mf">1e-3</span>
<span class="n">f</span> <span class="o">=</span> <span class="n">Constant</span><span class="p">((</span><span class="mi">0</span><span class="p">,</span><span class="o">-</span><span class="n">rho_g</span><span class="p">))</span>

<span class="n">V</span> <span class="o">=</span> <span class="n">VectorFunctionSpace</span><span class="p">(</span><span class="n">mesh</span><span class="p">,</span> <span class="s1">&#39;Lagrange&#39;</span><span class="p">,</span> <span class="n">degree</span><span class="o">=</span><span class="mi">2</span><span class="p">)</span>
<span class="n">du</span> <span class="o">=</span> <span class="n">TrialFunction</span><span class="p">(</span><span class="n">V</span><span class="p">)</span>
<span class="n">u_</span> <span class="o">=</span> <span class="n">TestFunction</span><span class="p">(</span><span class="n">V</span><span class="p">)</span>
<span class="n">a</span> <span class="o">=</span> <span class="n">inner</span><span class="p">(</span><span class="n">sigma</span><span class="p">(</span><span class="n">du</span><span class="p">),</span> <span class="n">eps</span><span class="p">(</span><span class="n">u_</span><span class="p">))</span><span class="o">*</span><span class="n">dx</span>
<span class="n">l</span> <span class="o">=</span> <span class="n">inner</span><span class="p">(</span><span class="n">f</span><span class="p">,</span> <span class="n">u_</span><span class="p">)</span><span class="o">*</span><span class="n">dx</span>
</pre></div>
</div>
</div>
<div class="section" id="resolution">
<h2>Resolution<a class="headerlink" href="#resolution" title="Permalink to this headline"></a></h2>
<p>Fixed displacements are imposed on the left part of the beam, the <code class="docutils literal"><span class="pre">solve</span></code>
function is then called and solution is plotted by deforming the mesh:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="k">def</span> <span class="nf">left</span><span class="p">(</span><span class="n">x</span><span class="p">,</span> <span class="n">on_boundary</span><span class="p">):</span>
    <span class="k">return</span> <span class="n">near</span><span class="p">(</span><span class="n">x</span><span class="p">[</span><span class="mi">0</span><span class="p">],</span><span class="mf">0.</span><span class="p">)</span>

<span class="n">bc</span> <span class="o">=</span> <span class="n">DirichletBC</span><span class="p">(</span><span class="n">V</span><span class="p">,</span> <span class="n">Constant</span><span class="p">((</span><span class="mf">0.</span><span class="p">,</span><span class="mf">0.</span><span class="p">)),</span> <span class="n">left</span><span class="p">)</span>

<span class="n">u</span> <span class="o">=</span> <span class="n">Function</span><span class="p">(</span><span class="n">V</span><span class="p">,</span> <span class="n">name</span><span class="o">=</span><span class="s2">&quot;Displacement&quot;</span><span class="p">)</span>
<span class="n">solve</span><span class="p">(</span><span class="n">a</span> <span class="o">==</span> <span class="n">l</span><span class="p">,</span> <span class="n">u</span><span class="p">,</span> <span class="n">bc</span><span class="p">)</span>

<span class="n">plot</span><span class="p">(</span><span class="mf">1e3</span><span class="o">*</span><span class="n">u</span><span class="p">,</span> <span class="n">mode</span><span class="o">=</span><span class="s2">&quot;displacement&quot;</span><span class="p">)</span>
</pre></div>
</div>
<p>The (amplified) solution should look like this:</p>
<a class="reference internal image-reference" href="../../_images/cantilever_deformed.png"><img alt="../../_images/cantilever_deformed.png" src="../../_images/cantilever_deformed.png" style="width: 754.8px; height: 247.5px;" /></a>
</div>
<div class="section" id="validation-and-post-processing">
<h2>Validation and post-processing<a class="headerlink" href="#validation-and-post-processing" title="Permalink to this headline"></a></h2>
<p>The maximal deflection is compared against the analytical solution from
Euler-Bernoulli beam theory which is here <span class="math">\(w_{beam} = \dfrac{qL^4}{8EI}\)</span>:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="nb">print</span><span class="p">(</span><span class="s2">&quot;Maximal deflection:&quot;</span><span class="p">,</span> <span class="o">-</span><span class="n">u</span><span class="p">(</span><span class="n">L</span><span class="p">,</span><span class="n">H</span><span class="o">/</span><span class="mf">2.</span><span class="p">)[</span><span class="mi">1</span><span class="p">])</span>
<span class="nb">print</span><span class="p">(</span><span class="s2">&quot;Beam theory deflection:&quot;</span><span class="p">,</span> <span class="nb">float</span><span class="p">(</span><span class="mi">3</span><span class="o">*</span><span class="n">rho_g</span><span class="o">*</span><span class="n">L</span><span class="o">**</span><span class="mi">4</span><span class="o">/</span><span class="mi">2</span><span class="o">/</span><span class="n">E</span><span class="o">/</span><span class="n">H</span><span class="o">**</span><span class="mi">3</span><span class="p">))</span>
</pre></div>
</div>
<p>One finds <span class="math">\(w_{FE} = 5.8638\text{e-3}\)</span> against <span class="math">\(w_{beam} = 5.8594\text{e-3}\)</span>
that is a 0.07% difference.</p>
<p>The stress tensor must be projected on an appropriate function space in order to
evaluate pointwise values or export it for Paraview vizualisation. Here we choose
to describe it as a (2D) tensor and project it onto a piecewise constant function
space:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="n">Vsig</span> <span class="o">=</span> <span class="n">TensorFunctionSpace</span><span class="p">(</span><span class="n">mesh</span><span class="p">,</span> <span class="s2">&quot;DG&quot;</span><span class="p">,</span> <span class="n">degree</span><span class="o">=</span><span class="mi">0</span><span class="p">)</span>
<span class="n">sig</span> <span class="o">=</span> <span class="n">Function</span><span class="p">(</span><span class="n">Vsig</span><span class="p">,</span> <span class="n">name</span><span class="o">=</span><span class="s2">&quot;Stress&quot;</span><span class="p">)</span>
<span class="n">sig</span><span class="o">.</span><span class="n">assign</span><span class="p">(</span><span class="n">project</span><span class="p">(</span><span class="n">sigma</span><span class="p">(</span><span class="n">u</span><span class="p">),</span> <span class="n">Vsig</span><span class="p">))</span>
<span class="nb">print</span><span class="p">(</span><span class="s2">&quot;Stress at (0,H):&quot;</span><span class="p">,</span> <span class="n">sig</span><span class="p">(</span><span class="mi">0</span><span class="p">,</span> <span class="n">H</span><span class="p">))</span>
</pre></div>
</div>
<p>Fields can be exported in a suitable format for vizualisation using Paraview.
VTK-based extensions (.pvd,.vtu) are not suited for multiple fields and parallel
writing/reading. Prefered output format is now .xdmf:</p>
<div class="highlight-default"><div class="highlight"><pre><span></span><span class="n">file_results</span> <span class="o">=</span> <span class="n">XDMFFile</span><span class="p">(</span><span class="s2">&quot;elasticity_results.xdmf&quot;</span><span class="p">)</span>
<span class="n">file_results</span><span class="o">.</span><span class="n">parameters</span><span class="p">[</span><span class="s2">&quot;flush_output&quot;</span><span class="p">]</span> <span class="o">=</span> <span class="kc">True</span>
<span class="n">file_results</span><span class="o">.</span><span class="n">parameters</span><span class="p">[</span><span class="s2">&quot;functions_share_mesh&quot;</span><span class="p">]</span> <span class="o">=</span> <span class="kc">True</span>
<span class="n">file_results</span><span class="o">.</span><span class="n">write</span><span class="p">(</span><span class="n">u</span><span class="p">,</span> <span class="mf">0.</span><span class="p">)</span>
<span class="n">file_results</span><span class="o">.</span><span class="n">write</span><span class="p">(</span><span class="n">sig</span><span class="p">,</span> <span class="mf">0.</span><span class="p">)</span>
</pre></div>
</div>
</div>
</div>


          </div>
        </div>
      </div>
        </div>
        <div class="sidebar">
          <h3>Table Of Contents</h3>
          <ul class="current">
<li class="toctree-l1"><a class="reference internal" href="../../intro.html">Introduction</a><ul>
<li class="toctree-l2"><a class="reference internal" href="../../intro.html#what-is-it-about">What is it about ?</a></li>
<li class="toctree-l2"><a class="reference internal" href="../../intro.html#how-do-i-get-started">How do I get started ?</a></li>
<li class="toctree-l2"><a class="reference internal" href="../../intro.html#about-the-author">About the author</a></li>
</ul>
</li>
<li class="toctree-l1 current"><a class="reference internal" href="../../linear_problems.html">Linear problems in solid mechanics</a><ul class="current">
<li class="toctree-l2 current"><a class="current reference internal" href="#">2D linear elasticity</a></li>
<li class="toctree-l2"><a class="reference internal" href="orthotropic_elasticity.py.html">Orthotropic linear elasticity</a></li>
<li class="toctree-l2"><a class="reference internal" href="axisymmetric_elasticity.html">Axisymmetric formulation for elastic structures of revolution</a></li>
<li class="toctree-l2"><a class="reference internal" href="../thermoelasticity/thermoelasticity.html">Linear thermoelasticity (weak coupling)</a></li>
<li class="toctree-l2"><a class="reference internal" href="../thermoelasticity/thermoelasticity_transient.html">Thermo-elastic evolution problem (full coupling)</a></li>
<li class="toctree-l2"><a class="reference internal" href="../modal_analysis_dynamics/cantilever_modal.py.html">Modal analysis of an elastic structure</a></li>
</ul>
</li>
<li class="toctree-l1"><a class="reference internal" href="../../homogenization.html">Homogenization of heterogeneous materials</a><ul>
<li class="toctree-l2"><a class="reference internal" href="../periodic_homog_elas/periodic_homog_elas.html">Periodic homogenization of linear elastic materials</a></li>
</ul>
</li>
<li class="toctree-l1"><a class="reference internal" href="../../nonlinear_problems.html">Nonlinear problems in solid mechanics</a><ul>
Jeremy BLEYER's avatar
Jeremy BLEYER committed
241
<li class="toctree-l2"><a class="reference internal" href="../viscoelasticity/linear_viscoelasticity.html">Linear viscoelasticity</a></li>
Jeremy BLEYER's avatar
Jeremy BLEYER committed
242 243 244
<li class="toctree-l2"><a class="reference internal" href="../2D_plasticity/vonMises_plasticity.py.html">Elasto-plastic analysis of a 2D von Mises material</a></li>
</ul>
</li>
Jeremy BLEYER's avatar
Jeremy BLEYER committed
245 246
<li class="toctree-l1"><a class="reference internal" href="../../beams_and_plates.html">Beams and plates</a><ul>
<li class="toctree-l2"><a class="reference internal" href="../timoshenko/beam_buckling.html">Eulerian buckling of a beam</a></li>
Jeremy BLEYER's avatar
Jeremy BLEYER committed
247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301
<li class="toctree-l2"><a class="reference internal" href="../reissner_mindlin/reissner_mindlin_quads.py.html">Reissner-Mindlin plate with Quadrilaterals</a></li>
<li class="toctree-l2"><a class="reference internal" href="../reissner_mindlin/reissner_mindlin_dg.py.html">Reissner-Mindlin plate with a Discontinuous-Galerkin approach</a></li>
</ul>
</li>
<li class="toctree-l1"><a class="reference internal" href="../../tips_and_tricks.html">Tips and Tricks</a><ul>
<li class="toctree-l2"><a class="reference internal" href="../../tips_and_tricks.html#efficient-projection-on-dg-or-quadrature-spaces">Efficient projection on DG or Quadrature spaces</a></li>
</ul>
</li>
</ul>

          <div role="search">
            <h3 style="margin-top: 1.5em;">Search</h3>
            <form class="search" action="../../search.html" method="get">
                <input type="text" name="q" />
                <input type="submit" value="Go" />
                <input type="hidden" name="check_keywords" value="yes" />
                <input type="hidden" name="area" value="default" />
            </form>
          </div>
        </div>
        <div class="clearer"></div>
      </div>
    </div>

    <div class="footer-wrapper">
      <div class="footer">
        <div class="left">
          <div role="navigation" aria-label="related navigaton">
            <a href="../../linear_problems.html" title="Linear problems in solid mechanics"
              >previous</a> |
            <a href="orthotropic_elasticity.py.html" title="Orthotropic linear elasticity"
              >next</a> |
            <a href="../../genindex.html" title="General Index"
              >index</a>
          </div>
          <div role="note" aria-label="source link">
              <br/>
              <a href="../../_sources/demo/elasticity/2D_elasticity.py.rst.txt"
                rel="nofollow">Show Source</a>
          </div>
        </div>

        <div class="right">
          
    <div class="footer" role="contentinfo">
        &#169; Copyright 2016, Jeremy Bleyer.
      Created using <a href="http://sphinx-doc.org/">Sphinx</a> 1.6.5.
    </div>
        </div>
        <div class="clearer"></div>
      </div>
    </div>

  </body>
</html>