turbulentTemperatureRadFixedGradientFvPatchScalarField.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-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
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 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 namespace incompressible
37 {
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
42 (
43  const fvPatch& p,
45 )
46 :
47  fixedGradientFvPatchScalarField(p, iF),
48  //QrName_("undefined"),
49  Qs(p.size(),0.0),
50  Qv(p.size(),0.0)
51 {
52 }
53 
54 
57 (
59  const fvPatch& p,
61  const fvPatchFieldMapper& mapper
62 )
63 :
64  fixedGradientFvPatchScalarField(psf, p, iF, mapper),
65  //QrName_("undefined"),
66  Qs(psf.Qs,mapper),
67  Qv(psf.Qv,mapper)
68 {}
69 
70 
73 (
74  const fvPatch& p,
76  const dictionary& dict
77 )
78 :
79  fixedGradientFvPatchScalarField(p, iF),
80  //QrName_(dict.lookupOrDefault<word>("Qr", "Qr")),
81  //Qs(dict.lookupOrDefault<scalarField>("Qs",scalarField(p.size(),0.0))),
82  //Qv(dict.lookupOrDefault<scalarField>("Qv",scalarField(p.size(),0.0)))
83  Qs("Qs",dict,p.size()),
84  Qv("Qv",dict,p.size())
85 {
86 
87 
88  if (debug)
89  {
90  Info << "initialeze the heat boundary conditions"<<endl;
91  }
92 
93  if (dict.found("value"))
94  {
95 
96  if (debug)
97  {
98  Info << "find value" <<endl;
99  }
100 
101  //操作符 = 来给 value 赋值
102 
104  (
105  scalarField("value", dict, p.size())
106  );
107  }
108  else
109  {
110  if (debug)
111  {
112  Info << "value can not found,"<<endl;
113  }
114 
116  (
117  0
118  );
119  }
120 
121  if (dict.found("gradient"))
122  {
123  if (debug)
124  {
125  Info << "find gradient" <<endl;
126  }
127 
128  gradient() = scalarField("gradient", dict, p.size());
129  fixedGradientFvPatchScalarField::updateCoeffs();
130  fixedGradientFvPatchScalarField::evaluate();
131  }
132  else
133  {
134 
135  if (debug)
136  {
137  Info << "gardient can not found,"<< endl;
138  Info << "the value is equal to the internal value near cell" << endl;
139  Info << "gradient default is 0" << endl;
140  }
141  fvPatchField<scalar>::operator=(patchInternalField());
142  gradient() = 0.0;
143 
144  }
145 
146 
147  fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
148 
149 }
150 
151 
154 (
157 )
158 :
159 
160  fixedGradientFvPatchScalarField(psf, iF),
161  //QrName_("undefined"),
162  Qs(psf.Qs),
163  Qv(psf.Qv)
164 {}
165 
166 
167 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
168 
170 {
171  // Since we're inside initEvaluate/evaluate there might be processor
172  // comms underway. Change the tag we use.
173  int oldTag = UPstream::msgType();
174  UPstream::msgType() = oldTag+1;
175 
176  if (debug)
177  {
178  Info << "update the boundary patch field Coeffs now ... " << endl;
179  }
180 
181  if (updated())
182  {
183  return;
184  }
185 
186  const dictionary& transportProperties
187  = db().lookupObject<IOdictionary>("transportProperties");
188 
189  dimensionedScalar Cp0(transportProperties.lookup("Cp"));
190 
191  dimensionedScalar rho0(transportProperties.lookup("rho"));
192 
193  const fvPatchField<scalar>& alphaEffW =
194  patch().lookupPatchField<volScalarField, scalar>("alphaEff");
195 
196  const fvPatchField<scalar>& QtotalFromRadiation =
197  patch().lookupPatchField<volScalarField, scalar>("QtotalRad_");
198 
199 // const fvPatchField<scalar>& alphaEffW =
200 // patch().lookupPatchField<volScalarField, scalar>("kappaEff");
201 
202 // Info << "cp= "<< Cp0 << endl;
203 // Info << "rho0= "<< rho0 << endl;
204  //Info << "alphaEffw=" << alphaEffW << endl;
205 // Info << "heatflux= " << heatFlux_ << endl;
206 
207  const scalar Ap = gSum(patch().magSf());
208  // gradient() = q_/(Ap*rhoCp0*alphaEffp);
209  // Info << "Area ="<< Ap;
210  // Info << "Cp ="<< Cp0.value();
211  // Info << "rho ="<< rho0.value();
212  // Info << "alpha ="<< alphaEffW;
213  //gradient() = heatFlux_ / (alphaEffW * Cp0.value() * rho0.value() );
214 
215 
216  //scalarField& p = *this;
217  // scalarField Qr(p.size(), 0.0);
218  // if (QrName_ != "none" || QrName_ != "undefined")
219  // {
220  // Info << "get the " << QrName_ << " value now....." << endl;
221 
222  // Qr = patch().lookupPatchField<volScalarField, scalar>(QrName_);
223  // }
224  //gradient() = (QtotalFromRadiation*Ap+Qs*Ap+Qv) / (Ap*alphaEffW * Cp0.value() * rho0.value() ); //热辐射W+单位面积传热量W/m2*面积+直接发热量W
225  gradient() = (Qs*Ap+Qv) / (Ap*alphaEffW * Cp0.value() * rho0.value() ); //热辐射W+单位面积传热量W/m2*面积+直接发热量W
226  //if (debug)
227  {
228  Info << "input radiation heat flux is : " << QtotalFromRadiation << endl;
229  Info << "input Qs heat flux is : " << Qs << "W/m2" << endl;
230  Info << "input Qv heat flux is : " << Qv << "W" << endl;
231  Info << "input Ap is : " << Ap << "m2" <<endl;
232  Info << "input alphaEffw is : " << alphaEffW << "m2/s"<<endl;
233 // Info << "input alpha is : " << alpha << "m2/s"<<endl;
234 // Info << "input alphat is : " << alpht << "m2/s"<<endl;
235  Info << "input Cp0 is : " << Cp0 << endl;
236  Info << "input rho0 is : " << rho0 << endl;
237  Info << "the gradient of temperature is : " << gradient() << endl;
238  }
239 
240  fixedGradientFvPatchScalarField::updateCoeffs();
241 //Info << "update coeffs" << endl;
242 
243 
244 
245 
246 
247 
248  if (updated())
249  {
250  return;
251  }
252 
253 
254  // Restore tag
255  UPstream::msgType() = oldTag;
256 }
257 
258 
260 (
261  Ostream& os
262 ) const
263 {
264  // os.writeKeyword("QrName")<< QrName_ << token::END_STATEMENT << nl;
266  Qs.writeEntry("Qs", os);
267  Qv.writeEntry("Qv", os);
268  gradient().writeEntry("gradient", os);
269  this->writeEntry("value", os);
270 }
271 
272 
273 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274 
276 (
279 );
280 
281 
282 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
283 } //End incompressible
284 } // End namespace Foam
285 
286 
Foam::fvPatchField< scalar >
volFields.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
Foam::fvPatchField::operator=
virtual void operator=(const UList< Type > &)
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::incompressible::turbulentTemperatureRadFixedGradientFvPatchScalarField::Qv
scalarField Qv
Name of Heat flux for W.
Definition: turbulentTemperatureRadFixedGradientFvPatchScalarField.H:152
p
p
Definition: pEqn.H:62
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::fvPatchField::operator
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvPatchFieldMapper.H
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::incompressible::turbulentTemperatureRadFixedGradientFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: turbulentTemperatureRadFixedGradientFvPatchScalarField.C:169
Foam::Info
messageStream Info
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::incompressible::turbulentTemperatureRadFixedGradientFvPatchScalarField
Definition: turbulentTemperatureRadFixedGradientFvPatchScalarField.H:137
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::incompressible::turbulentTemperatureRadFixedGradientFvPatchScalarField::Qs
scalarField Qs
Name of field for radiation heatFlux.
Definition: turbulentTemperatureRadFixedGradientFvPatchScalarField.H:148
turbulentTemperatureRadFixedGradientFvPatchScalarField.H
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
rho0
scalar rho0
Definition: readInitialConditions.H:97
Foam::incompressible::turbulentTemperatureRadFixedGradientFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: turbulentTemperatureRadFixedGradientFvPatchScalarField.C:260
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:452
Foam::incompressible::turbulentTemperatureRadFixedGradientFvPatchScalarField::turbulentTemperatureRadFixedGradientFvPatchScalarField
turbulentTemperatureRadFixedGradientFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: turbulentTemperatureRadFixedGradientFvPatchScalarField.C:42
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::incompressible::makePatchTypeField
makePatchTypeField(fvPatchScalarField, turbulentBoundaryCoupledFvPatchScalarField)