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

Go to the source code of this file.

Functions

void flux (const double *state, double *flux, double *eigenvalue)
 as well as a function which, given the state of each quantity, returns the fluxes and the minimum/maximum eigenvalues (i.e.
 
double update_conservation (scalar *conserved, scalar *updates, double dtmax)
 
void event_defaults (void)
 Event: defaults (i = 0)
 
event cleanup (i=end) free(evolving)
 At the end of the run we need to free the list (to avoid a memory leak).
 
event init (i=0)
 User initialisation happens here.
 
static double riemann (const double *right, const double *left, double Delta, double *f, int len, double dtmax)
 

Variables

scalarscalars
 
vectorvectors
 
scalarevolving
 The generic time-integration scheme in predictor-corrector.h needs to know which fields are updated i.e.
 

Function Documentation

◆ cleanup()

event cleanup ( i  = end)

At the end of the run we need to free the list (to avoid a memory leak).

◆ 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.

We switch to a pure minmod limiter by default for increased robustness.

With the MUSCL scheme we use the CFL depends on the dimension of the problem.

On trees we need to replace the default bilinear refinement/prolongation with linear so that reconstructed values also use slope limiting.

Definition at line 59 of file conservation.h.

References CFL, dimension, evolving, list_concat(), refine_linear(), restriction_volume_average(), s, scalars, set_prolongation(), set_restriction(), theta, update, update_conservation(), vectors, and x.

Here is the call graph for this function:

◆ flux()

void flux ( const double s,
double f,
double e 
)

as well as a function which, given the state of each quantity, returns the fluxes and the minimum/maximum eigenvalues (i.e.

eigenvalue[0]/eigenvalue[1]) of the corresponding Riemann problem.

as well as a function which, given the state of each quantity, returns the fluxes and the minimum/maximum eigenvalues (i.e.

The parameter passed to the function is the array s which contains the state variables for each conserved field, in the order of their definition above (i.e. scalars then vectors).

We first recover each value ( \(\rho\), \(E\), \(w_x\) and \(w_y\)) and then compute the corresponding fluxes (f[0], f[1], f[2] and f[3]).

The minimum and maximum eigenvalues for the Euler system are the characteristic speeds \(u \pm \sqrt(\gamma p / \rho)\).

Definition at line 63 of file compressible.h.

References c, dimension, E, f, gammao, i, p, rho, s, sq(), un, and x.

Here is the call graph for this function:

◆ init()

event init ( i  = 0)

User initialisation happens here.

◆ riemann()

static double riemann ( const double right,
const double left,
double  Delta,
double f,
int  len,
double  dtmax 
)
static

Computing fluxes

The core of the central-upwind scheme (see e.g. section 3.1 of Kurganov & Levy, 2002) is the approximate solution of the Riemann problem given by the left and right states to get the fluxes f.

Definition at line 110 of file conservation.h.

References a, CFL, dt, dtmax, f, flux, i, left, max, min, right, and x.

Referenced by update_conservation().

Here is the caller graph for this function:

◆ update_conservation()

double update_conservation ( scalar conserved,
scalar updates,
double  dtmax 
)

The gradients of each quantity are stored in a list of dynamically-allocated fields. First-order reconstruction is used for the gradient fields.

We allocated fields for storing fluxes for each scalar and vector quantity.

The predictor-corrector scheme treats all fields as scalars (stored in the conserved list). We need to recover vector and tensor quantities from these lists. To do so, knowing the number of scalar fields, we split the scalar list into a list of scalars and a list of vectors.

We then do the same for the gradients i.e. split the list of vectors into a list of vectors and a list of tensors.

And again for the fluxes.

We are ready to compute the fluxes through each face of the domain.

Left/right state reconstruction

We use the central values of each scalar/vector quantity and the pre-computed gradients to compute the left and right states.

Riemann problem

We then call the generic Riemann solver and store the resulting fluxes in the pre-allocated fields.

Update

The update for each scalar quantity is the divergence of the fluxes.

Cleanup

We finally deallocate the memory used to store lists and gradient fields.

Definition at line 134 of file conservation.h.

References _i, cm, dimension, dtmax, dx, f, fm, free(), fs, g, gradients(), scalar::i, i, l, list_copy(), list_len(), refine_linear(), riemann(), s, scalars, t, tensors_from_vectors(), v, vector, vectors, vectors_append(), vectors_copy(), vectors_from_scalars(), vector::x, tensor::x, x, and zero().

Referenced by event_defaults().

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

Variable Documentation

◆ evolving

scalar* evolving

The generic time-integration scheme in predictor-corrector.h needs to know which fields are updated i.e.

Time-integration

Setup

Time integration will be done with a generic predictor-corrector scheme. all the scalars + the components of all the vector fields. It also needs a function to compute the time-derivatives of the evolving variables.

Definition at line 55 of file conservation.h.

Referenced by event_defaults(), and vertical_fluxes().

◆ scalars

scalar* scalars
extern

A generic solver for systems of conservation laws

Using the ideas of Kurganov and Tadmor, 2000 it is possible to write a generic solver for systems of conservation laws of the form

\[ \partial_t\left(\begin{array}{c} s_i\\ \mathbf{v}_j\\ \end{array}\right) + \nabla\cdot\left(\begin{array}{c} \mathbf{F}_i\\ \mathbf{T}_j\\ \end{array}\right) = 0 \]

where \(s_i\) is a list of scalar fields, \(\mathbf{v}_j\) a list of vector fields and \(\mathbf{F}_i\), \(\mathbf{T}_j\) are the corresponding vector (resp. tensor) fluxes.

Note that the Saint-Venant solver is a particular case of this generic algorithm.

The user must provide the lists of conserved scalar and vector fields

Definition at line 52 of file compressible.h.

Referenced by event_defaults(), and update_conservation().

◆ vectors

vector* vectors
extern

Definition at line 53 of file compressible.h.

Referenced by event_defaults(), if(), and update_conservation().