Introduction
In this example, we check the correctness of SFEMaNS for a magnetic problem involving a conducting domain only. We consider Dirichlet boundary conditions. We use P2 finite elements for the magnetic field.
We solve the Maxwell equations:
\begin{align} \begin{cases} \partial_t (\mu^c \mathbf{H}) + \nabla \times \left(\frac{1}{\Rm \sigma}\nabla \times \mathbf{H} \right) = \nabla\times (\bu^\text{given} \times \mu^c \mathbf{H}) + \nabla \times \left(\frac{1}{\Rm \sigma} \nabla\times \mathbf{j} \right) & \text{in } \Omega_1, \\ \text{div} (\mu^c \mathbf {H}) = 0 & \text{in } \Omega_1 ,\\ \bH \times \bn = \bH_{\text{bdy}} \times \bn & \text{on } \Gamma_1,\\ \bH_{|t=0}= \bH_0, \end{cases} \end{align}
in the conducting domain \(\Omega_1= \{ (r,\theta,z) \in {R}^3 : (r,\theta,z) \in [0,0.5] \times [0,2\pi) \times [0,1]\} \) with \(\Gamma_1= \partial \Omega_1\). The data are the source term \(\mathbf{j}\), the given velocity \(\bu^\text{given}\), the boundary data \(\bH_{\text{bdy}}\), the initial datas \(\bH_0\). The parameter \(\Rm\) is the magnetic Reynolds number. The parameter \(\mu^c\) is the magnetic permeability of the conducting region \(\Omega_1\). The parameter \(\sigma\) is the conductivity in the conducting region \(\Omega_1\).
Manufactured solutions
We approximate the following analytical solutions:
\begin{align*} H_r(r,\theta,z,t) &= \frac{\alpha z}{\mu^c} \left( \cos(\theta) + \frac{r}{4}\cos(2\theta)+ \frac{r^2}{9}\cos(3\theta) \right)\cos(t) + \frac{\beta z}{\mu^c} \left(\sin(\theta) + \frac{r}{4}\sin(2\theta)+ \frac{r^2}{9}\sin(3\theta) \right)\cos(t), \\ H_{\theta}(r,\theta,z,t) &= \frac{\beta z}{\mu^c} \left( \cos(\theta) + \frac{r}{4}\cos(2\theta)+\frac{r^2}{9}\cos(3\theta) \right)\cos(t) -\frac{\alpha z}{\mu^c} \left(\sin(\theta) + \frac{r}{4}\sin(2\theta)+\frac{r^2}{9}\sin(3\theta) \right)\cos(t), \\ H_z(r,\theta,z,t) &= \frac{\alpha}{\mu^c} \left(r\cos(\theta) +\frac{r^2}{8}\cos(2\theta) + \frac{r^3}{27}\cos(3\theta) \right) \cos(t) + \frac{\beta}{\mu^c} \left(r\sin(\theta) +\frac{r^2}{8}\sin(2\theta) + \frac{r^3}{27}\sin(3\theta) \right) \cos(t), \end{align*}
where \(\alpha\) and \(\beta\) are parameters. We also set the given velocity field as follows:
\begin{align*} u^\text{given}_r(r,\theta,z,t) &= \alpha z \left(\cos(\theta) + \frac{r}{4}\cos(2\theta) + \frac{r^2}{9}\cos(3\theta) \right) + \beta z \left( \sin(\theta) + \frac{r}{4}\sin(2\theta) + \frac{r^2}{9}\sin(3\theta) \right), \\ u^\text{given}_{\theta}(r,\theta,z,t) &= \beta z \left( \cos(\theta) + \frac{r}{4}\cos(2\theta) + \frac{r^2}{9}\cos(3\theta)\right) -\alpha z \left( \sin(\theta) + \frac{r}{4}\sin(2\theta) + \frac{r^2}{9}\sin(3\theta) \right), \\ u^\text{given}_z(r,\theta,z,t) &= \alpha \left( r \cos(\theta) + \frac{r^2}{8}\cos(2\theta) + \frac{r^3}{27}\cos(3\theta)\right) +\beta \left( r \sin(\theta) + \frac{r^2}{8}\sin(2\theta) + \frac{r^3}{27}\sin(3\theta) \right). \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_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 images show the mesh for P1 finite elements.
Finite element mesh (P1).
|
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_10.f90. Here is a description of the subroutines and functions of interest.
-
First we define the parameters \(\alpha\) and \(\beta\) at the begining of the module so that every subroutines has access to these real numbers.
-
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: DO k=1,6
Hn1(:,k,i) = Hexact(H_mesh,k, H_mesh%rr,
list_mode(i), mu_H_field,
time)
IF (k<3) THEN
ENDIF
END IF
ENDDO
ENDDO
DO k=1,6
Hn(:,k,i) = Hexact(H_mesh,k, H_mesh%rr,
list_mode(i), mu_H_field,
time)
IF (k<3) THEN
ENDIF
END IF
ENDDO
ENDDO
The scalar potential \(\phi\) is initialized but we note it will not be used for this example.
-
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.
-
We check that the magnetic permeability is constant in the conducting domain.
IF (MAXVAL(mu_H_field) /= MINVAL(mu_H_field)) THEN
END IF
-
We define the magnetic permeability and the radial and vertical coordinates r, z.
muH=mu_H_field(1)
r = rr(1,:)
-
If the Fourier mode m is equal to 0, the magnetic field is set to zero.
IF (m==0) THEN
vv = 0
RETURN
END IF
-
For the Fourier modes \(m\in\{1,2,3\}\), 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
ELSEIF (TYPE == 2) THEN
vv =
beta *
z*(r**(m-1))*m
ELSEIF (TYPE ==3) THEN
vv =
beta *
z*(r**(m-1))*m
ELSEIF (TYPE == 4) THEN
ELSEIF (TYPE == 5) THEN
ELSEIF (TYPE == 6) THEN
ENDIF
vv = (vv/muH)*COS(
t)/m**3
RETURN
where \(t\) is the time.
-
The function
Vexact is used to define the given velocity field \(\bu^\text{given}\).
-
If the Fourier mode m is equal to 0,the velocity field is set to zero.
IF (m==0) THEN
vv = 0
RETURN
END IF
-
We define the radial and vertical coordinates r, z.
-
For the Fourier modes \(m\in\{1,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).
vv(:,1) =
alpha*
z*(r**(m-1))*m
vv(:,2) =
beta *
z*(r**(m-1))*m
vv(:,3) =
beta *
z*(r**(m-1))*m
vv(:,4) =-
alpha*
z*(r**(m-1))*m
VV = vv/m**3
RETURN
-
The function
Jexact_gauss is used to define the source term \(\textbf{j}\). Since \(\ROT \textbf{H}=0\) and \( \bu^\text{given} \times \textbf{H}=0 \), this source term has to satisfy the relation \( \ROT(\frac{1}{\sigma\Rm}\textbf{j})=\partial_t (\mu \textbf{H}) \). It is done by using the function Eexact_gauss as follows: vv = -sigma* Eexact_gauss(TYPE, rr, m, mu_phi, sigma, mu_H,
t)
RETURN
We note that we do not multiply by the magnetic Reynolds number \(\Rm\) because the sigma that is given to the function Jexact_gauss has already been multiplied by \(\Rm\).
-
The function
Eexact_gauss is used to define a vector, let's say \(\textbf{D}\), such that \( \ROT \textbf{D}=\partial_t (\mu^c \textbf{H})\).
-
We define the radial and vertical coordinates r, z.
-
We set Eexact_gauss to zero.
-
If the Fourier mode m is equal to 0, Eexact_gauss is already defined correcty.
-
For the Fourier modes \(m\in\{1,2,3\}\), we define it 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).
IF (TYPE == 1) THEN
vv = 0.d0
ELSEIF (TYPE == 2) THEN
vv = 0.d0
ELSEIF (TYPE ==3) THEN
vv =
alpha*(-1.d0/(m+2)*r**(m+1))
ELSEIF (TYPE == 4) THEN
vv =
beta *(-1.d0/(m+2)*r**(m+1))
ELSEIF (TYPE == 5) THEN
ELSEIF (TYPE == 6) THEN
ENDIF
RETURN
As \(\mu^c\) is set to one in the data file, we do not multiply by its value.
All the other subroutines present in the file condlim_test_10.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_10 and can be found in the following directory: ($SFEMaNS_DIR)/MHD_DATA_TEST_CONV_PETSC.
-
We use a formatted mesh by setting:
===Is mesh file formatted (true/false)?
.t.
-
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.
-
We use two processors in the meridian section. It means the finite element mesh is subdivised in two.
===Number of processors in meridian section
2
-
We solve the problem for \(3\) Fourier modes.
===Number of Fourier modes
3
-
We use \(3\) processors in Fourier space.
===Number of processors in Fourier space
3
It means that each processors is solving the problem for \(3/3=1\) Fourier modes.
-
We select specific Fourier modes to solve.
===Select Fourier modes? (true/false)
-
We give the list of the Fourier modes to solve.
===List of Fourier modes (if select_mode=.TRUE.)
1 2 3
-
We approximate the Maxwell equations by setting:
===Problem type: (nst, mxw, mhd, fhd)
'mxw'
-
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\). We note that we need to say if we use a restart for the velocity. If yes, the \(\bu^\text{given}\) comes from a suite file and not the function Vexact.
-
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
-
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
-
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.
-
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 4 5
-
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
-
We set the conductivity in each domains where the magnetic field is approximated.
===Conductivity in the conductive part (1:nb_dom_H)
1.d0
-
We set the type of finite element used to approximate the magnetic field.
===Type of finite element for magnetic field
2
-
We set the magnetic Reynolds number \(\Rm\).
===Magnetic Reynolds number
1.d0
-
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\).
-
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
-
We give information on how to solve the matrix associated to the time marching of the Maxwell equations.
-
===Maximum number of iterations
for Maxwell
solver100
-
===Relative tolerance
for Maxwell
solver1.d-6
===Absolute tolerance
for Maxwell
solver1.d-10
-
===Solver type for Maxwell (FGMRES, CG, ...)
GMRES
===Preconditionner type for Maxwell
solver (HYPRE, JACOBI, MUMPS...)
MUMPS
-
To get the total elapse time and the average time in loop minus initialization, we write
===Verbose divergence? (true/false)
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 three quantities:
-
The L2 norm of the error on the magnetic field \(\textbf{H}\).
-
The L2 norm of the error on the curl of the magnetic field \(\textbf{H}\).
-
The L2 norm of the divergence of the magnetic field \(\textbf{B}=\mu\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_10 in the directory ($SFEMaNS_DIR)/MHD_DATA_TEST_CONV_PETSC. They are equal to:
============================================
Mesh_10_form.FEM, P2
===Reference results
3.343870115552030E-006 L2 error on Hn
3.440137072254043E-005 L2 error on Curl(Hn)
1.051662042447933E-004 L2 norm of Div(mu Hn)
1.d0 No scalar potential
To conclude this test, we show the profile of the approximated magnetic field \(\textbf{H}\) at the final time. This figure is done in the plane \(y=0\) which is the union of the half plane \(\theta=0\) and \(\theta=\pi\).
Magnetic field magnitude in the plane plane y=0.
|