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

Go to the source code of this file.

Functions

trace mgstats diffusion (scalar f, double dt, vector D={{-1}}, scalar r={-1}, scalar beta={-1}, scalar theta={-1})
 The parameters of the diffusion() function are a scalar field f, scalar fields r and \(\beta\) defining the reactive term, the timestep dt and a vector field containing the diffusion coefficient D.
 

Function Documentation

◆ diffusion()

trace mgstats diffusion ( scalar  f,
double  dt,
vector  D = {{-1}},
scalar  r = {-1},
scalar  beta = {-1},
scalar  theta = {-1} 
)

The parameters of the diffusion() function are a scalar field f, scalar fields r and \(\beta\) defining the reactive term, the timestep dt and a vector field containing the diffusion coefficient D.

Time-implicit discretisation of reaction–diffusion equations

We want to discretise implicitly the reaction–diffusion equation

\[ \theta\partial_tf = \nabla\cdot(D\nabla f) + \beta f + r \]

where \(\beta f + r\) is a reactive term, \(D\) is the diffusion coefficient and \(\theta\) can be a density term.

Using a time-implicit backward Euler discretisation, this can be written

\[ \theta\frac{f^{n+1} - f^{n}}{dt} = \nabla\cdot(D\nabla f^{n+1}) + \beta f^{n+1} + r \]

Rearranging the terms we get

\[ \nabla\cdot(D\nabla f^{n+1}) + (\beta - \frac{\theta}{dt}) f^{n+1} = - \frac{\theta}{dt}f^{n} - r \]

This is a Poisson–Helmholtz problem which can be solved with a multigrid solver. If D or \(\theta\) are omitted they are set to one. If \(\beta\) is omitted it is set to zero. Both D and \(\beta\) may be constant fields.

Note that the r, \(\beta\) and \(\theta\) fields will be modified by the solver.

The function returns the statistics of the Poisson solver.

If dt is zero we don't do anything.

We define \(f\) and \(r\) for convenience.

We define a (possibly constant) field equal to \(\theta/dt\).

We use r to store the r.h.s. of the Poisson–Helmholtz solver.

If \(\beta\) is provided, we use it to store the diagonal term \(\lambda\).

Finally we solve the system.

Definition at line 42 of file diffusion.h.

Referenced by event_tracer_diffusion().

Here is the caller graph for this function: