mappedVelocityFluxFixedValueFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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 (
51  const fvPatch& p,
53  const dictionary& dict
54 )
55 :
56  fixedValueFvPatchVectorField(p, iF, dict),
57  phiName_(dict.lookupOrDefault<word>("phi", "phi"))
58 {
59  if (!isA<mappedPatchBase>(this->patch().patch()))
60  {
62  << "Patch type '" << p.type()
63  << "' not type '" << mappedPatchBase::typeName << "'"
64  << " for patch " << p.name()
65  << " of field " << internalField().name()
66  << " in file " << internalField().objectPath()
67  << exit(FatalError);
68  }
69 
70  const mappedPatchBase& mpp = refCast<const mappedPatchBase>
71  (
72  this->patch().patch()
73  );
74  if (mpp.mode() == mappedPolyPatch::NEARESTCELL)
75  {
77  << "Patch " << p.name()
78  << " of type '" << p.type()
79  << "' can not be used in 'nearestCell' mode"
80  << " of field " << internalField().name()
81  << " in file " << internalField().objectPath()
82  << exit(FatalError);
83  }
84 }
85 
86 
89 (
91  const fvPatch& p,
93  const fvPatchFieldMapper& mapper
94 )
95 :
96  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
97  phiName_(ptf.phiName_)
98 {
99  if (!isA<mappedPatchBase>(this->patch().patch()))
100  {
102  << "Patch type '" << p.type()
103  << "' not type '" << mappedPatchBase::typeName << "'"
104  << " for patch " << p.name()
105  << " of field " << internalField().name()
106  << " in file " << internalField().objectPath()
107  << exit(FatalError);
108  }
109 }
110 
111 
114 (
117 )
118 :
119  fixedValueFvPatchVectorField(ptf, iF),
120  phiName_(ptf.phiName_)
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
127 {
128  if (updated())
129  {
130  return;
131  }
132 
133  // Since we're inside initEvaluate/evaluate there might be processor
134  // comms underway. Change the tag we use.
135  int oldTag = UPstream::msgType();
136  UPstream::msgType() = oldTag+1;
137 
138  // Get the mappedPatchBase
139  const mappedPatchBase& mpp = refCast<const mappedPatchBase>
140  (
141  mappedVelocityFluxFixedValueFvPatchField::patch().patch()
142  );
143  const fvMesh& nbrMesh = refCast<const fvMesh>(mpp.sampleMesh());
144  const word& fieldName = internalField().name();
145  const volVectorField& UField =
146  nbrMesh.lookupObject<volVectorField>(fieldName);
147 
148  surfaceScalarField& phiField =
149  nbrMesh.lookupObjectRef<surfaceScalarField>(phiName_);
150 
151  vectorField newUValues;
152  scalarField newPhiValues;
153 
154  switch (mpp.mode())
155  {
157  {
158  vectorField allUValues(nbrMesh.nFaces(), Zero);
159  scalarField allPhiValues(nbrMesh.nFaces(), 0.0);
160 
161  forAll(UField.boundaryField(), patchi)
162  {
163  const fvPatchVectorField& Upf = UField.boundaryField()[patchi];
164  const scalarField& phipf = phiField.boundaryField()[patchi];
165 
166  label faceStart = Upf.patch().start();
167 
168  forAll(Upf, facei)
169  {
170  allUValues[faceStart + facei] = Upf[facei];
171  allPhiValues[faceStart + facei] = phipf[facei];
172  }
173  }
174 
175  mpp.distribute(allUValues);
176  newUValues.transfer(allUValues);
177 
178  mpp.distribute(allPhiValues);
179  newPhiValues.transfer(allPhiValues);
180 
181  break;
182  }
185  {
186  const label nbrPatchID =
187  nbrMesh.boundaryMesh().findPatchID(mpp.samplePatch());
188 
189  newUValues = UField.boundaryField()[nbrPatchID];
190  mpp.distribute(newUValues);
191 
192  newPhiValues = phiField.boundaryField()[nbrPatchID];
193  mpp.distribute(newPhiValues);
194 
195  break;
196  }
197  default:
198  {
200  << "patch can only be used in NEARESTPATCHFACE, "
201  << "NEARESTPATCHFACEAMI or NEARESTFACE mode" << nl
202  << abort(FatalError);
203  }
204  }
205 
206  operator==(newUValues);
207  phiField.boundaryFieldRef()[patch().index()] == newPhiValues;
208 
209  // Restore tag
210  UPstream::msgType() = oldTag;
211 
212  fixedValueFvPatchVectorField::updateCoeffs();
213 }
214 
215 
217 (
218  Ostream& os
219 ) const
220 {
222  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
223  writeEntry(os, "value", *this);
224 }
225 
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 namespace Foam
230 {
232  (
235  );
236 }
237 
238 
239 // ************************************************************************* //
Foam::surfaceFields.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:151
const polyMesh & sampleMesh() const
Get the region mesh.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
mappedVelocityFluxFixedValueFvPatchField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
This boundary condition maps the velocity and flux from a neighbour patch to this patch...
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
Macros for easy insertion into run-time selection tables.
const sampleMode & mode() const
What to sample.
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
static const zero Zero
Definition: zero.H:97
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:351
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const word & samplePatch() const
Patch (only if NEARESTPATCHFACE)
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
Namespace for OpenFOAM.