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

Go to the source code of this file.

Macros

#define rho(f)   (clamp(f,0.,1.)*(rho1 - rho2) + rho2)
 The density and viscosity are defined using arithmetic averages by default.
 
#define mu(f)   (clamp(f,0.,1.)*(mu1 - mu2) + mu2)
 
#define sf   f
 We have the option of using some "smearing" of the density/viscosity jump.
 

Functions

void event_defaults (void)
 Event: defaults (i = 0)
 
void event_tracer_advection (void)
 Event: tracer_advection (i++)
 
void event_properties (void)
 Event: properties (i++)
 

Variables

double rho1 = 1.
 
double mu1 = 0.
 
double rho2 = 1.
 
double mu2 = 0.
 
vector alphav []
 Auxilliary fields are necessary to define the (variable) specific volume \(\alpha=1/\rho\) as well as the cell-centered density.
 
scalar rhov []
 

Macro Definition Documentation

◆ mu

#define mu (   f)    (clamp(f,0.,1.)*(mu1 - mu2) + mu2)

Definition at line 42 of file two-phase-generic.h.

◆ rho

#define rho (   f)    (clamp(f,0.,1.)*(rho1 - rho2) + rho2)

The density and viscosity are defined using arithmetic averages by default.

The user can overload these definitions to use other types of averages (i.e. harmonic).

Definition at line 39 of file two-phase-generic.h.

◆ sf

#define sf   f

We have the option of using some "smearing" of the density/viscosity jump.

Definition at line 52 of file two-phase-generic.h.

Function Documentation

◆ event_defaults()

void event_defaults ( void  )

Event: defaults (i = 0)

Boundary conditions for VOF-advected tracers usually depend on boundary conditions for the VOF field.

Event: defaults (i = 0)

Event: defaults (i = 0)

Initialisation

We set the default values.

Electrohydrodynamic stresses

The EHD force density, \(\mathbf{f}_e\), can be computed as the divergence of the Maxwell stress tensor \(\mathbf{M}\),

\[ M_{ij} = \varepsilon (E_i E_j - \frac{E^2}{2}\delta_{ij}) \]

where \(E_i\) is the \(i\)-component of the electric field, \(\mathbf{E}=-\nabla \phi\) and \(\delta_{ij}\) is the Kronecker delta.

We need to add the corresponding acceleration to the Navier–Stokes solver.

If the acceleration vector a (defined by the Navier–Stokes solver) is constant, we make it variable.

Event: defaults (i = 0)

Event: defaults (i = 0 )

Defaults

On trees we need to ensure conservation of the tracer when refining/coarsening.

Event: defaults (i = 0)

it is an acceleration. If necessary, we allocate a new vector field to store it.

Event: defaults (i = 0)

Event: defaults (i = 0)

\[ \begin{aligned} 0 & = - \sum_k \nabla \cdot [\theta_H (hu)_k^{n + 1} + (1 - \theta_H) (hu)^n_k] \\ \frac{(hu)^{n + 1}_k - (hu)_k^n}{\Delta t} & = - \Delta tgh^{n + 1 / 2}_k (\theta_H \nabla \eta_r^{n + 1} + (1 - \theta_H) \nabla \eta_r^n) \end{aligned} \]

where \(\eta_r\) is the equivalent free-surface height (i.e. pressure) applied on the rigid lid.

Event: defaults (i = 0)

Setup

The \(w_k\) and \(\phi_k\) scalar fields are allocated and the \(w_k\) are added to the list of advected tracers.

Event: defaults (i = 0)

Initial conditions

Event: defaults (i = 0)

Event: defaults (i = 0)

The default convention in Basilisk is no-flow through the boundaries of the domain, i.e. they are a streamline i.e. \(\psi=\)constant on the boundary. We set the default value for the CFL (the default in utils.h is 0.5). This is done once at the beginning of the simulation.

Event: defaults (i = 0)

Event: defaults (i = 0)

We will need to add the acceleration term \(w^2/y\) in the evolution equation for \(u_y\). If the acceleration field is not allocated yet, we do so.

Event: defaults (i = 0)

Initialisation and cleanup

We use the main time loop (in the predictor-corrector scheme) to setup the initial defaults.

If the viscosity is non-zero, we need to allocate the face-centered viscosity field.

We add the interface to the default display.

Definition at line 13 of file two-phase-generic.h.

References alpha, alphav, display(), mu, mu1, mu2, reset, rho, and rhov.

Here is the call graph for this function:

◆ event_properties()

void event_properties ( void  )

Event: properties (i++)

Definition at line 89 of file two-phase-generic.h.

References _i, alphav, cm, fm, fraction_refine(), mu, mu1, mu2, muv, rho, rhov, set_prolongation(), sf, vector::x, and x.

Here is the call graph for this function:

◆ event_tracer_advection()

void event_tracer_advection ( void  )

Event: tracer_advection (i++)

When using smearing of the density jump, we initialise sf with the vertex-average of f.

Definition at line 56 of file two-phase-generic.h.

References _i, f, refine_bilinear(), set_prolongation(), sf, and x.

Here is the call graph for this function:

Variable Documentation

◆ alphav

vector alphav[]

Auxilliary fields are necessary to define the (variable) specific volume \(\alpha=1/\rho\) as well as the cell-centered density.

Definition at line 9 of file two-phase-generic.h.

Referenced by event_defaults(), and event_properties().

◆ mu1

double mu1 = 0.

Definition at line 3 of file two-phase-generic.h.

Referenced by event_defaults(), and event_properties().

◆ mu2

double mu2 = 0.

Definition at line 3 of file two-phase-generic.h.

Referenced by event_defaults(), and event_properties().

◆ rho1

double rho1 = 1.

Definition at line 3 of file two-phase-generic.h.

◆ rho2

double rho2 = 1.

Definition at line 3 of file two-phase-generic.h.

◆ rhov

scalar rhov[]

Definition at line 10 of file two-phase-generic.h.

Referenced by event_defaults(), and event_properties().