version 0.6.0
Energy equation

Namespaces

module  enum_energy_equation_formulation
 Enumerations associated to the energy equation.
 
module  fields_energy
 Contains the field arrays associated to the Energy equation.
 
module  mod_solve_cell_allencahn_equation
 Solve the Allen-Cahn equation (experimental module, not followed anymore)
 
module  variables_energy
 Declaration of scalar variables associated to the Energy equation.
 

Functions

subroutine solve_energy ()
 This routine solves the energy equation: More...
 
subroutine mod_solve_cell_allencahn_equation::solve_allencahn_equation (phase_field, phase_field_n, phase_field_nm1, time_step, time_step_n, time_order_discretization, has_ghost_boundary_cells, allencahn_solver, mobility, thickness, energy_wall, undercooling_energy)
 Solve the Allen-Cahn equation. More...
 

Detailed Description

This directory contains the routines necessary to prepare the code to solve the Energy equation, discretized in time as follow:

\[ \rho Cp \left( \frac {\alpha \mathbf{T}^{n+1} + \beta \mathbf{T}^n + \gamma \mathbf{T}^{n-1}} {\Delta t} + \mathbf{u^{n+1}} \cdot \nabla \tilde{T} \right) = \nabla \cdot \left( \lambda \nabla T^{n+1} \right) \]

where values of \( \alpha, \beta, \gamma \) helps to switch from Euler time discretization scheme of order 1 to the 2nd order backward differential one:

\( \alpha = 1, \beta = -1, \gamma = 0 \) for the Euler scheme

\( \alpha = \frac {3} {2}, \beta = -2, \gamma = \frac {1} {2} \) for 2nd order BDF

As regards the advection term, velocity at time \( t^{n+1} \) is know since the Navier-Stokes equations are solved before. This term can be treated explicitly ( \( \tilde{T}=T^{n} \)) or implicitly ( \( \tilde{T}=T^{n+1} \))

Description of the directory

Function Documentation

◆ solve_allencahn_equation()

subroutine mod_solve_cell_allencahn_equation::solve_allencahn_equation ( double precision, dimension(:,:,:), intent(inout)  phase_field,
double precision, dimension(:,:,:), intent(in)  phase_field_n,
double precision, dimension(:,:,:), intent(in), allocatable  phase_field_nm1,
double precision, intent(in)  time_step,
double precision, intent(in)  time_step_n,
integer, intent(in)  time_order_discretization,
logical, intent(in)  has_ghost_boundary_cells,
class(t_base_solver), intent(inout)  allencahn_solver,
double precision, intent(in)  mobility,
double precision, intent(in)  thickness,
double precision, intent(in)  energy_wall,
double precision, dimension(:,:,:), intent(in)  undercooling_energy 
)

The Allen-Cahn equation is

\begin{align} \frac1M \frac {\partial \phi} {\partial t} - \epsilon^2 \Delta \phi + 2 W \phi (1 - \phi) (1 - 2 \phi) + 30 δL\, \phi^2 (1 - \phi)^2 &= 0, \end{align}

where \(\phi\) is the phase field, \(M\) is the mobility, \(\epsilon\) is the interface thickness, \(W\) is the energy wall, and \(δL\) is the undercooling energy, The nonlinear part of the equation represent the local free energy \( f(\phi) = W\, g(\phi) + δL\, p(\phi) \), where \(g\) is the double-well function and \(p\) is the interpolation function

\begin{align} g'(\phi) &= 2 \phi (1 - \phi) (1 - 2 \phi), & p'(\phi) &= 30 \phi^2 (1 - \phi)^2. \end{align}

The undercooling develops as \(δL = S (T_M - T)\), where \(S\) is the melting entropy and \((T_M - T)\) is the temperature difference with respect to the melting temperature. In this implementation, \(W\) is a scalar while \(δL\) is a field that allow to account for inhomogeneities in the temperature field.

Warning
The remainder of this text is out-dated.

Discretization of the time-variation term is done through Backward Difference Formulas (BDF); see discretize_cell_temporal_term. Discretization of the diffusion term is done classically. Discretization of the double-well term is obtained explicitely

\begin{align} F_m'(\phi^{n+1}) &\approx F_m'(\phi^{n}) \end{align}

when BDF1 is in use, and

\begin{align} F_m'(\phi^{n+1}) &\approx 2F_m'(\phi^{n}) - F_m'(\phi^{n-1}) \end{align}

when BDF2 is in use.

Parameters
[in,out]phase_field\(\phi^{n+1}\)
[in]phase_field_n,phase_field_nm1\(\phi^{n}\), \(\phi^{n-1}\)
[in]time_step,time_step_n\(\Delta t^{n+1}\), \(\Delta t^{n}\)
[in]time_order_discretizationBDF scheme
[in,out]allencahn_solverLinear System solver
[in]mobility\(m\)
[in]thickness\(\epsilon\)
[in]energy_wall\(W\)
[in]undercooling_energy\(δL\)
[in]has_ghost_boundary_cells

◆ solve_energy()

subroutine solve_energy

This equation is solved thanks to the generic advection/diffusion equation. The coefficient variable of the temporal term of this equation is set to \( \rho C_p \).