Skip to content
Snippets Groups Projects
problem.hh 8.4 KiB
Newer Older
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
Bernd Flemisch's avatar
Bernd Flemisch committed
/*****************************************************************************
 *   See the file COPYING for full copying permissions.                      *
Bernd Flemisch's avatar
Bernd Flemisch committed
 *                                                                           *
 *   This program is free software: you can redistribute it and/or modify    *
Bernd Flemisch's avatar
Bernd Flemisch committed
 *   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.                                     *
Bernd Flemisch's avatar
Bernd Flemisch committed
 *                                                                           *
 *   This program 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 this program.  If not, see <http://www.gnu.org/licenses/>.   *
Bernd Flemisch's avatar
Bernd Flemisch committed
 *****************************************************************************/
 * \ingroup RichardsTests
 * \brief A water infiltration problem with a low-permeability lens
 *        embedded into a high-permeability domain which uses the
#ifndef DUMUX_RICHARDS_LENSPROBLEM_HH
#define DUMUX_RICHARDS_LENSPROBLEM_HH
#include <dumux/common/properties.hh>
#include <dumux/common/parameters.hh>
#include <dumux/common/boundarytypes.hh>
#include <dumux/common/numeqvector.hh>
Katharina Heck's avatar
Katharina Heck committed
#include <dumux/porousmediumflow/problem.hh>

Bernd Flemisch's avatar
Bernd Flemisch committed
/*!
 *
 * \brief A water infiltration problem with a low-permeability lens
 *        embedded into a high-permeability domain which uses the
Natalie Schröder's avatar
Natalie Schröder committed
 *        Richards model.
Bernd Flemisch's avatar
Bernd Flemisch committed
 *
 * The domain is box shaped. Left and right boundaries are Dirichlet
 * boundaries with fixed water pressure (fixed Saturation \f$S_w = 0\f$),
 * bottom boundary is closed (Neumann 0 boundary), the top boundary
 * (Neumann 0 boundary) is also closed except for infiltration
 * section, where water is infiltrating into an initially unsaturated
 * porous medium. This problem is very similar to the LensProblem
 * which uses the TwoPBoxModel, with the main difference being that
 * the domain is initially fully saturated by gas instead of water and
 * water instead of a %DNAPL infiltrates from the top.
 * This problem uses the \ref RichardsModel
Bernd Flemisch's avatar
Bernd Flemisch committed
 *
 * To run the simulation execute the following line in shell:
 * <tt>./test_boxrichards -parameterFile test_boxrichards.input -TimeManager.TEnd 10000000</tt>
 * <tt>./test_ccrichards -parameterFile test_ccrichards.input -TimeManager.TEnd 10000000</tt>
Bernd Flemisch's avatar
Bernd Flemisch committed
 *
 * where the initial time step is 100 seconds, and the end of the
 * simulation time is 10,000,000 seconds (115.7 days)
Bernd Flemisch's avatar
Bernd Flemisch committed
 */
Katharina Heck's avatar
Katharina Heck committed
class RichardsLensProblem : public PorousMediumFlowProblem<TypeTag>
Bernd Flemisch's avatar
Bernd Flemisch committed
{
Katharina Heck's avatar
Katharina Heck committed
    using ParentType = PorousMediumFlowProblem<TypeTag>;
    using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView;
    using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
    using BoundaryTypes = Dumux::BoundaryTypes<GetPropType<TypeTag, Properties::ModelTraits>::numEq()>;
    using NumEqVector = Dumux::NumEqVector<PrimaryVariables>;
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices;
    using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
Bernd Flemisch's avatar
Bernd Flemisch committed
    enum {
        // copy some indices for convenience
        pressureIdx = Indices::pressureIdx,
        conti0EqIdx = Indices::conti0EqIdx,
Timo Koch's avatar
Timo Koch committed
        bothPhases = Indices::bothPhases,
        dimWorld = GridView::dimensionworld
Bernd Flemisch's avatar
Bernd Flemisch committed
    };
    using Element = typename GridView::template Codim<0>::Entity;
    using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
Bernd Flemisch's avatar
Bernd Flemisch committed

public:
    RichardsLensProblem(std::shared_ptr<const GridGeometry> gridGeometry)
    : ParentType(gridGeometry)
Bernd Flemisch's avatar
Bernd Flemisch committed
    {
        name_ = getParam<std::string>("Problem.Name");
Bernd Flemisch's avatar
Bernd Flemisch committed
    }
Bernd Flemisch's avatar
Bernd Flemisch committed
    /*!
     * \name Problem parameters
     */
    // \{

    /*!
Bernd Flemisch's avatar
Bernd Flemisch committed
     *
     * This is used as a prefix for files generated by the simulation.
     */
    const std::string& name() const
Bernd Flemisch's avatar
Bernd Flemisch committed
    /*!
     * \brief Returns the temperature [K] within a finite volume
Bernd Flemisch's avatar
Bernd Flemisch committed
     *
     * This problem assumes a temperature of 10 degrees Celsius.
     */
    { return 273.15 + 10; }; // -> 10°C
     * \brief Returns the reference pressure [Pa] of the nonwetting
     *        fluid phase within a finite volume
     *
     * This problem assumes a constant reference pressure of 1 bar.
Bernd Flemisch's avatar
Bernd Flemisch committed
     */
    Scalar nonwettingReferencePressure() const
Bernd Flemisch's avatar
Bernd Flemisch committed
    // \}
Bernd Flemisch's avatar
Bernd Flemisch committed
    /*!
     * \name Boundary conditions
     */
    // \{

    /*!
     * \brief Specifies which kind of boundary condition should be
     *        used for which equation on a given boundary segment.
     * \param globalPos The position for which the boundary type is set
Bernd Flemisch's avatar
Bernd Flemisch committed
     */
    BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Bernd Flemisch's avatar
Bernd Flemisch committed
    {
        BoundaryTypes bcTypes;
        if (onLeftBoundary_(globalPos) || onRightBoundary_(globalPos))
            bcTypes.setAllDirichlet();
            bcTypes.setAllNeumann();
        return bcTypes;
Bernd Flemisch's avatar
Bernd Flemisch committed
    }
Bernd Flemisch's avatar
Bernd Flemisch committed
    /*!
     * \brief Evaluates the boundary conditions for a Dirichlet boundary segment.
Bernd Flemisch's avatar
Bernd Flemisch committed
     *
     * \param globalPos The position for which the Dirichlet value is set
Bernd Flemisch's avatar
Bernd Flemisch committed
     * For this method, the \a values parameter stores primary variables.
     */
    PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Bernd Flemisch's avatar
Bernd Flemisch committed
    {
Bernd Flemisch's avatar
Bernd Flemisch committed
    }
Bernd Flemisch's avatar
Bernd Flemisch committed
    /*!
     * \brief Evaluates the boundary conditions for a Neumann boundary segment.
Bernd Flemisch's avatar
Bernd Flemisch committed
     *
     * For this method, the \a values parameter stores the mass flux
     * in normal direction of each phase. Negative values mean influx.
     * \param globalPos The position for which the Neumann value is set
Bernd Flemisch's avatar
Bernd Flemisch committed
     */
    NumEqVector neumannAtPos(const GlobalPosition &globalPos) const
Bernd Flemisch's avatar
Bernd Flemisch committed
    {
        if (onInlet_(globalPos))
            values[conti0EqIdx] = -0.04; // kg/(m*s)
        return values;
Bernd Flemisch's avatar
Bernd Flemisch committed
    }

    /*!
     * \name Volume terms
     */
    // \{

    /*!
     * \brief Evaluates the initial values for a control volume.
Bernd Flemisch's avatar
Bernd Flemisch committed
     *
     * For this method, the \a values parameter stores primary
     * variables.
     * \param globalPos The position for which the boundary type is set
Bernd Flemisch's avatar
Bernd Flemisch committed
     */
    PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
    { return initial_(globalPos); };
    PrimaryVariables initial_(const GlobalPosition &globalPos) const
        PrimaryVariables values(0.0);
        const Scalar sw = 0.0;
        const Scalar pc = this->spatialParams().fluidMatrixInteractionAtPos(globalPos).pc(sw);
        values[pressureIdx] = nonwettingReferencePressure() - pc;
Timo Koch's avatar
Timo Koch committed
        values.setState(bothPhases);
    bool onLeftBoundary_(const GlobalPosition &globalPos) const
    {
        return globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_;
    }

    bool onRightBoundary_(const GlobalPosition &globalPos) const
    {
        return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_;
    }

    bool onLowerBoundary_(const GlobalPosition &globalPos) const
    {
        return globalPos[1] < this->gridGeometry().bBoxMin()[1] + eps_;
    }

    bool onUpperBoundary_(const GlobalPosition &globalPos) const
    {
        return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_;
    }

    bool onInlet_(const GlobalPosition &globalPos) const
    {
        Scalar width = this->gridGeometry().bBoxMax()[0] - this->gridGeometry().bBoxMin()[0];
        Scalar lambda = (this->gridGeometry().bBoxMax()[0] - globalPos[0])/width;
        return onUpperBoundary_(globalPos) && 0.5 < lambda + eps_ && lambda < 2.0/3.0 + eps_;
    static constexpr Scalar eps_ = 1.5e-7;
    GlobalPosition lensLowerLeft_;
    GlobalPosition lensUpperRight_;
Bernd Flemisch's avatar
Bernd Flemisch committed
};
Bernd Flemisch's avatar
Bernd Flemisch committed

#endif