From f7fb8166ac2d7a2615fd8196c02b784eb3317336 Mon Sep 17 00:00:00 2001 From: Benjamin Faigle <benjamin.faigle@posteo.de> Date: Thu, 4 Oct 2012 14:26:59 +0000 Subject: [PATCH] Removed headercheck errors by forward declaration of PropTags. Finished doxygen documentation of the pressure modules. git-svn-id: svn://svn.iws.uni-stuttgart.de/DUMUX/dumux/trunk@9195 2fb0f335-1f38-0410-981e-8018bf24f1b0 --- .../decoupled/2p2c/2p2cadaptiveproperties.hh | 1 + dumux/decoupled/2p2c/fvpressure2p2c.hh | 43 +++++------ .../decoupled/2p2c/fvpressure2p2cadaptive.hh | 52 +++++++++----- .../2p2c/fvpressure2p2cmultiphysics.hh | 71 +++++++++++++++++-- .../decoupled/2p2c/fvtransport2p2cadaptive.hh | 6 +- 5 files changed, 130 insertions(+), 43 deletions(-) diff --git a/dumux/decoupled/2p2c/2p2cadaptiveproperties.hh b/dumux/decoupled/2p2c/2p2cadaptiveproperties.hh index 7ee1a32a03..fc964db411 100644 --- a/dumux/decoupled/2p2c/2p2cadaptiveproperties.hh +++ b/dumux/decoupled/2p2c/2p2cadaptiveproperties.hh @@ -29,6 +29,7 @@ #define DUMUX_2P2CADAPTIVE_PROPERTIES_HH #include <dumux/decoupled/2p2c/2p2cproperties.hh> +#include <dumux/decoupled/common/fv/mpfa/fvmpfaproperties.hh> namespace Dumux { diff --git a/dumux/decoupled/2p2c/fvpressure2p2c.hh b/dumux/decoupled/2p2c/fvpressure2p2c.hh index 2f6360d2cf..fa87f2a14b 100644 --- a/dumux/decoupled/2p2c/fvpressure2p2c.hh +++ b/dumux/decoupled/2p2c/fvpressure2p2c.hh @@ -181,8 +181,8 @@ public: protected: Problem& problem_; - Scalar maxError_; - bool enableVolumeIntegral; + Scalar maxError_; //!> Maximum volume error of all cells + bool enableVolumeIntegral; //!> Enables the volume integral of the pressure equation Scalar ErrorTermFactor_; //!< Handling of error term: relaxation factor Scalar ErrorTermLowerBound_; //!< Handling of error term: lower bound for error dampening Scalar ErrorTermUpperBound_; //!< Handling of error term: upper bound for error dampening @@ -198,7 +198,7 @@ private: }; -//! function which assembles the system of equations to be solved +//! Assembles the source term /** for first == true, a source is implemented as in FVPressure2P. * for first == false, the source is translated into a volumentric source term: * \f[ V_i \sum_{\kappa} \frac{\partial v_{t}}{\partial C^{\kappa}} q^{\kappa}_i \f]. @@ -241,7 +241,8 @@ void FVPressure2P2C<TypeTag>::getSource(Dune::FieldVector<Scalar, 2>& sourceEntr return; } -//! function which assembles the system of equations to be solved + +//! Assembles the storage term /** for first == true, there is no storage contribution. * for first == false, the storage term comprises the compressibility (due to a change in * pressure from last timestep): @@ -328,18 +329,18 @@ void FVPressure2P2C<TypeTag>::getStorage(Dune::FieldVector<Scalar, 2>& storageEn //! Get flux at an interface between two cells /** for first == true, the flux is calculated in traditional fractional-flow forn as in FVPressure2P. - * for first == false, the flux thorugh \f$ \gamma{ij} \f$ is calculated via a volume balance formulation - * \f[ - A_{\gamma_{ij}} \cdot \mathbf{K} \cdot \mathbf{u} \cdot (\mathbf{n}_{\gamma_{ij}} \cdot \mathbf{u}) - \sum_{\alpha} \varrho_{\alpha} \lambda_{\alpha} \sum_{\kappa} \frac{\partial v_{t}}{\partial C^{\kappa}} X^{\kappa}_{\alpha} - \left( \frac{p_{\alpha,j}^t - p^{t}_{\alpha,i}}{\Delta x} + \varrho_{\alpha} \mathbf{g}\right) - + V_i \frac{A_{\gamma_{ij}}}{U_i} \mathbf{K} - \sum_{\alpha} \varrho_{\alpha} \lambda_{\alpha} \sum_{\kappa} -\frac{\frac{\partial v_{t,j}}{\partial C^{\kappa}_j}-\frac{\partial v_{t,i}}{\partial C^{\kappa}_i}}{\Delta x} X^{\kappa}_{\alpha} - \left( \frac{p_{\alpha,j}^t - p^{t}_{\alpha,i}}{\Delta x} + \varrho_{\alpha} \mathbf{g}\right) \f] + * for first == false, the flux thorugh \f$ \gamma \f$ is calculated via a volume balance formulation + * \f[ - A_{\gamma} \mathbf{n}^T_{\gamma} \mathbf{K} \sum_{\alpha} \varrho_{\alpha} \lambda_{\alpha} + \mathbf{d}_{ij} \left( \frac{p_{\alpha,j}^t - p^{t}_{\alpha,i}}{\Delta x} + \varrho_{\alpha} \mathbf{g}^T \mathbf{d}_{ij} \right) + \sum_{\kappa} X^{\kappa}_{\alpha} \frac{\partial v_{t}}{\partial C^{\kappa}} + + V_i \frac{A_{\gamma}}{U_i} \mathbf{d}^T \mathbf{K} \sum_{\alpha} \varrho_{\alpha} \lambda_{\alpha} + \mathbf{d}_{ij} \left( \frac{p_{\alpha,j}^t - p^{t}_{\alpha,i}}{\Delta x} + \varrho_{\alpha} \mathbf{g}^T \mathbf{d}_{ij} \right) + \sum_{\kappa} X^{\kappa}_{\alpha} + \frac{\frac{\partial v_{t,j}}{\partial C^{\kappa}_j}-\frac{\partial v_{t,i}}{\partial C^{\kappa}_i}}{\Delta x} \f] * This includes a boundary integral and a volume integral, because * \f$ \frac{\partial v_{t,i}}{\partial C^{\kappa}_i} \f$ is not constant. - * Here, \f$ \mathbf{u} \f$ is the normalized vector connecting the cell centers, and \f$ \mathbf{n}_{\gamma_{ij}} \f$ - * represents the normal of the face \f$ \gamma{ij} \f$. + * Here, \f$ \mathbf{d}_{ij} \f$ is the normalized vector connecting the cell centers, and \f$ \mathbf{n}_{\gamma} \f$ + * represents the normal of the face \f$ \gamma \f$. * \param entries The Matrix and RHS entries * \param intersection Intersection between cell I and J * \param cellDataI Data of cell I @@ -612,14 +613,16 @@ void FVPressure2P2C<TypeTag>::getFlux(Dune::FieldVector<Scalar, 2>& entries, //! Get flux on Boundary /** for first == true, the flux is calculated in traditional fractional-flow forn as in FVPressure2P. - * for first == false, the flux thorugh \f$ \gamma{ij} \f$ is calculated via a volume balance formulation - * \f[ - A_{\gamma_{ij}} \cdot \mathbf{K} \cdot \mathbf{u} \cdot (\mathbf{n}_{\gamma_{ij}} \cdot \mathbf{u}) - \sum_{\alpha} \varrho_{\alpha} \lambda_{\alpha} \sum_{\kappa} \frac{\partial v_{t}}{\partial C^{\kappa}} X^{\kappa}_{\alpha} - \left( \frac{p_{\alpha,j}^t - p^{t}_{\alpha,i}}{\Delta x} + \varrho_{\alpha} \mathbf{g}\right) \;, \f] + * for first == false, the flux thorugh \f$ \gamma \f$ is calculated via a volume balance formulation + * \f[ - A_{\gamma} \mathbf{n}^T_{\gamma} \mathbf{K} \sum_{\alpha} \varrho_{\alpha} \lambda_{\alpha} \mathbf{d}_{ij} + \left( \frac{p_{\alpha,j}^t - p^{t}_{\alpha,i}}{\Delta x} + \varrho_{\alpha} \mathbf{g}^T \mathbf{d}_{ij} \right) + \sum_{\kappa} \frac{\partial v_{t}}{\partial C^{\kappa}} X^{\kappa}_{\alpha} \;, \f] * where we skip the volume integral assuming \f$ \frac{\partial v_{t,i}}{\partial C^{\kappa}_i} \f$ * to be constant at the boundary. - * Here, \f$ \mathbf{u} \f$ is the normalized vector connecting the cell centers, and \f$ \mathbf{n}_{\gamma_{ij}} \f$ - * represents the normal of the face \f$ \gamma{ij} \f$. + * Here, \f$ \mathbf{d}_{ij} \f$ is the normalized vector connecting the cell centers, and \f$ \mathbf{n}_{\gamma} \f$ + * represents the normal of the face \f$ \gamma \f$. + * + * If a Neumann BC is set, the given (mass-)flux is directly multiplied by the volume derivative and inserted. * \param entries The Matrix and RHS entries * \param intersection Intersection between cell I and J * \param cellDataI Data of cell I diff --git a/dumux/decoupled/2p2c/fvpressure2p2cadaptive.hh b/dumux/decoupled/2p2c/fvpressure2p2cadaptive.hh index 98057772d8..db02b9fda3 100644 --- a/dumux/decoupled/2p2c/fvpressure2p2cadaptive.hh +++ b/dumux/decoupled/2p2c/fvpressure2p2cadaptive.hh @@ -31,14 +31,14 @@ #include <dune/istl/preconditioners.hh> // dumux environment +// include 2p mpfa pressure model +#include <dumux/decoupled/common/fv/mpfa/fvmpfaproperties.hh> +#include <dumux/decoupled/2p/diffusion/fvmpfa/lmethod/fvmpfal2pfaboundpressure2padaptive.hh> + #include <dumux/decoupled/2p2c/fvpressure2p2c.hh> #include <dumux/common/math.hh> #include <dumux/io/vtkmultiwriter.hh> -// include pressure model from Markus -#include <dumux/decoupled/common/fv/mpfa/fvmpfaproperties.hh> -#include <dumux/decoupled/2p/diffusion/fvmpfa/lmethod/fvmpfal2pfaboundpressure2padaptive.hh> - /** * @file * @brief Finite Volume Diffusion Model @@ -47,6 +47,13 @@ namespace Dumux { +namespace Properties +{ +// forward declaration of properties +NEW_PROP_TAG(GridAdaptEnableMultiPointFluxApproximation); +NEW_PROP_TAG(GridAdaptEnableSecondHalfEdge); +} + //! The finite volume model for the solution of the compositional pressure equation /*! \ingroup multiphase * Provides a Finite Volume implementation for the pressure equation of a compressible @@ -64,11 +71,12 @@ namespace Dumux * \f$ \rho_{\alpha} \f$ the phase density and \f$ \bf{g} \f$ the gravity constant and \f$ C^{\kappa} \f$ the total Component concentration. * See paper SPE 99619 or "Analysis of a Compositional Model for Fluid * Flow in Porous Media" by Chen, Qin and Ewing for derivation. - * - * The pressure base class FVPressure assembles the matrix and right-hand-side vector and solves for the pressure vector, - * whereas this class provides the actual entries for the matrix and RHS vector. * The partial derivatives of the actual fluid volume \f$ v_{total} \f$ are gained by using a secant method. * + * This adaptive implementation uses its own initialization and assembling methods for matrix and right-hand-side vector, + * but solution is done by the base class FVPressure. Fluxes near hanging nodes can be + * calculated using an \a mpfa method. + * * \tparam TypeTag The Type Tag */ template<class TypeTag> class FVPressure2P2CAdaptive @@ -151,6 +159,7 @@ public: void initializeMatrix(); //function which assembles the system of equations to be solved void assemble(bool first); + void getMpfaFlux(const IntersectionIterator&, const CellData&); // mpfa transmissibilities @@ -213,12 +222,12 @@ protected: GlobalPosition&, TransmissivityMatrix&); - // Matrix for vector rotation of mpfa - // introduce matrix R for vector rotation and R is initialized as zero matrix + //! Matrix for vector rotation used in mpfa DimMatrix R_; - bool enableVolumeIntegral; - bool enableMPFA; - bool enableSecondHalfEdge; + bool enableVolumeIntegral; //!> Enables the volume integral of the pressure equation + bool enableMPFA; //!> Enables mpfa method to calculate the fluxes near hanging nodes + bool enableSecondHalfEdge; //!> If possible, 2 interaction volumes are used for the mpfa method near hanging nodes + //! The 2p Mpfa pressure module, that is only used for the calulation of transmissibility of the second interaction volumes FVMPFAL2PFABoundPressure2PAdaptive<TypeTag>* pressureModelAdaptive2p_; }; @@ -468,7 +477,15 @@ void FVPressure2P2CAdaptive<TypeTag>::assemble(bool first) } //! Compute flux over an irregular interface using a \a mpfa method -/** TODO: Insert correct formulas! +/** A mpfa l-method is applied to calculate fluxes near hanging nodes, using: + * \f[ + - \sum_{\alpha} \varrho_{\alpha} \lambda_{\alpha} + \left( \sum_k \tau_{2k} p^t_{\alpha,k} + \varrho_{\alpha} \sum_k \tau_{2k} \mathbf{g}^T \mathbf{x}_{k} \right) + \sum_{\kappa} X^{\kappa}_{\alpha} \frac{\partial v_{t}}{\partial C^{\kappa}} + + \frac{ V_i}{U_i} \sum_{\alpha} \varrho_{\alpha} \lambda_{\alpha} + \left( \sum_k \tau_{2k} p^t_{\alpha,k} + \varrho_{\alpha} \sum_k \tau_{2k} \mathbf{g}^T \mathbf{x}_{k} \right) + \sum_{\kappa} X^{\kappa}_{\alpha} \frac{\frac{\partial v_{t,j}}{\partial C^{\kappa}_j}-\frac{\partial v_{t,i}}{\partial C^{\kappa}_i}}{\Delta x} + \f] * * We provide two options: Calculating the flux expressed by twice the flux * through the one unique interaction region on the hanging node if one @@ -1077,10 +1094,13 @@ int FVPressure2P2CAdaptive<TypeTag>::computeTransmissibilities(const Intersectio } } -// mpfa transmissibilities from mpfal2pfa - -/*! Indices used in a interaction volume of the MPFA-o method +//! An adapter to use the traditional 2p implementation for the second interaction region +/*! + * The second interaction region consists of 4 cells, hence the traditional non-adaptive + * implementation to calculate the transmissibility coefficients is used. This, however, + * uses its own naming conventions and local Indices used in the work of I. Aavatsmark. * + * Indices used in a interaction volume of the MPFA-l method * \verbatim | | | | 4-----------3-----------3 | diff --git a/dumux/decoupled/2p2c/fvpressure2p2cmultiphysics.hh b/dumux/decoupled/2p2c/fvpressure2p2cmultiphysics.hh index 02e90df986..76ed05bd32 100644 --- a/dumux/decoupled/2p2c/fvpressure2p2cmultiphysics.hh +++ b/dumux/decoupled/2p2c/fvpressure2p2cmultiphysics.hh @@ -189,7 +189,7 @@ protected: Dune::BlockVector<Dune::FieldVector<int,1> > nextSubdomain; //! vector holding next subdomain const GlobalPosition& gravity_; //!< vector including the gravity constant static constexpr int pressureType = GET_PROP_VALUE(TypeTag, PressureFormulation); //!< gives kind of pressure used (\f$ 0 = p_w \f$, \f$ 1 = p_n \f$, \f$ 2 = p_{global} \f$) - Dune::Timer timer_; + Dune::Timer timer_; //!> A timer for the time spent on the multiphysics framework. //! Indices of matrix and rhs entries /** @@ -298,6 +298,16 @@ void FVPressure2P2CMultiPhysics<TypeTag>::assemble(bool first) return; } +//! Assembles the source term +/** The source is translated into a volumentric source term: + * \f[ V_i \sum_{\kappa} \frac{1}{\varrho} q^{\kappa}_i \; , \f] + * because under singlephase conditions + * \f[ \frac{\partial v_{t}}{\partial C^{\kappa}} \approx \frac{1}{\varrho} \f]. + * \param sourceEntry The Matrix and RHS entries + * \param elementI The element I + * \param cellDataI Data of cell I + * \param first Flag if pressure field is unknown + */ template<class TypeTag> void FVPressure2P2CMultiPhysics<TypeTag>::get1pSource(Dune::FieldVector<Scalar, 2>& sourceEntry, const Element& elementI, const CellData& cellDataI) { @@ -317,7 +327,17 @@ void FVPressure2P2CMultiPhysics<TypeTag>::get1pSource(Dune::FieldVector<Scalar, return; } - +//! Assembles the storage term for a 1p cell in a multiphysics framework +/** The storage term comprises the (single-phase) compressibility (due to a change in + * pressure from last timestep): + * \f[ V_i c_{i} \frac{p^t_i - p^{t-\Delta t}_i}{\Delta t} \f] + * and the damped error introduced by the incorrect transport of the last timestep: + * \f[ V_i \alpha_r \frac{v_{t} - \phi}{\Delta t} \f]. + * The latter is damped according to Fritz 2011. + * \param storageEntry The Matrix and RHS entries + * \param elementI The element I + * \param cellDataI Data of cell I + */ template<class TypeTag> void FVPressure2P2CMultiPhysics<TypeTag>::get1pStorage(Dune::FieldVector<Scalar, 2>& storageEntry, const Element& elementI, @@ -422,8 +442,19 @@ void FVPressure2P2CMultiPhysics<TypeTag>::get1pStorage(Dune::FieldVector<Scalar, return; } - - +/*! \brief The compositional single-phase flux in the multiphysics framework + * + * If only single-phase conditions are encountered, the flux expression simplifies to (written for the + * case where the wetting phase is only present): + \f[ + A_{\gamma} \mathbf{n}_{\gamma}^T \mathbf{K} + \varrho_w \lambda_w \mathbf{d}_{ij} \left( \frac{p_{w,j}^t - p^{t}_{w,i}}{\Delta x} + \varrho_{w} \mathbf{g}^T \mathbf{d}_{ij} \right) . + \f] + * + * \param entries The Matrix and RHS entries + * \param intersection Intersection between cell I and J + * \param cellDataI Data of cell I + */ template<class TypeTag> void FVPressure2P2CMultiPhysics<TypeTag>::get1pFlux(Dune::FieldVector<Scalar, 2>& entries, const Intersection& intersection, const CellData& cellDataI) @@ -514,6 +545,22 @@ void FVPressure2P2CMultiPhysics<TypeTag>::get1pFlux(Dune::FieldVector<Scalar, 2> return; } +/*! \brief The compositional single-phase flux in the multiphysics framework + * + * If only single-phase conditions are encountered, the flux expression simplifies to (written for the + * case where the wetting phase is only present): + \f[ + A_{\gamma} \mathbf{n}_{\gamma}^T \mathbf{K} + \varrho_w \lambda_w \mathbf{d}_{i-Boundary} \left( \frac{p_{w,Boundary}^t - p^{t}_{w,i}}{\Delta x} + + \varrho_{w} \mathbf{g}^T \mathbf{d}_{i-Boundary} \right) . + \f] + * + * If a Neumann BC is set, the given (mass-)flux is directly multiplied by the volume derivative and inserted. + * + * \param entries The Matrix and RHS entries + * \param intersection Intersection between cell I and J + * \param cellDataI Data of cell I + */ template<class TypeTag> void FVPressure2P2CMultiPhysics<TypeTag>::get1pFluxOnBoundary(Dune::FieldVector<Scalar, 2>& entries, const Intersection& intersection, const CellData& cellDataI) @@ -661,9 +708,13 @@ void FVPressure2P2CMultiPhysics<TypeTag>::get1pFluxOnBoundary(Dune::FieldVector< //! constitutive functions are updated once if new concentrations are calculated and stored in the variables container /*! - * In contrast to the standard sequential 2p2c model, this method also holds routines - * to adapt the subdomain. The subdomain indicates weather we are in 1p domain (value = 1) + * In contrast to the standard sequential 2p2c model ( FVPressure2P2CMultiPhysics<TypeTag>::updateMaterialLaws() ), + * this method also holds routines to adapt the subdomain. The subdomain indicates weather we are in 1p domain (value = 1) * or in the two phase subdomain (value = 2). + * Note that the type of flash, i.e. the type of FluidState (FS), present in each cell does not have to + * coincide with the subdomain. If a cell will be simple and was complex, a complex FS is available, so next time step + * will use this complex FS, but updateMaterialLaw afterwards will finally transform that to simple FS. + * \param postTimeStep Flag indicating method is called from Problem::postTimeStep() */ template<class TypeTag> void FVPressure2P2CMultiPhysics<TypeTag>::updateMaterialLaws(bool postTimeStep) @@ -766,6 +817,14 @@ void FVPressure2P2CMultiPhysics<TypeTag>::updateMaterialLaws(bool postTimeStep) return; } + +//! updates secondary variables of one single phase cell +/*! For each element, the secondary variables are updated according to the + * primary variables. Only a simple flash calulation has to be carried out, + * as phase distribution is already known: single-phase. + * \param elementI The element + * \param postTimeStep Flag indicating if we have just completed a time step + */ template<class TypeTag> void FVPressure2P2CMultiPhysics<TypeTag>::update1pMaterialLawsInElement(const Element& elementI, CellData& cellData, bool postTimeStep) { diff --git a/dumux/decoupled/2p2c/fvtransport2p2cadaptive.hh b/dumux/decoupled/2p2c/fvtransport2p2cadaptive.hh index fa9dcd3dc6..cc5970f70a 100644 --- a/dumux/decoupled/2p2c/fvtransport2p2cadaptive.hh +++ b/dumux/decoupled/2p2c/fvtransport2p2cadaptive.hh @@ -23,7 +23,6 @@ #define DUMUX_FVTRANSPORT2P2C_ADAPTIVE_HH #include <dune/grid/common/gridenums.hh> -#include <dumux/decoupled/2p2c/2p2cadaptiveproperties.hh> #include <dumux/decoupled/2p2c/fvtransport2p2c.hh> #include <dumux/common/math.hh> @@ -35,6 +34,11 @@ namespace Dumux { +namespace Properties +{ +// forward declaration of properties +NEW_PROP_TAG(GridAdaptEnableMultiPointFluxApproximation); +} //! Miscible Transport step in a Finite Volume discretization /*! \ingroup multiphase * The finite volume model for the solution of the transport equation for compositional -- GitLab