HydrothermalHeatFluxFvPatchScalarField.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8License
9 This file is part of OpenFOAM.
10
11 OpenFOAM is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by
13 the Free Software Foundation, either version 3 of the License, or
14 (at your option) any later version.
15
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19 for more details.
20
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23
24\*---------------------------------------------------------------------------*/
25
27#include "addToRunTimeSelectionTable.H"
28#include "linear.H"
29// #include "fvPatchFieldMapper.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
34(
35 const fvPatch& p,
36 const DimensionedField<scalar, volMesh>& iF
37)
38 :
39 fixedGradientFvPatchScalarField(p, iF),
40 krName_("kr"),
41 q_(p.size(), 0.0),
42 shape_("fixed"),
43 qmax_(0),
44 qmin_(0),
45 x0_(0),
46 z0_(0),
47 c_(1)
48{}
49
51(
52 const HydrothermalHeatFlux& ptf,
53 const fvPatch& p,
54 const DimensionedField<scalar, volMesh>& iF,
55 const fvPatchFieldMapper& mapper
56)
57 :
58 fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
59 krName_(ptf.krName_),
60 q_(ptf.q_), //why q_(ptf.q_, mapper) is used in turbulentHeatFluxTemperatureFvPatchScalarField.C ?
61 shape_(ptf.shape_),
62 qmax_(ptf.qmax_),
63 qmin_(ptf.qmin_),
64 x0_(ptf.x0_),
65 z0_(ptf.z0_),
66 c_(ptf.c_),
67 x1_(ptf.x1_),
68 x2_(ptf.x2_),
69 q1_(ptf.q1_),
70 q2_(ptf.q2_)
71{
72 patchType() = ptf.patchType();
73
74 // Map gradient. Set unmapped values and overwrite with mapped ptf
75 gradient() = 0.0;
76 mapper(gradient(), ptf.gradient());
77
78 // Evaluate the value field from the gradient if the internal field is valid
79 if (notNull(iF) && iF.size())
80 {
81 scalarField::operator=
82 (
83 // patchInternalField() + gradient()/patch().deltaCoeffs()
84 // ***HGW Hack to avoid the construction of mesh.deltaCoeffs
85 // which fails for AMI patches for some mapping operations
86 patchInternalField() + gradient()*(patch().nf() & patch().delta())
87 );
88 }
89 else
90 {
91 // Enforce mapping of values so we have a valid starting value. This
92 // constructor is used when reconstructing fields
93 mapper(*this, ptf);
94 }
95}
96
98(
99 const fvPatch& p,
100 const DimensionedField<scalar, volMesh>& iF,
101 const dictionary& dict
102)
103 :
104 fixedGradientFvPatchScalarField(p, iF),
105 krName_(dict.lookupOrDefault<word>("kr", "kr")),
106 q_("q", dict, p.size()),
107 shape_(dict.lookupOrDefault<word>("shape", "fixed")),
108 qmax_(dict.lookupOrDefault<scalar>("qmax", 0)),
109 qmin_(dict.lookupOrDefault<scalar>("qmin", 0)),
110 x0_(dict.lookupOrDefault<scalar>("x0", 0)),
111 z0_(dict.lookupOrDefault<scalar>("z0_", 0)),
112 c_(dict.lookupOrDefault<scalar>("c", 1)),
113 x1_(dict.lookupOrDefault<scalar>("x1", 0)),
114 x2_(dict.lookupOrDefault<scalar>("x2", 1)),
115 q1_(dict.lookupOrDefault<scalar>("q1", 0)),
116 q2_(dict.lookupOrDefault<scalar>("q2", 0))
117{
118 if (dict.found("value") && dict.found("gradient"))
119 {
120 fvPatchField<scalar>::operator=
121 (
122 scalarField("value", dict, p.size())
123 );
124 gradient() = scalarField("gradient", dict, p.size());
125 }
126 else
127 {
128 fvPatchField<scalar>::operator=(patchInternalField());
129 gradient() = 0.0;
130 }
131}
132
134(
135 const HydrothermalHeatFlux& ptf
136)
137 :
138 fixedGradientFvPatchScalarField(ptf),
139 krName_(ptf.krName_),
140 q_(ptf.q_),
141 shape_(ptf.shape_),
142 qmax_(ptf.qmax_),
143 qmin_(ptf.qmin_),
144 x0_(ptf.x0_),
145 z0_(ptf.z0_),
146 c_(ptf.c_),
147 x1_(ptf.x1_),
148 x2_(ptf.x2_),
149 q1_(ptf.q1_),
150 q2_(ptf.q2_)
151{}
152
154(
155 const HydrothermalHeatFlux& ptf,
156 const DimensionedField<scalar, volMesh>& iF
157)
158 :
159 fixedGradientFvPatchScalarField(ptf, iF),
160 krName_(ptf.krName_),
161 q_(ptf.q_),
162 shape_(ptf.shape_),
163 qmax_(ptf.qmax_),
164 qmin_(ptf.qmin_),
165 x0_(ptf.x0_),
166 z0_(ptf.z0_),
167 c_(ptf.c_),
168 x1_(ptf.x1_),
169 x2_(ptf.x2_),
170 q1_(ptf.q1_),
171 q2_(ptf.q2_)
172{}
173
174// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
175
177{
178 if (updated())
179 {
180 return;
181 }
182
183 const IOdictionary& thermophysicalProperties = db().lookupObject<IOdictionary>("thermophysicalProperties");
184 // const scalar kr(readScalar(transportProperties.lookup("kr")));
185 dimensionedScalar kr(thermophysicalProperties.subDict("mixture").subDict("porousMedia").lookup(krName_));
186 gradient() = q_/kr.value();
187 // 1. gaussian shape
188 if(shape_=="gaussian2d")
189 {
190 scalarField x(this->patch().Cf().component(0));
191 gradient()=(qmin_ + (qmax_-qmin_)*exp(-((x-x0_)*(x-x0_))/(2*c_*c_)))/kr.value();
192 }else if(shape_=="gaussian3d")
193 {
194 scalarField x(this->patch().Cf().component(0));
195 scalarField y(this->patch().Cf().component(2));
196 gradient()=(qmin_ + (qmax_-qmin_)*exp(-((x-x0_)*(x-x0_) + (y-z0_)*(y-z0_))/(2*c_*c_)))/kr.value();
197 }else if(shape_=="linear2d")
198 {
199 scalar slope = (q2_ - q1_)/(x2_ - x1_);
200 scalarField x(this->patch().Cf().component(0));
201 gradient()=((x-x1_)*slope + q1_)/kr.value();
202 }
203 fixedGradientFvPatchScalarField::updateCoeffs();
204}
205
207{
208 fixedGradientFvPatchScalarField::write(os);
209 writeEntry(os, "q", q_);
210 writeEntry(os, "value", *this);
211 writeEntryIfDifferent<word>(os, "kr","kr",krName_);
212 writeEntry(os, "shape", shape_);
213 if(shape_=="gaussian2d")
214 {
215 writeEntry(os, "x0", x0_);
216 writeEntry(os, "c", c_);
217 writeEntry(os, "qmin", qmin_);
218 writeEntry(os, "qmax", qmax_);
219 }else if(shape_=="gaussian3d")
220 {
221 writeEntry(os, "x0", x0_);
222 writeEntry(os, "z0", z0_);
223 writeEntry(os, "c", c_);
224 writeEntry(os, "qmin", qmin_);
225 writeEntry(os, "qmax", qmax_);
226 }
227}
228
229// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230
231namespace Foam
232{
234(
235 fvPatchScalarField,
237);
238}
239
240// ************************************************************************* //
volScalarField & p
Definition: createFields.H:10
dimensionedScalar kr
Definition: createFields.H:48
HydrothermalHeatFlux(const fvPatch &, const DimensionedField< scalar, volMesh > &)
makePatchTypeField(fvPatchScalarField, HydrothermalHeatFlux)