|
Basilisk CFD
Adaptive Cartesian mesh PDE framework
|
#include <viscosity-embed.h>
Data Fields | |
| vector | mu |
| scalar | rho |
| double | dt |
| double(* | embed_flux )(Point, scalar, vector, double *) |
We consider the Stokes equations
\[ \rho \mathbf{u}_t = \rho \mathbf{g} + \nabla \cdot [\mu ( \nabla \mathbf{u} + \nabla^T \mathbf{u})] \]
with \(\mathbf{g}\) the acceleration and \(T\) the transpose.
In the case of incompressible flow and uniform viscosity, the viscous term can be further simplified. Since
\[ \nabla \cdot (\mu \nabla \mathbf{u}) = (\nabla \mu) \cdot \nabla \mathbf{u} + \mu \nabla (\nabla \cdot \mathbf{u}) = 0 \]
the viscous term reduces to
\[ \nabla \cdot [\mu ( \nabla \mathbf{u} + \nabla^T \mathbf{u})] = \nabla \cdot (\mu \nabla^T \mathbf{u}) = \nabla^2 (\mu \mathbf{u}) \]
and the equations for each velocity component are decoupled. For each component \(v_i\), the following discrete implicit equation is solved using the multigrid solver
\[ \frac{\Delta t} \nabla (\mu \nabla v_i^{n+1}) - (\rho + \lambda_i) v_i^{n+1} + \underbrace{(\rho v_i^{n} + g_i\, Delta t)}_{r_i} = 0 \]
\(\lambda_i\) is a possible extra term due to the metric.
This header file implicitly solves the viscous diffusion equation:
\[ \rho\frac{\partial\boldsymbol{u}}{\partial t} = \boldsymbol{\nabla}\cdot\left[\mu\left(\boldsymbol{\nabla u} + \boldsymbol{\nabla u}^T\right)\right]. \]
Temporally discretised, this equation reads,
\[ -\frac{\rho}{\Delta t}\boldsymbol{u}_{n+1} + \boldsymbol{\nabla}\cdot\left[\mu\left(\boldsymbol{\nabla u} + \boldsymbol{\nabla u}^T\right)\right]_{n+1} = -\frac{\rho}{\Delta t}\boldsymbol{u}_n, \]
and thus has the form,
\[ L(\boldsymbol{a}) = \boldsymbol{\nabla}\cdot\left[\mu\left(\boldsymbol{\nabla a} + \boldsymbol{\nabla a}^T\right)\right] = \boldsymbol{b}, \]
where \(L()\) is a linear operator, and \(\boldsymbol{a}\) and \(\boldsymbol{b}\) are vectors. This system of mutually coupled equations can therefore be solved efficiently using a multigrid solver, described for the SGN equations in Popinet, 2015. When solving time-dependent problems, a good initial guess \(\tilde{\boldsymbol{a}} = \boldsymbol{a} - d\boldsymbol{a}\) is available, where \(d\boldsymbol{a}\) is an unknown correction. Therefore, it is usually more efficient to solve for the equivalent problem,
\[ L(d\boldsymbol{a}) = \boldsymbol{b} - L(\tilde{\boldsymbol{a}}) = \boldsymbol{res}, \]
where \(\boldsymbol{res}\) is the residual. Owing to the linearity of the operator \(L()\), \(d\boldsymbol{a}\) can be added to the initial guess \(\tilde{\boldsymbol{a}}\), and the process is then repeated until the residual falls below a given tolerance. The procedure can be summarised by the following steps:
#. Compute the residual \(\boldsymbol{res} = \boldsymbol{b} - L(\tilde{\boldsymbol{a}})\). #. If \(\left\lVert\boldsymbol{res}\right\rVert < \epsilon\), \(\tilde{\boldsymbol{a}}\) is good enough, stop. #. Else, solve \(L(d\boldsymbol{a})\simeq\boldsymbol{res}\). #. Add \(d\boldsymbol{a}\) to \(\tilde{\boldsymbol{a}}\) and go back to step 1.
This generic strategy is implemented in the standard Poisson solver. We also define a data structure for the main parameters of the viscous problem.
Definition at line 35 of file viscosity-embed.h.
| double Viscosity::dt |
Definition at line 38 of file viscosity-embed.h.
Definition at line 39 of file viscosity-embed.h.
Referenced by viscosity().
| vector Viscosity::mu |
Definition at line 36 of file viscosity-embed.h.
| scalar Viscosity::rho |
Definition at line 37 of file viscosity-embed.h.