SFEMaNS  version 4.1 (work in progress)
Reference documentation for SFEMaNS
 All Classes Files Functions Variables Groups Pages
Test 18: Maxwell Equations with variable magnetic permeability in (r,z) with jumps

Introduction

In this example, we check the correctness of SFEMaNS for stationary magnetic problem involving a conducting. The magnetic permeability is variable in (r,z) and presents a jump in \(r=1\). The magnetic permeability is given. We consider Dirichlet boundary conditions. We use P2 finite elements. This set up is approximated over one time iteration with a large time step.

We solve the Maxwell equations with the magnetic field \(\bB=\mu \bH\) as dependent variable:

\begin{align} \begin{cases} \partial_t (\mathbf{B}) + \nabla \times \left(\frac{1}{\Rm \sigma}\nabla \times( \frac{1}{\mu^c}\mathbf{B} ) \right) = \nabla\times (\bu^\text{given} \times \mathbf{B}) + \nabla \times \left(\frac{1}{\Rm \sigma} \nabla\times \mathbf{j} \right) & \text{in } \Omega, \\ \text{div} (\mathbf{B}) = 0 & \text{in } \Omega ,\\ \bH_1 \times \bn_1 + \bH_2 \times \bn_2 = 0 & \text{on } \Sigma_\mu,\\ \mu^c_1\bH _1 \cdot \bn_1 + \mu^c_2 \bH _2 \cdot \bn_2 = 0 & \text{on } \Sigma_\mu,\\ {\bH \times \bn}_{|\Gamma} = {\bH_{\text{bdy}} \times \bn}_{|\Gamma}, &\\ \bH_{|t=0}= \bH_0, \\ \end{cases} \end{align}

in the domain \(\Omega= \{ (r,\theta,z) \in {R}^3 : (r,\theta,z) \in [0,2] \times [0,2\pi) \times [0.25,1]\} \). This domain is the union of two conduction domain \(\Omega_1=\{ (r,\theta,z) \in {R}^3 : (r,\theta,z) \in [0,1] \times [0,2\pi) \times [0.25,1]\}\) and \( \Omega_2= \{ (r,\theta,z) \in {R}^3 : (r,\theta,z) \in [1,2] \times [0,2\pi) \times [0.25,1]\}\). We denote by \(\Sigma_\mu\) the interface of \(\Omega_1\) and \(\Omega_2\). We also set \(\Gamma= \partial \Omega\). We define outward normals \(\bn_1, \bn_2\) to the surface \(\Sigma_\mu\) and the outward normals \(\bn\) to the surface \(\Gamma\). The data are the source term \(\mathbf{j}\), the given velocity \(\bu^\text{given}\), the boundary data \(\bH_{\text{bdy}}\) and the initial data \(\bH_0\). The parameter \(\Rm\) is the magnetic Reynolds number. The parameter \(\mu^c\) is the magnetic permeability of the conducting region \(\Omega\). The parameter \(\sigma\) is the conductivity in the conducting region \(\Omega\).

Manufactured solutions

We define the magnetic permeability in the conducting region as follows:

\begin{align*} \mu^c(r,z)= \begin{cases} 1 + r & \text{ in } \Omega_1, \\ 1+ r + 2 \lambda_\mu \frac{1+r}{z^2(3r+2)} &\text{ in } \Omega_2, \end{cases} \end{align*}

with \(\lambda_\mu= 10\).

We approximate the following analytical solutions:

\begin{align*} H_r(r,\theta,z,t) &= \begin{cases} r z & \text{ in } \Omega_1, \\ \frac{(3r+2)z^3r}{3z^2r+2z^2+ 2 \lambda_\mu} &\text{ in } \Omega_2, \end{cases} \\ H_{\theta}(r,\theta,z,t) &=0, \\ H_z(r,\theta,z,t) &= - \frac{z^2(3r+2)}{2(1+r)}, \end{align*}

and remind that \(\bB=\mu\bH\).

We also set the given velocity field to zero.

\begin{align*} u^\text{given}_r(r,\theta,z,t) &= 0, \\ u^\text{given}_{\theta}(r,\theta,z,t) &= 0, \\ u^\text{given}_z(r,\theta,z,t) &= 0. \end{align*}

The source term \( \mathbf{j}\) and the boundary data \(\bH_{\text{bdy}}\) are computed accordingly.

Generation of the mesh

The finite element mesh used for this test is named mesh_18_0.05.FEM. The mesh size for the P1 approximation is \(0.05\) in \(\Omega\). You can generate this mesh with the files in the following directory: ($SFEMaNS_MESH_GEN_DIR)/EXAMPLES/EXAMPLES_MANUFACTURED_SOLUTIONS/mesh_18_0.05.FEM. The following images show the mesh for P1 finite elements for the conducting and vacuum regions.

fig_mesh_18_0.05.png
Finite element mesh (P1) of the conducting region.

Information on the file condlim.f90

The initial conditions, boundary conditions and the forcing term \(\textbf{j}\) in the Maxwell equations are set in the file condlim_test_18.f90. This condlim also uses functions of the file test_18.f90 to define the magnetic permeability and its gradient. Here is a description of the subroutines and functions of interest of the file condlim_test_18.f90.

  1. The subroutine init_maxwell initializes the magnetic field and the scalar potential at the time \(-dt\) and \(0\) with \(dt\) being the time step. It is done by using the functions Hexact, Phiexact as follows:
    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
  2. 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.
    1. If the Fourier mode m is not equal to 0, we set the magnetic field to zero.
      IF (m/=0) THEN
      vv(:) = 0.d0
      RETURN
      END IF
    2. We define a tabular id so that each nodes n, id(n) is equal to 1 or 2 depending if the node is in \(\Omega_1\) or \(\Omega_2\).
      1. If rr contains the radial and vertical coordinates of each np nodes, we define id as follows.
        if (SIZE(rr,2)== H_mesh%np) THEN
        DO mm = 1, H_mesh%me !id is used to determine on which side of the interface we are
        id(H_mesh%jj(:,mm)) = H_mesh%i_d(mm) !used for initialization
        END DO
        This case is done during the initialization of the magnetic field.
      2. If rr does not contain the information on np nodes, it means Hexact is called to set boundary conditions. In that case, Hexact is called one time per nodes. So we define id, a tabular of dimension 1, as follows:
        ELSE
        IF (rr(1,1)<1) THEN !used for boundary condition
        id = 1
        ELSE
        id = 2
        END IF
        END IF
    3. For the Fourier mode \(m=0\), we define the magnetic 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) as follows:
      IF (TYPE==1) THEN
      DO n = 1, SIZE(rr,2)
      IF (id(n)==1) THEN
      vv(n) = rr(1,n)*rr(2,n)
      ELSE
      vv(n) = (3.d0*rr(1,n)+2.d0)*rr(2,n)**3*rr(1,n)/&
      (3.d0*rr(2,n)**2*rr(1,n) + 2.d0*rr(2,n)**2 + 2.d0*lambda_mu_T18 )
      END IF
      END DO
      ELSE IF (TYPE==5) THEN
      vv = - 0.5 * ( rr(2,:)**2*(3*rr(1,:)+2) )/( 1.d0 + rr(1,:) )
      ELSE
      vv = 0.d0
      END IF
      RETURN
      where lambda_mu_T18 is set to \(10\) in the file test_18.f90.
  3. The function Vexact is used to define the given velocity field \(\bu^\text{given}\). It is set to zero.
    vv = 0.d0
  4. The function Jexact_gauss is used to define the source term \(\textbf{j}\). We remind that \( \bu^\text{given}=0 \) and \(\partial_t \textbf{H}=0\) so it is defined such that \(\textbf{j}=\ROT \bH\).
    1. If the Fourier mode m is not equal to 0, it is set to zero.
      IF (m/=0) THEN
      vv = 0.d0
      RETURN
      END IF
    2. For the Fourier mode \(m=0\), we define the source term 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) as follows:
      IF (TYPE==3) THEN
      IF (mesh_id==1) THEN !Determine on which side of the interface we are
      vv = rr(1) - (1.d0/2.d0)*rr(2)**2*(3*rr(1)+2)/(1+rr(1))**2 &
      + (3.d0/2.d0)*rr(2)**2/(1+rr(1))
      ELSE
      vv = rr(1)*rr(2)**2*(3*rr(1)+2)*(3*rr(2)**2*rr(1)+2*rr(2)**2+6*lambda_mu_T18) &
      /(3*rr(2)**2*rr(1)+2*rr(2)**2+2*lambda_mu_T18)**2 &
      + (3.d0/2.d0)*rr(2)**2/(1+rr(1)) &
      - (1.d0/2.d0)*rr(2)**2*(3*rr(1)+2)/(1+rr(1))**2
      ENDIF
      ELSE
      vv = 0.d0
      END IF
      RETURN
  5. The function mu_bar_in_fourier_space defines the magnetic permeability that depends of the radial and vertical coordinates. It is done by using the function mu_bar_in_fourier_space_anal_T18 of the file test_18.f90.
    IF( PRESENT(pts) .AND. PRESENT(pts_ids) ) THEN
    vv=mu_bar_in_fourier_space_anal_T18(H_mesh,nb,ne,pts,pts_ids)
    ELSE
    vv=mu_bar_in_fourier_space_anal_T18(H_mesh,nb,ne)
    END IF
    RETURN
    When the magnetic permeability depends of the time or the azimuthal direction (set "===Is permeability variable in theta (true/false)?" to true in the data file), this function is used to define a stabilization term mu_bar. It needs to be smaller than mu everywhere in the domain. The magnetic permeability is then defined with the function mu_in_real_space.
  6. The function grad_mu_bar_in_fourier_space defines the gradient of the magnetic permeability. It is done by using the function grad_mu_bar_in_fourier_space_anal_T18 of the file test_18.f90.
    vv=grad_mu_bar_in_fourier_space_anal_T18(pt,pt_id)
    RETURN
    When the magnetic permeability depends of the time or the azimuthal direction, this function is used to define the gradient of a stabilization term mu_bar.

All the other subroutines present in the file condlim_test_18.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. The following describes the functions of interest of the file test_18.f90.

  1. First we define the real number \(\lambda_\mu\).
    REAL (KIND=8), PARAMETER, PUBLIC :: lambda_mu_T18=10.0d0
  2. The function mu_bar_in_fourier_space_anal_T18 defines the magnetic permeability.
    1. First we define the radial and vertical cylindrical coordinates of each nodes. We also define a tabular local_ids that is equal to 1 if the node is in \(\Omega_1\) and 2 if it's in \(\Omega_2\).
      IF( PRESENT(pts) .AND. PRESENT(pts_ids) ) THEN !Computing mu at pts
      r=pts(1,nb:ne)
      z=pts(2,nb:ne)
      local_ids=pts_ids
      ELSE
      r=H_mesh%rr(1,nb:ne) !Computing mu at nodes
      z=H_mesh%rr(2,nb:ne)
      DO m = 1, H_mesh%me
      global_ids(H_mesh%jj(:,m)) = H_mesh%i_d(m)
      END DO
      local_ids=global_ids(nb:ne)
      END IF
    2. We define the magnetic permeability.
      DO n = 1, ne - nb + 1
      IF(local_ids(n)==1) THEN
      vv(n) = 1.d0 + r(n) !mu1 , DCQ: If you change mu1_bar, you have to change
      !Jexact_gauss() as well
      ELSE
      vv(n) = 1.d0 + r(n) + 2*lambda_mu_T18*(1+r(n))/(z(n)**2*(3*r(n)+2)) !mu2
      END IF
      END DO
      RETURN
  3. The function grad_mu_bar_in_fourier_space_anal_T18 defines the gradient of the magnetic permeability on gauss points.
    1. We define the radial and vertical cylindrical coordinates.
      r=pt(1)
      z=pt(2)
    2. We define the gradient of the magnetic permeability.
      IF(pt_id(1)==1) THEN !grad_mu_1
      vv(1)= 1.d0
      vv(2)= 0.d0
      ELSE !grad_mu_2
      vv(1)=1.d0 + ( (3*r+2)*(2*lambda_mu_T18) - ( (2*lambda_mu_T18*(1+r)))*(3) ) /( z*(3*r+2) )**2
      vv(2)= (2*lambda_mu_T18*(1+r))/(3*r+2)*(-2/z**3)
      END IF
      RETURN

Setting in the data file

We describe the data file of this test. It is called debug_data_test_18 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_18_0.05.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 \(1\) 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 Fourier mode \(m=0\).
  7. We approximate the Maxwell equations by setting:
    ===Problem type: (nst, mxw, mhd, fhd)
    'mxw'
  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 \(10^{40}\) and solve the problem over \(1\) time iteration.
    ===Time step and number of time iterations
    1d40, 1
  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 Maxwell equations.
    ===Number of subdomains in magnetic field (H) mesh
    2
    ===List of subdomains for magnetic field (H) mesh
    1 2
  11. We set the number of interface in H_mesh and give their respective labels.
    ===Number of interfaces in H mesh
    1
    ===List of interfaces in H mesh
    5
    Such interfaces represent interfaces with discontinuities in magnetic permeability or interfaces between the magnetic field mesh and the temperature or the velocity field meshes.
  12. 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
    2 3 4
  13. The magnetic permeability is defined with a function of the file condlim_test_18.f90.
    ===Is permeability defined analytically (true/false)?
    .t.
  14. We construct the magnetic permeability on the gauss points without using its value on the nodes and a finite element interpolation.
    ===Use FEM Interpolation for magnetic permeability (true/false)?
    .f.
  15. We set the conductivity in each domains where the magnetic field is approximated.
    ===Conductivity in the conductive part (1:nb_dom_H)
    1.d0 1.d0
  16. We set the type of finite element used to approximate the magnetic field.
    ===Type of finite element for magnetic field
    2
  17. We set the magnetic Reynolds number \(\Rm\).
    ===Magnetic Reynolds number
    1.d0
  18. 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
  19. 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

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 divergence of the magnetic field \(\textbf{B}=\mu\textbf{H}\).
  2. The L2 norm of the error on the curl of the magnetic field \(\textbf{H}\).
  3. The L2 norm of the error on the magnetic field \(\textbf{H}\).
  4. The L2 norm of the magnetic filed \(\bH\).

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_18 in the directory ($SFEMaNS_DIR)/MHD_DATA_TEST_CONV_PETSC. They are equal to:

============================================
mesh_18_0.05.FEM, dt=1d40, it_max=1
===Reference results
5.053392038384935E-005 !L2 norm of div(mu (H-Hexact))
4.903853881658507E-005 !L2 norm of curl(H-Hexact)
6.193260400950516E-006 !L2 norm of H-Hexact
2.20415767167993 !L2 norm of H

To conclude this test, we show the profile of the approximated magnetic field \(\textbf{H}\) 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_18_mag_tfin.png
Magnetic field magnitude in the plane plane y=0.