Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
momentum.h File Reference
#include "all-mach.h"
#include "vof.h"
Include dependency graph for momentum.h:

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)
 

Functions

void event_defaults (void)
 Event: defaults (i = 0)
 
void event_properties (void)
 Event: properties (i++)
 
void event_vof (void)
 Event: vof (i++)
 
void event_tracer_advection (void)
 We set the list of interfaces back to its default value.
 

Variables

scalar f []
 
scalarinterfaces = {f}
 The height functions are stored in the vector field associated with each VOF tracer.
 
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\) and average viscosity \(\mu\) (on faces) as well as the cell-centered density.
 
vector muv []
 
scalar rhov []
 
static scalarinterfaces1 = NULL
 We overload the vof() event to transport consistently the volume fraction and the momentum of each phase.
 

Macro Definition Documentation

◆ mu

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

Definition at line 48 of file momentum.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 45 of file momentum.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.

We use (strict) minmod slope limiting for all components.

Definition at line 26 of file momentum.h.

References alpha, alphav, dimension, minmod2(), mu, muv, q, rho, rhov, theta, vector::x, and x.

Here is the call graph for this function:

◆ event_properties()

void event_properties ( void  )

Event: properties (i++)

Definition at line 52 of file momentum.h.

References _i, alphav, cm, f, fm, mu, muv, rho, rhov, vector::x, and x.

◆ event_tracer_advection()

void event_tracer_advection ( void  )

We set the list of interfaces back to its default value.

Event: tracer_advection (i++)

Definition at line 138 of file momentum.h.

References interfaces, and interfaces1.

◆ event_vof()

void event_vof ( void  )

Event: vof (i++)

We split the total momentum \(q\) into its two components \(q1\) and \(q2\) associated with \(f\) and \(1 - f\) respectively.

Momentum \(q2\) is associated with \(1 - f\), so we set the inverse attribute to true. We use the same limiting for q1 and q2.

The refinement function is modified by vof_advection(). To be able to restore it, we store its value.

We associate the transport of \(q1\) and \(q2\) with \(f\) and transport all fields consistently using the VOF scheme.

We recover the total momentum.

We restore the refinement function for the total momentum.

We set the list of interfaces to NULL so that the default vof() event does nothing (otherwise we would transport \(f\) twice).

Definition at line 69 of file momentum.h.

References _i, dimension, f, free(), i, interfaces, interfaces1, list_concat(), q, q1, q2, refine(), rho, rho1, rho2, s, set_prolongation(), tracers, u, vof_advection(), vector::x, 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\) and average viscosity \(\mu\) (on faces) as well as the cell-centered density.

Definition at line 22 of file momentum.h.

Referenced by event_defaults(), and event_properties().

◆ f

scalar f[]

Momentum-conserving formulation for two-phase interfacial flows

The interface between the fluids is tracked with a Volume-Of-Fluid method. The volume fraction in fluid 1 is \(f=1\) and \(f=0\) in fluid

  1. The densities and dynamic viscosities for fluid 1 and 2 are rho1, mu1*, rho2, mu2, respectively.

Definition at line 14 of file momentum.h.

Referenced by event_properties(), and event_vof().

◆ interfaces

scalar * interfaces = {f}

The height functions are stored in the vector field associated with each VOF tracer.

We will need basic functions for volume fraction computations.

They need to be updated every time the VOF field changes. For the centered Navier-Stokes solver, this means after initialisation and after VOF advection.

Note that strictly speaking this should be done for each sweep of the direction-split VOF advection, which we do not do here i.e. we use the normal at the beginning of the timestep and assume it is constant during each sweep. This seems to work fine.

Definition at line 14 of file momentum.h.

Referenced by event_tracer_advection(), and event_vof().

◆ interfaces1

scalar* interfaces1 = NULL
static

We overload the vof() event to transport consistently the volume fraction and the momentum of each phase.

Definition at line 66 of file momentum.h.

Referenced by event_tracer_advection(), and event_vof().

◆ mu1

double mu1 = 0.

Definition at line 15 of file momentum.h.

◆ mu2

double mu2 = 0.

Definition at line 15 of file momentum.h.

◆ muv

vector muv[]

Definition at line 22 of file momentum.h.

Referenced by event_defaults(), and event_properties().

◆ rho1

double rho1 = 1.

Definition at line 15 of file momentum.h.

Referenced by event_acceleration(), and event_vof().

◆ rho2

double rho2 = 1.

Definition at line 15 of file momentum.h.

Referenced by event_acceleration(), and event_vof().

◆ rhov

scalar rhov[]

Definition at line 23 of file momentum.h.

Referenced by event_defaults(), and event_properties().