SFEMaNS  version 4.1 (work in progress)
Reference documentation for SFEMaNS
 All Classes Files Functions Variables Groups Pages
Test 34: Ferrohydrodynamic periodic

Introduction

In this example, we check the correctness of SFEMaNS for a ferrohydrodynamic periodic in \(z\) problem. Maxwell, Navier-Stokes and temperature equations are solved. The domain presents a solid and a fluid subdomain. Maxwell equations are reduced to the magnetoststic equations due to the quasi-static approximation for the magnetic field. The fluid and solid regions do not have the same thermal properties. The volumetric heat capacity and the thermal conductivity are different. This case is usefull to study the cooling by convection and magnetoconvection of a solid by a ferrofluid, when a periodic assumption is used. The other difference with Test 33 is the non null pressure.

The domain of computation is \(\Omega= \{ (r,\theta,z) \in {R}^3 : (r,\theta,z) \in (0,1) \times [0,2\pi) \times (0,1)\} \). We note \(\Gamma= \partial \Omega\). It is composed of a solid and a fluid subdomain:

\begin{align*} \overline{\Omega} = \overline{\Omega_s} \cup \overline{\Omega_f}.\\ \end{align*}

The subdomains are defined followingly: \(\Omega_s = \{ (r,\theta,z) \in {R}^3 : (r,\theta,z) \in (0,1/2) \times [0,2\pi) \times (0,1)\} \) and \(\Omega_f= \{ (r,\theta,z) \in {R}^3 : (r,\theta,z) \in (1/2,1) \times [0,2\pi) \times (0,1)\} \). We note \(\Gamma_f = \partial \Omega_f\).

We solve the temperature equations:

\begin{align*} \begin{cases} c\partial_t T+ c\tilde{\bu} \cdot \GRAD T - \DIV (\lambda \GRAD T) &= f_T, \\ T_{|\Gamma} &= T_\text{bdy} , \\ T_{|t=0} &= T_0, \end{cases} \end{align*}

in the domain \(\Omega\). The extended velocity is defined by \(\tilde{\bu} = \bu\) in \(\Omega_f\) and \(0\) in \(\Omega_s\). The volumetric heat capacity \(c\) and the thermal conductivity \(\lambda\) are piecewise constant functions of space: \(c = c_f\) in \(\Omega_f\) and \(c_s\) in \(\Omega_s\), with \(c_f = 2 c_s\); \(\lambda = \lambda_f\) in \(\Omega_f\) and \(\lambda_s\) in \(\Omega_s\), with \(\lambda_s = 10 \lambda_f\). We recall that these parameters are dimensionless.

We solve the Navier-Stokes equations:

\begin{align*} \begin{cases} \partial_t\bu+\left(\ROT\bu\right)\CROSS\bu - \frac{1}{\Re}\LAP \bu +\GRAD p &= \alpha T \textbf{e}_z + \chi(T) \GRAD \left( \frac{H^2}{2} \right) + \bef, \\ \DIV \bu &= 0, \\ \bu_{|\Gamma_f} &= \bu_{\text{bdy}} , \\ \bu_{|t=0} &= \bu_0, \\ p_{|t=0} &= p_0, \end{cases} \end{align*}

in the domain \(\Omega_f\). We denote by \(\textbf{e}_z\) the unit vector in the vertical direction. The function of the temperature \(\chi\) in the Kelvin force is defined by

\begin{align*} \chi(T) = - T^2 \end{align*}

We solve the magnetostatics equations:

\begin{align*} \begin{cases} \ROT \bH &= \bJ,\\ \DIV (\mu \bH) &= 0, \\ \bH \CROSS \bn &= \bH_{\text{bdy}} \CROSS \bn, \text{ on } \Gamma,\\ \bH_{|t=0} &= \bH_0, \\ \end{cases} \end{align*}

in the domain \(\Omega\). The permeability \(\mu\) is constant in the whole domain. Solving magnetostatic equations is possible by using a very small electrical conductivity in the whole domain (see data file description). The other terms of the induction equation are thus neglectable in that case.

The data are the source terms \(f_T\), \(\bef\) and \(\bJ\), the boundary data \(T_\text{bdy}\), \(\bu_{\text{bdy}}\) and \(\bH_{\text{bdy}}\), the initial data \(T_0\), \(\bu_0\), \(p_0\) and \(\bH_0\). The parameters are the volumetric heat capacities \(c_s\) and \(c_f\), the thermal conductivities \(\lambda_s\) and \(\lambda_f\), the kinetic Reynolds number \(\Re\), the thermal gravity number \(\alpha\) and the permeability \(\mu\).

Manufactured solutions

We approximate the following analytical solutions:

\begin{align*} T(r,\theta,z,t) & = \lambda^{-1}r^2(r-r_0)\sin(2\pi z)(1+\cos(\theta))\cos(t), \\ u_r(r,\theta,z,t) &= -2\pi(r-r_0)^2\cos(2\pi z)(1+\cos(\theta))\cos(t), \\ u_{\theta}(r,\theta,z,t) &= 2\pi(r-r_0)^2\cos(2\pi z)(1+\cos(\theta))\cos(t), \\ u_z(r,\theta,z,t) &= \frac{r-r_0}{r} \sin(2\pi z) ((3r-r_0)(1+\cos(\theta))+(r-r_0)\sin(\theta))\cos(t), \\ p(r,\theta,z,t) &= r^3 \sin(2\pi z) \cos(\theta) \cos(t), \\ H_r(r,\theta,z,t) &= 2\pi r^3 \sin(2\pi z) (1 + \sin(\theta)) \cos(t), \\ H_\theta(r,\theta,z,t) &= - 2\pi r^3 \sin(2\pi z) (1 + \sin(\theta)) \cos(t), \\ H_z(r,\theta,z,t) &= - r^2 \cos(2\pi z) (4 - \cos(\theta) + 4 \sin(\theta)) \cos(t) \end{align*}

where \(r_0 = 1/2\) is the limit between solid and fluid regions. The velocity is the curl of a vector field, it is thus divergence free. Same thing for the magnetic field.

The source terms \(f_T, \bef\), \(\bJ\) and the boundary data \(T_\text{bdy}, \bu_{\text{bdy}}\), \(\bH_{\text{bdy}}\) are computed accordingly.

Generation of the mesh

The finite element mesh used for this test is named SOLID_FLUID_10.FEM and has a mesh size of \(0.1\) for the P1 approximation. You can generate this mesh with the files in the following directory: ($SFEMaNS_MESH_GEN_DIR)/EXAMPLES/EXAMPLES_MANUFACTURED_SOLUTIONS/SOLID_FLUID_10. The following image shows the mesh for P1 finite elements.

fig_SOLID_FLUID_10.png
Finite element mesh (P1).

Information on the file condlim.f90

The initial conditions, boundary conditions and the forcing terms are set in the file condlim_test_34.f90. Here is a description of the subroutines and functions of interest.

  1. The subroutine init_velocity_pressure initializes the velocity field and the pressure at the times \(-dt\) and \(0\) with \(dt\) being the time step. This is done by using the functions vv_exact and pp_exact as follows:
    time = 0.d0
    DO i= 1, SIZE(list_mode)
    mode = list_mode(i)
    DO j = 1, 6
    !===velocity
    un_m1(:,j,i) = vv_exact(j,mesh_f%rr,mode,time-dt)
    un (:,j,i) = vv_exact(j,mesh_f%rr,mode,time)
    END DO
    DO j = 1, 2
    !===pressure
    pn_m2(:) = pp_exact(j,mesh_c%rr,mode,time-2*dt)
    pn_m1 (:,j,i) = pp_exact(j,mesh_c%rr,mode,time-dt)
    pn (:,j,i) = pp_exact(j,mesh_c%rr,mode,time)
    phin_m1(:,j,i) = pn_m1(:,j,i) - pn_m2(:)
    phin (:,j,i) = Pn (:,j,i) - pn_m1(:,j,i)
    ENDDO
    ENDDO
  2. The subroutine init_temperature initializes the temperature at the times \(-dt\) and \(0\) with \(dt\) the time step. This is done by using the function temperature_exact as follows:
    time = 0.d0
    DO i= 1, SIZE(list_mode)
    mode = list_mode(i)
    DO j = 1, 2
    tempn_m1(:,j,i) = temperature_exact(j, mesh%rr, mode, time-dt)
    tempn (:,j,i) = temperature_exact(j, mesh%rr, mode, time)
    ENDDO
    ENDDO
  3. The function vv_exact contains the analytical velocity field. It is used to initialize the velocity field and to impose Dirichlet boundary conditions on the velocity field.
    1. The limit between fluid and solid region and \(\pi\) are defined:
      REAL(KIND=8) :: r0 = 0.5d0, pi = acos(-1.d0)
    2. We construct the radial and vertical coordinates r, z.
      r = rr(1,:)
      z = rr(2,:)
    3. We define the velocity field depending of its TYPE (1 and 2 for the component radial cosine and sine, 3 and 4 for the component azimuthal cosine and sine, 5 and 6 for the component vertical cosine and sine) and of its mode m as follows:
      IF (TYPE==1) THEN
      IF ((m==0).OR.(m==1)) THEN
      vv = -2*pi*(r-r0)**2*cos(t)*cos(2*pi*z)
      ELSE
      vv = 0.d0
      END IF
      ELSE IF (TYPE==3) THEN
      IF ((m==0).OR.(m==1)) THEN
      vv = 2*pi*(r-r0)**2*cos(t)*cos(2*pi*z)
      ELSE
      vv = 0.d0
      END IF
      ELSE IF (TYPE==5) THEN
      IF ((m==0).OR.(m==1)) THEN
      vv = (r-r0)*cos(t)*sin(2*pi*z)/r * (3*r-r0)
      ELSE
      vv = 0.d0
      END IF
      ELSE IF (TYPE==6) THEN
      IF (m==1) THEN
      vv = (r-r0)*cos(t)*sin(2*pi*z)/r * (r-r0)
      ELSE
      vv = 0.d0
      END IF
      ELSE
      vv = 0.d0
      END IF
      where \(t\) is the time. It is important to specify the null types or modes to avoid nonsense results.
  4. The function pp_exact contains the analytical pressure. It is used to initialize the pressure.
    1. \(\pi\) is defined:
      REAL(KIND=8) :: pi = acos(-1.d0)
    2. We construct the radial and vertical coordinates r, z.
      r = rr(1,:)
      z = rr(2,:)
    3. We define the pressure depending on its TYPE (1 and 2 for cosine and sine) and on its mode as follows:
      IF ((TYPE==1).AND.(m==1)) THEN
      vv = r**3*sin(2*pi*z)*cos(t)
      ELSE
      vv = 0.d0
      END IF
  5. The function temperature_exact contains the analytical temperature. It is used to initialize the temperature and to impose Dirichlet boundary condition on the temperature.
    1. An array for the conductivity \(\lambda\) at each node is declared:
      REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv, r, z, lambda
    2. The limit between fluid and solid region and \(\pi\) are defined:
      REAL(KIND=8) :: r0=0.5d0, pi = acos(-1.d0)
    3. We construct the radial and vertical coordinates r, z.
      r = rr(1,:)
      z = rr(2,:)
    4. We construct the conductivity array.
      DO i=1,SIZE(rr,2)
      IF (r(i).LE.r0) THEN
      lambda(i) = inputs%temperature_diffusivity(1)
      ELSE
      lambda(i) = inputs%temperature_diffusivity(2)
      END IF
      END DO
    5. We define the temperature depending on its TYPE (1 and 2 for cosine and sine) and on its mode as follows:
      IF ((TYPE==1).AND.((m==0).OR.(m==1))) THEN
      vv = r**2*(r-r0)*sin(2*pi*z)*cos(t) / lambda
      ELSE
      vv = 0.d0
      END IF
  6. The function Hexact contains the analytical magnetic field. It is used to initialize the magnetic field and enforce the boundary conditions.
    1. The limit between fluid and solid region and \(\pi\) are defined:
      REAL(KIND=8) :: r0 = 0.5d0, pi = acos(-1.d0)
    2. We construct the radial and vertical coordinates r, z.
      r = rr(1,:)
      z = rr(2,:)
    3. We define the magnetic field depending on its TYPE (1 and 2 for the component radial cosine and sine, 3 and 4 for the component azimuthal cosine and sine, 5 and 6 for the component vertical cosine and sine) and on its mode as follows:
      IF (TYPE == 1) THEN
      IF (m == 0) THEN
      vv = 2 * pi * r**3 * sin(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      ELSE IF (TYPE == 2) THEN
      IF (m == 1) THEN
      vv = 2 * pi * r**3 * sin(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      ELSE IF (TYPE == 3) THEN
      IF (m == 0) THEN
      vv = - 2 * pi * r**3 * sin(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      ELSE IF (TYPE == 4) THEN
      IF (m == 1) THEN
      vv = - 2 * pi * r**3 * sin(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      ELSE IF (TYPE == 5) THEN
      IF (m == 0) THEN
      vv = 4 * r**2 * cos(2*pi*z) * cos(t)
      ELSE IF (m == 1) THEN
      vv = - r**2 * cos(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      ELSE
      IF (m == 1) THEN
      vv = 4 * r**2 * cos(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      END IF
  7. The function source_in_temperature computes the source term \(f_T\) of the temperature equations.
    1. Arrays c and lambda for the value of the volumetric heat capacity and thermal conductivity at each node must be declared.
      REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv, r, z, c, lambda
    2. The limit between fluid and solid region and \(\pi\) are defined:
      REAL(KIND=8) :: r0 = 0.5d0, pi = acos(-1.d0)
    3. We construct the radial and vertical coordinates r, z.
      r = rr(1,:)
      z = rr(2,:)
    4. The c array is filled based on the data file. The solid volumetric heat capacity is used for the region \(r \le r_0\) and the fluid volumetric heat capacity is used in the fluid region \(r > r_0\).
      DO i=1,SIZE(rr,2)
      IF (r(i).LE.r0) THEN
      c(i) = inputs%vol_heat_capacity(1)
      ELSE
      c(i) = inputs%vol_heat_capacity(2)
      END IF
      END DO
    5. The lambda array is filled based on the data file. The solid conductivity is used for the region \(r \le r_0\) and the fluid conductivity is used in the fluid region \(r > r_0\). The temperature_diffusivity input contains the thermal conductivity when a volumetric heat capacity is used.
      DO i=1,SIZE(rr,2)
      IF (r(i).LE.r0) THEN
      lambda(i) = inputs%temperature_diffusivity(1)
      ELSE
      lambda(i) = inputs%temperature_diffusivity(2)
      END IF
      END DO
    6. The source term \( f_T = c\partial_t T+ c\tilde{\bu} \cdot \GRAD T - \DIV(\lambda\GRAD T) \) is defined in two parts. Firstly, we define the part \( c\partial_t T - \DIV (\lambda\GRAD T) \):
      IF (TYPE==1) THEN
      IF (m==0) THEN
      vv = ((-9*r + 4*pi**2*r**3 + 4*r0 - 4*pi**2*r**2*r0) * lambda * Cos(t) + &
      c * r**2 * (-r + r0) * Sin(t)) * Sin(2*pi*z) / lambda
      ELSE IF (m==1) THEN
      vv = - ((-3*r0 + 4*r*(2.d0 + pi**2*r*(-r + r0))) * lambda * Cos(t) + &
      c * r**2 * (r - r0) * Sin(t)) * Sin(2*pi*z) / lambda
      ELSE
      vv = 0.d0
      END IF
      ELSE
      vv = 0.d0
      END IF
      Secondly, we add the part \( c\tilde{\bu} \cdot \GRAD T \), which is different from 0 only in the fluid:
      IF (TYPE==1) THEN
      IF (m==0) THEN
      DO i=1,size(rr,2)
      IF (r(i)>r0) THEN
      vv(i) = vv(i) + (3*c(i) * pi * r(i) * (r(i) - r0)**2 * r0 &
      * Cos(t)**2 * Sin(4*pi*z(i))) / (2 * lambda(i))
      END IF
      END DO
      ELSE IF (m==1) THEN
      DO i=1,size(rr,2)
      IF (r(i)>r0) THEN
      vv(i) = vv(i) + (2*c(i) * pi * r(i) * (r(i) - r0)**2 * r0 &
      * Cos(t)**2 * Sin(4*pi*z(i))) / lambda(i)
      END IF
      END DO
      ELSE IF (m==2) THEN
      DO i=1,size(rr,2)
      IF (r(i)>r0) THEN
      vv(i) = vv(i) + (c(i) * pi * r(i) * (r(i) - r0)**2 * r0 &
      * Cos(t)**2 * Sin(4*pi*z(i))) / (2 *lambda(i))
      END IF
      END DO
      END IF
      END IF
      In the second part, the array is filled cell by cell because we have to test if the node is in the fluid region.
  8. The function source_in_NS_momentum computes the source term \(\alpha T \textbf{e}_z+\bef\) of the Navier-Stokes equations.
    1. The coordinates \(r\) and \(z\), and the conductivity \(\lambda\) are declared. The conductivity is required to compute the source term corresponding to the Kelvin force.
      REAL(KIND=8), DIMENSION(SIZE(rr,2)) :: vv, r, z, lambda
    2. The coefficient \(\alpha\) of the Boussinesq force is declared and the limit between fluid and solid region and \(\pi\) are defined:
      REAL(KIND=8) :: alpha, r0 = 0.5d0, pi = acos(-1.d0)
    3. The coefficient \(\alpha\) is defined based on the data file.
      alpha = inputs%gravity_coefficient
    4. We construct the radial and vertical coordinates r, z.
      r = rr(1,:)
      z = rr(2,:)
    5. We construct the conductivity array.
      DO j=1,SIZE(rr,2)
      IF (r(j).LE.r0) THEN
      lambda(j) = inputs%temperature_diffusivity(1)
      ELSE
      lambda(j) = inputs%temperature_diffusivity(2)
      END IF
      END DO
    6. We construct the first part of the source term containing the Boussinesq force \(\alpha T \be_z\) and the piece of \(\bef\) that cancels it.
      IF (TYPE==5) THEN
      vv = alpha*(opt_tempn(:,1,i) - temperature_exact(1,rr,mode,time))
      ELSE IF (TYPE==6) THEN
      vv = alpha*(opt_tempn(:,2,i) - temperature_exact(2,rr,mode,time))
      ELSE
      vv = 0.d0
      END IF
      The array opt_temp is an input argument and produces the Boussinesq force. The temperature_exact function is defined in the same file and leads to the cancelling term of \(\bef\).
    7. The part of the source term \( \partial_t\bu+\left(\ROT\bu\right)\CROSS\bu - \frac{1}{\Re}\LAP \bu +\GRAD p\) is then added. It depends on the TYPE (1-6) and the mode (0-2).
      IF (TYPE==1) THEN
      IF (mode==0) THEN
      vv = vv -(4*pi*r*(4*pi**2*r**4 - 8*pi**2*r**3*r0 + r0**2 + r**2*(-3 + 4*pi**2*r0**2))*Cos(time)*Cos(2*pi*z) + &
      (r - r0)*Re*Cos(time)**2*(14*r**3 - 5*r**2*r0 - 5*r*r0**2 + 2*r0**3 + &
      (36*pi**2*r**5 - 84*pi**2*r**4*r0 + 5*r*r0**2 - 2*r0**3 + 2*r**3*(-7 + 30*pi**2*r0**2) + &
      r**2*(5*r0 - 12*pi**2*r0**3))*Cos(4*pi*z)) - &
      4*pi*r**3*(r - r0)**2*Re*Cos(2*pi*z)*Sin(time))/(2.*r**3*Re)
      ELSE IF (mode==1) THEN
      vv = vv -2*(2*pi*r*(2*pi**2*r**4 - r*r0 - 4*pi**2*r**3*r0 + r0**2 + r**2*(-1 + 2*pi**2*r0**2))*Cos(time)*Cos(2*pi*z) + &
      ((3*r**2 - 4*r*r0 + r0**2)*Re*Cos(time)**2*(3*r**2 - r0**2 + &
      (8*pi**2*r**4 - 16*pi**2*r**3*r0 + r0**2 + r**2*(-3 + 8*pi**2*r0**2))*Cos(4*pi*z)))/2. - &
      pi*r**3*(r - r0)**2*Re*Cos(2*pi*z)*Sin(time))/(r**3*Re)
      ELSE IF (mode==2) THEN
      vv = vv -((r - r0)*Cos(time)**2*(4*r**2 - r*r0 - r0**2 + (12*pi**2*r**4 - &
      28*pi**2*r**3*r0 + r0**2 + 4*r**2*(-1.d0 + 5*pi**2*r0**2) + &
      r*(r0 - 4*pi**2*r0**3))*Cos(4*pi*z)))/(2.*r**2)
      END IF
      ELSE IF (TYPE==2) THEN
      IF (mode==1) THEN
      vv = vv + ((r - r0)**2*Cos(time)*(-4*pi*r*Cos(2*pi*z) + Re*Cos(time)*(4*pi**2*r**4 - r*r0 - 8*pi**2*r**3*r0 + r0**2 + &
      r**2*(-3.d0 + 4*pi**2*r0**2) + (3*r**2 + r*r0 - r0**2)*Cos(4*pi*z))))/(r**3*Re)
      ELSE IF (mode==2) THEN
      vv = vv + ((r - r0)**2*Cos(time)**2*(4*pi**2*r**4 - r*r0 - 8*pi**2*r**3*r0 + r0**2 + &
      r**2*(-3.d0 + 4*pi**2*r0**2) + (3*r**2 + r*r0 - r0**2)*Cos(4*pi*z)))/(2.*r**3)
      END IF
      ELSE IF (TYPE==3) THEN
      IF (mode==0) THEN
      vv = vv + (2*pi*(-3*pi*r*(r - r0)**3*(3*r - r0)*Re*Cos(time)**2 + &
      (4*pi**2*r**4 - 8*pi**2*r**3*r0 + r0**2 + r**2*(-3.d0 + 4*pi**2*r0**2))*Cos(time)*Cos(2*pi*z) - &
      r**2*(r - r0)**2*Re*Cos(2*pi*z)*Sin(time)))/(r**2*Re)
      ELSE IF (mode==1) THEN
      vv = vv + (4*pi*r*(2*pi**2*r**4 - r*r0 - 4*pi**2*r**3*r0 + r0**2 + &
      r**2*(-1 + 2*pi**2*r0**2))*Cos(time)*Cos(2*pi*z) + &
      ((r - r0)**3*(3*r - r0)*Re*Cos(time)**2*(-1.d0 - 16*pi**2*r**2 + Cos(4*pi*z)))/2. - &
      2*pi*r**3*(r - r0)**2*Re*Cos(2*pi*z)*Sin(time))/(r**3*Re)
      ELSE IF (mode==2) THEN
      vv = vv + ((r - r0)**3*(3*r - r0)*Cos(time)**2*(-1.d0 - 4*pi**2*r**2 + Cos(4*pi*z)))/(2.*r**3)
      END IF
      ELSE IF (TYPE==4) THEN
      IF (mode==1) THEN
      vv = vv + ((r - r0)**2*Cos(time)*(-8*pi*r*Cos(2*pi*z) + Re*Cos(time)*((-3*r + r0)**2 + &
      (8*pi**2*r**4 + 6*r*r0 - 16*pi**2*r**3*r0 - r0**2 + r**2*(-9.d0 + 8*pi**2*r0**2))*Cos(4*pi*z))))/(2.*r**3*Re)
      ELSE IF (mode==2) THEN
      vv = vv + ((r - r0)**2*Cos(time)**2*(2*r - r0 + (2*r*(-1.d0 + pi*(r - r0))*(1 + pi*(r - r0)) + r0)*Cos(4*pi*z)))/r**2
      END IF
      ELSE IF (TYPE==5) THEN
      IF (mode==0) THEN
      vv = vv + ((4*(12*pi**2*r**4 - 16*pi**2*r**3*r0 - r0**2 + &
      r**2*(-3.d0 + 4*pi**2*r0**2))*Cos(time)*Sin(2*pi*z) - &
      4*r**2*(3*r**2 - 4*r*r0 + r0**2)*Re*Sin(time)*Sin(2*pi*z) + &
      pi*r*(r - r0)**2*(12*pi**2*r**4 - r*r0 - 24*pi**2*r**3*r0 + 2*r0**2 + &
      4*r**2*(-1.d0 + 3*pi**2*r0**2))*Re*4*Cos(time)**2*Sin(4*pi*z)))/(4.*r**3*Re)
      ELSE IF (mode==1) THEN
      vv = vv + (4*(-r0 + pi**2*r*(3*r**2 - 4*r*r0 + r0**2))*Cos(time)*Sin(2*pi*z) - &
      r*(3*r**2 - 4*r*r0 + r0**2)*Re*Sin(time)*Sin(2*pi*z) + &
      pi*(r - r0)**2*(16*pi**2*r**4 - 2*r*r0 - 32*pi**2*r**3*r0 + 3*r0**2 + &
      r**2*(-5.d0 + 16*pi**2*r0**2))*Re*Cos(time)**2*Sin(4*pi*z))/(r**2*Re)
      ELSE IF (mode==2) THEN
      vv = vv + (pi*(r - r0)**2*(4*pi**2*r**4 - r*r0 - 8*pi**2*r**3*r0 + r0**2 + &
      r**2*(-1.d0 + 4*pi**2*r0**2))*Cos(time)**2*Sin(4*pi*z))/r**2
      END IF
      ELSE IF (TYPE==6) THEN
      IF (mode==1) THEN
      vv = vv + (2*(2*pi**2*r*(r - r0)**2 - r0)*Cos(time)*Sin(2*pi*z) - r*(r - r0)**2*Re*Sin(time)*Sin(2*pi*z) -&
      4*pi*r*(r - r0)**3*Re*Cos(time)**2*Sin(4*pi*z))/(r**2*Re)
      ELSE IF (mode==2) THEN
      vv = vv -2*pi*(r - r0)**3*Cos(time)**2*Sin(4*pi*z)/r
      END IF
      END IF
      IF ((TYPE==1).AND.(mode==1)) THEN
      vv = vv + 3*r**2*cos(time)*sin(2*pi*z)
      ELSE IF ((TYPE==4).AND.(mode==1)) THEN
      vv = vv - r**2*cos(time)*sin(2*pi*z)
      ELSE IF ((TYPE==5).AND.(mode==1)) THEN
      vv = vv + 2*pi*r**3*cos(time)*cos(2*pi*z)
      END IF
    8. The last part of the source term \(- \chi(T) \GRAD \left( \frac{H^2}{2} \right)\) is then added. It depends on the TYPE (1-6) and the mode (0-4).
      IF (TYPE == 1) THEN
      IF (mode == 0) THEN
      vv = vv + r**7 * (r - r0)**2 * cos(time)**4 * (215.d0 + 204 * pi**2 * r**2 + &
      (215.d0 - 204 * pi**2 * r**2) * cos(4*pi * z)) * sin(2*pi*z)**2 / &
      (8 * lambda**2)
      ELSE IF (mode == 1) THEN
      vv = vv - 5 * r**7 * (r - r0)**2 * cos(time)**4 * (-11.d0 - 12 * pi**2 * r**2 + &
      (-11.d0 + 12 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      (2 * lambda**2)
      ELSE IF (mode == 2) THEN
      vv = vv - 7 * r**7 * (r - r0)**2 * cos(time)**4 * sin(4*pi*z)**2 / &
      (2 * lambda**2)
      ELSE IF (mode == 3) THEN
      vv = vv + r**7 * (r - r0)**2 * cos(time)**4 * (-19.d0 - 12 * pi**2 *r**2 + &
      (-19.d0 + 12 * pi**2 *r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      (2 * lambda**2)
      ELSE IF (mode == 4) THEN
      vv = vv + 3 * r**7 * (r - r0)**2 * cos(time)**4 * (-5.d0 - 4 * pi**2 * r**2 + &
      (-5.d0 + 4 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      (8 * lambda**2)
      END IF
      ELSE IF (TYPE == 2) THEN
      IF (mode == 1) THEN
      vv = vv - 6 * r**7 * (r - r0)**2 * cos(time)**4 * (-6.d0 - 5 * pi**2 * r**2 + &
      (-6.d0 + 5.d0 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      lambda**2
      ELSE IF (mode == 2) THEN
      vv = vv + 2 * r**7 * (r - r0)**2 * cos(time)**4 * (13.d0 + 12 * pi**2 * r**2 + &
      (13.d0 - 12 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      lambda**2
      ELSE IF (mode == 3) THEN
      vv = vv + 2 * r**7 * (r - r0)**2 * cos(time)**4 * (2.d0 + 3 * pi**2 * r**2 + &
      (2.d0 - 3 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      lambda**2
      ELSE IF (mode == 4) THEN
      vv = vv - r**7 * (r - r0)**2 * cos(time)**4 * sin(4*pi*z)**2 / &
      (2 * lambda**2)
      END IF
      ELSE IF (TYPE == 3) THEN
      IF (mode == 0) THEN
      vv = vv + r**7 * (r - r0)**2 * cos(time)**4 * (15.d0 + 8 * pi**2 * r**2 + &
      (15.d0 - 8 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      (2 * lambda**2)
      ELSE IF (mode == 1) THEN
      vv = vv + r**7 * (r - r0)**2 * cos(time)**4 * (12.d0 + 7 * pi**2 * r**2 + &
      (12.d0 - 7 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      lambda**2
      ELSE IF (mode == 2) THEN
      vv = vv + r**7 * (r - r0)**2 * cos(time)**4 * (5.d0 + 4 * pi**2 * r**2 + &
      (5.d0 - 4 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      lambda**2
      ELSE IF (mode == 3) THEN
      vv = vv + 2 * pi**2 * r**9 * (r - r0)**2 * cos(time)**4 * sin(2*pi*z)**4 / &
      lambda**2
      ELSE IF (mode == 4) THEN
      vv = vv - r**7 * (r - r0)**2 * cos(time)**4 * sin(4*pi*z)**2 / &
      (4 * lambda**2)
      END IF
      ELSE IF (TYPE == 4) THEN
      IF (mode == 1) THEN
      vv = vv + r**7 * (r - r0)**2 * cos(time)**4 * (25.d0 + 8 * pi**2 * r**2 + &
      (25.d0 - 8 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      (4 * lambda**2)
      ELSE IF (mode == 2) THEN
      vv = vv + r**7 * (r - r0)**2 * cos(time)**4 * (61.d0 + 24 * pi**2 * r**2 + &
      (61.d0 - 24 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      (8 * lambda**2)
      ELSE IF (mode == 3) THEN
      vv = vv + r**7 * (r - r0)**2 * cos(time)**4 * (17.d0 + 8 * pi**2 * r**2 + &
      (17.d0 - 8 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      (4 * lambda**2)
      ELSE IF (mode == 4) THEN
      vv = vv + r**7 * (r - r0)**2 * cos(time)**4 * (15.d0 + 8 * pi**2 * r**2 + &
      (15.d0 - 8 * pi**2 * r**2) * cos(4*pi*z)) * sin(2*pi*z)**2 / &
      (16 * lambda**2)
      END IF
      ELSE IF (TYPE == 5) THEN
      IF (mode == 0) THEN
      vv = vv + pi * r**8 * (-215.d0 + 136 * pi**2 * r**2) * (r - r0)**2 * cos(time)**4 * &
      cos(2*pi*z) * sin(2*pi*z)**3 / (4 * lambda**2)
      ELSE IF (mode == 1) THEN
      vv = vv + 5 * pi * r**8 * (-11.d0 + 8 * pi**2 * r**2) * (r - r0)**2 * cos(time)**4 * &
      cos(2*pi*z) * sin(2*pi*z)**3 / lambda**2
      ELSE IF (mode == 2) THEN
      vv = vv + 14 * pi * r**8 * (r - r0)**2 * cos(time)**4 * &
      cos(2*pi*z) * sin(2*pi*z)**3 / lambda**2
      ELSE IF (mode == 3) THEN
      vv = vv - pi * r**8 * (-19.d0 + 8 * pi**2 * r**2) * (r - r0)**2 * cos(time)**4 * &
      cos(2*pi*z) * sin(2*pi*z)**3 / lambda**2
      ELSE IF (mode == 4) THEN
      vv = vv - pi * r**8 * (-15.d0 + 8 * pi**2 * r**2) * (r - r0)**2 * cos(time)**4 * &
      cos(2*pi*z) * sin(2*pi*z)**3 / (4 * lambda**2)
      END IF
      ELSE IF (TYPE == 6) THEN
      IF (mode == 1) THEN
      vv = vv + 8 * pi * r**8 * (-9.d0 + 5 * pi**2 * r**2) * (r - r0)**2 * cos(time)**4 * &
      cos(2*pi*z) * sin(2*pi*z)**3 / lambda**2
      ELSE IF (mode == 2) THEN
      vv = vv + 4 * pi * r**8 * (-13.d0 + 8 * pi**2 * r**2) * (r - r0)**2 * cos(time)**4 * &
      cos(2*pi*z) * sin(2*pi*z)**3 / lambda**2
      ELSE IF (mode == 3) THEN
      vv = vv + 8 * pi * r**8 * (-1.d0 + pi**2 * r**2) * (r - r0)**2 * cos(time)**4 * &
      cos(2*pi*z) * sin(2*pi*z)**3 / lambda**2
      ELSE IF (mode == 4) THEN
      vv = vv + 2 * pi * r**8 * (r - r0)**2 * cos(time)**4 * &
      cos(2*pi*z) * sin(2*pi*z)**3 / lambda**2
      END IF
      END IF
  9. The function Jexact_gauss computes the source term \(\bJ\) of the magnetostatic eaquations.
    1. The limit between fluid and solid region and \(\pi\) are defined:
      REAL(KIND=8) :: r0 = 0.5d0, pi = acos(-1.d0)
    2. We construct the radial and vertical coordinates r, z.
      r = rr(1)
      z = rr(2)
    3. The current density is defined by type and mode.
      IF (TYPE == 1) THEN
      IF (m == 0) THEN
      vv = 4 * pi**2 * r**3 * cos(2*pi*z) * cos(t)
      ELSE IF (m == 1) THEN
      vv = 4 * r * cos(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      ELSE IF (TYPE == 2) THEN
      IF (m == 1) THEN
      vv = r * (1.d0 + 4 * pi**2 * r**2) * cos(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      ELSE IF (TYPE == 3) THEN
      IF (m == 0) THEN
      vv = 4 * r * (-2.d0 + pi**2 * r**2) * cos(2*pi*z) * cos(t)
      ELSE IF (m == 1) THEN
      vv = 2 * r * cos(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      ELSE IF (TYPE == 4) THEN
      IF (m == 1) THEN
      vv = 4 * r * (-2.d0 + pi**2 * r**2) * cos(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      ELSE IF (TYPE == 5) THEN
      IF (m == 0) THEN
      vv = - 8 * pi * r**2 * sin(2*pi*z) * cos(t)
      ELSE IF (m == 1) THEN
      vv = - 2 * pi * r**2 * sin(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      ELSE
      IF (m == 1) THEN
      vv = - 8 * pi * r**2 * sin(2*pi*z) * cos(t)
      ELSE
      vv = 0.d0
      END IF
      END IF
  10. The function kelvin_force_coeff computes the coefficient \(\chi(T)\) of the Kelvin force. It takes a real \(T\) and gives the real \(\chi(T)\). In this case \(\chi(T) = - T^2\) so we implement:
    vv = - temp**2

All the other subroutines present in the file condlim_test_34.f90 are not used in this test. We refer to the section Fortran file condlim.f90 for a description of all the subroutines of the condlim file.

Setting in the data file

We describe the data file of this test. It is called debug_data_test_34 and can be found in the directory ($SFEMaNS_DIR)/MHD_DATA_TEST_CONV_PETSC.

  1. We use a formatted mesh by setting:
    ===Is mesh file formatted (true/false)?
    .t.
  2. The path and the name of the mesh are specified with the two following lines:
    ===Directory and name of mesh file
    '.' 'SOLID_FLUID_10.FEM'
    where '.' refers to the directory where the data file is, meaning ($SFEMaNS_DIR)/MHD_DATA_TEST_CONV_PETSC.
  3. We use one processor in the meridian section. It means the finite element mesh is not subdivised.
    ===Number of processors in meridian section
    1
  4. We solved the problem for \(5\) Fourier modes.
    ===Number of Fourier modes
    5
    The Fourier modes are not detailed so the first three modes \(0,1,2,3,4\) are solved.
  5. We use \(5\) processors in Fourier space.
    ===Number of processors in Fourier space
    5
    It means that each processors is solving the problem for \(5/5=1\) Fourier modes.
  6. We approximate the ferrohydrodynamic equations by setting:
    ===Problem type: (nst, mxw, mhd, fhd)
    'fhd'
  7. We do not restart the computations from previous results.
    ===Restart on velocity (true/false)
    .f.
    ===Restart on magnetic field (true/false)
    .f.
    It means the computation starts from the time \(t=0\).
  8. We use a time step of \(0.01\) and solve the problem over \(100\) time iterations.
    ===Time step and number of time iterations
    1.d-2 100
  9. We set the periodic boundaries.
    1. We indicate the number of couples of periodic interfaces. In our case, there is only one, the couple of interfaces 4 (bottom) and 2 (top)
      ===How many pieces of periodic boundary?
      1
      There could have been a second couple if we had given a different number to the top (and bottom) boundaries of solid and fluid subdomains.
    2. We indicate the interfaces of each couple, and the coordinates in the plan \((r,z)\) of the vector that joins them. The vector that joins interfaces 4 and 2 is the vector \(\be_z\).
      ===Indices of periodic boundaries and corresponding vectors
      4 2 0.d0 1.d0
      If there was additional couples, corresponding lines should be written below.
  10. We set the number of domains and their label, see the files associated to the generation of the mesh, where the code approximates the Navier-Stokes equations.
    ===Number of subdomains in Navier-Stokes mesh
    1
    ===List of subdomains for Navier-Stokes mesh
    2
  11. We set the number of boundaries with Dirichlet conditions on the velocity field and give their respective labels.
    ===How many boundary pieces for full Dirichlet BCs on velocity?
    2
    ===List of boundary pieces for full Dirichlet BCs on velocity
    3 5
  12. We set the kinetic Reynolds number \(\Re\).
    ===Reynolds number
    1.d0
  13. We give information on how to solve the matrix associated to the time marching of the velocity.
    1. ===Maximum number of iterations for velocity solver
      100
    2. ===Relative tolerance for velocity solver
      1.d-6
      ===Absolute tolerance for velocity solver
      1.d-10
    3. ===Solver type for velocity (FGMRES, CG, ...)
      GMRES
      ===Preconditionner type for velocity solver (HYPRE, JACOBI, MUMPS...)
      MUMPS
  14. We give information on how to solve the matrix associated to the time marching of the pressure.
    1. ===Maximum number of iterations for pressure solver
      100
    2. ===Relative tolerance for pressure solver
      1.d-6
      ===Absolute tolerance for pressure solver
      1.d-10
    3. ===Solver type for pressure (FGMRES, CG, ...)
      GMRES
      ===Preconditionner type for pressure solver (HYPRE, JACOBI, MUMPS...)
      MUMPS
  15. We give information on how to solve the mass matrix.
    1. ===Maximum number of iterations for mass matrix solver
      100
    2. ===Relative tolerance for mass matrix solver
      1.d-6
      ===Absolute tolerance for mass matrix solver
      1.d-10
    3. ===Solver type for mass matrix (FGMRES, CG, ...)
      CG
      ===Preconditionner type for mass matrix solver (HYPRE, JACOBI, MUMPS...)
      MUMPS
  16. We solve the temperature equation.
    ===Is there a temperature field?
    .t.
  17. We set the coefficient \(\alpha\) of Boussinesq force.
    ===Non-dimensional gravity coefficient
    1.d0
  18. We set the number of domains and their label, see the files associated to the generation of the mesh, where the code approximates the temperature equation.
    ===Number of subdomains in temperature mesh
    2
    ===List of subdomains for temperature mesh
    1 2
  19. We set the volumetric heat capacity \(c\) and the thermal conductivity \(\lambda\).
    ===Volumetric heat capacity (1:nb_dom_temp)
    1.d0 2.d0
    ===Thermal conductivity (1:nb_dom_temp)
    10.d0 1.d0
  20. We set the number of boundaries with Dirichlet conditions on the velocity and give their respective labels.
    ===How many boundary pieces for Dirichlet BCs on temperature?
    1
    ===List of boundary pieces for Dirichlet BCs on temperature
    5
  21. Setting the interfaces between regions where only the temperature is solved and regions where velocity and temperature are solved is not required for 'fhd' problems.
  22. We give information on how to solve the matrix associated to the time marching of the temperature.
    1. ===Maximum number of iterations for temperature solver
      100
    2. ===Relative tolerance for temperature solver
      1.d-6
      ===Absolute tolerance for temperature solver
      1.d-10
    3. ===Solver type for temperature (FGMRES, CG, ...)
      GMRES
      ===Preconditionner type for temperature solver (HYPRE, JACOBI, MUMPS...)
      MUMPS
  23. We set the number of domains and their label, see the files associated to the generation of the mesh, where the code approximates the magnetic field.
    ===Number of subdomains in magnetic field (H) mesh
    2
    ===List of subdomains for magnetic field (H) mesh
    1 2
  24. We set the number of interfaces in H_mesh and give their respective labels.
    ===Number of interfaces in H mesh
    1
    ===List of interfaces in H mesh
    3
    There is no permeability jump but it is necessary to indicate that there is an interface in H_mesh in order to impose boundary conditions on the velocity on this interface.
  25. We set the number of boundaries with Dirichlet conditions on the magnetic field and give their respective labels.
    ===Number of Dirichlet sides for Hxn
    1
    ===List of Dirichlet sides for Hxn
    5
  26. We set the magnetic permeability in each domains where the magnetic field is approximated.
    ===Permeability in the conductive part (1:nb_dom_H)
    1.d0 1.d0
  27. We set the conductivity in each domains where the magnetic field is approximated.
    ===Conductivity in the conductive part (1:nb_dom_H)
    1.d-20 1.d-20
  28. We set the type of finite element used to approximate the magnetic field.
    ===Type of finite element for magnetic field
    2
  29. We set the magnetic Reynolds number \(\Rm\).
    ===Magnetic Reynolds number
    1.d0
  30. We set the stabilization coefficients for the divergence of the magnetic field and the penalty of the Dirichlet and interface terms.
    ===Stabilization coefficient (divergence)
    1.d0
    ===Stabilization coefficient for Dirichlet H and/or interface H/H
    1.d0
    We note that these coefficients are usually set to \(1\).
  31. We set the number of domains and their label, see the files associated to the generation of the mesh, where the code approximates the potential scalar: none, only an H formulation is used.
    ===Number of subdomains in magnetic potential (phi) mesh
    0
  32. We give information on how to solve the matrix associated to the time marching of the Maxwell equations.
    1. ===Maximum number of iterations for Maxwell solver
      100
    2. ===Relative tolerance for Maxwell solver
      1.d-6
      ===Absolute tolerance for Maxwell solver
      1.d-10
    3. ===Solver type for Maxwell (FGMRES, CG, ...)
      GMRES
      ===Preconditionner type for Maxwell solver (HYPRE, JACOBI, MUMPS...)
      MUMPS
  33. To get the total elapse time and the average time in loop minus initialization, we write:
    ===Verbose timing? (true/false)
    .t.
    These informations are written in the file lis when you run the shell debug_SFEMaNS_template.

Outputs and value of reference

The outputs of this test are computed with the file post_processing_debug.f90 that can be found in: ($SFEMaNS_DIR)/MHD_DATA_TEST_CONV_PETSC.

To check the well behavior of the code, we compute four quantities:

  1. The L2-norm of error on u divided by the L2-norm of u exact.
  2. The L2-norm of error on p divided by L2-norm of p exact.
  3. The L2-norm of error on T divided by L2-norm of T exact.
  4. The L2-norm of error on H divided by L2-norm of H exact.

These quantities are computed at the final time \(t=1\). They are compared to reference values to attest of the correct behavior of the code. These values of reference are in the last lines of the file debug_data_test_34 in ($SFEMaNS_DIR)/MHD_DATA_TEST_CONV_PETSC. They are equal to:

============================================
(SOLID_FLUID_10.FEM)
===Reference results
3.6382183420384366E-004 L2-norm of error on u / L2-norm of u exact
6.1696264206533347E-002 L2-norm of error on p / L2-norm of p exact
2.2000110457751656E-004 L2-norm of error on T / L2-norm of T exact
2.8695295848810841E-003 L2-norm of error on H / L2-norm of H exact

To conclude this test, we display the profile of the approximated pressure, velocity magnitude, temperature and magnetic field magnitude at the final time. These figures are done in the plane \(y=0\) which is the union of the half plane \(\theta=0\) and \(\theta=\pi\).

fig_test_34_pre_tfin.png
Pressure in the plane y=0.
fig_test_34_vel_tfin.png
Velocity magnitude in the plane y=0.
fig_test_34_temp_tfin.png
Temperature in the plane y=0.
fig_test_34_mag_tfin.png
Magnetic field magnitude in the plane y=0.