integrand_ladvection Derived Type

type, public, extends(integrand_tester_object) :: integrand_ladvection

type~~integrand_ladvection~~InheritsGraph type~integrand_ladvection integrand_ladvection interpolator_object interpolator_object interpolator_object->type~integrand_ladvection interpolator type~integrand_tester_object integrand_tester_object type~integrand_tester_object->type~integrand_ladvection type~integrand_object integrand_object type~integrand_object->type~integrand_tester_object
Help


1D linear advection field.

It is a FOODIE integrand class concrete extension.

1D linear advection field

The 1D linear advection equation is a conservation law that reads as $$ \begin{matrix} u_t = R(u) \Leftrightarrow u_t = F(u)_x \\ F(u) = a * u $$ where a is scalar constant coefficient. The PDE must completed with the proper initial and boundary conditions.

Numerical grid organization

The finite volume, Godunov's like approach is employed. The conservative variables (and the primitive ones) are co-located at the cell center. The cell and (inter)faces numeration is as follow.

                cell            (inter)faces
                 |                   |
                 v                   v
     |-------|-------|-.....-|-------|-------|-------|-------|-.....-|-------|-------|-------|-.....-|-------|-------|
     | 1-Ng  | 2-Ng  | ..... |  -1   |   0   |   1   |  2    | ..... |  Ni   | Ni+1  | Ni+2  | ..... |Ni+Ng-1| Ni+Ng |
     |-------|-------|-.....-|-------|-------|-------|-------|-.....-|-------|-------|-------|-.....-|-------|-------|
    0-Ng                             -1      0       1       2      Ni-1     Ni                                    Ni+Ng

Where Ni are the finite volumes (cells) used for discretizing the domain and Ng are the ghost cells used for imposing the left and right boundary conditions (for a total of 2Ng cells).

Inherited By

type~~integrand_ladvection~~InheritedByGraph type~integrand_ladvection integrand_ladvection type~test_object test_object type~integrand_ladvection->type~test_object ladvection_0
Help

Source Code


Components

TypeVisibility AttributesNameInitial
real(kind=R_P), public :: CFL =0._R_P

CFL value.

real(kind=R_P), public :: Dx =0._R_P

Space step.

integer(kind=I_P), public :: Ng =0

Ghost cells number.

integer(kind=I_P), public :: Ni =0

Space dimension.

real(kind=R_P), public :: a =0._R_P

Advection coefficient.

character(len=99), public :: initial_state

Initial state.

class(interpolator_object), public, allocatable:: interpolator

WENO interpolator.

real(kind=R_P), public :: length =0._R_P

Domain length.

real(kind=R_P), public, allocatable:: u(:)

Integrand (state) variable.

character(len=99), public :: w_scheme =''

WENO Scheme used.

real(kind=R_P), public :: weno_eps =0._R_P

WENO epsilon to avoid division by zero, default value.

integer(kind=I_P), public :: weno_order =0

WENO reconstruction order.


Type-Bound Procedures

generic, public :: add_fast => integrand_add_integrand_fast

Overloading add_fast method.

procedure, public, pass(lhs) :: assign_integrand

= operator.

procedure, public, pass(lhs) :: assign_real

= real operator.

  • private pure subroutine assign_real(lhs, rhs)

    Assign one real to an advection field.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(inout) :: lhs

    Left hand side.

    real(kind=R_P), intent(in) :: rhs(1:)

    Right hand side.

generic, public :: assignment(=) => assign_integrand, assign_real

Overloading = assignament.

  • private pure subroutine assign_integrand(lhs, rhs)

    = operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_oscillation), intent(inout) :: lhs

    Left hand side.

    class(integrand_object), intent(in) :: rhs

    Right hand side.

  • private pure subroutine assign_real(lhs, rhs)

    = real operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_oscillation), intent(inout) :: lhs

    Left hand side.

    real(kind=R_P), intent(in) :: rhs(1:)

    Right hand side.

procedure, public, pass(self) :: compute_dx

Compute the space step by means of CFL condition.

  • private pure function compute_dx(self, Dt) result(Dx)

    Compute the space step step by means of CFL condition.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: self

    Advection field.

    real(kind=R_P), intent(in) :: Dt

    Time step.

    Return Value real(kind=R_P)

    Space step.

procedure, public, pass(self) :: description

Return an informative description of the test.

  • private pure function description(self, prefix) result(desc)

    Return informative integrator description.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: self

    Integrand.

    character(len=*), intent(in), optional :: prefix

    Prefixing string.

    Return Value character(len=:), allocatable

    Description.

procedure, public, pass(self) :: destroy

Destroy field.

  • private subroutine destroy(self)

    Destroy field.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(inout) :: self

    Advection field.

procedure, public, pass(self) :: dt => compute_dt

Compute the current time step by means of CFL condition.

  • private pure function compute_dt(self, final_time, t) result(Dt)

    Compute the current time step by means of CFL condition.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: self

    Advection field.

    real(kind=R_P), intent(in) :: final_time

    Maximum integration time.

    real(kind=R_P), intent(in), optional :: t

    Time.

    Return Value real(kind=R_P)

    Time step.

procedure, public, pass(self) :: error

Return error.

  • private pure function error(self, t, t0, U0)

    Return error.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: self

    Integrand.

    real(kind=R_P), intent(in) :: t

    Time.

    real(kind=R_P), intent(in), optional :: t0

    Initial time.

    class(integrand_object), intent(in), optional :: U0

    Initial conditions.

    Return Value real(kind=R_P), allocatable, (:)

    Error.

procedure, public, pass(self) :: exact_solution

Return exact solution.

  • private pure function exact_solution(self, t, t0, U0) result(exact)

    Return exact solution.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: self

    Integrand.

    real(kind=R_P), intent(in) :: t

    Time.

    real(kind=R_P), intent(in), optional :: t0

    Initial time.

    class(integrand_object), intent(in), optional :: U0

    Initial conditions.

    Return Value real(kind=R_P), allocatable, (:)

    Exact solution.

procedure, public, pass(self) :: export_tecplot

Export integrand to Tecplot file.

  • private subroutine export_tecplot(self, file_name, t, scheme, close_file, with_exact_solution, U0)

    Export integrand to Tecplot file.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: self

    Advection field.

    character(len=*), intent(in), optional :: file_name

    File name.

    real(kind=R_P), intent(in), optional :: t

    Time.

    character(len=*), intent(in), optional :: scheme

    Scheme used to integrate integrand.

    logical, intent(in), optional :: close_file

    Flag for closing file.

    logical, intent(in), optional :: with_exact_solution

    Flag for export also exact solution.

    class(integrand_object), intent(in), optional :: U0

    Initial conditions.

procedure, private, pass(self) :: impose_boundary_conditions

Impose boundary conditions.

  • private pure subroutine impose_boundary_conditions(self, u)

    Impose boundary conditions.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: self

    Advection field.

    real(kind=R_P), intent(inout) :: u(1-self%Ng:)

    Conservative variables.

procedure, public, pass(self) :: initialize

Initialize integrand.

  • private subroutine initialize(self, Dt)

    Initialize integrand.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(inout) :: self

    Integrand.

    real(kind=R_P), intent(in) :: Dt

    Time step.

procedure, public, pass(lhs) :: integrand_add_integrand

+ operator.

  • private pure function integrand_add_integrand(lhs, rhs) result(opr)

    + operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: lhs

    Left hand side.

    class(integrand_object), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

procedure, public, pass(opr) :: integrand_add_integrand_fast

+ fast operator.

procedure, public, pass(lhs) :: integrand_add_real

+ real operator.

  • private pure function integrand_add_real(lhs, rhs) result(opr)

    + real operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: lhs

    Left hand side.

    real(kind=R_P), intent(in) :: rhs(1:)

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

procedure, public, pass(self) :: integrand_dimension

Return integrand dimension.

  • private pure function integrand_dimension(self)

    return integrand dimension.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: self

    integrand.

    Return Value integer(kind=I_P)

    integrand dimension.

procedure, public, pass(lhs) :: integrand_multiply_integrand

* operator.

procedure, public, pass(opr) :: integrand_multiply_integrand_fast

* fast operator.

procedure, public, pass(lhs) :: integrand_multiply_real

* real operator.

  • private pure function integrand_multiply_real(lhs, rhs) result(opr)

    * real_scalar operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: lhs

    Left hand side.

    real(kind=R_P), intent(in) :: rhs(1:)

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

procedure, public, pass(lhs) :: integrand_multiply_real_scalar

* real_scalar operator.

  • private pure function integrand_multiply_real_scalar(lhs, rhs) result(opr)

    * real_scalar operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: lhs

    Left hand side.

    real(kind=R_P), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

procedure, public, pass(opr) :: integrand_multiply_real_scalar_fast

* real_scalar fast operator.

procedure, public, pass(lhs) :: integrand_sub_integrand

- operator.

  • private pure function integrand_sub_integrand(lhs, rhs) result(opr)

    - operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: lhs

    Left hand side.

    class(integrand_object), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

procedure, public, pass(lhs) :: integrand_sub_real

- real operator.

  • private pure function integrand_sub_real(lhs, rhs) result(opr)

    - real operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: lhs

    Left hand side.

    real(kind=R_P), intent(in) :: rhs(1:)

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

procedure, public, pass(opr) :: integrand_subtract_integrand_fast

- fast operator.

procedure, public, pass(lhs) :: local_error

||integrand_ladvection - integrand_ladvection|| operator.

  • private function local_error(lhs, rhs) result(error)

    Estimate local truncation error between 2 advection approximations.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: lhs

    Left hand side.

    class(integrand_object), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P)

    Error estimation.

Overloading multiply_fast method.

Overloading * operator.

  • private pure function integrand_multiply_integrand(lhs, rhs) result(opr)

    * operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_oscillation), intent(in) :: lhs

    Left hand side.

    class(integrand_object), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

  • private pure function integrand_multiply_real(lhs, rhs) result(opr)

    * real_scalar operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_oscillation), intent(in) :: lhs

    Left hand side.

    real(kind=R_P), intent(in) :: rhs(1:)

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

  • private pure function real_multiply_integrand(lhs, rhs) result(opr)

    real_scalar * operator.

    Arguments

    Type IntentOptional AttributesName
    real(kind=R_P), intent(in) :: lhs(1:)

    Left hand side.

    class(integrand_oscillation), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

  • private pure function integrand_multiply_real_scalar(lhs, rhs) result(opr)

    * real_scalar operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_oscillation), intent(in) :: lhs

    Left hand side.

    real(kind=R_P), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

  • private pure function real_scalar_multiply_integrand(lhs, rhs) result(opr)

    real_scalar * operator.

    Arguments

    Type IntentOptional AttributesName
    real(kind=R_P), intent(in) :: lhs

    Left hand side.

    class(integrand_oscillation), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

Overloading + operator.

  • private pure function integrand_add_integrand(lhs, rhs) result(opr)

    + operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_oscillation), intent(in) :: lhs

    Left hand side.

    class(integrand_object), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

  • private pure function integrand_add_real(lhs, rhs) result(opr)

    + real operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_oscillation), intent(in) :: lhs

    Left hand side.

    real(kind=R_P), intent(in) :: rhs(1:)

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

  • private pure function real_add_integrand(lhs, rhs) result(opr)

    real + operator.

    Arguments

    Type IntentOptional AttributesName
    real(kind=R_P), intent(in) :: lhs(1:)

    Left hand side.

    class(integrand_oscillation), intent(in) :: rhs

    Left hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

Overloading - operator.

  • private pure function integrand_sub_integrand(lhs, rhs) result(opr)

    - operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_oscillation), intent(in) :: lhs

    Left hand side.

    class(integrand_object), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

  • private pure function integrand_sub_real(lhs, rhs) result(opr)

    - real operator.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_oscillation), intent(in) :: lhs

    Left hand side.

    real(kind=R_P), intent(in) :: rhs(1:)

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

  • private pure function real_sub_integrand(lhs, rhs) result(opr)

    real - operator.

    Arguments

    Type IntentOptional AttributesName
    real(kind=R_P), intent(in) :: lhs(1:)

    Left hand side.

    class(integrand_oscillation), intent(in) :: rhs

    Left hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

generic, public :: operator(.lterror.) => local_error

Estimate local truncation error.

  • private pure function local_error(lhs, rhs) result(error)

    Estimate local truncation error between 2 oscillation approximations.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_oscillation), intent(in) :: lhs

    Left hand side.

    class(integrand_object), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P)

    Error estimation.

procedure, public, pass(self) :: output

Extract integrand state field.

  • private pure function output(self) result(state)

    Output the advection field state.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: self

    Advection field.

    Return Value real(kind=R_P), allocatable, (:)

    Advection state

procedure, public, pass(self) :: parse_cli

Initialize from command line interface.

  • private subroutine parse_cli(self, cli)

    Initialize from command line interface.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(inout) :: self

    Advection field.

    type(command_line_interface), intent(inout) :: cli

    Command line interface handler.

procedure, public, pass(rhs) :: real_add_integrand

real + operator.

  • private pure function real_add_integrand(lhs, rhs) result(opr)

    real + operator.

    Arguments

    Type IntentOptional AttributesName
    real(kind=R_P), intent(in) :: lhs(1:)

    Left hand side.

    class(integrand_ladvection), intent(in) :: rhs

    Left hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

procedure, public, pass(rhs) :: real_multiply_integrand

real * operator.

  • private pure function real_multiply_integrand(lhs, rhs) result(opr)

    real_scalar * operator.

    Arguments

    Type IntentOptional AttributesName
    real(kind=R_P), intent(in) :: lhs(1:)

    Left hand side.

    class(integrand_ladvection), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

procedure, public, pass(rhs) :: real_scalar_multiply_integrand

real_scalar * operator.

  • private pure function real_scalar_multiply_integrand(lhs, rhs) result(opr)

    real_scalar * operator.

    Arguments

    Type IntentOptional AttributesName
    real(kind=R_P), intent(in) :: lhs

    Left hand side.

    class(integrand_ladvection), intent(in) :: rhs

    Right hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

procedure, public, pass(rhs) :: real_sub_integrand

real - operator.

  • private pure function real_sub_integrand(lhs, rhs) result(opr)

    real - operator.

    Arguments

    Type IntentOptional AttributesName
    real(kind=R_P), intent(in) :: lhs(1:)

    Left hand side.

    class(integrand_ladvection), intent(in) :: rhs

    Left hand side.

    Return Value real(kind=R_P), allocatable, (:)

    Operator result.

procedure, private, pass(self) :: reconstruct_interfaces

Reconstruct interface states.

  • private subroutine reconstruct_interfaces(self, conservative, r_conservative)

    Reconstruct interfaces states.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: self

    Advection field.

    real(kind=R_P), intent(in) :: conservative(1-self%Ng:)

    Conservative variables.

    real(kind=R_P), intent(inout) :: r_conservative(1:,0:)

    Reconstructed conservative vars.

procedure, public, nopass :: set_cli

Set command line interface.

  • private subroutine set_cli(cli)

    Set command line interface.

    Arguments

    Type IntentOptional AttributesName
    type(command_line_interface), intent(inout) :: cli

    Command line interface handler.

procedure, private, pass(self) :: set_sin_wave_initial_state

Set initial state as a sin wave.

procedure, private, pass(self) :: set_square_wave_initial_state

Set initial state as a square wave.

generic, public :: subtract_fast => integrand_subtract_integrand_fast

Overloading subtract_fast method.

procedure, public, pass(self) :: t => dU_dt

Time derivative, residuals.

  • private function dU_dt(self, t) result(dState_dt)

    Time derivative of advection field, the residuals function.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_ladvection), intent(in) :: self

    Advection field.

    real(kind=R_P), intent(in), optional :: t

    Time.

    Return Value real(kind=R_P), allocatable, (:)

    Advection field time derivative.

procedure, public, pass(self) :: t_fast

Time derivative, residuals, fast mode.

  • private subroutine t_fast(self, t)

    Time derivative function of integrand class, i.e. the residuals function. Fast mode acting directly on self.

    Arguments

    Type IntentOptional AttributesName
    class(integrand_object), intent(inout) :: self

    Integrand field.

    real(kind=R_P), intent(in), optional :: t

    Time.

Source Code

type, extends(integrand_tester_object) :: integrand_ladvection
   !< 1D linear advection field.
   !<
   !< It is a FOODIE integrand class concrete extension.
   !<
   !<### 1D linear advection field
   !< The 1D linear advection equation is a conservation law that reads as
   !<$$
   !<\begin{matrix}
   !<u_t = R(u)  \Leftrightarrow u_t = F(u)_x \\
   !<F(u) = a * u
   !<$$
   !< where `a` is scalar constant coefficient. The PDE must completed with the proper initial and boundary conditions.
   !<
   !<#### Numerical grid organization
   !< The finite volume, Godunov's like approach is employed. The conservative variables (and the primitive ones) are co-located at
   !< the cell center. The cell and (inter)faces numeration is as follow.
   !<```
   !<                cell            (inter)faces
   !<                 |                   |
   !<                 v                   v
   !<     |-------|-------|-.....-|-------|-------|-------|-------|-.....-|-------|-------|-------|-.....-|-------|-------|
   !<     | 1-Ng  | 2-Ng  | ..... |  -1   |   0   |   1   |  2    | ..... |  Ni   | Ni+1  | Ni+2  | ..... |Ni+Ng-1| Ni+Ng |
   !<     |-------|-------|-.....-|-------|-------|-------|-------|-.....-|-------|-------|-------|-.....-|-------|-------|
   !<    0-Ng                             -1      0       1       2      Ni-1     Ni                                    Ni+Ng
   !<```
   !< Where *Ni* are the finite volumes (cells) used for discretizing the domain and *Ng* are the ghost cells used for imposing the
   !< left and right boundary conditions (for a total of *2Ng* cells).
   character(99)                           :: w_scheme=''     !< WENO Scheme used.
   integer(I_P)                            :: weno_order=0    !< WENO reconstruction order.
   real(R_P)                               :: weno_eps=0._R_P !< WENO epsilon to avoid division by zero, default value.
   real(R_P)                               :: CFL=0._R_P      !< CFL value.
   integer(I_P)                            :: Ni=0            !< Space dimension.
   integer(I_P)                            :: Ng=0            !< Ghost cells number.
   real(R_P)                               :: length=0._R_P   !< Domain length.
   real(R_P)                               :: Dx=0._R_P       !< Space step.
   real(R_P)                               :: a=0._R_P        !< Advection coefficient.
   real(R_P), allocatable                  :: u(:)            !< Integrand (state) variable.
   class(interpolator_object), allocatable :: interpolator    !< WENO interpolator.
   character(99)                           :: initial_state   !< Initial state.
   contains
      ! auxiliary methods
      procedure, pass(self), public :: destroy          !< Destroy field.
      procedure, pass(self), public :: dt => compute_dt !< Compute the current time step by means of CFL condition.
      procedure, pass(self), public ::       compute_dx !< Compute the space step by means of CFL condition.
      procedure, pass(self), public :: output           !< Extract integrand state field.
      ! integrand_tester_object deferred methods
      procedure, pass(self), public :: description    !< Return an informative description of the test.
      procedure, pass(self), public :: error          !< Return error.
      procedure, pass(self), public :: exact_solution !< Return exact solution.
      procedure, pass(self), public :: export_tecplot !< Export integrand to Tecplot file.
      procedure, pass(self), public :: initialize     !< Initialize integrand.
      procedure, pass(self), public :: parse_cli      !< Initialize from command line interface.
      procedure, nopass,     public :: set_cli        !< Set command line interface.
      ! integrand_object deferred methods
      procedure, pass(self), public :: integrand_dimension !< Return integrand dimension.
      procedure, pass(self), public :: t => dU_dt          !< Time derivative, residuals.
      ! operators
      procedure, pass(lhs), public :: local_error !<`||integrand_ladvection - integrand_ladvection||` operator.
      ! +
      procedure, pass(lhs), public :: integrand_add_integrand !< `+` operator.
      procedure, pass(lhs), public :: integrand_add_real      !< `+ real` operator.
      procedure, pass(rhs), public :: real_add_integrand      !< `real +` operator.
      ! *
      procedure, pass(lhs), public :: integrand_multiply_integrand   !< `*` operator.
      procedure, pass(lhs), public :: integrand_multiply_real        !< `* real` operator.
      procedure, pass(rhs), public :: real_multiply_integrand        !< `real *` operator.
      procedure, pass(lhs), public :: integrand_multiply_real_scalar !< `* real_scalar` operator.
      procedure, pass(rhs), public :: real_scalar_multiply_integrand !< `real_scalar *` operator.
      ! -
      procedure, pass(lhs), public :: integrand_sub_integrand !< `-` operator.
      procedure, pass(lhs), public :: integrand_sub_real      !< `- real` operator.
      procedure, pass(rhs), public :: real_sub_integrand      !< `real -` operator.
      ! =
      procedure, pass(lhs), public :: assign_integrand !< `=` operator.
      procedure, pass(lhs), public :: assign_real      !< `= real` operator.
      ! override fast operators
      ! procedure, pass(self), public :: t_fast                              !< Time derivative, residuals, fast mode.
      ! procedure, pass(opr),  public :: integrand_add_integrand_fast        !< `+` fast operator.
      ! procedure, pass(opr),  public :: integrand_multiply_integrand_fast   !< `*` fast operator.
      ! procedure, pass(opr),  public :: integrand_multiply_real_scalar_fast !< `* real_scalar` fast operator.
      ! procedure, pass(opr),  public :: integrand_subtract_integrand_fast   !< `-` fast operator.
      ! private methods
      procedure, pass(self), private :: impose_boundary_conditions    !< Impose boundary conditions.
      procedure, pass(self), private :: reconstruct_interfaces        !< Reconstruct interface states.
      procedure, pass(self), private :: set_sin_wave_initial_state    !< Set initial state as a sin wave.
      procedure, pass(self), private :: set_square_wave_initial_state !< Set initial state as a square wave.
endtype integrand_ladvection