// -*- 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::GroundWaterProblem */ #ifndef EWOMS_GROUND_WATER_PROBLEM_HH #define EWOMS_GROUND_WATER_PROBLEM_HH #include #include #include #include #include #include #include #include #include #include #include #include namespace Opm { template class GroundWaterProblem; } namespace Opm::Properties { namespace TTag { struct GroundWaterBaseProblem {}; } template struct LensLowerLeftX { using type = UndefinedProperty; }; template struct LensLowerLeftY { using type = UndefinedProperty; }; template struct LensLowerLeftZ { using type = UndefinedProperty; }; template struct LensUpperRightX { using type = UndefinedProperty; }; template struct LensUpperRightY { using type = UndefinedProperty; }; template struct LensUpperRightZ { using type = UndefinedProperty; }; template struct Permeability { using type = UndefinedProperty; }; template struct PermeabilityLens { using type = UndefinedProperty; }; template struct Fluid { private: using Scalar = GetPropType; public: using type = Opm::LiquidPhase >; }; // Set the grid type template struct Grid { using type = Dune::YaspGrid<2>; }; // struct Grid { using type = Dune::SGrid<2, 2>; }; template struct Problem { using type = Opm::GroundWaterProblem; }; template struct LensLowerLeftX { using type = GetPropType; static constexpr type value = 0.25; }; template struct LensLowerLeftY { using type = GetPropType; static constexpr type value = 0.25; }; template struct LensLowerLeftZ { using type = GetPropType; static constexpr type value = 0.25; }; template struct LensUpperRightX { using type = GetPropType; static constexpr type value = 0.75; }; template struct LensUpperRightY { using type = GetPropType; static constexpr type value = 0.75; }; template struct LensUpperRightZ { using type = GetPropType; static constexpr type value = 0.75; }; template struct Permeability { using type = GetPropType; static constexpr type value = 1e-10; }; template struct PermeabilityLens { using type = GetPropType; static constexpr type value = 1e-12; }; // Use the conjugated gradient linear solver with the default preconditioner (i.e., // ILU-0) from dune-istl template struct LinearSolverSplice { using type = TTag::ParallelIstlLinearSolver; }; template struct LinearSolverWrapper { using type = Opm::Linear::SolverWrapperConjugatedGradients; }; } // namespace Opm::Properties namespace Opm::Parameters { // Enable gravity template struct EnableGravity { static constexpr bool value = true; }; // The default for the end time of the simulation template struct EndTime { using type = GetPropType; static constexpr type value = 1; }; // The default DGF file to load template struct GridFile { static constexpr auto value = "./data/groundwater_2d.dgf"; }; // The default for the initial time step size of the simulation template struct InitialTimeStepSize { using type = GetPropType; static constexpr type value = 1; }; } // namespace Opm::Parameters namespace Opm { /*! * \ingroup TestProblems * * \brief Test for the immisicible VCVF discretization with only a single phase * * This problem is inspired by groundwater flow. Don't expect it to be * realistic, though: For two dimensions, the domain size is 1m times * 1m. On the left and right of the domain, no-flow boundaries are * used, while at the top and bottom free flow boundaries with a * pressure of 2 bar and 1 bar are used. The center of the domain is * occupied by a rectangular lens of lower permeability. */ template class GroundWaterProblem : public GetPropType { using ParentType = GetPropType; using GridView = GetPropType; using Scalar = GetPropType; using FluidSystem = GetPropType; // copy some indices for convenience using Indices = GetPropType; enum { numPhases = FluidSystem::numPhases, // Grid and world dimension dim = GridView::dimension, dimWorld = GridView::dimensionworld, // indices of the primary variables pressure0Idx = Indices::pressure0Idx }; using Simulator = GetPropType; using EqVector = GetPropType; using RateVector = GetPropType; using BoundaryRateVector = GetPropType; using PrimaryVariables = GetPropType; using Model = GetPropType; using CoordScalar = typename GridView::ctype; using GlobalPosition = Dune::FieldVector; using DimMatrix = Dune::FieldMatrix; public: /*! * \copydoc Doxygen::defaultProblemConstructor */ GroundWaterProblem(Simulator& simulator) : ParentType(simulator) { } /*! * \copydoc FvBaseProblem::finishInit */ void finishInit() { ParentType::finishInit(); eps_ = 1.0e-3; lensLowerLeft_[0] = Parameters::get(); if (dim > 1) lensLowerLeft_[1] = Parameters::get(); if (dim > 2) lensLowerLeft_[2] = Parameters::get(); lensUpperRight_[0] = Parameters::get(); if (dim > 1) lensUpperRight_[1] = Parameters::get(); if (dim > 2) lensUpperRight_[2] = Parameters::get(); intrinsicPerm_ = this->toDimMatrix_(Parameters::get()); intrinsicPermLens_ = this->toDimMatrix_(Parameters::get()); } /*! * \copydoc FvBaseMultiPhaseProblem::registerParameters */ static void registerParameters() { ParentType::registerParameters(); Parameters::registerParam ("The x-coordinate of the lens' lower-left corner [m]."); Parameters::registerParam ("The x-coordinate of the lens' upper-right corner [m]."); if (dimWorld > 1) { Parameters::registerParam ("The y-coordinate of the lens' lower-left corner [m]."); Parameters::registerParam ("The y-coordinate of the lens' upper-right corner [m]."); } if (dimWorld > 2) { Parameters::registerParam ("The z-coordinate of the lens' lower-left corner [m]."); Parameters::registerParam ("The z-coordinate of the lens' upper-right corner [m]."); } Parameters::registerParam ("The intrinsic permeability [m^2] of the ambient material."); Parameters::registerParam ("The intrinsic permeability [m^2] of the lens."); } /*! * \name Problem parameters */ // \{ /*! * \copydoc FvBaseProblem::name */ std::string name() const { std::ostringstream oss; oss << "groundwater_" << Model::name(); return oss.str(); } /*! * \copydoc FvBaseProblem::endTimeStep */ void endTimeStep() { #ifndef NDEBUG this->model().checkConservativeness(); // Calculate storage terms EqVector storage; this->model().globalStorage(storage); // Write mass balance information for rank 0 if (this->gridView().comm().rank() == 0) { std::cout << "Storage: " << storage << std::endl << std::flush; } #endif // NDEBUG } /*! * \copydoc FvBaseMultiPhaseProblem::temperature */ template Scalar temperature(const Context& /*context*/, unsigned /*spaceIdx*/, unsigned /*timeIdx*/) const { return 273.15 + 10; } // 10C /*! * \copydoc FvBaseMultiPhaseProblem::porosity */ template Scalar porosity(const Context& /*context*/, unsigned /*spaceIdx*/, unsigned /*timeIdx*/) const { return 0.4; } /*! * \copydoc FvBaseMultiPhaseProblem::intrinsicPermeability */ template const DimMatrix& intrinsicPermeability(const Context& context, unsigned spaceIdx, unsigned timeIdx) const { if (isInLens_(context.pos(spaceIdx, timeIdx))) return intrinsicPermLens_; else return intrinsicPerm_; } //! \} /*! * \name Boundary conditions */ //! \{ /*! * \copydoc FvBaseProblem::boundary */ template void boundary(BoundaryRateVector& values, const Context& context, unsigned spaceIdx, unsigned timeIdx) const { const GlobalPosition& globalPos = context.pos(spaceIdx, timeIdx); if (onLowerBoundary_(globalPos) || onUpperBoundary_(globalPos)) { Scalar pressure; Scalar T = temperature(context, spaceIdx, timeIdx); if (onLowerBoundary_(globalPos)) pressure = 2e5; else // on upper boundary pressure = 1e5; Opm::ImmiscibleFluidState fs; fs.setSaturation(/*phaseIdx=*/0, 1.0); fs.setPressure(/*phaseIdx=*/0, pressure); fs.setTemperature(T); typename FluidSystem::template ParameterCache paramCache; paramCache.updateAll(fs); for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) { fs.setDensity(phaseIdx, FluidSystem::density(fs, paramCache, phaseIdx)); fs.setViscosity(phaseIdx, FluidSystem::viscosity(fs, paramCache, phaseIdx)); } // impose an freeflow boundary condition values.setFreeFlow(context, spaceIdx, timeIdx, fs); } else { // no flow boundary values.setNoFlow(); } } //! \} /*! * \name Volumetric terms */ //! \{ /*! * \copydoc FvBaseProblem::initial */ template void initial(PrimaryVariables& values, const Context& /*context*/, unsigned /*spaceIdx*/, unsigned /*timeIdx*/) const { // const GlobalPosition& globalPos = context.pos(spaceIdx, timeIdx); values[pressure0Idx] = 1.0e+5; // + 9.81*1.23*(20-globalPos[dim-1]); } /*! * \copydoc FvBaseProblem::source */ template void source(RateVector& rate, const Context& /*context*/, unsigned /*spaceIdx*/, unsigned /*timeIdx*/) const { rate = Scalar(0.0); } //! \} private: bool onLowerBoundary_(const GlobalPosition& pos) const { return pos[dim - 1] < eps_; } bool onUpperBoundary_(const GlobalPosition& pos) const { return pos[dim - 1] > this->boundingBoxMax()[dim - 1] - eps_; } bool isInLens_(const GlobalPosition& pos) const { return lensLowerLeft_[0] <= pos[0] && pos[0] <= lensUpperRight_[0] && lensLowerLeft_[1] <= pos[1] && pos[1] <= lensUpperRight_[1]; } GlobalPosition lensLowerLeft_; GlobalPosition lensUpperRight_; DimMatrix intrinsicPerm_; DimMatrix intrinsicPermLens_; Scalar eps_; }; } // namespace Opm #endif