Newer
Older
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* See the file COPYING for full copying permissions. *
* This program 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. *
* 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/>. *
*****************************************************************************/

Andreas Lauser
committed
/*!
* \file
*
* \brief A water infiltration problem with a low-permeability lens
* embedded into a high-permeability domain which uses the

Andreas Lauser
committed
* Richards box model.
*/
#ifndef DUMUX_RICHARDS_LENSPROBLEM_HH
#define DUMUX_RICHARDS_LENSPROBLEM_HH
#include <dune/grid/io/file/dgfparser.hh>
#include <dumux/boxmodels/richards/richardsmodel.hh>

Bernd Flemisch
committed
#include <dumux/material/components/simpleh2o.hh>
#include <dumux/material/fluidsystems/liquidphase.hh>
#include "richardslensspatialparams.hh"
namespace Dumux
{
template <class TypeTag>
//////////
// Specify the properties for the lens problem
//////////
namespace Properties
{
NEW_TYPE_TAG(RichardsLensProblem, INHERITS_FROM(BoxRichards, RichardsLensSpatialParams));
// Use 2d YaspGrid
SET_TYPE_PROP(RichardsLensProblem, Grid, Dune::YaspGrid<2>);
// Set the physical problem to be solved

Andreas Lauser
committed
{ typedef Dumux::RichardsLensProblem<TypeTag> type; };
SET_PROP(RichardsLensProblem, WettingPhase)
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;

Bernd Flemisch
committed
typedef Dumux::LiquidPhase<Scalar, Dumux::SimpleH2O<Scalar> > type;
};
// Enable gravity
SET_BOOL_PROP(RichardsLensProblem, EnableGravity, true);
// Enable partial reassembly of the Jacobian matrix
SET_BOOL_PROP(RichardsLensProblem, EnablePartialReassemble, true);
// Enable re-use of the Jacobian matrix for the first iteration of a time step
SET_BOOL_PROP(RichardsLensProblem, EnableJacobianRecycling, true);
// Use forward differences to approximate the Jacobian matrix
SET_INT_PROP(RichardsLensProblem, NumericDifferenceMethod, +1);
// Set the maximum number of newton iterations of a time step
SET_INT_PROP(RichardsLensProblem, NewtonMaxSteps, 28);
// Set the "desireable" number of newton iterations of a time step
SET_INT_PROP(RichardsLensProblem, NewtonTargetSteps, 18);
// Do not write the intermediate results of the newton method
SET_BOOL_PROP(RichardsLensProblem, NewtonWriteConvergence, false);

Andreas Lauser
committed
* \ingroup RichardsModel
* \ingroup BoxTestProblems

Andreas Lauser
committed
*
* \brief A water infiltration problem with a low-permeability lens
* embedded into a high-permeability domain which uses the

Andreas Lauser
committed
* Richards box model.

Andreas Lauser
committed
* The domain is box shaped. Left and right boundaries are Dirichlet

Andreas Lauser
committed
* boundaries with fixed water pressure (fixed Saturation \f$S_w = 0\f$),

Andreas Lauser
committed
* 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 the the LensProblem
* which uses the TwoPBoxModel, with the main difference being that
* the domain is initally fully saturated by gas instead of water and

Andreas Lauser
committed
* water instead of a %DNAPL infiltrates from the top.
*
* To run the simulation execute the following line in shell:
* <tt>./test_richards -parameterFile test_richards.input -TimeManager.TEnd 10000000</tt>

Andreas Lauser
committed
* where the initial time step is 100 seconds, and the end of the
* simulation time is 10,000,000 seconds (115.7 days)

Andreas Lauser
committed
template <class TypeTag>
class RichardsLensProblem : public RichardsBoxProblem<TypeTag>
typedef RichardsBoxProblem<TypeTag> ParentType;
typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
typedef typename GET_PROP_TYPE(TypeTag, BoundaryTypes) BoundaryTypes;
typedef typename GET_PROP_TYPE(TypeTag, TimeManager) TimeManager;
typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry;
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
pwIdx = Indices::pwIdx,
contiEqIdx = Indices::contiEqIdx,
dimWorld = GridView::dimensionworld
typedef typename GridView::template Codim<0>::Entity Element;
typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition;

Andreas Lauser
committed
/*!
* \brief Constructor
*
* \param timeManager The Dumux TimeManager for simulation management.
* \param gridView The grid view on the spatial domain of the problem
*/
RichardsLensProblem(TimeManager &timeManager,
const GridView &gridView)
: ParentType(timeManager, gridView)
eps_ = 3e-6;
pnRef_ = 1e5;
lensLowerLeft_[0] = 1.0;
lensLowerLeft_[1] = 2.0;
lensUpperRight_[0] = 4.0;
lensUpperRight_[1] = 3.0;
this->spatialParams().setLensCoords(lensLowerLeft_, lensUpperRight_);
}
/*!
* \name Problem parameters
*/
// \{
/*!

Andreas Lauser
committed
* \brief The problem name
*
* This is used as a prefix for files generated by the simulation.
*/
const char *name() const

Andreas Lauser
committed
* \brief Returns the temperature [K] within a finite volume
*
* This problem assumes a temperature of 10 degrees Celsius.
*/

Andreas Lauser
committed
Scalar temperature() const

Andreas Lauser
committed
{ return 273.15 + 10; }; // -> 10°C

Andreas Lauser
committed
* \brief Returns the reference pressure [Pa] of the non-wetting
* fluid phase within a finite volume
*
* This problem assumes a constant reference pressure of 1 bar.

Andreas Lauser
committed
* \param element The DUNE Codim<0> entity which intersects with
* the finite volume in question
* \param fvElemGeom The finite volume geometry of the element
* \param scvIdx The sub control volume index inside the finite

Andreas Lauser
committed
* volume geometry

Andreas Lauser
committed
Scalar referencePressure(const Element &element,
const FVElementGeometry &fvElemGeom,
const int scvIdx) const

Andreas Lauser
committed
{ return pnRef_; };
void sourceAtPos(PrimaryVariables &values,
const GlobalPosition &globalPos) const
{
values = 0;
}
// \}
/*!
* \name Boundary conditions
*/
// \{
/*!
* \brief Specifies which kind of boundary condition should be
* used for which equation on a given boundary segment.
*
* \param values The boundary types for the conservation equations
* \param globalPos The position for which the boundary type is set
void boundaryTypesAtPos(BoundaryTypes &values,

Andreas Lauser
committed
const GlobalPosition &globalPos) const
onRightBoundary_(globalPos))
{
values.setAllNeumann();
}
/*!
* \brief Evaluate the boundary conditions for a dirichlet
* boundary segment.
*
* \param values The dirichlet values for the primary variables
* \param globalPos The position for which the Dirichlet value is set
* For this method, the \a values parameter stores primary variables.
*/
void dirichletAtPos(PrimaryVariables &values,

Andreas Lauser
committed
const GlobalPosition &globalPos) const
// use initial values as boundary conditions
}
/*!
* \brief Evaluate the boundary conditions for a neumann
* boundary segment.
*
* For this method, the \a values parameter stores the mass flux
* in normal direction of each phase. Negative values mean influx.

Andreas Lauser
committed
*
* \param values The neumann values for the conservation equations
* \param globalPos The position for which the Neumann value is set
void neumannAtPos(PrimaryVariables &values,

Andreas Lauser
committed
const GlobalPosition &globalPos) const
values = 0.0;
// inflow of water
values[contiEqIdx] = -0.04; // kg / (m * s)
}
/*!
* \name Volume terms
*/
// \{
/*!

Andreas Lauser
committed
* \brief Evaluate the initial values for a control volume.
*
* For this method, the \a values parameter stores primary
* variables.

Andreas Lauser
committed
*
* \param values Storage for all primary variables of the initial condition
* \param pos The position for which the boundary type is set
void initialAtPos(PrimaryVariables &values,
const GlobalPosition &pos) const

Andreas Lauser
committed
{ initial_(values, pos); };
// \}
private:
void initial_(PrimaryVariables &values, const GlobalPosition &pos) const
{
Scalar Sw = 0.0;
Scalar pc =
MaterialLaw::pC(this->spatialParams().materialLawParams(pos),

Andreas Lauser
committed
values[pwIdx] = pnRef_ - pc;
bool onLeftBoundary_(const GlobalPosition &globalPos) const
{
return globalPos[0] < this->bboxMin()[0] + eps_;
}
bool onRightBoundary_(const GlobalPosition &globalPos) const
{
return globalPos[0] > this->bboxMax()[0] - eps_;
}
bool onLowerBoundary_(const GlobalPosition &globalPos) const
{
return globalPos[1] < this->bboxMin()[1] + eps_;
}
bool onUpperBoundary_(const GlobalPosition &globalPos) const
{
return globalPos[1] > this->bboxMax()[1] - eps_;
}
bool onInlet_(const GlobalPosition &globalPos) const
{
Scalar width = this->bboxMax()[0] - this->bboxMin()[0];
Scalar lambda = (this->bboxMax()[0] - globalPos[0])/width;
return onUpperBoundary_(globalPos) && 0.5 < lambda && lambda < 2.0/3.0;

Andreas Lauser
committed
Scalar pnRef_;
GlobalPosition lensLowerLeft_;
GlobalPosition lensUpperRight_;