celldata2p2cadaptive.hh 13 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
// -*- 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/>.   *
 *****************************************************************************/
#ifndef DUMUX_ELEMENTDATA2P2C_ADAPTIVE_HH
#define DUMUX_ELEMENTDATA2P2C_ADAPTIVE_HH

Markus Wolff's avatar
Markus Wolff committed
22
23
#include <dumux/decoupled/2p2c/celldata2p2c.hh>
#include <dumux/decoupled/2p2c/celldata2p2cmultiphysics.hh>
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
/**
 * @file
 * @brief  Class including the variables and data of discretized data of the constitutive relations for one element
 * @author Benjamin Faigle, Markus Wolff
 */
namespace Dumux
{
/*!
 * \ingroup Adaptive2p2c
 */
//! Class including the data of a grid cell needed if an adaptive grid is used.
/*! The class provides model-specific functions needed to adapt the stored cell data to a new (adapted) grid.
 * Additionally, it provides the storage-infrastructure for explicit front tracking.
 *
 * @tparam TypeTag The Type Tag
 */
template<class TypeTag>
41
class CellData2P2CAdaptive: public CellData2P2CMultiPhysics<TypeTag>
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
{
private:
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GridView::Grid Grid;
    typedef typename GET_PROP_TYPE(TypeTag, Problem) Problem;
    typedef typename GET_PROP_TYPE(TypeTag, FluidState) FluidState;

    enum
    {
        dim = GridView::dimension
    };

    typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;

    enum
    {
        wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx,
        wCompIdx = Indices::wCompIdx, nCompIdx = Indices::nCompIdx
    };
    enum
    {
        numPhases = GET_PROP_VALUE(TypeTag, NumPhases)
    };
    typedef typename GridView::Traits::template Codim<0>::Entity Element;

68
    //! gives kind of pressure used (\f$ 0 = p_w \f$, \f$ 1 = p_n \f$, \f$ 2 = p_{global} \f$)
69
    static constexpr int pressureType = GET_PROP_VALUE(TypeTag, PressureFormulation);
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
    int upwindError_[numPhases];

public:
    //! A container for all necessary variables to map an old solution to a new grid
    /**
     * If the primary variables (pressure, total concentrations) are mapped to a new grid,
     * the secondary variables can be calulated. For the mapping between sons and father, it
     * is in addition necessary to know about how many sons live in each father ("count").
     * While only one phase pressure is PV, the method updateMaterialLaws() that
     * recalculates the secondary variables needs both phase pressures (initiated via the
     * capillary pressure of the last time step) to start the iteration to determine
     * appropriate phase composition and pc. Hence, both phase pressures are mapped
     * to the new solution.
     */
    struct AdaptedValues
    {
        Dune::FieldVector<Scalar,2> totalConcentration_; //!< Transport primary variables
        Dune::FieldVector<Scalar,2> pressure_; //!< Pressure primary variables
        //! Storage for volume derivatives, as transport estimate on old pressure field is not correct after refinement
        Dune::FieldVector<Scalar,3> volumeDerivatives_;
        Scalar cellVolume; //!< Cell volume for transformation of volume-specific primary variables
        FluxData2P2C<TypeTag> fluxData_; //!< Information on old flux direction
        int subdomain; //!< subdomain
        int count; //!< counts the number of cells averaged
        int isRefined; //!< Indicates if cell is refined
        AdaptedValues()
        {
            totalConcentration_=0.;
            pressure_ = 0.;
            volumeDerivatives_ =0.;
            cellVolume = 0.;
            subdomain=-1;
            count = 0;
            isRefined = false;
        }
    };


    //! Constructs an adaptive CellData object
109
    CellData2P2CAdaptive() : CellData2P2CMultiPhysics<TypeTag>()
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
    {
        for (int i = 0; i < numPhases;i++)
            upwindError_[i] = 0;
    }

    //! stores leaf cell primary variables to transfer to new indexing
    /**
     * Stores values to be adapted from the current CellData objects into
     * the adaptation container in order to be mapped on a new grid.
     *
     * @param adaptedValues Container for model-specific values to be adapted
     * @param element The element to be stored
     */
    void storeAdaptionValues(AdaptedValues& adaptedValues, const Element& element)
    {
        adaptedValues.totalConcentration_[wCompIdx]= this->massConcentration(wCompIdx);
        adaptedValues.totalConcentration_[nCompIdx]= this->massConcentration(nCompIdx);
        adaptedValues.pressure_[wPhaseIdx] = this->pressure(wPhaseIdx);
        adaptedValues.pressure_[nPhaseIdx] = this->pressure(nPhaseIdx);
        adaptedValues.volumeDerivatives_[wCompIdx] = this->dv(wPhaseIdx);
        adaptedValues.volumeDerivatives_[nCompIdx] = this->dv(nPhaseIdx);
        adaptedValues.volumeDerivatives_[2] = this->dv_dp();
        adaptedValues.cellVolume = element.geometry().volume();
        adaptedValues.subdomain = this->subdomain();
        adaptedValues.fluxData_=this->fluxData();
    }

Nicolas Schwenck's avatar
Nicolas Schwenck committed
137
    //! adds cell information to father element for possible averaging / coarsening
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
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
    /**
     * Sum up the adaptedValues (sons values) into father element. We store from leaf
     * upwards, so sons are stored first, then cells on the next leaf (=fathers)
     * can be averaged.
     *
     * @param adaptedValues Container for model-specific values to be adapted
     * @param adaptedValuesFather Values to be adapted of father cell
     * @param fatherElement The element of the father
     */
    static void storeAdaptionValues(AdaptedValues& adaptedValues,
                                    AdaptedValues& adaptedValuesFather,
                                    const Element& fatherElement)
    {
        adaptedValuesFather.totalConcentration_[wCompIdx]
                += adaptedValues.cellVolume* adaptedValues.totalConcentration_[wCompIdx];
        adaptedValuesFather.totalConcentration_[nCompIdx]
                += adaptedValues.cellVolume* adaptedValues.totalConcentration_[nCompIdx];
        // if all cells are summed up, re-convert mass into total concentrations
        Scalar fatherVolume = fatherElement.geometry().volume();
        if(adaptedValuesFather.count == (pow(2,dim)))
        {
            adaptedValuesFather.totalConcentration_[wCompIdx] /= fatherVolume;
            adaptedValuesFather.totalConcentration_[nCompIdx] /= fatherVolume;
        }
        adaptedValuesFather.cellVolume = fatherVolume;

        adaptedValuesFather.pressure_[wPhaseIdx] += adaptedValues.pressure_[wPhaseIdx] / adaptedValues.count;
        adaptedValuesFather.pressure_[nPhaseIdx] += adaptedValues.pressure_[nPhaseIdx] / adaptedValues.count;
        // apply maximum complexity for new cell
        adaptedValuesFather.subdomain = std::max(adaptedValuesFather.subdomain, adaptedValues.subdomain);
    }
    //! Set adapted values in CellData
    /**
     * This methods stores reconstructed values into the cellData object, by
     * this setting a newly mapped solution to the storage container of the
     * decoupled models.
     * In new cells, update estimate does not give meaningful results. We therefore
     * copy volume derivatives from old time step, and indicate that those are already availabe.
     *
     * @param adaptedValues Container for model-specific values to be adapted
     * @param element The element where things are stored.
     */
    void setAdaptionValues(AdaptedValues& adaptedValues, const Element& element)
    {
        // in new cells, there is a cellData object, but not yet a fluidstate.
        // To write the adapted values in this fluidstate, we have to enshure that
        // it was created. We therefore specify the type of FluidState that should be stored.
        this->setSubdomainAndFluidStateType(adaptedValues.subdomain);
        this->setMassConcentration(wCompIdx,
                adaptedValues.totalConcentration_[wCompIdx]);
        this->setMassConcentration(nCompIdx,
                adaptedValues.totalConcentration_[nCompIdx]);
        this->setPressure(wPhaseIdx,
                adaptedValues.pressure_[wPhaseIdx] / adaptedValues.count);
        this->setPressure(nPhaseIdx,
                adaptedValues.pressure_[nPhaseIdx] / adaptedValues.count);

        //copy flux directions
        this->fluxData()=adaptedValues.fluxData_;

        //indicate if cell was just refined in this TS
        this->wasRefined()= adaptedValues.isRefined;
        if(adaptedValues.isRefined)
        {
            this->dv(wPhaseIdx) = adaptedValues.volumeDerivatives_[wCompIdx];
            this->dv(nPhaseIdx) = adaptedValues.volumeDerivatives_[nCompIdx];
            this->dv_dp() = adaptedValues.volumeDerivatives_[2];
            this->volumeDerivativesAvailable(true);
        }
        else
            this->volumeDerivativesAvailable(false);  // recalculate volume derivatives
    }
    //! Reconstructs sons entries from data of father cell
    /**
     * Reconstructs an new solution from a father cell into for a newly
     * generated son cell. The new cell is stored into the global
     * adaptationMap.
     *
     * @param adaptationMap Global map storing all values to be adapted
     * @param father Entity Pointer to the father cell
     * @param son Entity Pointer to the newly created son cell
     * @param problem The problem
     */
    static void reconstructAdaptionValues(Dune::PersistentContainer<Grid, AdaptedValues>& adaptationMap,
            const Element& father, const Element& son, const Problem& problem)
    {
        // acess father and son
        AdaptedValues& adaptedValues = adaptationMap[son];
        AdaptedValues& adaptedValuesFather = adaptationMap[father];

        adaptedValues.subdomain = adaptedValuesFather.subdomain;
        adaptedValues.totalConcentration_[wCompIdx]
                   = adaptedValuesFather.totalConcentration_[wCompIdx] / adaptedValuesFather.count;
        adaptedValues.totalConcentration_[nCompIdx]
                   = adaptedValuesFather.totalConcentration_[nCompIdx] / adaptedValuesFather.count;

        // Introduce a hydrostatic pressure distribution inside the father cell
        Scalar gTimesHeight = problem.gravity()
                    * (son.geometry().center() - father.geometry().center());
237
238
        gTimesHeight *= (adaptedValuesFather.totalConcentration_[wCompIdx]+adaptedValuesFather.totalConcentration_[nCompIdx])/
                         problem.spatialParams().porosity(son);
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
270
271
272
273
274
275
276
//        int globalIdxSon = problem.variables().index(son);


        // recalculate new Primary variable and store pc (i.e. other pressure)
        if(pressureType == wPhaseIdx)
        {
            // recompute pressure and pc
            Scalar pressure = adaptedValuesFather.pressure_[wPhaseIdx] / adaptedValuesFather.count;
            Scalar pc = adaptedValuesFather.pressure_[nPhaseIdx] / adaptedValuesFather.count
                        - pressure;

            adaptedValues.pressure_[wPhaseIdx]
                    = pressure + gTimesHeight  ;
            adaptedValues.pressure_[nPhaseIdx]
                    = pc + adaptedValues.pressure_[wPhaseIdx];
        }
        else
        {
            // recompute pressure and pc
            Scalar pressure = adaptedValuesFather.pressure_[nPhaseIdx] / adaptedValuesFather.count;
            Scalar pc = pressure
                        - adaptedValuesFather.pressure_[wPhaseIdx] / adaptedValuesFather.count;

            adaptedValues.pressure_[nPhaseIdx]
                    = pressure + gTimesHeight  ;
            adaptedValues.pressure_[wPhaseIdx]
                    = adaptedValues.pressure_[nPhaseIdx] - pc;
        }
        // copy volume derivatives and compressibility, and indicate that cell was refined
        adaptedValues.volumeDerivatives_[wCompIdx] = adaptedValuesFather.volumeDerivatives_[wCompIdx];
        adaptedValues.volumeDerivatives_[nCompIdx] = adaptedValuesFather.volumeDerivatives_[nCompIdx];
        adaptedValues.volumeDerivatives_[2] = adaptedValuesFather.volumeDerivatives_[2];
        adaptedValues.isRefined = true;
    }
};

}
#endif