Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
mac.h File Reference
#include "run.h"
#include "timestep.h"
#include "poisson.h"
Include dependency graph for mac.h:

Go to the source code of this file.

Functions

event set_dtmax (i++, last) dtmax
 
void event_stability (void)
 Event: stability (i++,last)
 
void event_advance (void)
 Event: advance (i++,last)
 
void event_projection (void)
 Event: projection (i++,last)
 

Variables

scalar p []
 The Markers-And-Cells (MAC) formulation was first described in the pioneering paper of Harlow and Welch, 1965.
 
vector u []
 
double nu = 0.
 The only parameter is the viscosity coefficient \(\nu\).
 
mgstats mgp
 
u n [left] = 0
 We need to explicitly set the zero normal velocity condition.
 
double dtmax
 

Function Documentation

◆ event_advance()

void event_advance ( void  )

Event: advance (i++,last)

We allocate a local symmetric tensor field. To be able to compute the divergence of the tensor at the face locations, we need to compute the diagonal components at the center of cells and the off-diagonal component at the vertices.

Staggering of \f$\mathbf{u}\f$ and \f$\mathbf{S}\f$

We average the velocity components at the center to compute the diagonal components.

We average horizontally and vertically to compute the off-diagonal component at the vertices.

Finally we compute

\[ \mathbf{u}_* = \mathbf{u}_n + dt\nabla\cdot\mathbf{S} \]

Definition at line 83 of file mac.h.

References _i, dimension, dt, nu, S, sq(), u, vector::x, x, and vector::y.

Here is the call graph for this function:

◆ event_projection()

void event_projection ( void  )

Event: projection (i++,last)

Projection

In a second step we compute

\[ \mathbf{u}_{n+1} = \mathbf{u}_* - \Delta t\nabla p \]

with the condition

\[ \nabla\cdot\mathbf{u}_{n+1} = 0 \]

This gives the Poisson equation for the pressure

\[ \nabla\cdot(\nabla p) = \frac{\nabla\cdot\mathbf{u}_*}{\Delta t} \]

Definition at line 140 of file mac.h.

References dt, mgp, mgstats::nrelax, p, project(), u, and unity.

Here is the call graph for this function:

◆ event_stability()

void event_stability ( void  )

Event: stability (i++,last)

We need to overload the stability event so that the CFL is taken into account (because we set stokes to true).

We first compute the minimum and maximum values of \(\alpha/f_m = 1/\rho\), as well as \(\Delta_{min}\).

The maximum timestep is set using the sum of surface tension coefficients.

Definition at line 78 of file mac.h.

References dt, dtmax, dtnext(), timestep(), and u.

Here is the call graph for this function:

◆ set_dtmax()

event set_dtmax ( i++  ,
last   
)

Variable Documentation

◆ dtmax

double dtmax

Time integration

Advection–Diffusion

In a first step, we compute \(\mathbf{u}_*\)

\[ \frac{\mathbf{u}_* - \mathbf{u}_n}{dt} = \nabla\cdot\mathbf{S} \]

with \(\mathbf{S}\) the symmetric tensor

\[ \mathbf{S} = - \mathbf{u}\otimes\mathbf{u} + \nu\nabla\mathbf{u} = \left(\begin{array}{cc} - u_x^2 + 2\nu\partial_xu_x & - u_xu_y + \nu(\partial_yu_x + \partial_xu_y)\\ \ldots & - u_y^2 + 2\nu\partial_yu_y \end{array}\right) \]

The timestep for this iteration is controlled by the CFL condition (and the timing of upcoming events).

Definition at line 73 of file mac.h.

Referenced by event_stability().

◆ mgp

mgstats mgp

Definition at line 42 of file mac.h.

Referenced by event_projection().

◆ n

u n[bottom] = 0

We need to explicitly set the zero normal velocity condition.

Definition at line 47 of file mac.h.

◆ nu

double nu = 0.

The only parameter is the viscosity coefficient \(\nu\).

The statistics for the (multigrid) solution of the Poisson problem are stored in mgp.

Definition at line 41 of file mac.h.

Referenced by event_advance().

◆ p

scalar p[]

The Markers-And-Cells (MAC) formulation was first described in the pioneering paper of Harlow and Welch, 1965.

Incompressible Navier–Stokes solver (MAC formulation)

We wish to approximate numerically the incompressible Navier–Stokes equations

\[ \partial_t\mathbf{u}+\nabla\cdot(\mathbf{u}\otimes\mathbf{u}) = -\nabla p + \nabla\cdot(\nu\nabla\mathbf{u}) \]

\[ \nabla\cdot\mathbf{u} = 0 \]

We will use the generic time loop, a CFL-limited timestep and we will need to solve a Poisson problem. It relies on a face discretisation of the velocity components u.x and u.y, relative to the (centered) pressure p. This guarantees the consistency of the discrete gradient, divergence and Laplacian operators and leads to a stable (mode-free) integration.

Definition at line 32 of file mac.h.

Referenced by event_projection().

◆ u

vector u[]

Definition at line 33 of file mac.h.

Referenced by event_advance(), event_projection(), and event_stability().