mappedVelocityFluxFixedValueFvPatchField.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 
27 #include "fvPatchFieldMapper.H"
28 #include "mappedPatchBase.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
32 
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
43  fixedValueFvPatchVectorField(p, iF),
44  phiName_("phi")
45 {}
46 
47 
50 (
52  const fvPatch& p,
54  const fvPatchFieldMapper& mapper
55 )
56 :
57  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
58  phiName_(ptf.phiName_)
59 {
60  if (!isA<mappedPatchBase>(this->patch().patch()))
61  {
63  << "Patch type '" << p.type()
64  << "' not type '" << mappedPatchBase::typeName << "'"
65  << " for patch " << p.name()
66  << " of field " << dimensionedInternalField().name()
67  << " in file " << dimensionedInternalField().objectPath()
68  << exit(FatalError);
69  }
70 }
71 
72 
75 (
76  const fvPatch& p,
78  const dictionary& dict
79 )
80 :
81  fixedValueFvPatchVectorField(p, iF, dict),
82  phiName_(dict.lookupOrDefault<word>("phi", "phi"))
83 {
84  if (!isA<mappedPatchBase>(this->patch().patch()))
85  {
87  << "Patch type '" << p.type()
88  << "' not type '" << mappedPatchBase::typeName << "'"
89  << " for patch " << p.name()
90  << " of field " << dimensionedInternalField().name()
91  << " in file " << dimensionedInternalField().objectPath()
92  << exit(FatalError);
93  }
94 
95  const mappedPatchBase& mpp = refCast<const mappedPatchBase>
96  (
97  this->patch().patch()
98  );
99  if (mpp.mode() == mappedPolyPatch::NEARESTCELL)
100  {
102  << "Patch " << p.name()
103  << " of type '" << p.type()
104  << "' can not be used in 'nearestCell' mode"
105  << " of field " << dimensionedInternalField().name()
106  << " in file " << dimensionedInternalField().objectPath()
107  << exit(FatalError);
108  }
109 }
110 
111 
114 (
116 )
117 :
118  fixedValueFvPatchVectorField(ptf),
119  phiName_(ptf.phiName_)
120 {}
121 
122 
125 (
128 )
129 :
130  fixedValueFvPatchVectorField(ptf, iF),
131  phiName_(ptf.phiName_)
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
138 {
139  if (updated())
140  {
141  return;
142  }
143 
144  // Since we're inside initEvaluate/evaluate there might be processor
145  // comms underway. Change the tag we use.
146  int oldTag = UPstream::msgType();
147  UPstream::msgType() = oldTag+1;
148 
149  // Get the mappedPatchBase
150  const mappedPatchBase& mpp = refCast<const mappedPatchBase>
151  (
152  mappedVelocityFluxFixedValueFvPatchField::patch().patch()
153  );
154  const fvMesh& nbrMesh = refCast<const fvMesh>(mpp.sampleMesh());
155  const word& fieldName = dimensionedInternalField().name();
156  const volVectorField& UField =
157  nbrMesh.lookupObject<volVectorField>(fieldName);
158 
159  surfaceScalarField& phiField = const_cast<surfaceScalarField&>
160  (
162  );
163 
164  vectorField newUValues;
165  scalarField newPhiValues;
166 
167  switch (mpp.mode())
168  {
170  {
171  vectorField allUValues(nbrMesh.nFaces(), vector::zero);
172  scalarField allPhiValues(nbrMesh.nFaces(), 0.0);
173 
174  forAll(UField.boundaryField(), patchI)
175  {
176  const fvPatchVectorField& Upf = UField.boundaryField()[patchI];
177  const scalarField& phipf = phiField.boundaryField()[patchI];
178 
179  label faceStart = Upf.patch().start();
180 
181  forAll(Upf, faceI)
182  {
183  allUValues[faceStart + faceI] = Upf[faceI];
184  allPhiValues[faceStart + faceI] = phipf[faceI];
185  }
186  }
187 
188  mpp.distribute(allUValues);
189  newUValues.transfer(allUValues);
190 
191  mpp.distribute(allPhiValues);
192  newPhiValues.transfer(allPhiValues);
193 
194  break;
195  }
198  {
199  const label nbrPatchID =
200  nbrMesh.boundaryMesh().findPatchID(mpp.samplePatch());
201 
202  newUValues = UField.boundaryField()[nbrPatchID];
203  mpp.distribute(newUValues);
204 
205  newPhiValues = phiField.boundaryField()[nbrPatchID];
206  mpp.distribute(newPhiValues);
207 
208  break;
209  }
210  default:
211  {
213  << "patch can only be used in NEARESTPATCHFACE, "
214  << "NEARESTPATCHFACEAMI or NEARESTFACE mode" << nl
215  << abort(FatalError);
216  }
217  }
218 
219  operator==(newUValues);
220  phiField.boundaryField()[patch().index()] == newPhiValues;
221 
222  // Restore tag
223  UPstream::msgType() = oldTag;
224 
225  fixedValueFvPatchVectorField::updateCoeffs();
226 }
227 
228 
230 (
231  Ostream& os
232 ) const
233 {
235  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
236  this->writeEntry("value", os);
237 }
238 
239 
240 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
241 
242 namespace Foam
243 {
245  (
248  );
249 }
250 
251 
252 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
mappedVelocityFluxFixedValueFvPatchField.H
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: polyBoundaryMesh.C:678
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::mappedPatchBase::NEARESTPATCHFACEAMI
@ NEARESTPATCHFACEAMI
Definition: mappedPatchBase.H:113
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::mappedPatchBase::NEARESTFACE
@ NEARESTFACE
Definition: mappedPatchBase.H:115
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::operator==
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
surfaceFields.H
Foam::surfaceFields.
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::mappedVelocityFluxFixedValueFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: mappedVelocityFluxFixedValueFvPatchField.C:137
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::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::mappedVelocityFluxFixedValueFvPatchField
This boundary condition maps the velocity and flux from a neighbour patch to this patch.
Definition: mappedVelocityFluxFixedValueFvPatchField.H:93
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::mappedVelocityFluxFixedValueFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: mappedVelocityFluxFixedValueFvPatchField.C:230
Foam::mappedPatchBase::mode
const sampleMode & mode() const
What to sample.
Definition: mappedPatchBaseI.H:27
dict
dictionary dict
Definition: searchingEngine.H:14
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::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::mappedVelocityFluxFixedValueFvPatchField::phiName_
word phiName_
Name of flux field.
Definition: mappedVelocityFluxFixedValueFvPatchField.H:100
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::mappedVelocityFluxFixedValueFvPatchField::mappedVelocityFluxFixedValueFvPatchField
mappedVelocityFluxFixedValueFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: mappedVelocityFluxFixedValueFvPatchField.C:38
Foam::mappedPatchBase::NEARESTPATCHFACE
@ NEARESTPATCHFACE
Definition: mappedPatchBase.H:112
Foam::fvPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
mappedPatchBase.H
Foam::mappedPatchBase::samplePatch
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE)
Definition: mappedPatchBaseI.H:59