Opm::TofDiscGalReorder Class Reference

#include <TofDiscGalReorder.hpp>

Inheritance diagram for Opm::TofDiscGalReorder:
Inheritance graph

Public Member Functions

 TofDiscGalReorder (const UnstructuredGrid &grid, const parameter::ParameterGroup &param)
 
void solveTof (const double *darcyflux, const double *porevolume, const double *source, std::vector< double > &tof_coeff)
 
void solveTofTracer (const double *darcyflux, const double *porevolume, const double *source, const SparseTable< int > &tracerheads, std::vector< double > &tof_coeff, std::vector< double > &tracer_coeff)
 

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 discontinuous Galerkin 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. The user may specify the polynomial degree of the basis function space used, but only degrees 0 and 1 are supported so far.

Constructor & Destructor Documentation

Opm::TofDiscGalReorder::TofDiscGalReorder ( const UnstructuredGrid grid,
const parameter::ParameterGroup param 
)

Construct solver.

Parameters
[in]gridA 2d or 3d grid.
[in]paramParameters for the solver. The following parameters are accepted (defaults):
  • dg_degree (0) – Polynomial degree of basis functions.
  • use_tensorial_basis (false) – Use tensor-product basis, interpreting dg_degree as bi/tri-degree not total degree.
  • use_cvi (false) – Use ECVI velocity interpolation.
  • use_limiter (false) – Use a slope limiter. If true, the next three parameters are used.
  • limiter_relative_flux_threshold (1e-3) – Ignore upstream fluxes below this threshold, relative to total cell flux.
  • limiter_method ("MinUpwindFace") – Limiter method used. Accepted methods are:
    • MinUpwindFace – Limit cell tof to >= inflow face tofs.
    • MinUpwindAverage – Limit cell tof to >= inflow cell average tofs.
  • limiter_usage ("DuringComputations") – Usage pattern for limiter. Accepted choices are:
    • DuringComputations – Apply limiter to cells as they are computed, so downstream cells' solutions may be affected by limiting in upstream cells.
    • AsPostProcess – Apply in dependency order, but only after computing (unlimited) solution.
    • AsSimultaneousPostProcess – Apply to each cell independently, using un- limited solution in neighbouring cells.

Member Function Documentation

const std::vector<int>& Opm::ReorderSolverInterface::components ( ) const
protectedinherited
void Opm::ReorderSolverInterface::reorderAndTransport ( const UnstructuredGrid grid,
const double *  darcyflux 
)
protectedinherited
const std::vector<int>& Opm::ReorderSolverInterface::sequence ( ) const
protectedinherited
void Opm::TofDiscGalReorder::solveTof ( const double *  darcyflux,
const double *  porevolume,
const double *  source,
std::vector< double > &  tof_coeff 
)

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]tof_coeffArray of time-of-flight solution coefficients. The values are ordered by cell, meaning that the K coefficients corresponding to the first cell come before the K coefficients corresponding to the second cell etc. K depends on degree and grid dimension.
void Opm::TofDiscGalReorder::solveTofTracer ( const double *  darcyflux,
const double *  porevolume,
const double *  source,
const SparseTable< int > &  tracerheads,
std::vector< double > &  tof_coeff,
std::vector< double > &  tracer_coeff 
)

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]tof_coeffArray of time-of-flight solution coefficients. The values are ordered by cell, meaning that the K coefficients corresponding to the first cell comes before the K coefficients corresponding to the second cell etc. K depends on degree and grid dimension.
[out]tracer_coeffArray of tracer solution coefficients. N*K per cell, where N is equal to tracerheads.size(). All K coefs for a tracer are consecutive, and all tracers' coefs for a cell come before those for the next cell.

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