2012-09-26 01:58:03 -05:00
|
|
|
/*
|
|
|
|
Copyright 2012 SINTEF ICT, Applied Mathematics.
|
|
|
|
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
*/
|
|
|
|
|
2012-11-16 09:06:01 -06:00
|
|
|
#ifndef OPM_TOFDISCGALREORDER_HEADER_INCLUDED
|
|
|
|
#define OPM_TOFDISCGALREORDER_HEADER_INCLUDED
|
2012-09-26 01:58:03 -05:00
|
|
|
|
2012-11-16 09:24:54 -06:00
|
|
|
#include <opm/core/transport/reorder/ReorderSolverInterface.hpp>
|
2013-08-08 05:32:01 -05:00
|
|
|
#include <memory>
|
2012-09-26 01:58:03 -05:00
|
|
|
#include <vector>
|
|
|
|
#include <map>
|
|
|
|
#include <ostream>
|
2012-10-16 04:11:33 -05:00
|
|
|
|
2012-09-26 01:58:03 -05:00
|
|
|
struct UnstructuredGrid;
|
|
|
|
|
|
|
|
namespace Opm
|
|
|
|
{
|
|
|
|
|
|
|
|
class IncompPropertiesInterface;
|
2012-10-16 04:11:33 -05:00
|
|
|
class VelocityInterpolationInterface;
|
2013-01-16 08:13:45 -06:00
|
|
|
class DGBasisInterface;
|
2013-01-08 06:14:26 -06:00
|
|
|
namespace parameter { class ParameterGroup; }
|
2013-04-24 03:39:50 -05:00
|
|
|
template <typename T> class SparseTable;
|
2012-09-26 01:58:03 -05:00
|
|
|
|
2012-10-09 02:54:54 -05:00
|
|
|
/// Implements a discontinuous Galerkin solver for
|
|
|
|
/// (single-phase) time-of-flight using reordering.
|
|
|
|
/// The equation solved is:
|
2013-03-21 08:51:49 -05:00
|
|
|
/// \f[v \cdot \nabla\tau = \phi\f]
|
2013-03-22 03:45:00 -05:00
|
|
|
/// in which \f$ v \f$ is the fluid velocity, \f$ \tau \f$ is time-of-flight and
|
|
|
|
/// \f$ \phi \f$ is the porosity. This is a boundary value problem, and
|
2013-03-21 08:51:49 -05:00
|
|
|
/// \f$ \tau \f$ is specified to be zero on all inflow boundaries.
|
2012-10-09 02:54:54 -05:00
|
|
|
/// The user may specify the polynomial degree of the basis function space
|
|
|
|
/// used, but only degrees 0 and 1 are supported so far.
|
2012-11-16 09:06:01 -06:00
|
|
|
class TofDiscGalReorder : public ReorderSolverInterface
|
2012-09-26 01:58:03 -05:00
|
|
|
{
|
|
|
|
public:
|
|
|
|
/// Construct solver.
|
|
|
|
/// \param[in] grid A 2d or 3d grid.
|
2013-01-08 06:34:50 -06:00
|
|
|
/// \param[in] param Parameters for the solver.
|
2013-03-21 08:51:49 -05:00
|
|
|
/// The following parameters are accepted (defaults):\n
|
2013-03-22 03:55:12 -05:00
|
|
|
/// - \c dg_degree (0) -- Polynomial degree of basis functions.
|
|
|
|
/// - \c use_tensorial_basis (false) -- Use tensor-product basis, interpreting dg_degree as
|
2013-03-21 08:51:49 -05:00
|
|
|
/// bi/tri-degree not total degree.
|
2013-03-22 03:55:12 -05:00
|
|
|
/// - \c use_cvi (false) -- Use ECVI velocity interpolation.
|
|
|
|
/// - \c use_limiter (false) -- Use a slope limiter. If true, the next three parameters are used.
|
|
|
|
/// - \c limiter_relative_flux_threshold (1e-3) -- Ignore upstream fluxes below this threshold,
|
|
|
|
/// relative to total cell flux.
|
|
|
|
/// - \c 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.
|
|
|
|
/// - \c limiter_usage ("DuringComputations") -- Usage pattern for limiter. Accepted choices are:
|
|
|
|
/// - DuringComputations -- Apply limiter to cells as they are computed,
|
2013-03-21 08:51:49 -05:00
|
|
|
/// so downstream cells' solutions may be affected
|
|
|
|
/// by limiting in upstream cells.
|
2013-03-22 03:55:12 -05:00
|
|
|
/// - AsPostProcess -- Apply in dependency order, but only after
|
2013-03-21 08:51:49 -05:00
|
|
|
/// computing (unlimited) solution.
|
2013-03-22 03:55:12 -05:00
|
|
|
/// - AsSimultaneousPostProcess -- Apply to each cell independently, using un-
|
2013-03-21 08:51:49 -05:00
|
|
|
/// limited solution in neighbouring cells.
|
2012-11-16 09:06:01 -06:00
|
|
|
TofDiscGalReorder(const UnstructuredGrid& grid,
|
2013-03-14 10:18:39 -05:00
|
|
|
const parameter::ParameterGroup& param);
|
2012-09-26 01:58:03 -05:00
|
|
|
|
2012-10-09 02:54:54 -05:00
|
|
|
|
|
|
|
/// Solve for time-of-flight.
|
2012-09-26 01:58:03 -05:00
|
|
|
/// \param[in] darcyflux Array of signed face fluxes.
|
|
|
|
/// \param[in] porevolume Array of pore volumes.
|
2012-10-09 02:54:54 -05:00
|
|
|
/// \param[in] source Source term. Sign convention is:
|
|
|
|
/// (+) inflow flux,
|
|
|
|
/// (-) outflow flux.
|
2012-09-26 01:58:03 -05:00
|
|
|
/// \param[out] tof_coeff Array of time-of-flight solution coefficients.
|
|
|
|
/// The values are ordered by cell, meaning that
|
|
|
|
/// the K coefficients corresponding to the first
|
2013-04-24 03:39:50 -05:00
|
|
|
/// cell come before the K coefficients corresponding
|
2012-09-26 01:58:03 -05:00
|
|
|
/// to the second cell etc.
|
|
|
|
/// K depends on degree and grid dimension.
|
|
|
|
void solveTof(const double* darcyflux,
|
|
|
|
const double* porevolume,
|
|
|
|
const double* source,
|
|
|
|
std::vector<double>& tof_coeff);
|
|
|
|
|
2013-04-24 03:39:50 -05:00
|
|
|
/// Solve for time-of-flight and a number of tracers.
|
|
|
|
/// \param[in] darcyflux Array of signed face fluxes.
|
|
|
|
/// \param[in] porevolume Array of pore volumes.
|
|
|
|
/// \param[in] source Source term. Sign convention is:
|
|
|
|
/// (+) inflow flux,
|
|
|
|
/// (-) outflow flux.
|
|
|
|
/// \param[in] tracerheads Table containing one row per tracer, and each
|
|
|
|
/// row contains the source cells for that tracer.
|
|
|
|
/// \param[out] tof_coeff Array 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.
|
|
|
|
/// \param[out] tracer_coeff Array 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.
|
|
|
|
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);
|
|
|
|
|
2012-09-26 01:58:03 -05:00
|
|
|
private:
|
|
|
|
virtual void solveSingleCell(const int cell);
|
|
|
|
virtual void solveMultiCell(const int num_cells, const int* cells);
|
|
|
|
|
|
|
|
private:
|
2012-10-16 04:11:33 -05:00
|
|
|
// Disable copying and assignment.
|
2012-11-16 09:06:01 -06:00
|
|
|
TofDiscGalReorder(const TofDiscGalReorder&);
|
|
|
|
TofDiscGalReorder& operator=(const TofDiscGalReorder&);
|
2012-10-16 04:11:33 -05:00
|
|
|
|
2012-12-18 07:15:31 -06:00
|
|
|
// Data members
|
2012-09-26 01:58:03 -05:00
|
|
|
const UnstructuredGrid& grid_;
|
2013-08-08 05:32:01 -05:00
|
|
|
std::shared_ptr<VelocityInterpolationInterface> velocity_interpolation_;
|
2012-12-06 03:18:57 -06:00
|
|
|
bool use_cvi_;
|
2012-12-18 07:15:31 -06:00
|
|
|
bool use_limiter_;
|
2013-01-08 04:04:43 -06:00
|
|
|
double limiter_relative_flux_threshold_;
|
|
|
|
enum LimiterMethod { MinUpwindFace, MinUpwindAverage };
|
|
|
|
LimiterMethod limiter_method_;
|
|
|
|
enum LimiterUsage { DuringComputations, AsPostProcess, AsSimultaneousPostProcess };
|
|
|
|
LimiterUsage limiter_usage_;
|
2012-09-26 01:58:03 -05:00
|
|
|
const double* darcyflux_; // one flux per grid face
|
|
|
|
const double* porevolume_; // one volume per cell
|
|
|
|
const double* source_; // one volumetric source term per cell
|
2013-08-08 05:32:01 -05:00
|
|
|
std::shared_ptr<DGBasisInterface> basis_func_;
|
2012-09-26 01:58:03 -05:00
|
|
|
double* tof_coeff_;
|
2013-04-24 03:39:50 -05:00
|
|
|
// For tracers.
|
|
|
|
double* tracer_coeff_;
|
|
|
|
int num_tracers_;
|
|
|
|
enum { NoTracerHead = -1 };
|
|
|
|
std::vector<int> tracerhead_by_cell_;
|
|
|
|
// Used by solveSingleCell().
|
|
|
|
std::vector<double> rhs_; // single-cell right-hand-sides
|
2012-09-26 06:30:54 -05:00
|
|
|
std::vector<double> jac_; // single-cell jacobian
|
2013-04-24 03:39:50 -05:00
|
|
|
std::vector<double> orig_rhs_; // single-cell right-hand-sides (copy)
|
2012-10-17 05:40:43 -05:00
|
|
|
std::vector<double> orig_jac_; // single-cell jacobian (copy)
|
2012-09-26 06:30:54 -05:00
|
|
|
std::vector<double> coord_;
|
2013-01-21 07:55:27 -06:00
|
|
|
mutable std::vector<double> basis_;
|
|
|
|
mutable std::vector<double> basis_nb_;
|
2012-09-26 06:30:54 -05:00
|
|
|
std::vector<double> grad_basis_;
|
|
|
|
std::vector<double> velocity_;
|
2013-04-23 06:38:47 -05:00
|
|
|
int num_singlesolves_;
|
|
|
|
// Used by solveMultiCell():
|
|
|
|
double gauss_seidel_tol_;
|
|
|
|
int num_multicell_;
|
|
|
|
int max_size_multicell_;
|
|
|
|
int max_iter_multicell_;
|
2012-12-18 07:15:31 -06:00
|
|
|
|
|
|
|
// Private methods
|
2013-01-08 04:04:43 -06:00
|
|
|
|
|
|
|
// Apply some limiter, writing to array tof
|
|
|
|
// (will read data from tof_coeff_, it is ok to call
|
|
|
|
// with tof_coeff as tof argument.
|
|
|
|
void applyLimiter(const int cell, double* tof);
|
2013-01-21 03:46:06 -06:00
|
|
|
void applyMinUpwindLimiter(const int cell, const bool face_min, double* tof);
|
2013-01-08 04:04:43 -06:00
|
|
|
void applyLimiterAsPostProcess();
|
|
|
|
void applyLimiterAsSimultaneousPostProcess();
|
2013-01-21 07:55:27 -06:00
|
|
|
double totalFlux(const int cell) const;
|
|
|
|
double minCornerVal(const int cell, const int face) const;
|
2013-04-24 04:27:04 -05:00
|
|
|
|
|
|
|
// Apply a simple (restrict to [0,1]) limiter.
|
|
|
|
// Intended for tracers.
|
|
|
|
void applyTracerLimiter(const int cell, double* local_coeff);
|
2012-09-26 01:58:03 -05:00
|
|
|
};
|
|
|
|
|
|
|
|
} // namespace Opm
|
|
|
|
|
|
|
|
#endif // OPM_TRANSPORTMODELTRACERTOFDISCGAL_HEADER_INCLUDED
|