stokesncnimodel.hh 8.33 KB
Newer Older
Timo Koch's avatar
Timo Koch committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
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
// -*- 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/>.   *
 *****************************************************************************/
/*!
 * \file
 *
 * \brief Adaption of the BOX scheme to the non-isothermal
 *        compositional stokes model (with n components).
 */
#ifndef DUMUX_STOKESNCNI_MODEL_HH
#define DUMUX_STOKESNCNI_MODEL_HH

#include <dumux/freeflow/stokesnc/stokesncmodel.hh>

#include "stokesncnilocalresidual.hh"
#include "stokesncniproperties.hh"

namespace Dumux {
/*!
 * \ingroup BoxStokesncniModel
 * \brief Adaption of the BOX scheme to the non-isothermal compositional n-component Stokes model.
 *
 * This model implements a non-isothermal n-component Stokes flow of a fluid
 * solving a momentum balance, a mass balance, a conservation equation for one component,
 * and one balance equation for the energy.
 *
 * Momentum Balance:
 * \f[
\frac{\partial \left(\varrho_g {\boldsymbol{v}}_g\right)}{\partial t}
+ \boldsymbol{\nabla} \boldsymbol{\cdot} \left(p_g {\bf {I}}
- \mu_g \left(\boldsymbol{\nabla} \boldsymbol{v}_g
+ \boldsymbol{\nabla} \boldsymbol{v}_g^T\right)\right)
- \varrho_g {\bf g} = 0,
 * \f]
 *
 * Mass balance equation:
 * \f[
\frac{\partial \varrho_g}{\partial t} + \boldsymbol{\nabla}\boldsymbol{\cdot}\left(\varrho_g {\boldsymbol{v}}_g\right) - q_g = 0
 * \f]
 *
 * Component mass balance equation:
 * \f[
 \frac{\partial \left(\varrho_g X_g^\kappa\right)}{\partial t}
 + \boldsymbol{\nabla} \boldsymbol{\cdot} \left( \varrho_g {\boldsymbol{v}}_g X_g^\kappa
 - D^\kappa_g \varrho_g \frac{M^\kappa}{M_g} \boldsymbol{\nabla} x_g^\kappa \right)
 - q_g^\kappa = 0
 * \f]
 *
 * Energy balance equation:
 * \f[
\frac{\partial (\varrho_g  u_g)}{\partial t}
67
+ \boldsymbol{\nabla} \boldsymbol{\cdot} \left( \varrho_g h_g {\boldsymbol{v}}_g
Timo Koch's avatar
Timo Koch committed
68
69
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
109
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
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
- \lambda_g \boldsymbol{\nabla} T \right) - q_T = 0
 * \f]
 *
 * This is discretized using a fully-coupled vertex
 * centered finite volume (box) scheme as spatial and
 * the implicit Euler method as temporal discretization.
 *
 */
template<class TypeTag>
class StokesncniModel : public StokesncModel<TypeTag>
{
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;

    enum { dim = GridView::dimension };
    enum { transportCompIdx = Indices::transportCompIdx };
    enum { phaseIdx = GET_PROP_VALUE(TypeTag, PhaseIdx) };
    enum { useMoles = GET_PROP_VALUE(TypeTag, UseMoles) };
	enum { numComponents = Indices::numComponents };

    typedef typename GridView::template Codim<0>::Iterator ElementIterator;

    typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry;
    typedef typename GET_PROP_TYPE(TypeTag, ElementBoundaryTypes) ElementBoundaryTypes;
    typedef typename GET_PROP_TYPE(TypeTag, SolutionVector) SolutionVector;

    typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
    typedef typename GET_PROP_TYPE(TypeTag, VolumeVariables) VolumeVariables;

public:
    //! \copydoc ImplicitModel::addOutputVtkFields
    template <class MultiWriter>
    void addOutputVtkFields(const SolutionVector &sol,
                            MultiWriter &writer)
    {
        typedef Dune::BlockVector<Dune::FieldVector<Scalar, 1> > ScalarField;
        typedef Dune::BlockVector<Dune::FieldVector<Scalar, dim> > VelocityField;
		
        const Scalar scale_ = GET_PROP_VALUE(TypeTag, Scaling);
		
        // create the required scalar fields
        unsigned numVertices = this->gridView_().size(dim);
        ScalarField &pn = *writer.allocateManagedBuffer(numVertices);
        ScalarField &delP = *writer.allocateManagedBuffer(numVertices);
		ScalarField &T = *writer.allocateManagedBuffer(numVertices);
        ScalarField &h = *writer.allocateManagedBuffer(numVertices);
		ScalarField *moleFraction[numComponents];
		for (int i = 0; i < numComponents; ++i)
        moleFraction[i] = writer.template allocateManagedBuffer<Scalar, 1>(numVertices);
		
		ScalarField *massFraction[numComponents];
		for (int i = 0; i < numComponents; ++i)
        massFraction[i] = writer.template allocateManagedBuffer<Scalar, 1>(numVertices);
        
		ScalarField &rho = *writer.allocateManagedBuffer(numVertices);
        ScalarField &mu = *writer.allocateManagedBuffer(numVertices);
        VelocityField &velocity = *writer.template allocateManagedBuffer<Scalar, dim> (numVertices);
		
        unsigned numElements = this->gridView_().size(0);
        ScalarField &rank = *writer.allocateManagedBuffer(numElements);
		
        FVElementGeometry fvGeometry;
        VolumeVariables volVars;
        ElementBoundaryTypes elemBcTypes;
		
        ElementIterator elemIt = this->gridView_().template begin<0>();
        ElementIterator endit = this->gridView_().template end<0>();
        for (; elemIt != endit; ++elemIt)
        {
            int idx = this->elementMapper().map(*elemIt);
            rank[idx] = this->gridView_().comm().rank();
			
            fvGeometry.update(this->gridView_(), *elemIt);
            elemBcTypes.update(this->problem_(), *elemIt, fvGeometry);
			
            int numLocalVerts = elemIt->template count<dim>();
            for (int i = 0; i < numLocalVerts; ++i)
            {
                int globalIdx = this->vertexMapper().map(*elemIt, i, dim);
                volVars.update(sol[globalIdx],
                               this->problem_(),
                               *elemIt,
                               fvGeometry,
                               i,
                               false);
				
                pn[globalIdx] = volVars.pressure()*scale_;
                delP[globalIdx] = volVars.pressure()*scale_ - 1e5;
				for (int compIdx = 0; compIdx < numComponents; ++compIdx)
				{
159
160
					(*moleFraction[compIdx])[globalIdx]= volVars.moleFraction(compIdx);
					(*massFraction[compIdx])[globalIdx]= volVars.massFraction(compIdx);
Timo Koch's avatar
Timo Koch committed
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
					Valgrind::CheckDefined((*moleFraction[compIdx])[globalIdx]);
					Valgrind::CheckDefined((*massFraction[compIdx])[globalIdx]);
				}
				
				T   [globalIdx] = volVars.temperature();
                
				rho[globalIdx] = volVars.density()*scale_*scale_*scale_;
                mu[globalIdx] = volVars.dynamicViscosity()*scale_;
                h[globalIdx] = volVars.enthalpy();
                velocity[globalIdx] = volVars.velocity();
                velocity[globalIdx] *= 1/scale_;
            }
        }
		writer.attachVertexData(T, "temperature");
        writer.attachVertexData(pn, "pg");
        writer.attachVertexData(delP, "delP");
		
		
		if(useMoles)
		{
			for (int j = 0; j < numComponents; ++j)
			{
				std::ostringstream oss;
				oss << "x"
				<< FluidSystem::componentName(j)
				<< "g";
				writer.attachVertexData(*moleFraction[j], oss.str().c_str());
			}
		} else {
			for (int j = 0; j < numComponents; ++j)
			{
				std::ostringstream oss;
				oss << "X"
				<< FluidSystem::componentName(j)
				<< "g";
				writer.attachVertexData(*massFraction[j], oss.str().c_str());
			}
		}
        
        writer.attachVertexData(h, "h");
        writer.attachVertexData(rho, "rho");
        writer.attachVertexData(mu, "mu");
        writer.attachVertexData(velocity, "v", dim);

    }
};

}

#include "stokesncnipropertydefaults.hh"

#endif