Rodin::Advection namespace

Module which provides models for the resolution of unsteady advection equations.

The Advection module implements semi-Lagrangian (characteristic-based) transport of scalar fields on unstructured meshes. The primary class is Lagrangian, which solves the advection equation:

\[ \frac{\partial \phi}{\partial t} + \mathbf{v} \cdot \nabla \phi = 0 \]

Algorithm

At each time step, the Lagrangian method:

  1. Traces characteristics backward from each mesh vertex along the velocity field $ \mathbf{v} $ for a duration $ \Delta t $ using the Flow operator
  2. Interpolates the solution at the departure points
  3. Projects the interpolated values back onto the finite element space via an $ L^2 $ projection (mass-matrix solve using CG)

The variational formulation at each step is:

\[ \int_\Omega u^{n+1} v \, dx = \int_\Omega u^n(\mathbf{X}) \, v \, dx \]

where $ \mathbf{X} $ is the departure point traced backward along $ \mathbf{v} $ .

Time Stepping and Boundary Policies

The characteristic tracing uses configurable time-stepping schemes. The default is Runge-Kutta 4 (Math::RungeKutta::RK4); custom steppers can be passed as the fifth constructor argument.

Two boundary policies handle characteristics that hit the domain boundary:

PolicyBehavior
StopInsideBoundaryPolicyClamps the departure point to an interior point with a small inward offset (eps_in_phys, default $ 10^{-12} $ ). Used automatically for the first time step (from initial condition).
TaylorBoundaryShiftPolicy<Field>Applies a first-order Taylor correction at the boundary hit: shifts the point inward and accumulates the gradient correction $ \nabla\phi \cdot \delta\mathbf{x} $ . Used automatically for subsequent steps (from the current solution).

The Lagrangian class tracks time internally (m_t) and switches between the initial condition (first step) and the current solution (later steps).

Typical Usage

P1 Vh(mesh);
TrialFunction u(Vh);
TestFunction  v(Vh);

GridFunction phi(Vh);
// ... initialize phi ...

VectorFunction velocity{1.0, 0.0};  // uniform horizontal advection

// Default stepper: RK4
Advection::Lagrangian advection(u, v, phi, velocity);

Real dt = 0.01;
for (int step = 0; step < 100; ++step)
  advection.step(dt);

Classes

template<class ... Params>
class Lagrangian
Lagrangian variational advection for scalar fields.
template<class FES, class Data, class Initial, class VectorField, class Step>
class Lagrangian<Variational::TrialFunction<Variational::GridFunction<FES, Data>, FES>, Variational::TestFunction<FES>, Initial, VectorField, Step>
Lagrangian advection specialization for trial/test function formulation.
class StopInsideBoundaryPolicy
Boundary policy for characteristic tracing: stop-at-boundary with inward offset.
template<class Field>
class TaylorBoundaryShiftPolicy
First-order boundary-hit handler for semi-Lagrangian tracing with an inward shift.

Functions

template<class FES, class Data, class Initial, class VVel>
Lagrangian(Variational::TrialFunction<Variational::GridFunction<FES, Data>, FES>&, Variational::TestFunction<FES>&, Initial&&, VVel&&) -> Lagrangian< Variational::TrialFunction< Variational::GridFunction< FES, Data >, FES >, Variational::TestFunction< FES >, Initial, VVel, Math::RungeKutta::RK4 >
Deduction guide for Lagrangian with default RK4 stepper.
template<class FES, class Data, class Initial, class VVel, class SStep>
Lagrangian(Variational::TrialFunction<Variational::GridFunction<FES, Data>, FES>&, Variational::TestFunction<FES>&, Initial&&, VVel&&, SStep&&) -> Lagrangian< Variational::TrialFunction< Variational::GridFunction< FES, Data >, FES >, Variational::TestFunction< FES >, Initial, VVel, SStep >
Deduction guide for Lagrangian with custom stepper.

Function documentation

template<class FES, class Data, class Initial, class VVel>
Rodin::Advection::Lagrangian(Variational::TrialFunction<Variational::GridFunction<FES, Data>, FES>&, Variational::TestFunction<FES>&, Initial&&, VVel&&) -> Lagrangian< Variational::TrialFunction< Variational::GridFunction< FES, Data >, FES >, Variational::TestFunction< FES >, Initial, VVel, Math::RungeKutta::RK4 >

Deduction guide for Lagrangian with default RK4 stepper.

Allows construction without explicitly specifying the Step template parameter.

template<class FES, class Data, class Initial, class VVel, class SStep>
Rodin::Advection::Lagrangian(Variational::TrialFunction<Variational::GridFunction<FES, Data>, FES>&, Variational::TestFunction<FES>&, Initial&&, VVel&&, SStep&&) -> Lagrangian< Variational::TrialFunction< Variational::GridFunction< FES, Data >, FES >, Variational::TestFunction< FES >, Initial, VVel, SStep >

Deduction guide for Lagrangian with custom stepper.

Allows construction with explicit time-stepping scheme specification.