opm-simulators/examples/tutorialproblem_coupled.hh

217 lines
8.8 KiB
C++
Raw Normal View History

// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* Copyright (C) 2008-2009 by Melanie Darcis, Klaus Mosthaf *
* Copyright (C) 2009-2012 by Andreas Lauser *
* Institute for Modelling Hydraulic and Environmental Systems *
* University of Stuttgart, Germany *
* email: <givenname>.<name>@iws.uni-stuttgart.de *
* *
* 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/>. *
*****************************************************************************/
2010-11-11 09:01:37 -06:00
/*!
* \file
*
* \brief Tutorial problem for a fully coupled twophase box model.
*/
#ifndef DUMUX_TUTORIAL_PROBLEM_COUPLED_HH // guardian macro /*@\label{tutorial-coupled:guardian1}@*/
#define DUMUX_TUTORIAL_PROBLEM_COUPLED_HH // guardian macro /*@\label{tutorial-coupled:guardian2}@*/
// The numerical model
2010-08-03 12:11:34 -05:00
#include <dumux/boxmodels/2p/2pmodel.hh>
// The DUNE grid used
#include <dune/grid/yaspgrid.hh>
#include "tutorialspatialparameters_coupled.hh"
// The components that are used
#include <dumux/material/components/h2o.hh>
#include <dumux/material/components/oil.hh>
#include <dumux/common/cubegridcreator.hh>
namespace Dumux{
// Forward declaration of the problem class
template <class TypeTag>
class TutorialProblemCoupled;
namespace Properties {
// Create a new type tag for the problem
NEW_TYPE_TAG(TutorialProblemCoupled, INHERITS_FROM(BoxTwoP, TutorialSpatialParametersCoupled)); /*@\label{tutorial-coupled:create-type-tag}@*/
// Set the "Problem" property
SET_PROP(TutorialProblemCoupled, Problem) /*@\label{tutorial-coupled:set-problem}@*/
{ typedef Dumux::TutorialProblemCoupled<TypeTag> type;};
// Set grid and the grid creator to be used
SET_TYPE_PROP(TutorialProblemCoupled, Grid, Dune::YaspGrid</*dim=*/2>); /*@\label{tutorial-coupled:set-grid}@*/
SET_TYPE_PROP(TutorialProblemCoupled, GridCreator, Dumux::CubeGridCreator<TypeTag>); /*@\label{tutorial-coupled:set-gridcreator}@*/
// Set the wetting phase
SET_PROP(TutorialProblemCoupled, WettingPhase) /*@\label{tutorial-coupled:2p-system-start}@*/
{
private: typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
public: typedef Dumux::LiquidPhase<Scalar, Dumux::H2O<Scalar> > type; /*@\label{tutorial-coupled:wettingPhase}@*/
};
// Set the non-wetting phase
SET_PROP(TutorialProblemCoupled, NonwettingPhase)
{
private: typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
public: typedef Dumux::LiquidPhase<Scalar, Dumux::Oil<Scalar> > type; /*@\label{tutorial-coupled:nonwettingPhase}@*/
}; /*@\label{tutorial-coupled:2p-system-end}@*/
SET_TYPE_PROP(TutorialProblemCoupled, FluidSystem, Dumux::TwoPImmiscibleFluidSystem<TypeTag>);
// Disable gravity
SET_BOOL_PROP(TutorialProblemCoupled, EnableGravity, false); /*@\label{tutorial-coupled:gravity}@*/
}
2010-11-11 09:01:37 -06:00
/*!
* \ingroup TwoPBoxModel
*
* \brief Tutorial problem for a fully coupled twophase box model.
*/
template <class TypeTag>
2010-08-03 12:11:34 -05:00
class TutorialProblemCoupled : public TwoPProblem<TypeTag> /*@\label{tutorial-coupled:def-problem}@*/
{
2010-08-03 12:11:34 -05:00
typedef TwoPProblem<TypeTag> ParentType;
typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
// Grid dimension
enum { dim = GridView::dimension };
// 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;
typedef Dune::FieldVector<Scalar, dim> GlobalPosition;
2010-08-03 12:11:34 -05:00
// Dumux specific types
typedef typename GET_PROP_TYPE(TypeTag, TimeManager) TimeManager;
typedef typename GET_PROP_TYPE(TypeTag, TwoPIndices) Indices;
typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
typedef typename GET_PROP_TYPE(TypeTag, BoundaryTypes) BoundaryTypes;
typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry;
public:
2010-08-03 12:11:34 -05:00
TutorialProblemCoupled(TimeManager &timeManager,
const GridView &gridView)
: ParentType(timeManager, gridView)
, eps_(3e-6)
2010-10-15 10:11:37 -05:00
{
}
//! Specifies the problem name. This is used as a prefix for files
//! generated by the simulation.
2010-10-15 10:11:37 -05:00
const char *name() const
{ return "tutorial_coupled"; }
//! Returns true if a restart file should be written.
2010-11-11 08:49:20 -06:00
bool shouldWriteRestartFile() const /*@\label{tutorial-coupled:restart}@*/
{ return false; }
2010-11-11 08:49:20 -06:00
//! Returns true if the current solution should be written to disk
//! as a VTK file
2010-11-11 08:49:20 -06:00
bool shouldWriteOutput() const /*@\label{tutorial-coupled:output}@*/
{
2011-06-27 11:41:17 -05:00
return
this->timeManager().timeStepIndex() > 0 &&
(this->timeManager().timeStepIndex() % 1 == 0);
2010-11-11 08:49:20 -06:00
}
//! Returns the temperature within a finite volume. We use constant
//! 10 degrees Celsius.
Scalar temperature() const
{ return 283.15; };
//! Specifies which kind of boundary condition should be used for
//! which equation for a finite volume on the boundary.
void boundaryTypes(BoundaryTypes &BCtypes, const Vertex &vertex) const
{
const GlobalPosition &pos = vertex.geometry().center();
if (pos[0] < eps_) // Dirichlet conditions on left boundary
BCtypes.setAllDirichlet();
else // neuman for the remaining boundaries
BCtypes.setAllNeumann();
}
//! Evaluates the Dirichlet boundary conditions for a finite volume
//! on the grid boundary. Here, the 'values' parameter stores
//! primary variables.
void dirichlet(PrimaryVariables &values, const Vertex &vertex) const
{
values[Indices::pwIdx] = 200.0e3; // 200 kPa = 2 bar
values[Indices::SnIdx] = 0.0; // 0 % oil saturation on left boundary
}
//! Evaluates the boundary conditions for a Neumann boundary
//! segment. Here, the 'values' parameter stores the mass flux in
//! [kg/(m^2 * s)] in normal direction of each phase. Negative
//! values mean influx.
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
{
const GlobalPosition &pos =
fvElemGeom.boundaryFace[boundaryFaceIdx].ipGlobal;
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
if (pos[0] > right - eps_) {
// oil outflux of 30 g/(m * s) on the right boundary.
values[Indices::contiWEqIdx] = 0;
values[Indices::contiNEqIdx] = 3e-2;
} else {
// no-flow on the remaining Neumann-boundaries.
values[Indices::contiWEqIdx] = 0;
values[Indices::contiNEqIdx] = 0;
}
}
//! Evaluates the initial value for a control volume. For this
//! method, the 'values' parameter stores primary variables.
void initial(PrimaryVariables &values,
2010-08-03 12:11:34 -05:00
const Element &element,
const FVElementGeometry &fvElemGeom,
2010-08-03 12:11:34 -05:00
int scvIdx) const
{
values[Indices::pwIdx] = 200.0e3; // 200 kPa = 2 bar
values[Indices::SnIdx] = 1.0;
}
//! Evaluates the source term for all phases within a given
//! sub-control-volume. In this case, the 'values' parameter
//! stores the rate mass generated or annihilated per volume unit
//! in [kg / (m^3 * s)]. Positive values mean that mass is created.
void source(PrimaryVariables &values,
2010-08-03 12:11:34 -05:00
const Element &element,
const FVElementGeometry &fvElemGeom,
2010-08-03 12:11:34 -05:00
int scvIdx) const
{
values[Indices::contiWEqIdx] = 0.0;
values[Indices::contiNEqIdx]= 0.0;
}
private:
// small epsilon value
Scalar eps_;
};
}
#endif