Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
viscosity.h File Reference
#include "poisson.h"
Include dependency graph for viscosity.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Data Structures

struct  Viscosity
 

Macros

#define lambda   ((coord){1.,1.,1.})
 

Functions

static void relax_viscosity (scalar *a, scalar *b, int l, void *data)
 
static double residual_viscosity (scalar *a, scalar *b, scalar *resl, void *data)
 
trace mgstats viscosity (vector u, vector mu, scalar rho, double dt, int nrelax=4, scalar *res=NULL)
 
trace mgstats viscosity_explicit (vector u, vector mu, scalar rho, double dt)
 

Macro Definition Documentation

◆ lambda

#define lambda   ((coord){1.,1.,1.})

Axisymmetry

In Basilisk, axisymmetric simulations are treated in a 2D Cartesian formulation for code generalisation purposes. The differential element \(dV\) is then accounted for in the metric (axi.h) incorporated in both \(\rho\) and \(\mu\). The strain rate tensor \(\boldsymbol{E}\) is of 3D nature:

\[ 2\boldsymbol{E} = \boldsymbol{\nabla u} + \boldsymbol{\nabla u}^T = \begin{bmatrix} 2\frac{\partial u_r}{\partial r} & 0 & \frac{\partial u_r}{\partial z} + \frac{\partial u_z}{\partial r}\\[.1cm] 0 & 2\frac{u_r}{r} & 0\\[.1cm] \frac{\partial u_r}{\partial z} + \frac{\partial u_z}{\partial r} & 0 & 2\frac{\partial u_z}{\partial z} \end{bmatrix}, \]

so that in cylindrical coordinates, the viscous diffusion equation reads,

\[ \rho\frac{\partial u_r}{\partial t} = \frac{\partial}{\partial r}\left(2\mu\frac{\partial u_r}{\partial r}\right) + \frac{\partial}{\partial z}\left[\mu\left(\frac{\partial u_r}{\partial z} + \frac{\partial u_z}{\partial r}\right)\right] + \frac{2\mu}{r}\left(\frac{\partial u_r}{\partial r} - \frac{u_r}{r}\right), \]

\[ \rho\frac{\partial u_z}{\partial t} = \frac{\partial}{\partial r}\left[\mu\left(\frac{\partial u_r}{\partial z} + \frac{\partial u_z}{\partial r}\right)\right] + \frac{\partial}{\partial z}\left(2\mu\frac{\partial u_z}{\partial z}\right) + \frac{\mu}{r}\left(\frac{\partial u_r}{\partial z} + \frac{\partial u_z}{\partial r}\right). \]

Conventionally, in basilisk, \(\boldsymbol{e}_y = \boldsymbol{e}_r\) and \(\boldsymbol{e}_x = \boldsymbol{e}_z\). For consistency reasons with 2D Cartesian simulations, the axisymmetric viscous diffusion equation in basilisk reads,

\[ \rho'\frac{\partial u}{\partial t} = \frac{\partial}{\partial x}\left(2\mu'\frac{\partial u}{\partial x}\right) + \frac{\partial}{\partial y}\left[\mu'\left(\frac{\partial u}{\partial y} + \frac{\partial v}{\partial x}\right)\right], \]

\[ \rho'\frac{\partial v}{\partial t} = \frac{\partial}{\partial x}\left[\mu'\left(\frac{\partial u}{\partial y} + \frac{\partial v}{\partial x}\right)\right] + \frac{\partial}{\partial y}\left(2\mu'\frac{\partial v}{\partial y}\right) - \lambda^*, \]

where \(\rho' = \rho y\) and \(\mu' = \mu y\). If one performs the straightforward derivation, one finds that the equation along \(x\) is identical to the equation along \(z\) in cylindrical coordinates. This is not the case for the equations along \(y\) and \(r\). That is why the variable \(\lambda^*\) is added to the equation along \(y\) in the 2D Cartesian formulation. Performing the derivations and equating both equations yields

\[ \lambda^* = \frac{2\mu'}{y^2}v. \]

The expression under "lambda.y" below stems from \(\lambda^*\) after discretisation. In non AXI simulations, \(\rho' = \rho\), \(\mu' = \mu\) and \(\lambda^* = 0\).

Definition at line 115 of file viscosity.h.

Function Documentation

◆ relax_viscosity()

static void relax_viscosity ( scalar a,
scalar b,
int  l,
void data 
)
static

Relaxation function

This function solves for the correction \(d\boldsymbol{a}\) in step 3 of the previously described algorithm. It is analogous to the relaxation function written for the Poisson-Helmholtz equation, with the only difference that the viscous diffusion equation is vectorial and yields a system of coupled scalar equations, the number of which is the dimension of the numerical simulation. This function is passed as an argument to the multigrid cycle.

We have the option of using red/black Gauss-Seidel relaxation or "re-use as soon as computed" relaxation. On GPUs (and probably also with OpenMP) red/black Gauss-Seidel converges much better (but requires two for (int _i = 0; _i < _N; _i++) /* foreach

Definition at line 131 of file viscosity.h.

References _GPU, _i, a, b, boundary_level, data, dimension, dt, Point::i, l, lambda, level, mu, p, point, rho, sq(), u, vector, w, vector::x, and x.

Referenced by viscosity().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ residual_viscosity()

static double residual_viscosity ( scalar a,
scalar b,
scalar resl,
void data 
)
static

Residual computation

This function computes the residual \(\boldsymbol{res}\) in step 1 of the previously described algorithm. It is passed as an argument to the multigrid solver.

We manually apply boundary conditions, so that all components are treated simultaneously. Otherwise (automatic) BCs would be applied component by component before each for (int _i = 0; _i < _N; _i++) /* foreach_face

Definition at line 222 of file viscosity.h.

References _i, a, b, boundary, d, data, dimension, dt, lambda, loop, max, mu, p, rho, sq(), u, vector, vector::x, and x.

Referenced by viscosity(), and viscosity_explicit().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ viscosity()

trace mgstats viscosity ( vector  u,
vector  mu,
scalar  rho,
double  dt,
int  nrelax = 4,
scalar res = NULL 
)

User interface

A user interface is provided for the solution of the viscous diffusion equation.

Implicit treatment

The velocity field \(\boldsymbol{u}_n\) is provided as an initial guess \(\tilde{\boldsymbol{a}}\).

We need \(\mu\) and \(\rho\) on all levels of the grid.

Definition at line 307 of file viscosity.h.

References _i, dimension, dt, mg_solve(), mu, p, relax_viscosity(), residual_viscosity(), restriction, rho, u, vector::x, and x.

Here is the call graph for this function:

◆ viscosity_explicit()

trace mgstats viscosity_explicit ( vector  u,
vector  mu,
scalar  rho,
double  dt 
)

Explicit treatment

This function does not make use of the multigrid solver. Instead, it explicitly advances the velocity field in time, in only one iteration.

Definition at line 337 of file viscosity.h.

References _i, dimension, dt, mu, p, residual_viscosity(), rho, u, vector::x, and x.

Here is the call graph for this function: