wallHeatFluxFvPatchScalarField.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) 1991-2008 OpenCFD Ltd.
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 the
13  Free Software Foundation; either version 2 of the License, or (at your
14  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, write to the Free Software Foundation,
23  Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
24 
25 \*---------------------------------------------------------------------------*/
26 
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 (
41  const fvPatch& p,
43 )
44 :
45  fixedGradientFvPatchScalarField(p, iF),
46  heatFlux_(p.size(), 0.0)
47 {
48 }
49 
50 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
60  heatFlux_(ptf.heatFlux_, mapper)
61 {}
62 
63 
65 (
66  const fvPatch& p,
68  const dictionary& dict
69 )
70 :
71  fixedGradientFvPatchScalarField(p, iF),
72  heatFlux_("heatFlux", dict, p.size())
73 {
74 
75  if (debug)
76  {
77  Info << "initialeze the heatFlux"<<endl;
78  }
79 
80  if (dict.found("value"))
81  {
82 
83  if (debug)
84  {
85  Info << "find value" <<endl;
86  }
88  (
89  scalarField("value", dict, p.size())
90  );
91  }
92  else
93  {
94  if (debug)
95  {
96  Info << "value can not found,"<<endl;
97  }
99  (
100  0
101  );
102  }
103 
104  if (dict.found("gradient"))
105  {
106  if (debug)
107  {
108  Info << "find gradient" <<endl;
109  }
110 
111  gradient() = scalarField("gradient", dict, p.size());
112  fixedGradientFvPatchScalarField::updateCoeffs();
113  fixedGradientFvPatchScalarField::evaluate();
114  }
115  else
116  {
117 
118  if (debug)
119  {
120  Info << "gardient can not found,"<< endl;
121  Info << "the value is equal to the internal value near cell" << endl;
122  Info << "gradient default is 0" << endl;
123  }
124  fvPatchField<scalar>::operator=(patchInternalField());
125  gradient() = 0.0;
126 
127  }
128 }
129 
131 (
132  const wallHeatFluxFvPatchScalarField& wbppsf
133 )
134 :
135  fixedGradientFvPatchScalarField(wbppsf),
136  heatFlux_(wbppsf.heatFlux_)
137 {}
138 
139 
141 (
142  const wallHeatFluxFvPatchScalarField& wbppsf,
144 )
145 :
146  fixedGradientFvPatchScalarField(wbppsf, iF),
147  heatFlux_(wbppsf.heatFlux_)
148 {
149 }
150 
151 
152 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153 
154 // Update the coefficients associated with the patch field
156 {
157  if (debug)
158  {
159  Info << "update the boundary patch field Coeffs now ... " << endl;
160  }
161 
162  if (updated())
163  {
164  return;
165  }
166 
167  const dictionary& transportProperties
168  = db().lookupObject<IOdictionary>("transportProperties");
169 
170  dimensionedScalar Cp0(transportProperties.lookup("Cp"));
171 
172  dimensionedScalar rho0(transportProperties.lookup("rho"));
173 
174  const fvPatchField<scalar>& alphaEffW =
175  patch().lookupPatchField<volScalarField, scalar>("alphaEff");
176 
177 // const fvPatchField<scalar>& alphaEffW =
178 // patch().lookupPatchField<volScalarField, scalar>("kappaEff");
179 
180 // Info << "cp= "<< Cp0 << endl;
181 // Info << "rho0= "<< rho0 << endl;
182  //Info << "alphaEffw=" << alphaEffW << endl;
183 // Info << "heatflux= " << heatFlux_ << endl;
184 
185  const scalar Ap = gSum(patch().magSf());
186  // gradient() = q_/(Ap*rhoCp0*alphaEffp);
187  // Info << "Area ="<< Ap;
188  // Info << "Cp ="<< Cp0.value();
189  // Info << "rho ="<< rho0.value();
190  // Info << "alpha ="<< alphaEffW;
191  //gradient() = heatFlux_ / (alphaEffW * Cp0.value() * rho0.value() );
192  gradient() = heatFlux_ / (Ap*alphaEffW * Cp0.value() * rho0.value() );
193 // if (debug)
194  {
195  Info << "input heat flux is : " << heatFlux_ << endl;
196  Info << "input Ap is : " << Ap << endl;
197  Info << "input alphaEffw is : " << alphaEffW << endl;
198  Info << "input Cp0 is : " << Cp0 << endl;
199  Info << "input rho0 is : " << rho0 << endl;
200  Info << "the gradient of temperature is : " << gradient() << endl;
201  }
202 
203  fixedGradientFvPatchScalarField::updateCoeffs();
204 //Info << "update coeffs" << endl;
205 }
206 
207 // Write
209 {
211  heatFlux_.writeEntry("heatFlux", os);
212  gradient().writeEntry("gradient", os);
213  this->writeEntry("value", os);
214 }
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
219 
220 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
221 
222 } // End namespace Foam
223 
224 // ************************************************************************* //
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::wallHeatFluxFvPatchScalarField
Definition: wallHeatFluxFvPatchScalarField.H:48
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::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::wallHeatFluxFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: wallHeatFluxFvPatchScalarField.C:207
Foam::Info
messageStream Info
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
wallHeatFluxFvPatchScalarField.H
Foam::wallHeatFluxFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: wallHeatFluxFvPatchScalarField.C:154
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
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
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::Field::writeEntry
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:700
rho0
scalar rho0
Definition: readInitialConditions.H:97
Foam::wallHeatFluxFvPatchScalarField::wallHeatFluxFvPatchScalarField
wallHeatFluxFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: wallHeatFluxFvPatchScalarField.C:39
Foam::wallHeatFluxFvPatchScalarField::heatFlux_
scalarField heatFlux_
Definition: wallHeatFluxFvPatchScalarField.H:55
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