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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "fvPatchFieldMapper.H"
30 #include "mappedPatchBase.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
34 
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
40 (
41  const fvPatch& p,
43 )
44 :
45  fixedValueFvPatchVectorField(p, iF),
46  phiName_("phi")
47 {}
48 
49 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
60  phiName_(ptf.phiName_)
61 {
62  if (!isA<mappedPatchBase>(this->patch().patch()))
63  {
65  << "Patch type '" << p.type()
66  << "' not type '" << mappedPatchBase::typeName << "'"
67  << " for patch " << p.name()
68  << " of field " << internalField().name()
69  << " in file " << internalField().objectPath()
70  << exit(FatalError);
71  }
72 }
73 
74 
77 (
78  const fvPatch& p,
80  const dictionary& dict
81 )
82 :
83  fixedValueFvPatchVectorField(p, iF, dict),
84  phiName_(dict.lookupOrDefault<word>("phi", "phi"))
85 {
86  if (!isA<mappedPatchBase>(this->patch().patch()))
87  {
89  << "Patch type '" << p.type()
90  << "' not type '" << mappedPatchBase::typeName << "'"
91  << " for patch " << p.name()
92  << " of field " << internalField().name()
93  << " in file " << internalField().objectPath()
94  << exit(FatalError);
95  }
96 
97  const mappedPatchBase& mpp = refCast<const mappedPatchBase>
98  (
99  this->patch().patch(),
100  dict
101  );
102  if (mpp.mode() == mappedPolyPatch::NEARESTCELL)
103  {
105  << "Patch " << p.name()
106  << " of type '" << p.type()
107  << "' can not be used in 'nearestCell' mode"
108  << " of field " << internalField().name()
109  << " in file " << internalField().objectPath()
110  << exit(FatalError);
111  }
112 }
113 
114 
117 (
119 )
120 :
121  fixedValueFvPatchVectorField(ptf),
122  phiName_(ptf.phiName_)
123 {}
124 
125 
128 (
131 )
132 :
133  fixedValueFvPatchVectorField(ptf, iF),
134  phiName_(ptf.phiName_)
135 {}
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
141 {
142  if (updated())
143  {
144  return;
145  }
146 
147  // Since we're inside initEvaluate/evaluate there might be processor
148  // comms underway. Change the tag we use.
149  int oldTag = UPstream::msgType();
150  UPstream::msgType() = oldTag+1;
151 
152  // Get the mappedPatchBase
153  const mappedPatchBase& mpp = refCast<const mappedPatchBase>
154  (
156  );
157  const fvMesh& nbrMesh = refCast<const fvMesh>(mpp.sampleMesh());
158  const word& fieldName = internalField().name();
159  const volVectorField& UField =
160  nbrMesh.lookupObject<volVectorField>(fieldName);
161 
162  const surfaceScalarField& phiField =
163  nbrMesh.lookupObject<surfaceScalarField>(phiName_);
164 
165  vectorField newUValues;
166  scalarField newPhiValues;
167 
168  switch (mpp.mode())
169  {
171  {
172  vectorField allUValues(nbrMesh.nFaces(), Zero);
173  scalarField allPhiValues(nbrMesh.nFaces(), Zero);
174 
175  forAll(UField.boundaryField(), patchi)
176  {
177  const fvPatchVectorField& Upf = UField.boundaryField()[patchi];
178  const scalarField& phipf = phiField.boundaryField()[patchi];
179 
180  label faceStart = Upf.patch().start();
181 
182  forAll(Upf, facei)
183  {
184  allUValues[faceStart + facei] = Upf[facei];
185  allPhiValues[faceStart + facei] = phipf[facei];
186  }
187  }
188 
189  mpp.distribute(allUValues);
190  newUValues.transfer(allUValues);
191 
192  mpp.distribute(allPhiValues);
193  newPhiValues.transfer(allPhiValues);
194 
195  break;
196  }
199  {
200  const label nbrPatchID =
201  nbrMesh.boundaryMesh().findPatchID(mpp.samplePatch());
202 
203  newUValues = UField.boundaryField()[nbrPatchID];
204  mpp.distribute(newUValues);
205 
206  newPhiValues = phiField.boundaryField()[nbrPatchID];
207  mpp.distribute(newPhiValues);
208 
209  break;
210  }
211  default:
212  {
214  << "patch can only be used in NEARESTPATCHFACE, "
215  << "NEARESTPATCHFACEAMI or NEARESTFACE mode" << nl
216  << abort(FatalError);
217  }
218  }
219 
220  operator==(newUValues);
221  const_cast<surfaceScalarField&>
222  (
223  phiField
224  ).boundaryFieldRef()[patch().index()] == newPhiValues;
225 
226  // Restore tag
227  UPstream::msgType() = oldTag;
228 
229  fixedValueFvPatchVectorField::updateCoeffs();
230 }
231 
232 
234 (
235  Ostream& os
236 ) const
237 {
239  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
240  this->writeEntry("value", os);
241 }
242 
243 
244 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
245 
246 namespace Foam
247 {
249  (
252  );
253 }
254 
255 
256 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:364
Foam::Ostream::writeEntryIfDifferent
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:231
p
volScalarField & p
Definition: createFieldRefs.H:8
mappedVelocityFluxFixedValueFvPatchField.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::mappedPatchBase::NEARESTPATCHFACEAMI
nearest patch face + AMI interpolation
Definition: mappedPatchBase.H:117
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::fvPatch::start
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:157
Foam::mappedPatchBase::NEARESTFACE
nearest face
Definition: mappedPatchBase.H:119
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
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:105
Foam::mappedVelocityFluxFixedValueFvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: mappedVelocityFluxFixedValueFvPatchField.C:140
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::mappedPatchBase::sampleMesh
const polyMesh & sampleMesh() const
Get the region mesh.
Definition: mappedPatchBase.C:1210
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:62
Foam::Field< vector >
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::mappedVelocityFluxFixedValueFvPatchField
This boundary condition maps the velocity and flux from a neighbour patch to this patch.
Definition: mappedVelocityFluxFixedValueFvPatchField.H:94
Foam::mappedPatchBase::distribute
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Definition: mappedPatchBaseTemplates.C:29
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Foam::mappedVelocityFluxFixedValueFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: mappedVelocityFluxFixedValueFvPatchField.C:234
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Definition: polyBoundaryMesh.C:766
Foam::mappedPatchBase::mode
const sampleMode & mode() const
What to sample.
Definition: mappedPatchBaseI.H:29
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
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:84
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::UPstream::msgType
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:491
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::mappedVelocityFluxFixedValueFvPatchField::mappedVelocityFluxFixedValueFvPatchField
mappedVelocityFluxFixedValueFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: mappedVelocityFluxFixedValueFvPatchField.C:40
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::mappedPatchBase::NEARESTPATCHFACE
nearest face on selected patch
Definition: mappedPatchBase.H:116
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::fvPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:343
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
mappedPatchBase.H
Foam::mappedPatchBase::samplePatch
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE)
Definition: mappedPatchBaseI.H:61