// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
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 2 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 .
Consult the COPYING file in the top-level source directory of this
module for the precise wording of the license and the list of
copyright holders.
*/
/*!
* \file
*
* \copydoc Opm::VcfvDiscretization
*/
#ifndef EWOMS_VCFV_DISCRETIZATION_HH
#define EWOMS_VCFV_DISCRETIZATION_HH
#include
#include "vcfvproperties.hh"
#include "vcfvstencil.hh"
#include "p1fegradientcalculator.hh"
#include "vcfvgridcommhandlefactory.hh"
#include "vcfvbaseoutputmodule.hh"
#include
#include
#if HAVE_DUNE_FEM
#include
#include
#endif
namespace Opm {
template
class VcfvDiscretization;
} // namespace Opm
namespace Opm::Properties {
//! Set the stencil
template
struct Stencil
{
private:
typedef GetPropType GridView;
typedef typename GridView::ctype CoordScalar;
public:
typedef Opm::VcfvStencil type;
};
//! Mapper for the degrees of freedoms.
template
struct DofMapper { using type = GetPropType; };
//! The concrete class which manages the spatial discretization
template
struct Discretization { using type = Opm::VcfvDiscretization; };
//! The base class for the output modules (decides whether to write
//! element or vertex based fields)
template
struct DiscBaseOutputModule
{ using type = Opm::VcfvBaseOutputModule; };
//! Calculates the gradient of any quantity given the index of a flux approximation point
template
struct GradientCalculator
{ using type = Opm::P1FeGradientCalculator; };
//! The class to create grid communication handles
template
struct GridCommHandleFactory
{ using type = Opm::VcfvGridCommHandleFactory; };
//! Use two-point gradients by default for the vertex centered finite volume scheme.
template
struct UseP1FiniteElementGradients { static constexpr bool value = false; };
#if HAVE_DUNE_FEM
//! Set the DiscreteFunctionSpace
template
struct DiscreteFunctionSpace
{
private:
typedef GetPropType Scalar;
typedef GetPropType GridPart;
enum { numEq = getPropValue() };
typedef Dune::Fem::FunctionSpace FunctionSpace;
public:
// Lagrange discrete function space with unknowns at the cell vertices
typedef Dune::Fem::LagrangeDiscreteFunctionSpace< FunctionSpace, GridPart, 1 > type;
};
#endif
//! Set the border list creator for vertices
template
struct BorderListCreator
{ private:
typedef GetPropType VertexMapper;
typedef GetPropType GridView;
public:
typedef Opm::Linear::VertexBorderListFromGrid type;
};
//! For the vertex centered finite volume method, ghost and overlap elements must _not_
//! be assembled to avoid accounting twice for the fluxes over the process boundary faces
//! of the local process' grid partition
template
struct LinearizeNonLocalElements { static constexpr bool value = false; };
} // namespace Opm::Properties
namespace Opm {
/*!
* \ingroup VcfvDiscretization
*
* \brief The base class for the vertex centered finite volume discretization scheme.
*/
template
class VcfvDiscretization : public FvBaseDiscretization
{
typedef FvBaseDiscretization ParentType;
typedef GetPropType Implementation;
typedef GetPropType DofMapper;
typedef GetPropType GridView;
typedef GetPropType Simulator;
enum { dim = GridView::dimension };
public:
VcfvDiscretization(Simulator& simulator)
: ParentType(simulator)
{ }
/*!
* \brief Returns a string of discretization's human-readable name
*/
static std::string discretizationName()
{ return "vcfv"; }
/*!
* \brief Returns the number of global degrees of freedom (DOFs) due to the grid
*/
size_t numGridDof() const
{ return static_cast(this->gridView_.size(/*codim=*/dim)); }
/*!
* \brief Mapper to convert the Dune entities of the
* discretization's degrees of freedoms are to indices.
*/
const DofMapper& dofMapper() const
{ return this->vertexMapper(); }
/*!
* \brief Serializes the current state of the model.
*
* \tparam Restarter The type of the serializer class
*
* \param res The serializer object
*/
template
void serialize(Restarter& res)
{ res.template serializeEntities*codim=*/dim>(asImp_(), this->gridView_); }
/*!
* \brief Deserializes the state of the model.
*
* \tparam Restarter The type of the serializer class
*
* \param res The serializer object
*/
template
void deserialize(Restarter& res)
{
res.template deserializeEntities*codim=*/dim>(asImp_(), this->gridView_);
this->solution(/*timeIdx=*/1) = this->solution(/*timeIdx=*/0);
}
private:
Implementation& asImp_()
{ return *static_cast(this); }
const Implementation& asImp_() const
{ return *static_cast(this); }
};
} // namespace Opm
#endif