Test Case of Elasto-Magnetism in Stationary and Axisymmetrical Case

1. Introduction

This page presents the simulation and result of page Elastic Equations with Electromagnetism and Thermic in Axisymmetrical case : Elastic, Thermic and Electromagnetism problem coupled by Thermal Dilatation in stationnary and axisymmetrical case.

2. Run the calculation

The command line to run this case is :

    mpirun -np 16 feelpp_toolbox_coefficientformpdes --config-file=elasto-thermo-mag.cfg --cfpdes.gmsh.hsize=5e-3

To compute with only the Thermal Dilatation, please, put :

        "bool_laplace":0,
        "bool_dilatation":1,

On Parameter function of .json file elasto-thermo-mag.json.

3. Data Files

The case data files are available in Github here :

4. Equation

In this subsection, we couple the equation (Static Elasticity Axis) of Elastic equation and (AV Axis) AV-Formulation in axisymmetrical coordinates.

The domain of resolution of electromagnetism part is Ωaxis with bounds Γaxis, ΓaxisD the bound of Dirichlet conditions and ΓaxisN the bound of Neumann conditions such that Γaxis=ΓaxisNΓaxisD.

The domain of resolution of elastic part is ΩaxiscΩaxis (and the domain of definition of electrical potential V and electrical conductivity σ) with bounds Γaxisc, ΓaxisDelas the bound of Dirichlet conditions and ΓaxisNelas the bound of Neumann conditions such that Γaxisc=ΓaxisDelasΓaxisNelas. The domain Ωaxisc correspond to the conductor.

Differential Formulation of Coupled Elasticity Equation and AV-Formulation in Stationary and Axisymmetrical case
(Elasticity-AV Axis Static){1μΔΦ+2μrΦr+σVθ=0 on Ωaxis(AV-1 Axis)2Vθ2=0 on Ωaxisc(AV-2 Axis)ˉˉσ=0on Ωaxisc(Elas-1)Aθ=0 on ΓaxisD(D Axis)Aθnaxis=0 on ΓaxisN(N Axis)u=0on ΓaxisDelas(D elas)ˉˉσn=0on ΓaxisNelas(N elas)

With :

  • Lamé’s coefficients : λ=Ev(12v)(1+v), μ=E2(1+v) with

    • Young modulus E

    • Poisson’s coefficient v

  • Magnetic potential field : A thus :

    • B=×A

    • Φ=rAθ

  • Displacement u=(uxuyuz)

  • Electric Potential Scalar : V

  • Stress Tensor with Hooke Law : ˉˉσ=ˉˉσE+ˉˉσT=λTr(ˉˉϵ)+2μˉˉϵ(3λ+2μ)αT(TT0)

  • Electrical conductivity : σ

And notations :

  • Divergence of tensor : ˉˉσ=(ˉˉσr,:+ˉˉσrrˉˉσθθrˉˉσz,:+ˉˉσrzr)

  • Scalar produce of tensor : ˉˉσn=(ˉˉσr,:naxisˉˉσz,:naxis)=(ˉˉσrrnr+ˉˉσrznzˉˉσzrnr+ˉˉσzznz)

The equation become coupled in one senses, the second equation (Elas-1 Axis) depends of potentials.
We don’t care of geometrical deformation which can change the mesh or the density. We suppose the geometry is fixed.

5. Geometry

The geometry is a torus of the conductor in cartesian coordinates (x,y,z) or rectangle in axisymmetric coordinates (r,z), surrounded by air.

1torus axis(1)
Geometry in axisymmetrical cut loop on Conductor
1torus axis(2)
Geometry in axisymmetrical cut
1torus 3d(1)
Geometry in three dimensions
1torus 3d(2)
Geometry in three dimensions cut loop on Conductor

The geometrical domains are :

  • Conductor : the torus is composed by conductor material

    • Interior : interior surface of conductor ring

    • Exterior : exterior surface of conductor ring

    • Upper : upper of conductor ring

    • Bottom : bottom of conductor ring

  • Air : the air surround Conductor

    • zAxis : a bound of Air, correspond to Oz axis ({(z,r),z=0})

    • infty : the rest of bound of Air

Symbol

Description

value

unit

rint

interior radius of torus

75e3

m

rext

exterior radius of torus

100.2e3

m

z1

half-height of torus

300e3

m

rinfty

radius of infty border

5rext

m

6. Boundary Conditions

We impose the boundary conditions :

  • Magnetic Equation :

    • Strong Dirichlet :

      • On zAxis : Φ=0 (Aθ=0 by symetric argument)

      • On infty : Φ=0 (Aθ=0 we consider the bound of resolution like infty for magnetic field)

  • Heat Equation :

    • Neumann : Tn=0 on Interior and Exterior

    • Robin : kTn=h(TTc) on Upper and Bottom. The Robin condition represents the cooling by water.

  • Elastic Equation :

    • Strong Dirichlet : u=0 on Upper and Bottom. The Dirichlet condition represents the embedding of mechanical part.

    • Neumann : ˉˉσn=0 on Interior and Exterior. The Neumann condition represents the freedom of displacement.

On JSON file, the boundary conditions are writed :

Boundary conditions on JSON file
    "BoundaryConditions":
    {
        "magnetic":
        {
            "Dirichlet":
            {
                "magdir":
                {
                    "markers":["ZAxis","Infty"],
                    "expr":"0"
                }
            }
        },
        "heat":
        {
            "Robin":
            {
                "heatdir":
                {
                    "markers":["Interior","Exterior"],
                    "expr1":"h*x:h:x",
                    "expr2":"h*T_c*x:h:T_c:x"
                }
            }
        },
        "elastic":
        {
            "Dirichlet":
            {
                "elasdir":
                {
                    "markers":["Upper","Bottom"],
                    "expr":"{0,0}"
                }
            }
        }
    }

7. Weak Formulation

We obtain :

Weak Formulation of Elastic Thermic Electromagnetism equations in Axisymmetric and Stationnary Case Vectorial
\tiny{ \text{(Weak Thermo-Mag Axis Static)} \\ \left\{ \begin{eqnarray*} \int_{\Omega^{axis}}{ \frac{1}{\mu} \tilde{\nabla} \Phi \cdot \tilde{\nabla} \mathbf{ϕ} \ r \ drdz } + \int_{\Omega^{axis}}{\frac{2}{\mu \, r} \frac{\partial \Phi}{\partial r} \, \mathbf{ϕ} \ r \ drdz} = - \int_{\Omega_c^{axis}}{ \sigma \frac{\partial V}{\partial \theta} \mathbf{ϕ} \ r \ drdz } & \text{(Weak AV Axis)} \\ \int_{\Omega_c^{axis}}{ k \, \tilde{\nabla} T \cdot \tilde{\nabla} \mathbf{ϕ} \, r \, drdz } + \int_{\Gamma_R^{axis}}{ h \, T \, \mathbf{ϕ} \, r \, d\Gamma } = \int_{\Omega_c^{axis}}{ Q \, \mathbf{ϕ} \, r \, drdz } + \int_{\Gamma_R^{axis}}{ h \, T_c \, \mathbf{ϕ} \, r \, d\Gamma } & \text{(Weak Heat Axis)} \\ \int_{\Omega_c^{axis}} \left( r \left( \lambda \left( \frac{\partial u_r}{\partial r} + \frac{1}{r} u_r + \frac{\partial u_z}{\partial z} \right) - (3\lambda + 2 \mu) \, \alpha_T \, \left( T-T_0 \right) \right) \left( \frac{\partial \xi_r}{\partial r} + \frac{1}{r} \xi_r \right) + r \mu \left[ 2 \left( \frac{\partial u_r}{\partial r} \frac{\partial \xi_r}{\partial r} + \frac{1}{r^2} u_r \xi_r \right) + \left( \frac{\partial u_r}{\partial z} \frac{\partial \xi_r}{\partial z} + \frac{\partial u_z}{\partial r} \frac{\partial \xi_r}{\partial z} \right) \right] \\ + r \mu \left( \frac{\partial u_r}{\partial z} + \frac{\partial u_z}{\partial r} \right) \frac{\partial \xi_z}{\partial r} + r \left[ \lambda \left( \frac{\partial u_r}{\partial r} + \frac{1}{r} u_r + \frac{\partial u_z}{\partial z} \right) + 2\mu \frac{\partial u_z}{\partial z} - (3\lambda + 2 \mu) \, \alpha_T \, \left( T-T_0 \right) \right] \frac{\partial \xi_z}{\partial z} \right) dr dz \\ = \int_{\Gamma_{N \hspace{0.05cm} elas}^{axis}}{ r \, \mathbf{g} \cdot \mathbf{ξ} \ dr dz } + \int_{\Gamma_{R \hspace{0.05cm} elas}^{axis}}{ r \left( \mathbf{g} - k \mathbf{u} \right) \cdot \mathbf{ξ} \ dr dz } & \hspace{0.5cm} \text{(Weak Static Elasticity Axis)} \\ \text{for } \mathbf{ϕ} \in H^1(\Omega)^2, \ \psi \in H^1(\Omega_c) , \ \eta \in H^1(\Omega_c), \ \mathbf{ξ} = \begin{pmatrix} xi_x \\ \xi_y \\ \xi_z \end{pmatrix} \in H^1(\Omega_c)^2 \ \text{and} \ \mathbf{F}_{laplace} = \begin{pmatrix} F_{laplace \, r} \\ 0 \\ F_{laplace \, z} \end{pmatrix} = \begin{pmatrix} - \frac{\sigma}{r^2} \frac{\partial V_{\theta}}{\partial \theta} \frac{\partial \Phi}{\partial r} \\ 0 \\ - \frac{\sigma}{r^2} \frac{\partial V_{\theta}}{\partial \theta} \frac{\partial \Phi}{\partial z} \end{pmatrix}_{cyl} \end{eqnarray*} \right. }

8. Parameters

The parameters of problem are :

  • On Conductor :

Symbol

Description

Value

Unit

V

scalar electrical potential

U \, \frac{\theta}{2\pi}

Volt

U

electrical potential

0.2

Volt

\sigma

electrical conductivity at T=T_0

58e6

S/m

j_{\theta}

orthoradial current density

- \sigma \frac{U}{2\pi} \, \frac{1}{r}

A/m^2

E

Young Modulus

2.1e6

Pa

v

Poisson’s coefficient

0.33

dimensionless

\lambda

Lame’s coefficient

\frac{E \, v}{(1-2v)(1+v)}

Pa

\mu

Lame’s coefficient

\frac{E}{2 (1+v)}

Pa

  • On Air :

Symbol

Description

Value

Unit

\mu=\mu_0

magnetic permeability of vacuum

4\pi.10^{-7}

kg \, m / A^2 / S^2

On JSON file, the parmeters are writed :

Parameters on JSON file
    "Parameters":
    {
        "bool_laplace":0,
        "bool_dilatation":1,

        "h":80000,      // W/m2/K
        "T_c":293,      // K
        "T0":293,       // K

        // Constants of analytical solve
        "a":77.32,      // K
        "b":0.40041,    // K
        "rmax":0.0861910719118454,  // m
        "Tmax":295.85   // K
    }

9. Coefficient Form PDEs

We use the application Coefficient Form PDEs. The coefficient associate to Weak Formulation are :

  • For MQS equation (Weak MQS Axis) :

    • On Conductor :

    Coefficient

    Description

    Expression

    d

    damping or mass coefficient

    \sigma

    c

    diffusion coefficient

    \frac{r}{\mu}

    \beta

    convection coefficient

    \begin{pmatrix} \frac{2}{\mu} \\ 0 \end{pmatrix}

    f

    source term

    j_{\theta} \, r^2

    • On Air :

    Coefficient

    Description

    Expression

    c

    diffusion coefficient

    \frac{r}{\mu}

  • For heat equation, on Conductor (the temperature isn’t computed on Air)

Coefficient

Description

Expression

d

damping or mass coefficient

\rho \, C_p

c

diffusion coefficient

k

f

source term

- \sigma \left( \frac{U}{2\pi \, r} \right)^2

  • For elastic equation, on Conductor (the displacement isn’t computed on Air) :

    Coefficient

    Description

    Expression

    c

    diffusion coefficient

    r \mu

    \gamma

    conservative flux source term

    \begin{pmatrix} - \lambda \left( r \frac{\partial u_r}{\partial r} + u_r + r \frac{\partial u_z}{\partial z} \right) - \mu \, r \, \frac{\partial u_r}{\partial r} - r \sigma_T & - \mu \, r \, \frac{\partial u_z}{\partial r} \\ - \mu \, r \, \frac{\partial u_r}{\partial z} & - \lambda \left( r \, \frac{\partial u_r}{\partial r} + u_r + r \, \frac{\partial u_z}{\partial z} \right) - \mu \, r \, \frac{\partial u_z}{\partial z} - r \sigma_T \end{pmatrix}

    f

    source term

    \begin{pmatrix} r F_r - \lambda \frac{\partial u_r}{\partial r} - (\lambda+2\mu) \frac{u_r}{r} - \lambda \frac{\partial u_z}{\partial z} - \sigma_T \\ r F_z \end{pmatrix}

On JSON file, the coefficients are writed :

CFPDEs coefficients on JSON file
   "Materials":
    {
        "Conductor":
        {
            // Magnetic Part
            "U":0.2,            // Volt
            "sigma":58e+6,      // S.m-1
            "mu":"4*pi*1e-7",   // kg.m/A2/S2

            "j_th":"-sigma*U/2/pi/x:sigma:U:x",

            "magnetic_c":"x/mu:x:mu",
            "magnetic_beta":"{2/mu,0}:mu",
            "magnetic_f":"-sigma*U/2/pi*x:sigma:U:x",

            // Heat Part
            "k":380,            // W/m/K
            "heat_c":"k*x:k:x",
            "heat_f":"sigma*U/(2*pi)*U/(2*pi)/x:sigma:U:x",

            // Elastic Part
            "E":2.1e6,          // Pa
            "v":0.33,           // dimensionless
            "alpha_T":"17e-6",  // K-1

            "Lame_lambda":"E*v/(1-2*v)/(1+v):E:v",
            "Lame_mu":"E/(2*(1+v)):E:v",

            "F_laplace":"{1/x*j_th*magnetic_grad_phi_0, 1/x*j_th*magnetic_grad_phi_1}::x:j_th:magnetic_grad_phi_0:magnetic_grad_phi_0:magnetic_grad_phi_1",
            "sigma_T":"-(3*Lame_lambda+2*Lame_mu)*alpha_T*(heat_T-T0):Lame_lambda:Lame_mu:alpha_T:heat_T:T0",

            "elastic_c":"x*Lame_mu:x:Lame_mu",
            "elastic_gamma":"{-Lame_lambda*(x*elastic_grad_u_00+elastic_u_0+x*elastic_grad_u_11)-x*Lame_mu*elastic_grad_u_00-bool_dilatation*x*sigma_T, -x*Lame_mu*elastic_grad_u_10, -x*Lame_mu*elastic_grad_u_01, -Lame_lambda*(x*elastic_grad_u_00+elastic_u_0+x*elastic_grad_u_11)-x*Lame_mu*elastic_grad_u_11-bool_dilatation*x*sigma_T}:bool_dilatation:x:Lame_lambda:Lame_mu:elastic_u_0:elastic_grad_u_00:elastic_grad_u_01:elastic_grad_u_10:elastic_grad_u_11:sigma_T",
            "elastic_f":"{bool_laplace*x*F_laplace_0-Lame_lambda*elastic_grad_u_00-(Lame_lambda+2*Lame_mu)*elastic_u_0/x-Lame_lambda*elastic_grad_u_11-bool_dilatation*sigma_T, bool_laplace*x*F_laplace_1}:bool_laplace:bool_dilatation:x:F_laplace_0:F_laplace_1:Lame_mu:Lame_lambda:elastic_u_0:elastic_grad_u_00:elastic_grad_u_01:elastic_grad_u_11:sigma_T",

            "sigma_rr":"(Lame_lambda+2*Lame_mu)*elastic_grad_u_00+Lame_lambda*elastic_u_0/x+Lame_lambda*elastic_grad_u_11+bool_dilatation*sigma_T:bool_dilatation:Lame_lambda:Lame_mu:x:elastic_u_0:elastic_grad_u_00:elastic_grad_u_11:elastic_grad_u_22:sigma_T",
            "sigma_thth":"Lame_lambda*elastic_grad_u_00+(Lame_lambda+2*Lame_mu)*Lame_lambda*elastic_u_0/x+Lame_lambda*elastic_grad_u_11+bool_dilatation*sigma_T:bool_dilatation:Lame_lambda:Lame_mu:x:elastic_u_0:elastic_grad_u_00:elastic_grad_u_11:elastic_grad_u_22:sigma_T",
            "sigma_zz":"Lame_lambda*elastic_grad_u_00+Lame_lambda*Lame_lambda*elastic_u_0/x+(Lame_lambda+2*Lame_mu)*elastic_grad_u_11+bool_dilatation*sigma_T:bool_dilatation:Lame_lambda:Lame_mu:x:elastic_u_0:elastic_grad_u_00:elastic_grad_u_11:elastic_grad_u_22:sigma_T",
            "sigma_rz":"Lame_mu*(elastic_grad_u_01+elastic_grad_u_10):Lame_mu:elastic_grad_u_01:elastic_grad_u_10",

            "sigma_1":"Lame_lambda/x*elastic_u_0+(Lame_lambda+Lame_mu)*(elastic_grad_u_00+elastic_grad_u_11)+bool_dilatation*sigma_T+Lame_mu*sqrt((elastic_grad_u_00-elastic_grad_u_11)*(elastic_grad_u_00-elastic_grad_u_11)+4*(elastic_grad_u_01+elastic_grad_u_10)*(elastic_grad_u_01+elastic_grad_u_10)):bool_dilatation:x:Lame_lambda:Lame_mu:elastic_u_0:elastic_grad_u_00:elastic_grad_u_11:elastic_grad_u_01:elastic_grad_u_10:sigma_T",
            "sigma_2":"Lame_lambda*elastic_grad_u_00+(Lame_lambda+2*Lame_mu)/x*elastic_u_0+Lame_lambda*elastic_grad_u_11+bool_dilatation*sigma_T:bool_dilatation:x:Lame_lambda:Lame_mu:elastic_u_0:elastic_grad_u_00:elastic_grad_u_11:sigma_T",
            "sigma_3":"Lame_lambda/x*elastic_u_0+(Lame_lambda+Lame_mu)*(elastic_grad_u_00+elastic_grad_u_11)+bool_dilatation*sigma_T-Lame_mu*sqrt((elastic_grad_u_00-elastic_grad_u_11)*(elastic_grad_u_00-elastic_grad_u_11)+4*(elastic_grad_u_01+elastic_grad_u_10)*(elastic_grad_u_01+elastic_grad_u_10)):bool_dilatation:x:Lame_lambda:Lame_mu:elastic_u_0:elastic_grad_u_00:elastic_grad_u_11:elastic_grad_u_01:elastic_grad_u_10:sigma_T"
        },
        "Air":
        {
            "physics":"magnetic",
            "mu":"4*pi*1e-7",   // kg.m/A2/s2
            "magnetic_c":"x/mu:x:mu",
            "magnetic_beta":"{2/mu,0}:mu"
        }
    }

10. Numeric Parameters

  • Mesh size :

    • Interior of torus : 0.01 m

    • Far from the torus : 0.2 m

  • Element Type : I run the simulation in P_2 elements.

1torus axis mesh
Mesh of Geometry

11. Result

11.1. Displacements

deformation
Deformation due to Thermal Dilatation Force : The white solid is initial geometry and the coloured solid is the final geometry (the colours represent Magnitude of displacement u (m))

The Thermal Dilatation deform the cylinder like a barrel filled with water. The cylinder deforms outwards.

We can note that near of embeded (Upper and Bottom), the cylinder deforms inwards.

ur
u_r(m)
uz
u_z(m)

The radial displacement u_r is maximum on Or axis near of the interior and decreases to zeros near of the embedded (on Upper and Bottom).

The z-displacement is symmetric by Or axis. Near of the interior, it positive in bottom part and negativ to upper part, the cylinder compresses toward the Or axis in z-displacement.

11.2. Tensor of Deformation \bar{\bar{\epsilon}}

epsilon rr
\bar{\bar{\epsilon}}_{rr} = \frac{\partial u_r}{\partial z}
epsilon thth
\bar{\bar{\epsilon}}_{\theta \theta} = \frac{1}{r} \, u_r
epsilon zz
\bar{\bar{\epsilon}}_{zz} = \frac{\partial u_z}{\partial z}
epsilon rz
\bar{\bar{\epsilon}}_{rz} = \frac{1}{2} \left( \frac{\partial u_r}{\partial z} + \frac{\partial u_z}{\partial r} \right)

The radial deformation is uniform by z. It maximum on middle radius (\frac{r_{int}+r_{ext}}{2} \approx 87 mm)

The orthoradial deformation is the greatest contribution of deformation. It maximum in Or axis near of the interior and decrease to negative value near of embedding (Upper and Bottom).

The z-deformation near of the center is symetric by Oz axis. The cylinder compresses near of the interior and it expands near of the exterior.

The rz-shear is great around the embededding and null in the rest of cylinder.

11.3. Stress Tensor \bar{\bar{\sigma}}

sigma rr
\scriptsize{\bar{\bar{\sigma}}_{rr} = \frac{\lambda}{r} \, u_r + \left( \lambda + 2\mu \right) \, \frac{\partial u_r}{\partial r} + \lambda \, \frac{\partial u_z}{\partial z} - (3\lambda + 2\mu) \, \alpha_T \, \left(T-T_0\right) (Pa)}
sigma thth
\scriptsize{\bar{\bar{\sigma}}_{\theta \theta} = \frac{\lambda + 2\mu}{r} \, u_r + \lambda \, \frac{\partial u_r}{\partial r} + \lambda \, \frac{\partial u_z}{\partial z}- (3\lambda + 2\mu) \, \alpha_T \, \left(T-T_0\right) (Pa)}
sigma zz
\scriptsize{\bar{\bar{\sigma}}_{zz} = \frac{\lambda}{r} \, u_r + \lambda \, \frac{\partial u_r}{\partial r} + \left( \lambda + 2\mu \right) \, \frac{\partial u_z}{\partial z} - (3\lambda + 2\mu) \, \alpha_T \, \left( T-T_0 \right) (Pa)}
sigma rz
\bar{\bar{\sigma}}_{rz} = \mu \, \left( \frac{\partial u_r}{\partial z} + \frac{\partial u_z}{\partial r} \right) (Pa)

The radial stress is uniform by the Oz axis and have a form. The minimum is reached around the middle.

The orthoradial stress is the greatest contribution of deformation. The orthoradial stress is uniform by the Oz axis and have a form of an upturned bell. The minimum is reached around the middle.

The z-deformation near of the center is symetric by Oz axis. The cylinder compresses near of the interior and it expands near of the exterior.

The rz-shear is great around the embededding and null in the rest of cylinder.

11.4. Other Tensor

11.4.1. Criterion of Von Mises \bar{\bar{\sigma}}_{vm}

The tensor of Von Mises permits to say if the matters is in elastic state or not. If the tensor exceeds an matter’s parameter, the matter pass to elastic state to plastic state.

\begin{eqnarray*} \bar{\bar{\sigma}}_{Von \, Misses} &=& \sqrt{3 \, J_2} \\ &=& \sqrt{ \frac{13}{2} \, \left( \bar{\bar{\sigma}}_1^2+\bar{\bar{\sigma}}_2^2+\bar{\bar{\sigma}}_3^2-\bar{\bar{\sigma}}_1 \bar{\bar{\sigma}}_2 - \bar{\bar{\sigma}}_1 \bar{\bar{\sigma}}_3 - \bar{\bar{\sigma}}_2 \bar{\bar{\sigma}}_3 \right) } \end{eqnarray*}

The \sigma_i are the principles tensors, the eigenvalues of stress tensor \bar{\bar{\sigma}}.

sigma vm Or
\bar{\bar{\sigma}}_{vm} (Pa)
sigma vm
\bar{\bar{\sigma}}_{vm} (Pa)

11.4.2. Criterion of Tresca \bar{\bar{\sigma}}_{tr}

The tensor of Tresca permits to say if the matters is in elastic state or not. If the tensor exceeds an matter’s parameter, the matter pass to elastic state to plastic state.

\bar{\bar{\sigma}}_{tr} = max \left( |\bar{\bar{\sigma}}_1-\bar{\bar{\sigma}}_2|, \, |\bar{\bar{\sigma}}_1-\bar{\bar{\sigma}}_3|, \, |\bar{\bar{\sigma}}_2-\bar{\bar{\sigma}}_3| \right)

The \sigma_i are the principles tensors, the eigenvalues of stress tensor \bar{\bar{\sigma}}.

sigma tr Or
\bar{\bar{\sigma}}_{tr} (Pa)
sigma tr
\bar{\bar{\sigma}}_{tr} (Pa)