henry1p2cproblem.hh 8.5 KB
Newer Older
Lena Walter's avatar
Lena Walter committed
1
2
3
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
4
5
 *   See the file COPYING for full copying permissions.                      *
 *                                                                           *
Lena Walter's avatar
Lena Walter committed
6
7
8
9
10
11
12
 *   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          *
13
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the            *
Lena Walter's avatar
Lena Walter committed
14
15
16
17
18
19
20
 *   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/>.   *
 *****************************************************************************/
/**
 * \file
Lena Walter's avatar
Lena Walter committed
21
22
 * \brief Definition of the 1p2c Henry problem
  *
Lena Walter's avatar
Lena Walter committed
23
24
25
26
 */
#ifndef DUMUX_HENRY1P2C_PROBLEM_HH
#define DUMUX_HENRY1P2C_PROBLEM_HH

27
28
#include <dumux/porousmediumflow/1p2c/implicit/model.hh>
#include <dumux/porousmediumflow/implicit/problem.hh>
Lena Walter's avatar
Lena Walter committed
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43

#include "watersaltfluidsystem.hh"
#include "henry1p2cspatialparameters.hh"

namespace Dumux
{

template <class TypeTag>
class Henry1p2cProblem;

namespace Properties
{
NEW_TYPE_TAG(Henry1p2cProblem, INHERITS_FROM(BoxOnePTwoC));

// Set the grid type
Timo Koch's avatar
Timo Koch committed
44
SET_TYPE_PROP(Henry1p2cProblem, Grid, Dune::YaspGrid<2>);
Lena Walter's avatar
Lena Walter committed
45
46

// Set the problem property
Katharina Heck's avatar
Katharina Heck committed
47
SET_TYPE_PROP(Henry1p2cProblem, Problem, Henry1p2cProblem<TypeTag>);
Lena Walter's avatar
Lena Walter committed
48

Lena Walter's avatar
Lena Walter committed
49
50
51
52
53
SET_TYPE_PROP(Henry1p2cProblem, FluidSystem, WaterSaltFluidSystem<TypeTag>);

// Set the spatial parameters
SET_TYPE_PROP(Henry1p2cProblem,
              SpatialParams,
Katharina Heck's avatar
Katharina Heck committed
54
              Henry1p2cSpatialParams<TypeTag>);
Lena Walter's avatar
Lena Walter committed
55
56
57
58
59
60
61
62

//Define whether mole(true) or mass (false) fractions are used
SET_BOOL_PROP(Henry1p2cProblem, UseMoles, false);
SET_BOOL_PROP(Henry1p2cProblem, ProblemEnableGravity, true);
}


template <class TypeTag>
Bernd Flemisch's avatar
Bernd Flemisch committed
63
class Henry1p2cProblem : public  ImplicitPorousMediaProblem<TypeTag>
Lena Walter's avatar
Lena Walter committed
64
{
Bernd Flemisch's avatar
Bernd Flemisch committed
65
    typedef ImplicitPorousMediaProblem<TypeTag> ParentType;
Lena Walter's avatar
Lena Walter committed
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83

    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry;
    typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
    typedef typename GET_PROP_TYPE(TypeTag, BoundaryTypes) BoundaryTypes;
    typedef typename GET_PROP_TYPE(TypeTag, TimeManager) TimeManager;
    typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;

    // copy some indices for convenience
    typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
      enum {
        // Grid and world dimension
        dim = GridView::dimension,
        dimWorld = GridView::dimensionworld,

        // indices of the primary variables
        pressureIdx = Indices::pressureIdx,
Bernd Flemisch's avatar
Bernd Flemisch committed
84
        massOrMoleFracIdx = Indices::massOrMoleFracIdx,
Lena Walter's avatar
Lena Walter committed
85
86

        // indices of the equations
Bernd Flemisch's avatar
Bernd Flemisch committed
87
88
        conti0EqIdx = Indices::conti0EqIdx,
        transportEqIdx = Indices::transportEqIdx
Lena Walter's avatar
Lena Walter committed
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
    };


    typedef typename GridView::template Codim<0>::Entity Element;
    typedef typename GridView::template Codim<0>::Iterator ElementIterator;
    typedef typename GridView::template Codim<dim>::Entity Vertex;
    typedef typename GridView::Intersection Intersection;

    typedef Dune::FieldVector<Scalar, dimWorld> GlobalPosition;

public:
    Henry1p2cProblem(TimeManager &timeManager,
            const GridView &gridView)
    : ParentType(timeManager, gridView)
      {

105
106
107
//         FVElementGeometry fvGeom;
//         ElementIterator elemIt = gridView.template begin<0>();
//         const ElementIterator endIt = gridView.template end<0>();
Lena Walter's avatar
Lena Walter committed
108
109
110
111
112

        freshWaterFluxRate_= GET_RUNTIME_PARAM(TypeTag, Scalar, Problem.freshWaterFluxRate);
      }


Lena Walter's avatar
Lena Walter committed
113
    bool shouldWriteRestartFile() const
Lena Walter's avatar
Lena Walter committed
114
           {
Lena Walter's avatar
Lena Walter committed
115
               return false;
Lena Walter's avatar
Lena Walter committed
116
117
           }

Lena Walter's avatar
Lena Walter committed
118
119
120
121
122
123
124
125
//    bool shouldWriteOutput() const
//           {
//               return
//                   this->timeManager().timeStepIndex() == 0 ||
//                   this->timeManager().timeStepIndex() %5 ==0 ||
//                   this->timeManager().willBeFinished();
//           }

Lena Walter's avatar
Lena Walter committed
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186



    /*!
     * \name Problem parameters
     */
    // \{

    /*!
     * \brief The problem name.
     *
     * This is used as a prefix for files generated by the simulation.
     */
    const char *name() const
    { return "Henry1p2c"; }

    /*!
     * \brief Returns the temperature within the domain.
     *
     * This problem assumes a temperature of 36 degrees Celsius.
     */
    Scalar temperature() const
    { return 273.15 + 20; }; // in [K]

    // \}

    /*!
     * \name Boundary conditions
     */
    // \{

    /*!
     * \brief Specifies which kind of boundary condition should be
     *        used for which equation on a given boundary segment.
     */
    void boundaryTypes(BoundaryTypes &values, const Vertex &vertex) const
    {
        const GlobalPosition globalPos = vertex.geometry().center();

        values.setAllDirichlet();

       if(globalPos[0]<eps_  || globalPos[1]<eps_ || globalPos[1]>1-eps_)
        {
            values.setAllNeumann();
        }

    }

    /*!
     * \brief Evaluate the boundary conditions for a dirichlet
     *        boundary segment.
     *
     * For this method, the \a values parameter stores primary variables.
     */
    void dirichlet(PrimaryVariables &values, const Vertex &vertex) const
    {
        const GlobalPosition globalPos = vertex.geometry().center();

        initial_(values, globalPos);
        if(globalPos[0]>2-eps_ && globalPos[1]<0.8+eps_)
        {
Bernd Flemisch's avatar
Bernd Flemisch committed
187
            values[massOrMoleFracIdx] = 0.03922;
Lena Walter's avatar
Lena Walter committed
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
        }

    }

    /*!
     * \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 component. Negative values mean
     * influx.
     */
    void neumann(PrimaryVariables &values,
                 const Element &element,
                 const FVElementGeometry &fvElemGeom,
                 const Intersection &is,
                 int scvIdx,
                 int boundaryFaceIdx) const
    {
        const GlobalPosition &globalPos = element.geometry().corner(scvIdx);
        values = 0;

        if(globalPos[0]<eps_)
              {
Bernd Flemisch's avatar
Bernd Flemisch committed
212
            values[conti0EqIdx] = -freshWaterFluxRate_; //-6.6E-2; //6.6e-5*m2/s*1000kg/m^3
Lena Walter's avatar
Lena Walter committed
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
              }

    }
    // \}

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

    /*!
     * \brief Evaluate the source term for all phases within a given
     *        sub-control-volume.
     *
     * For this method, the \a values parameter stores the rate mass
     * of a component is generated or annihilate per volume
     * unit. Positive values mean that mass is created, negative ones
     * mean that it vanishes.
     */
    void sourceAtPos(PrimaryVariables &values,
                     const GlobalPosition &globalPos) const
    {
        values = Scalar(0.0);

    }

    /*!
     * \brief Evaluate the initial value for a control volume.
     *
     * For this method, the \a values parameter stores primary
     * variables.
     */
    void initial(PrimaryVariables &values,
                 const Element &element,
                 const FVElementGeometry &fvElemGeom,
                 int scvIdx) const
    {
        const GlobalPosition &globalPos
            = element.geometry().corner(scvIdx);

        initial_(values, globalPos);
    }



    // \}

private:



    // the internal method for the initial condition
    void initial_(PrimaryVariables &values,
                  const GlobalPosition &globalPos) const
    {
        Scalar densityW = 1025;
        values[pressureIdx] = 1.0133e5+(depthBOR_-globalPos[1])*densityW*9.81; //initial condition for the pressure
Bernd Flemisch's avatar
Bernd Flemisch committed
270
        values[massOrMoleFracIdx] = 0.0; //initial condition for the salt molefraction
Lena Walter's avatar
Lena Walter committed
271
272
273
274
275

    }

    Scalar freshWaterFluxRate_;
    Scalar totalInjectionVolume_;
276
277
    static constexpr Scalar eps_ = 1e-6;
    static constexpr Scalar depthBOR_ = 1.0; // [m]
Lena Walter's avatar
Lena Walter committed
278
279
280
}; //end namespace
}
#endif