turbulentBoundaryCoupledFvPatchScalarField.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  \*---------------------------------------------------------------------------*/
27 #include "fvPatchFieldMapper.H"
28 #include "volFields.H"
29 #include "mappedPatchBase.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace incompressible
36  {
37 
38  // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
42  (
43  const fvPatch& p,
45  )
46  :
47  mixedFvPatchScalarField(p, iF),
48  TnbrName_("undefined-Tnbr")
49  {
50  this->refValue() = 0.0;
51  this->refGrad() = 0.0;
52  this->valueFraction() = 1.0;
53  }
54 
55 
58  (
60  const fvPatch& p,
62  const fvPatchFieldMapper& mapper
63  )
64  :
65  mixedFvPatchScalarField(ptf, p, iF, mapper),
66  TnbrName_(ptf.TnbrName_)
67  {}
68 
69 
72  (
73  const fvPatch& p,
75  const dictionary& dict
76  )
77  :
78  mixedFvPatchScalarField(p, iF),
79  TnbrName_(dict.lookup("Tnbr"))
80  {
81  if (!isA<mappedPatchBase>(this->patch().patch()))
82  {
84  << "' not type '" << mappedPatchBase::typeName << "'"
85  << "\n for patch " << p.name()
86  << " of field " << dimensionedInternalField().name()
87  << " in file " << dimensionedInternalField().objectPath()
88  << exit(FatalError);
89  }
90 
92 
93  if (dict.found("refValue"))
94  {
95  // Full restart
96  refValue() = scalarField("refValue", dict, p.size());
97  refGrad() = scalarField("refGradient", dict, p.size());
98  valueFraction() = scalarField("valueFraction", dict, p.size());
99  }
100  else
101  {
102  // Start from user entered data. Assume fixedValue.
103  refValue() = *this;
104  refGrad() = 0.0;
105  valueFraction() = 1.0;
106  }
107  }
108 
109 
112  (
115  )
116  :
117  mixedFvPatchScalarField(wtcsf, iF),
118  TnbrName_(wtcsf.TnbrName_)
119  {}
120 
121 
122  // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
125  {
126  if (updated())
127  {
128  return;
129  }
130 
131  // Since we're inside initEvaluate/evaluate there might be processor
132  // comms underway. Change the tag we use.
133  int oldTag = UPstream::msgType();
134  UPstream::msgType() = oldTag+1;
135 
136  // Get the coupling information from the mappedPatchBase
137  const mappedPatchBase& mpp =
138  refCast<const mappedPatchBase>(patch().patch());
139  const polyMesh& nbrMesh = mpp.sampleMesh();
140  const label samplePatchI = mpp.samplePolyPatch().index();
141  const fvPatch& nbrPatch =
142  refCast<const fvMesh>(nbrMesh).boundary()[samplePatchI];
143 
144  // Calculate the temperature by harmonic averaging
145  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
146 
148  refCast
149  <
151  >
152  (
153  nbrPatch.lookupPatchField<volScalarField, scalar>
154  (
155  TnbrName_
156  )
157  );
158 
159  // Swap to obtain full local values of neighbour internal field
160  tmp<scalarField> nbrIntFld(new scalarField(nbrField.size(), 0.0));
161  tmp<scalarField> nbrKDelta(new scalarField(nbrField.size(), 0.0));
162 
163  mpp.distribute(nbrIntFld());
164  mpp.distribute(nbrKDelta());
165 
166  // Both sides agree on
167  // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
168  // - gradient : (temperature-fld)*delta
169  // We've got a degree of freedom in how to implement this in a mixed bc.
170  // (what gradient, what fixedValue and mixing coefficient)
171  // Two reasonable choices:
172  // 1. specify above temperature on one side (preferentially the high side)
173  // and above gradient on the other. So this will switch between pure
174  // fixedvalue and pure fixedgradient
175  // 2. specify gradient and temperature such that the equations are the
176  // same on both sides. This leads to the choice of
177  // - refGradient = zero gradient
178  // - refValue = neighbour value
179  // - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
180 
181  this->refValue() = nbrIntFld();
182  this->refGrad() = 0.0;
183  this->valueFraction() = nbrKDelta()/(nbrKDelta());
184 
185  //这个updateCoeffs方法并不属于这个类或者这个命名空间,这个方法是属于
186  //mixedFvPatchScalarField的。
187  //#include "mixedFvPatchFields.H"
188  //makePatchTypeFieldTypedefs(mixed);
189  //这里定义了这个 mixedFvPatchScalarField.
190  //updateCoeffs这个函数的实现是在fvPatchField.H
191  //Update the coefficients associated with the patch field
192  //Sets Updated to true
193  // virtual void updateCoeffs()
194  // {
195  // updated_ = true;
196  // }
197  // updateCoeffs()的作用仅仅就是把一布尔类型的变量变成true
198  // turbulentTemperatureCoupledBaffleMixedFvPatchScalarField,这个边界条件类
199  // 的名字就可以看出,这是一个湍流温度耦合的边界,继承自 mixedFvPatchScalarField
200  //既然是继承的,如果updateCoefs是public的,那它自己也可以重载,但是没有
201  //所以就只能用::来调用这个updateCoefs()函数了。
202  mixedFvPatchScalarField::updateCoeffs();
203 
204  if (debug)
205  {
206  Info<< patch().boundaryMesh().mesh().name() << ':'
207  << patch().name() << ':'
208  << this->dimensionedInternalField().name() << " <- "
209  << nbrMesh.name() << ':'
210  << nbrPatch.name() << ':'
211  << this->dimensionedInternalField().name() << " :"
212  << endl;
213  }
214 
215  // Restore tag
216  UPstream::msgType() = oldTag;
217  }
218 
219 
221  (
222  Ostream& os
223  ) const
224  {
226  os.writeKeyword("Tnbr")<< TnbrName_
227  << token::END_STATEMENT << nl;
228  }
229 
230 
231  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
232 
234  (
237  );
238 
239 
240  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242  } // End namespace compressible
243 } // End namespace Foam
244 
245 
246 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
turbulentBoundaryCoupledFvPatchScalarField.H
volFields.H
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::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::incompressible::turbulentBoundaryCoupledFvPatchScalarField::TnbrName_
const word TnbrName_
Name of field on the neighbour region.
Definition: turbulentBoundaryCoupledFvPatchScalarField.H:144
p
p
Definition: pEqn.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::mappedPatchBase::samplePolyPatch
const polyPatch & samplePolyPatch() const
Get the patch on the region.
Definition: mappedPatchBase.C:1246
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvPatchFieldMapper.H
Foam::mappedPatchBase
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
Definition: mappedPatchBase.H:101
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::fvPatch::name
const word & name() const
Return name.
Definition: fvPatch.H:149
Foam::mappedPatchBase::sampleMesh
const polyMesh & sampleMesh() const
Get the region mesh.
Definition: mappedPatchBase.C:1237
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::mappedPatchBase::distribute
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Definition: mappedPatchBaseTemplates.C:27
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::incompressible::turbulentBoundaryCoupledFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: turbulentBoundaryCoupledFvPatchScalarField.C:124
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::incompressible::turbulentBoundaryCoupledFvPatchScalarField::turbulentBoundaryCoupledFvPatchScalarField
turbulentBoundaryCoupledFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: turbulentBoundaryCoupledFvPatchScalarField.C:42
Foam::FatalError
error FatalError
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
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:452
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
Foam::incompressible::turbulentBoundaryCoupledFvPatchScalarField
Definition: turbulentBoundaryCoupledFvPatchScalarField.H:136
Foam::fvPatch::lookupPatchField
const GeometricField::PatchFieldType & lookupPatchField(const word &name, const GeometricField *=NULL, const Type *=NULL) const
Lookup and return the patchField of the named field from the.
Definition: fvPatchFvMeshTemplates.C:32
Foam::incompressible::turbulentBoundaryCoupledFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: turbulentBoundaryCoupledFvPatchScalarField.C:221
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::patchIdentifier::index
label index() const
Return the index of this patch in the boundaryMesh.
Definition: patchIdentifier.H:133
write
Tcoeff write()
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)
mappedPatchBase.H