SFEMaNS  version 4.1 (work in progress)
Reference documentation for SFEMaNS
 All Classes Files Functions Variables Groups Pages
Test 21: MHD with variable density and electrical conductivity

Introduction

In this example, we check the correctness of SFEMaNS for a magnetohydrodynamic problem with variable electrical conductivity involving Dirichlet boundary conditions. The mass equation is not approximated. We consider a level set, solution of the same advection equation, that is used to reconstruct the density and the other fluid's properties.

We solve the level set equation:

\begin{align*} \partial_t \varphi + \bu \cdot \GRAD \varphi = f_\varphi. \end{align*}

We recontruct the density, dynamical viscosity and electrical conductivity as follows:

\begin{align*} \rho=\rho_1+(\rho_2-\rho_1) F(\varphi), \\ \eta=\eta_1+(\eta_2-\rho_1) F(\varphi), \end{align*}

with \(\eta_i\) and \(\rho_i\) data to define. The function \(F(\varphi)\) is either the identity (called linear reconstruction) or a piece-wise polynomial function (called reg reconstruction).

We solve the Navier-Stokes equations with the momentum \(\bm=\rho\bu\) as dependent variable:

\begin{align*} \partial_t\bm + \DIV(\bm\otimes\bu) - \frac{2}{\Re} \DIV(\eta \epsilon(\bu)) +\GRAD p &=\bef + (\ROT \textbf{H}) \times (\mu^c \textbf{H}), \\ \DIV \bu &= 0, \\ \bu_{|\Gamma} &= \bu_{\text{bdy}} , \\ \bu_{|t=0} &= \bu_0, \\ p_{|t=0} &= p_0, \end{align*}

where \(\epsilon(\bu)=\frac{1}{2}(\GRAD \bu + \GRAD^s \bu)\) is the strain rate tensor.

We solve the Maxwell equations:

\begin{align*} \partial_t (\mu^c \mathbf{H}) + \nabla \times \left(\frac{1}{\Rm \sigma}\nabla \times \mathbf{H} \right) & = \nabla\times (\bu \times \mu^c \mathbf{H}) + \nabla \times \left(\frac{1}{\Rm \sigma} \nabla\times \mathbf{j} \right), \\ \text{div} (\mu^c \mathbf {H}) &= 0 ,\\ {\bH \times \bn}_{|\Gamma} &= {\bH_{\text{bdy}} \times \bn}_{|\Gamma},\\ \bH_{|t=0}&= \bH_0. \end{align*}

These equations are solved in the domain \(\Omega= \{ (r,\theta,z) \in {R}^3 : (r,\theta,z) \in [0,1/2] \times [0,2\pi) \times [0,1]\} \) with \(\Gamma= \partial \Omega \). The data are the source terms \(f_\varphi\), \(\bef\) and \(\textbf{j}\), the boundary data \(\bu_{\text{bdy}}\), the initial datas \(\bu_0\) and \(p_0\). The parameter \(\Re\) is the kinetic Reynolds number$, the magnetic Reynolds number \(\Rm\), the magnetic permeability \(\mu^c\), the densities \(\rho_i\), the dynamical viscosities \(\eta_i\) and the electrical conductivities \(\sigma_i\).

Remarks:

  1. The level set and the momentum equations can be stabilized with the entropy viscosity method used in test 15 (called LES).
  2. For physical problem, the level set has to take values in [0,1] such that the interface is represented by \(\varphi^{-1}(\{1/2\})\). The fluids area are respectively represented by \(\varphi^{-1}(\{0\})\) and \(\varphi^{-1}(\{1\})\). This test does not consider immiscible fluids. It involves manufactured solution with a smooth variable density. As a consequence, the level set does not take values in [0,1].
  3. A compression term can also be added in the level set equation. This term allows the level set to remains sharp near the fluids interface represented by \( \varphi^{-1}(\{1/2\}) \).

We refer to the section Extension to multiphase flow problem for more details on the algorithms implemented in SFEMaNS for multiphase problem.

Manufactured solutions

We approximate the following analytical solutions:

\begin{align*} u_r(r,\theta,z,t) &= 0, \\ u_{\theta}(r,\theta,z,t) &= 0 , \\ u_z(r,\theta,z,t) &=0 , \\ p(r,\theta,z,t) &= 0, \\H_r(r,\theta,z,t) &= r^3 , \\H_{\theta}(r,\theta,z,t) &= 0 , \\H_z(r,\theta,z,t) &= - 4 r^2 z , \\ \varphi(r,\theta,z,t) &=1 + z + r \cos(\theta), \\ \rho(r,\theta,z,t) &=2 + z + r \cos(\theta), \\ \eta(r,\theta,z,t) &=1, \\ \sigma(r,\theta,z,t) &=2 + z + r \cos(\theta). \end{align*}

where the source terms \(f_\varphi\), \(\bef\), \(\textbf{j}\) and the boundary data \( \bu_{\text{bdy}}\) are computed accordingly.

Generation of the mesh

The finite element mesh used for this test is named Mesh_10_form.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/Mesh_10_form. The following image shows the mesh for P1 finite elements.

fig_Mesh_10_form.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_21.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 time \(-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 function init_level_set initializes the level set at the time \(-dt\) and \(0\) with \(dt\) being the time step. This is done by using the function level_set_exact as follows:
    time = 0.d0
    DO i= 1, SIZE(list_mode)
    mode = list_mode(i)
    DO j = 1, 2
    !===level_set
    DO n = 1, inputs%nb_fluid -1
    level_set_m1(n,:,j,i) = level_set_exact(n,j,vv_mesh%rr,mode,time-dt)
    level_set (n,:,j,i) = level_set_exact(n,j,vv_mesh%rr,mode,time)
    END DO
    END DO
    END DO
    We note there is one level set per interface and the different phases are stratified. So for nb_fluid fluids, we have inputs nb_fluid-1 interfaces.
  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. It is set to zero.
    vv = 0.d0
    RETURN
  4. The function pp_exact contains the analytical pressure. It is used to initialize the pressure to zero.
    vv=0.d0
    RETURN
  5. The function level_set_exact is used to initialize the level set.
    1. We define the level set of the mode and its TYPE (1 for cosine and 2 for sine) as follows:
      IF (interface_nb==1) THEN
      IF (m==0 .AND. TYPE ==1) THEN
      vv = 1.d0 + rr(2,:)
      ELSE IF (m==1 .AND. TYPE==1) THEN
      vv = rr(1,:)
      ELSE
      vv = 0.d0
      END IF
    2. If more than one level set is considered, the computation is stopped.
      ELSE
      CALL error_petsc(' BUG in level_set_exact, we should compute only 1 level set')
      END IF
      RETURN
      Indeed with two fluids there is one interface. So we should compute only one level set.
  6. The function source_in_NS_momentum computes the source term \(\bef\) of the Navier-Stokes equations.
  7. The function source_in_level_set computes the source term \(f_\varphi\) of the level set equations. It is equal to zero.
  8. The subroutine init_maxwell initializes the magnetic field at the time \(-dt\) and \(0\) with \(dt\) being the time step. It is done by using the function Hexact.
    time = -dt
    DO k=1,6
    DO i=1, SIZE(list_mode)
    Hn1(:,k,i) = Hexact(H_mesh,k, H_mesh%rr, list_mode(i), mu_H_field, time)
    IF (inputs%nb_dom_phi>0) THEN
    IF (k<3) THEN
    phin1(:,k,i) = Phiexact(k, phi_mesh%rr, list_mode(i) , mu_phi, time)
    ENDIF
    END IF
    ENDDO
    ENDDO
    time = time + dt
    DO k=1,6
    DO i=1, SIZE(list_mode)
    Hn(:,k,i) = Hexact(H_mesh,k, H_mesh%rr, list_mode(i), mu_H_field, time)
    IF (inputs%nb_dom_phi>0) THEN
    IF (k<3) THEN
    phin(:,k,i) = Phiexact(k, phi_mesh%rr, list_mode(i), mu_phi, time)
    ENDIF
    END IF
    ENDDO
    ENDDO
    We note the scalar potential is also initialized with the function Phiexact in the case that nb_dom_phi>0. As there is no scalar potential in this test, we have nb_dom_phi=0.
  9. The function Hexact contains the analytical magnetic field. It is used to initialize the magnetic field and to impose Dirichlet boundary conditions on \(\textbf{H}\times\textbf{n}\) with \(\textbf{n}\) the outter normal vector.
    IF (m==0.AND. TYPE==1) THEN
    vv = rr(1,:)**3
    ELSE IF (m==0.AND.TYPE==5) THEN
    vv = -4.d0*rr(1,:)**2*rr(2,:)
    ELSE
    vv = 0.d0
    END IF
    RETURN
  10. The function Phiexact contains the analytical scalar potential. It is not used in this test so we stop the code if this subroutine is called as follows:
    vv=0.d0
    CALL error_petsc('Phiexact: should not be called for this test')
    RETURN
  11. The function Jexact_gauss is used to define the source term \(\textbf{j}\). Since \(\partial_t (\mu^c \textbf{H})=0 \) and \( \bu=0 \), this source term has to satisfy the relation \(\textbf{j}= \ROT(\textbf{H}) \).
    IF(m==0.AND.TYPE==3) THEN
    vv = 8.d0*rr(1)*rr(2)
    ELSE
    vv = 0.d0
    END IF
    RETURN

All the other subroutines present in the file condlim_test_21.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_21 and can be found in the following 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
    '.' 'Mesh_10_form.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 solve the problem for \(4\) Fourier modes.
    ===Number of Fourier modes
    1
  5. We use \(1\) processors in Fourier space.
    ===Number of processors in Fourier space
    1
    It means that each processors is solving the problem for \(1/1=1\) Fourier modes.
  6. We do not select specific Fourier modes to solve.
    ===Select Fourier modes? (true/false)
    .f.
    As a consequence, the code approximates the problem on the first \(4\) Fourier modes.
  7. We approximate the MHD equations by setting:
    ===Problem type: (nst, mxw, mhd, fhd)
    'mhd'
  8. 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\).
  9. We use a time step of \(0.01\) and solve the problem over \(100\) time iterations.
    ===Time step and number of time iterations
    .01d0, 100
  10. We do not apply mass correction on the level set.
    ===Do we apply mass correction? (true/false)
    .f.
    The default value is true.
  11. We don't use a level set \(\varphi\) with values in \([0,1]\). So we do not kill the overshoot of the level set with respect of the interval \([0,1]\).
    ===Do we kill level set overshoot? (true/false)
    .f.
    The default value if false so these two lines are not required.
  12. 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
    1
  13. 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?
    3
    ===List of boundary pieces for full Dirichlet BCs on velocity
    2 4 5
  14. We use the momentum as dependent variable for the Navier-Stokes equations.
    ===Solve Navier-Stokes with u (true) or m (false)?
    .f.
    If the density or the viscosity are variable, this parameter needs to be false. The default value is true (constant density and viscosity).
  15. We use a BDF1 approximation of the time derivatives in the level set and momentum equations.
    ===Do we solve momentum with bdf2 (true/false)?
    .f.
    The default value is false.
  16. We set the kinetic Reynolds number \(\Re\).
    ===Reynolds number
    250.d0
  17. We use the entropy viscosity method to stabilize the level set equation.
    ===Use LES? (true/false)
    .t.
    This parameter needs to be true for multiphase problem.
  18. We define the coefficient \(c_\text{e}\) of the entropy viscosity.
    ===Coefficient multiplying residual
    0.d0
    It means the entropy viscosity is zero.
  19. We give information on how to solve the matrix associated to the time marching of the velocity (or momentum in this case).
    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
  20. 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
  21. 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
  22. We solve the level set equation.
    ===Is there a level set?
    .t.
  23. We want to consider one level set \(\varphi\), so we set:
    ===How many fluids?
    2
    We note this test does not consider two immiscible fluids.
  24. We do not use compression tools.
    ===Compression factor for level set
    0.d0
    This parameters is only relevant when we want to get sharp interface for immiscible fluids.
  25. We define the parameters \((\rho_1,\rho_2)\) used to reconstruct the density.
    ===Density of fluid 0, fluid 1, ...
    1.d0 2.d0
  26. We define the parameters \((\eta_1,\eta_2)\) used to reconstruct the dynamical viscosity.
    ===Dynamic viscosity of fluid 0, fluid 1, ...
    1.d0 1.d0
  27. We define a multiplier coefficient.
    ===multiplier for h_min for level set
    1.d0
    This multiplier times the smallest mesh size is stored in the variable inputs%h_min_distance. It can be used in the condlim file to set the wideness of the initial interface. It is not used in this case as the level set does not represent an interface between two immiscible fluids.
  28. We define the parameters \((\sigma_1,\sigma_2)\) used to reconstruct the electrical conductivity.
    ===Conductivity of fluid 0, fluid 1, ...
    1.d0 2.d0
  29. We use a linear reconstruction, meaning \(F(\varphi)=\varphi\).
    ===How are the variables reconstructed from the level set function? (lin, reg)
    'lin'
  30. We do not impose Dirichlet conditions on the level set.
    ===How many boundary pieces for Dirichlet BCs on level set?
    0
  31. We give information on how to solve the matrix associated to the time marching of the level set.
    1. ===Maximum number of iterations for level set solver
      100
    2. ===Relative tolerance for level set solver
      1.d-6
      ===Absolute tolerance for level set solver
      1.d-10
    3. ===Solver type for level set (FGMRES, CG, ...)
      GMRES
      ===Preconditionner type for level set solver (HYPRE, JACOBI, MUMPS...)
      MUMPS
  32. We use the magnetic field B as dependent variable for the Maxwell equations.
    ===Solve Maxwell with H (true) or B (false)?
    .f.
  33. 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
    1
    ===List of subdomains for magnetic field (H) mesh
    1
  34. We set the number of interface in H_mesh.
    ===Number of interfaces in H mesh
    0
    Such interfaces represent interfaces with discontinuities in magnetic permeability or interfaces between the magnetic field mesh and the temperature or the velocity field meshes.
  35. We set the number of boundaries with Dirichlet conditions on the magnetic field and give their respective labels.
    ===Number of Dirichlet sides for Hxn
    3
    ===List of Dirichlet sides for Hxn
    5 2 4
  36. 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
  37. We set the conductivity in each domains where the magnetic field is approximated.
    ===Conductivity in the conductive part (1:nb_dom_H)
    1.d0
    This value is not used by the code since the conductivity is reconstructed with the level set. These lines are required. They allow to approximate the Maxwell equations is another domain than the fluid one with a conductivity sets in the data.
  38. We set the type of finite element used to approximate the magnetic field.
    ===Type of finite element for magnetic field
    2
  39. We set the magnetic Reynolds number \(\Rm\).
    ===Magnetic Reynolds number
    2.d0
  40. We set stabilization coefficient 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 these coefficients are usually set to \(1\).
  41. We set the number of domains and their label, see the files associated to the generation of the mesh, where the code approximates the scalar potential.
    ===Number of subdomains in magnetic potential (phi) mesh
    0
  42. 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
  43. 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 the following directory: ($SFEMaNS_DIR)/MHD_DATA_TEST_CONV_PETSC.

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

  1. The L2 norm of the error on the velocity field.
  2. The H1 norm of the error on the velocity field.
  3. The L2 norm of the error on the level set.
  4. The L2 norm of the error on the the magnetic field.

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

============================================
(Mesh_10_form.FEM)
===Reference results
1.3906183032940818E-005 L2 error on velocity
3.0086983840053983E-006 L2 error on pressure
6.7632134413749881E-006 L2 error on level set
1.8335168044284966E-005 L2 error on Magnetic field

To conclude this test, we show the profile of the approximated level set, density, pressure and velocity 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_21_level_set_tfin.png
Level set in the plane plane y=0.
fig_test_21_pre_tfin.png
Presure in the plane plane y=0.
fig_test_21_vel_tfin.png
Velocity magnitude in the plane plane y=0.
fig_test_21_mag_tfin.png
Magnetic field magnitude in the plane plane y=0.