Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
redistance.h File Reference
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Macros

#define BGHOSTS   2
 

Functions

static double minmod3 (double a, double b)
 
static void dphidt (scalar phi, scalar dphi, scalar phi0, const double cfl, const double phixxmin)
 
trace int redistance (scalar phi, int imax=1, double cfl=0.5, int order=3, double eps=1e-6, double band=HUGE, scalar resf={-1}, const double phixxmin=1./HUGE)
 

Macro Definition Documentation

◆ BGHOSTS

#define BGHOSTS   2

Definition at line 73 of file redistance.h.

Function Documentation

◆ dphidt()

static void dphidt ( scalar  phi,
scalar  dphi,
scalar  phi0,
const double  cfl,
const double  phixxmin 
)
static

Time derivative

This function fills dphi with \(\partial_t\phi\) .

We first calculate the inputs of the Hamiltonian

\[ D^+\phi_{ij} = \dfrac{\phi_{i+1,j} - \phi_{ij}}{\Delta x} - \dfrac{\Delta} {2}\text{minmod}(D_{xx}\phi_{ij},D_{xx}\phi_{i+1,j}) \]

\[ D^-\phi_{ij} = \dfrac{\phi_{i,j} - \phi_{i-1,j}}{\Delta x} - \dfrac{\Delta} {2}\text{minmod}(D_{xx}\phi_{ij},D_{xx}\phi{i-1,j}) \]

where

\[ D_{xx}\phi_{ij} = \dfrac{\phi_{i-1,j} - 2\phi_{ij} + \phi_{i+1,j}}{\Delta x^2} \]

We check for interfacial cells.

Near the interface, i.e. for cells where

\[ \phi^0_i\phi^0_{i+1} \leq 0 \text{ or } \phi^0_i\phi^0_{i-1} \leq 0 \]

The scheme must stay truly upwind, meaning that the movement of the 0 level-set of the function must be as small as possible. Therefore the upwind numerical scheme is modified to

\[ D_x^+ = \dfrac{0-\phi_{ij}}{\Delta x^+} - \dfrac{\Delta x^+}{2} \text{minmod}(D_ {xx}\phi_{ij},D_{xx}\phi_{i+1,j}) \text{ if } \phi_{ij}\phi_{i+1,j} < 0 \]

\[ D_x^- = \dfrac{\phi_{ij}-0}{\Delta x^-} + \dfrac{\Delta x^-}{2} \text{minmod}(D_ {xx}\phi_{ij},D_{xx}\phi_{i-1,j}) \text{ if } \phi_{ij}\phi_{i+1,j} < 0 \]

which is the correction by Min & Gibou 2007.

We compute the subcell fix near the interface (see section 2.2 in Min, 2010).

\[ \Delta x^+ = \left\{ \begin{array}{ll} \Delta x \cdot \left( \frac{1}{2} + \dfrac{\phi^0_{i,j}-\phi^0_{i+1,j} - \text{sgn}(\phi^0_{i,j}-\phi^0_{i+1,j})\sqrt{D}}{}\right) \text{ if } \left| \phi^0_{xx}\right| >\epsilon \\ \Delta x \cdot \dfrac{\phi^0_{ij}}{\phi^0_{i,j}-\phi^0_{i+1,j}} \text{ else.} \\ \end{array} \right. \]

with

\[ \phi_{xx}^0 = \text{minmod}(\phi^0_{i-1,j}-2\phi^0_{ij}+\phi^0_{i+1,j}, \phi^0_{i,j}-2\phi^0_{i+1j}+\phi^0_{i+2,j}) \\ D = \left( \phi^0_{xx}/2 - \phi_{ij}^0 - \phi_{i+1,j} \right)^2 - 4\phi_{ij}^0\phi_{i+1,j}^0 \]

For the \(\Delta x^-\) calculation, replace all the \(+\) subscript by \(-\), this is dealt with properly with the j parameter below.

The Godunov Hamiltonian is

\[ H_G(a,b,c,d) = \left\{ \begin{array}{ll} \sqrt{\text{max}((a^-)^2,(b^+)^2 + (c^-)^2,(d^+)^2)} \text { when } \text{sign}(\phi^0_{ij}) \geq 0\\ \sqrt{\text{max}((a^+)^2,(b^-)^2 + (c^+)^2,(d^-)^2)} \text { when } \text{sign}(\phi^0_{ij}) < 0 \end{array} \right. \]

Definition at line 81 of file redistance.h.

References _i, a, b, cfl, D, dimension, dt, dx, HUGE, i, interfacial(), j, max, min, minmod3(), phi, sign2(), size, sq(), x, and coord::x.

Here is the call graph for this function:

◆ minmod3()

static double minmod3 ( double  a,
double  b 
)
inlinestatic

Redistancing of a distance field

The original implementation was by Alexandre Limare used in particular in Limare et al., 2022.

This file implements redistancing of a distance field with subcell correction, see the work of Russo & Smereka, 2000 with corrections by Min & Gibou, 2007 and by Min, 2010.

Let \(\phi\) be a function close to a signed function that has been perturbed by numerical diffusion (more precisely, a non-zero tangential velocity). By iterating on the eikonal equation

\[ \left\{\begin{array}{ll} \phi_t + \text{sign}(\phi^{0}) \left(\left| \nabla \phi\right| - 1 \right) = 0\\ \phi(x,0) = \phi^0(x) \end{array} \right. \]

we can correct or redistance \(\phi\) to make it a signed function.

We use a Godunov Hamiltonian approximation for \(\left| \nabla \phi\right|\)

\[ \left| \nabla \phi \right|_{ij} = H_G(D_x^+\phi_{ij}, D_x^-\phi_{ij}, D_y^+\phi_{ij}, D_y^-\phi_{ij}) \]

where \(D^\pm\phi_{ij}\) denotes the one-sided ENO difference finite difference in the x- direction

\[ D_x^+ = \dfrac{\phi_{i+1,j}-\phi_{i,j}}{\Delta} - \dfrac{\Delta}{2}\text{minmod}(D_{xx}\phi_{ij}, D_{xx}\phi_{i+1,j}) \]

\[ D_x^- = \dfrac{\phi_{i,j}-\phi_{i-1,j}}{\Delta} + \dfrac{\Delta}{2}\text{minmod}(D_{xx}\phi_{ij}, D_{xx}\phi_{i+1,j}) \]

here \(D_{xx}\phi_{ij} = (\phi_{i-1,j} - 2\phi{ij} + \phi_{i+1,j})/\Delta^2\).

The minmod function is zero when the two arguments have different signs, and takes the argument with smaller absolute value when the two have the same sign.

The Godunov Hamiltonian \(H_G\) is given as

\[ H_G(a,b,c,d) = \left\{ \begin{array}{ll} \sqrt{\text{max}((a^-)^2,(b^+)^2 + (c^-)^2,(d^+)^2)} \text { when } \text{sign}(\phi^0_{ij}) \geq 0\\ \sqrt{\text{max}((a^+)^2,(b^-)^2 + (c^+)^2,(d^-)^2)} \text { when } \text{sign}(\phi^0_{ij}) < 0 \end{array} \right. \]

with

\[ x^+ = \text{max}(0, x)\\ x^- = \text{min}(0, x)\\ \]

We use a minmod limiter.

Definition at line 66 of file redistance.h.

References a, b, and x.

Referenced by dphidt().

Here is the caller graph for this function:

◆ redistance()

trace int redistance ( scalar  phi,
int  imax = 1,
double  cfl = 0.5,
int  order = 3,
double  eps = 1e-6,
double  band = HUGE,
scalar  resf = {-1},
const double  phixxmin = 1./HUGE 
)

The redistance() function

We create phi0[], a copy of the initial level-set function before the iterations.

Time integration iteration loop.

We use either a RK2 scheme...

... or a RK3 compact scheme from Shu and Osher, 1988.

Definition at line 224 of file redistance.h.

Referenced by event_properties().

Here is the caller graph for this function: