2010-08-24 06:11:18 -05:00
|
|
|
// $Id$
|
2009-06-10 07:08:28 -05:00
|
|
|
/*****************************************************************************
|
2010-10-18 03:59:51 -05:00
|
|
|
* Copyright (C) 2008-2009 by Melanie Darcis, Klaus Mosthaf *
|
2009-06-10 07:08:28 -05:00
|
|
|
* Copyright (C) 2009 by Andreas Lauser *
|
|
|
|
* Institute of Hydraulic Engineering *
|
|
|
|
* University of Stuttgart, Germany *
|
|
|
|
* email: <givenname>.<name>@iws.uni-stuttgart.de *
|
|
|
|
* *
|
2011-01-04 09:49:28 -06:00
|
|
|
* This program is free software: you can redistribute it and/or modify *
|
2009-06-10 07:08:28 -05:00
|
|
|
* it under the terms of the GNU General Public License as published by *
|
2011-01-04 09:49:28 -06:00
|
|
|
* the Free Software Foundation, either version 2 of the License, or *
|
|
|
|
* (at your option) any later version. *
|
2009-06-10 07:08:28 -05:00
|
|
|
* *
|
2011-01-04 09:49:28 -06:00
|
|
|
* 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/>. *
|
2009-06-10 07:08:28 -05:00
|
|
|
*****************************************************************************/
|
2010-11-11 09:01:37 -06:00
|
|
|
/*!
|
|
|
|
* \file
|
|
|
|
*
|
|
|
|
* \brief Tutorial problem for a fully coupled twophase box model.
|
|
|
|
*/
|
2011-02-04 09:11:02 -06:00
|
|
|
#ifndef DUMUX_TUTORIALPROBLEM_COUPLED_HH // guardian macro /*@\label{tutorial-coupled:guardian1}@*/
|
|
|
|
#define DUMUX_TUTORIALPROBLEM_COUPLED_HH // guardian macro /*@\label{tutorial-coupled:guardian2}@*/
|
2009-06-10 07:08:28 -05:00
|
|
|
|
|
|
|
// the numerical model
|
2010-08-03 12:11:34 -05:00
|
|
|
#include <dumux/boxmodels/2p/2pmodel.hh>
|
2009-06-10 07:08:28 -05:00
|
|
|
|
2010-11-08 12:02:44 -06:00
|
|
|
// the DUNE grid used
|
|
|
|
#include <dune/grid/sgrid.hh>
|
2009-06-10 07:08:28 -05:00
|
|
|
|
2010-11-08 12:02:44 -06:00
|
|
|
// spatialy dependent parameters
|
2010-04-07 04:49:17 -05:00
|
|
|
#include "tutorialspatialparameters_coupled.hh"
|
2009-06-10 07:08:28 -05:00
|
|
|
|
2010-03-25 08:02:05 -05:00
|
|
|
namespace Dumux
|
2009-06-10 07:08:28 -05:00
|
|
|
{
|
|
|
|
|
|
|
|
// forward declaration of the problem class
|
|
|
|
template <class TypeTag>
|
|
|
|
class TutorialProblemCoupled;
|
|
|
|
|
|
|
|
namespace Properties
|
|
|
|
{
|
|
|
|
// create a new type tag for the problem
|
2009-06-18 11:46:37 -05:00
|
|
|
NEW_TYPE_TAG(TutorialProblemCoupled, INHERITS_FROM(BoxTwoP)); /*@\label{tutorial-coupled:create-type-tag}@*/
|
2009-06-10 07:08:28 -05:00
|
|
|
|
|
|
|
// Set the "Problem" property
|
2009-06-18 11:46:37 -05:00
|
|
|
SET_PROP(TutorialProblemCoupled, Problem) /*@\label{tutorial-coupled:set-problem}@*/
|
2010-11-08 12:02:44 -06:00
|
|
|
{ typedef Dumux::TutorialProblemCoupled<TypeTag> type;};
|
2009-06-10 07:08:28 -05:00
|
|
|
|
2009-06-18 11:46:37 -05:00
|
|
|
// Set the grid
|
|
|
|
SET_PROP(TutorialProblemCoupled, Grid) /*@\label{tutorial-coupled:set-grid}@*/
|
2010-06-15 04:52:10 -05:00
|
|
|
{
|
2009-06-18 11:46:37 -05:00
|
|
|
typedef Dune::SGrid<2,2> type;
|
|
|
|
static type *create() /*@\label{tutorial-coupled:create-grid-method}@*/
|
|
|
|
{
|
2010-04-06 04:47:25 -05:00
|
|
|
typedef typename type::ctype ctype;
|
2011-02-04 09:11:02 -06:00
|
|
|
Dune::FieldVector<int, 2> cellRes; // vector holding resolution of the grid
|
|
|
|
Dune::FieldVector<ctype, 2> lowerLeft(0.0); // Coordinate of lower left corner of the grid
|
|
|
|
Dune::FieldVector<ctype, 2> upperRight; // Coordinate of upper right corner of the grid
|
2010-12-02 05:14:57 -06:00
|
|
|
cellRes[0] = 100;
|
|
|
|
cellRes[1] = 1;
|
2009-06-18 11:46:37 -05:00
|
|
|
upperRight[0] = 300;
|
2009-08-05 09:11:46 -05:00
|
|
|
upperRight[1] = 60;
|
2009-06-18 11:46:37 -05:00
|
|
|
return new Dune::SGrid<2,2>(cellRes,
|
|
|
|
lowerLeft,
|
|
|
|
upperRight);
|
|
|
|
}
|
|
|
|
};
|
2009-06-10 07:08:28 -05:00
|
|
|
|
2010-10-18 03:59:51 -05:00
|
|
|
// Set the wetting phase
|
|
|
|
SET_PROP(TutorialProblemCoupled, WettingPhase) /*@\label{tutorial-coupled:2p-system-start}@*/
|
2010-04-07 04:49:17 -05:00
|
|
|
{
|
2010-11-08 12:02:44 -06:00
|
|
|
private: typedef typename GET_PROP_TYPE(TypeTag, PTAG(Scalar)) Scalar;
|
|
|
|
public: typedef Dumux::LiquidPhase<Scalar, Dumux::H2O<Scalar> > type; /*@\label{tutorial-coupled:wettingPhase}@*/
|
2010-04-07 04:49:17 -05:00
|
|
|
};
|
2010-04-06 04:47:25 -05:00
|
|
|
|
2010-10-18 03:59:51 -05:00
|
|
|
// Set the non-wetting phase
|
|
|
|
SET_PROP(TutorialProblemCoupled, NonwettingPhase)
|
|
|
|
{
|
2010-11-08 12:02:44 -06:00
|
|
|
private: typedef typename GET_PROP_TYPE(TypeTag, PTAG(Scalar)) Scalar;
|
|
|
|
public: typedef Dumux::LiquidPhase<Scalar, Dumux::Oil<Scalar> > type; /*@\label{tutorial-coupled:nonwettingPhase}@*/
|
2010-10-18 03:59:51 -05:00
|
|
|
}; /*@\label{tutorial-coupled:2p-system-end}@*/
|
|
|
|
|
|
|
|
|
2010-04-12 10:00:00 -05:00
|
|
|
// Set the spatial parameters
|
|
|
|
SET_PROP(TutorialProblemCoupled, SpatialParameters) /*@\label{tutorial-coupled:set-spatialparameters}@*/
|
2010-11-08 12:02:44 -06:00
|
|
|
{ typedef Dumux::TutorialSpatialParametersCoupled<TypeTag> type; };
|
2009-06-10 07:08:28 -05:00
|
|
|
|
2009-06-18 11:46:37 -05:00
|
|
|
// Disable gravity
|
|
|
|
SET_BOOL_PROP(TutorialProblemCoupled, EnableGravity, false); /*@\label{tutorial-coupled:gravity}@*/
|
2009-06-10 07:08:28 -05:00
|
|
|
}
|
|
|
|
|
2010-11-11 09:01:37 -06:00
|
|
|
/*!
|
|
|
|
* \ingroup TwoPBoxModel
|
|
|
|
*
|
|
|
|
* \brief Tutorial problem for a fully coupled twophase box model.
|
|
|
|
*/
|
|
|
|
|
2009-06-18 11:46:37 -05:00
|
|
|
// Definition of the actual problem
|
2010-11-08 12:02:44 -06:00
|
|
|
template <class TypeTag>
|
2010-08-03 12:11:34 -05:00
|
|
|
class TutorialProblemCoupled : public TwoPProblem<TypeTag> /*@\label{tutorial-coupled:def-problem}@*/
|
2009-06-10 07:08:28 -05:00
|
|
|
{
|
2010-08-03 12:11:34 -05:00
|
|
|
typedef TwoPProblem<TypeTag> ParentType;
|
2010-11-08 12:02:44 -06:00
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, PTAG(Scalar)) Scalar;
|
2010-08-03 12:11:34 -05:00
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, PTAG(GridView)) GridView;
|
2009-06-10 07:08:28 -05:00
|
|
|
|
2010-11-08 12:02:44 -06:00
|
|
|
// Grid dimension
|
|
|
|
enum { dim = GridView::dimension };
|
2009-06-10 07:08:28 -05:00
|
|
|
|
2010-11-08 12:02:44 -06:00
|
|
|
// Types from DUNE-Grid
|
2010-08-03 12:11:34 -05:00
|
|
|
typedef typename GridView::template Codim<0>::Entity Element;
|
|
|
|
typedef typename GridView::template Codim<dim>::Entity Vertex;
|
|
|
|
typedef typename GridView::Intersection Intersection;
|
2010-11-08 12:02:44 -06:00
|
|
|
typedef Dune::FieldVector<Scalar, dim> GlobalPosition;
|
2010-08-03 12:11:34 -05:00
|
|
|
|
2010-11-08 12:02:44 -06:00
|
|
|
// Dumux specific types
|
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, PTAG(TimeManager)) TimeManager;
|
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, PTAG(TwoPIndices)) Indices;
|
2010-08-05 07:11:43 -05:00
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, PTAG(PrimaryVariables)) PrimaryVariables;
|
2010-08-03 12:11:34 -05:00
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, PTAG(BoundaryTypes)) BoundaryTypes;
|
2009-06-10 07:08:28 -05:00
|
|
|
typedef typename GET_PROP_TYPE(TypeTag, PTAG(FVElementGeometry)) FVElementGeometry;
|
2010-08-05 07:11:43 -05:00
|
|
|
|
2009-06-10 07:08:28 -05:00
|
|
|
public:
|
2010-08-03 12:11:34 -05:00
|
|
|
TutorialProblemCoupled(TimeManager &timeManager,
|
|
|
|
const GridView &gridView)
|
|
|
|
: ParentType(timeManager, gridView)
|
2010-10-15 10:11:37 -05:00
|
|
|
{
|
|
|
|
}
|
|
|
|
|
2010-12-10 05:05:31 -06:00
|
|
|
// Specifies the problem name. This is used as a prefix for files
|
2010-11-08 12:02:44 -06:00
|
|
|
// generated by the simulation.
|
2010-10-15 10:11:37 -05:00
|
|
|
const char *name() const
|
|
|
|
{ return "tutorial_coupled"; }
|
2009-06-10 07:08:28 -05:00
|
|
|
|
2010-11-11 08:49:20 -06:00
|
|
|
//! Returns true if a restart file should be written.
|
|
|
|
/* The default behaviour is to write no restart file.
|
|
|
|
*/
|
|
|
|
bool shouldWriteRestartFile() const /*@\label{tutorial-coupled:restart}@*/
|
|
|
|
{
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
|
|
|
|
//! Returns true if the current solution should be written to disk (i.e. as a VTK file)
|
|
|
|
/*! The default behaviour is to write out the solution for
|
2010-12-10 05:05:31 -06:00
|
|
|
* every time step. Else, the user has to change the divisor in this function.
|
2010-11-11 08:49:20 -06:00
|
|
|
*/
|
|
|
|
bool shouldWriteOutput() const /*@\label{tutorial-coupled:output}@*/
|
|
|
|
{
|
|
|
|
return this->timeManager().timeStepIndex() > 0 &&
|
|
|
|
(this->timeManager().timeStepIndex() % 1 == 0);
|
|
|
|
}
|
|
|
|
|
2010-11-08 12:02:44 -06:00
|
|
|
// Return the temperature within a finite volume. We use constant
|
|
|
|
// 10 degrees Celsius.
|
2010-08-03 12:11:34 -05:00
|
|
|
Scalar temperature(const Element &element,
|
2009-11-02 08:02:06 -06:00
|
|
|
const FVElementGeometry &fvElemGeom,
|
2010-08-03 12:11:34 -05:00
|
|
|
int scvIdx) const
|
2009-06-10 07:08:28 -05:00
|
|
|
{ return 283.15; };
|
|
|
|
|
2009-06-18 11:46:37 -05:00
|
|
|
// Specifies which kind of boundary condition should be used for
|
2010-11-08 12:02:44 -06:00
|
|
|
// which equation for a finite volume on the boundary.
|
|
|
|
void boundaryTypes(BoundaryTypes &BCtypes, const Vertex &vertex) const
|
2010-06-15 04:52:10 -05:00
|
|
|
{
|
2010-10-15 08:37:28 -05:00
|
|
|
const GlobalPosition &pos = vertex.geometry().center();
|
2010-11-08 12:02:44 -06:00
|
|
|
if (pos[0] < eps_) // Dirichlet conditions on left boundary
|
|
|
|
BCtypes.setAllDirichlet();
|
2009-06-18 11:46:37 -05:00
|
|
|
else // neuman for the remaining boundaries
|
2010-11-08 12:02:44 -06:00
|
|
|
BCtypes.setAllNeumann();
|
2009-06-10 07:08:28 -05:00
|
|
|
|
|
|
|
}
|
|
|
|
|
2010-11-08 12:02:44 -06:00
|
|
|
// Evaluate the Dirichlet boundary conditions for a finite volume
|
|
|
|
// on the grid boundary. Here, the 'values' parameter stores
|
2009-06-18 11:46:37 -05:00
|
|
|
// primary variables.
|
2010-10-15 08:37:28 -05:00
|
|
|
void dirichlet(PrimaryVariables &values, const Vertex &vertex) const
|
2009-06-10 07:08:28 -05:00
|
|
|
{
|
2010-04-06 04:47:25 -05:00
|
|
|
values[Indices::pwIdx] = 200.0e3; // 200 kPa = 2 bar
|
|
|
|
values[Indices::SnIdx] = 0.0; // 0 % oil saturation on left boundary
|
2009-06-10 07:08:28 -05:00
|
|
|
}
|
|
|
|
|
2010-11-08 12:02:44 -06:00
|
|
|
// Evaluate the boundary conditions for a Neumann boundary
|
|
|
|
// segment. Here, the 'values' parameter stores the mass flux in
|
|
|
|
// normal direction of each phase. Negative values mean influx.
|
2010-08-05 07:11:43 -05:00
|
|
|
void neumann(PrimaryVariables &values,
|
2010-08-03 12:11:34 -05:00
|
|
|
const Element &element,
|
|
|
|
const FVElementGeometry &fvElemGeom,
|
|
|
|
const Intersection &isIt,
|
|
|
|
int scvIdx,
|
|
|
|
int boundaryFaceIdx) const
|
2009-06-10 07:08:28 -05:00
|
|
|
{
|
2010-06-15 04:52:10 -05:00
|
|
|
const GlobalPosition &pos =
|
2009-08-05 09:11:46 -05:00
|
|
|
fvElemGeom.boundaryFace[boundaryFaceIdx].ipGlobal;
|
2009-06-10 07:36:33 -05:00
|
|
|
Scalar right = this->bboxMax()[0];
|
2010-11-11 08:49:20 -06:00
|
|
|
// extraction of oil on the right boundary for approx. 1.e6 seconds
|
2010-12-02 05:14:57 -06:00
|
|
|
if (pos[0] > right - eps_) {
|
|
|
|
// oil outflux of 30 g/(m * s) on the right boundary.
|
2010-04-06 04:47:25 -05:00
|
|
|
values[Indices::contiWEqIdx] = 0;
|
2010-12-02 05:14:57 -06:00
|
|
|
values[Indices::contiNEqIdx] = 3e-2;
|
2009-06-10 07:08:28 -05:00
|
|
|
} else {
|
2010-11-08 12:02:44 -06:00
|
|
|
// no-flow on the remaining Neumann-boundaries.
|
2010-04-06 04:47:25 -05:00
|
|
|
values[Indices::contiWEqIdx] = 0;
|
|
|
|
values[Indices::contiNEqIdx] = 0;
|
2009-06-10 07:08:28 -05:00
|
|
|
}
|
|
|
|
}
|
2009-06-18 11:46:37 -05:00
|
|
|
|
|
|
|
// Evaluate the initial value for a control volume. For this
|
|
|
|
// method, the 'values' parameter stores primary variables.
|
2010-08-05 07:11:43 -05:00
|
|
|
void initial(PrimaryVariables &values,
|
2010-08-03 12:11:34 -05:00
|
|
|
const Element &element,
|
2009-06-10 07:08:28 -05:00
|
|
|
const FVElementGeometry &fvElemGeom,
|
2010-08-03 12:11:34 -05:00
|
|
|
int scvIdx) const
|
2009-06-10 07:08:28 -05:00
|
|
|
{
|
2010-04-06 04:47:25 -05:00
|
|
|
values[Indices::pwIdx] = 200.0e3; // 200 kPa = 2 bar
|
|
|
|
values[Indices::SnIdx] = 1.0;
|
2009-06-10 07:08:28 -05:00
|
|
|
}
|
|
|
|
|
2009-06-18 11:46:37 -05:00
|
|
|
// Evaluate the source term for all phases within a given
|
2010-11-08 12:02:44 -06:00
|
|
|
// sub-control-volume. In this case, the 'values' parameter stores
|
|
|
|
// the rate mass generated or annihilate per volume unit. Positive
|
|
|
|
// values mean that mass is created.
|
2010-08-05 07:11:43 -05:00
|
|
|
void source(PrimaryVariables &values,
|
2010-08-03 12:11:34 -05:00
|
|
|
const Element &element,
|
2009-06-18 11:46:37 -05:00
|
|
|
const FVElementGeometry &fvElemGeom,
|
2010-08-03 12:11:34 -05:00
|
|
|
int scvIdx) const
|
2009-06-10 07:08:28 -05:00
|
|
|
{
|
2010-04-06 04:47:25 -05:00
|
|
|
values[Indices::contiWEqIdx] = 0.0;
|
|
|
|
values[Indices::contiNEqIdx]= 0.0;
|
2009-06-10 07:08:28 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
private:
|
2010-11-08 12:02:44 -06:00
|
|
|
// small epsilon value
|
2009-06-10 07:08:28 -05:00
|
|
|
static const Scalar eps_ = 3e-6;
|
|
|
|
};
|
2009-06-18 11:46:37 -05:00
|
|
|
}
|
2009-06-10 07:08:28 -05:00
|
|
|
|
|
|
|
#endif
|