#include <TofReorder.hpp>

Inheritance diagram for Opm::TofReorder:
Inheritance graph

Public Member Functions

 TofReorder (const UnstructuredGrid &grid, const bool use_multidim_upwind=false)
 
void solveTof (const double *darcyflux, const double *porevolume, const double *source, std::vector< double > &tof)
 
void solveTofTracer (const double *darcyflux, const double *porevolume, const double *source, const SparseTable< int > &tracerheads, std::vector< double > &tof, std::vector< double > &tracer)
 

Protected Member Functions

void reorderAndTransport (const UnstructuredGrid &grid, const double *darcyflux)
 
const std::vector< int > & sequence () const
 
const std::vector< int > & components () const
 

Detailed Description

Implements a first-order finite volume solver for (single-phase) time-of-flight using reordering. The equation solved is:

\[v \cdot \nabla\tau = \phi\]

in which $ v $ is the fluid velocity, $ \tau $ is time-of-flight and $ \phi $ is the porosity. This is a boundary value problem, and $ \tau $ is specified to be zero on all inflow boundaries.

Constructor & Destructor Documentation

◆ TofReorder()

Opm::TofReorder::TofReorder ( const UnstructuredGrid &  grid,
const bool  use_multidim_upwind = false 
)

Construct solver.

Parameters
[in]gridA 2d or 3d grid.
[in]use_multidim_upwindIf true, use multidimensional tof upwinding.

Member Function Documentation

◆ components()

const std::vector< int > & Opm::ReorderSolverInterface::components ( ) const
protectedinherited

◆ reorderAndTransport()

void Opm::ReorderSolverInterface::reorderAndTransport ( const UnstructuredGrid &  grid,
const double *  darcyflux 
)
protectedinherited

◆ sequence()

const std::vector< int > & Opm::ReorderSolverInterface::sequence ( ) const
protectedinherited

◆ solveTof()

void Opm::TofReorder::solveTof ( const double *  darcyflux,
const double *  porevolume,
const double *  source,
std::vector< double > &  tof 
)

Solve for time-of-flight.

Parameters
[in]darcyfluxArray of signed face fluxes.
[in]porevolumeArray of pore volumes.
[in]sourceSource term. Sign convention is: (+) inflow flux, (-) outflow flux.
[out]tofArray of time-of-flight values.

◆ solveTofTracer()

void Opm::TofReorder::solveTofTracer ( const double *  darcyflux,
const double *  porevolume,
const double *  source,
const SparseTable< int > &  tracerheads,
std::vector< double > &  tof,
std::vector< double > &  tracer 
)

Solve for time-of-flight and a number of tracers.

Parameters
[in]darcyfluxArray of signed face fluxes.
[in]porevolumeArray of pore volumes.
[in]sourceSource term. Sign convention is: (+) inflow flux, (-) outflow flux.
[in]tracerheadsTable containing one row per tracer, and each row contains the source cells for that tracer.
[out]tofArray of time-of-flight values (1 per cell).
[out]tracerArray of tracer values. N per cell, where N is equalt to tracerheads.size().

The documentation for this class was generated from the following file: