mappedFlowRateFvPatchVectorField.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  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "volFields.H"
33 #include "mappedPatchBase.H"
34 #include "surfaceFields.H"
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvPatch& p,
42 )
43 :
45  nbrPhiName_("none"),
46  phiName_("phi"),
47  rhoName_("rho")
48 {}
49 
50 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
60  nbrPhiName_(ptf.nbrPhiName_),
61  phiName_(ptf.phiName_),
62  rhoName_(ptf.rhoName_)
63 {}
64 
65 
67 (
68  const fvPatch& p,
70  const dictionary& dict
71 )
72 :
74  nbrPhiName_(dict.getOrDefault<word>("nbrPhi", "phi")),
75  phiName_(dict.getOrDefault<word>("phi", "phi")),
76  rhoName_(dict.getOrDefault<word>("rho", "rho"))
77 {}
78 
79 
81 (
83 )
84 :
86  nbrPhiName_(ptf.nbrPhiName_),
87  phiName_(ptf.phiName_),
88  rhoName_(ptf.rhoName_)
89 {}
90 
91 
93 (
96 )
97 :
99  nbrPhiName_(ptf.nbrPhiName_),
100  phiName_(ptf.phiName_),
101  rhoName_(ptf.rhoName_)
102 {}
103 
104 
105 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
106 
108 {
109  if (updated())
110  {
111  return;
112  }
113 
114  // Since we're inside initEvaluate/evaluate there might be processor
115  // comms underway. Change the tag we use.
116  const int oldTag = UPstream::incrMsgType();
117 
118  // Get the coupling information from the mappedPatchBase
119  const mappedPatchBase& mpp = refCast<const mappedPatchBase>
120  (
121  patch().patch()
122  );
123  const polyMesh& nbrMesh = mpp.sampleMesh();
124  const fvPatch& nbrPatch = refCast<const fvMesh>
125  (
126  nbrMesh
127  ).boundary()[mpp.samplePolyPatch().index()];
128 
130  (
131  nbrPatch.lookupPatchField<surfaceScalarField>(nbrPhiName_)
132  );
133 
134  mpp.distribute(phi);
135 
136 
137  const surfaceScalarField& phiName =
138  db().lookupObject<surfaceScalarField>(phiName_);
139 
140  scalarField U(-phi/patch().magSf());
141 
142  vectorField n(patch().nf());
143 
144  if (phiName.dimensions() == dimVolume/dimTime)
145  {
146  // volumetric flow-rate
147  operator==(n*U);
148  }
149  else if (phiName.dimensions() == dimMass/dimTime)
150  {
151  const auto& rhop = patch().lookupPatchField<volScalarField>(rhoName_);
152 
153  // mass flow-rate
154  operator==(n*U/rhop);
155 
156  if (debug)
157  {
158  scalar phi = gSum(rhop*(*this) & patch().Sf());
159  Info<< patch().boundaryMesh().mesh().name() << ':'
160  << patch().name() << ':'
161  << this->internalField().name() << " <- "
162  << nbrMesh.name() << ':'
163  << nbrPatch.name() << ':'
164  << this->internalField().name() << " :"
165  << " mass flux[Kg/s]:" << -phi
166  << endl;
167  }
168  }
169  else
170  {
172  << "dimensions of " << phiName_ << " are incorrect" << nl
173  << " on patch " << this->patch().name()
174  << " of field " << this->internalField().name()
175  << " in file " << this->internalField().objectPath()
176  << nl << exit(FatalError);
177  }
178 
179  UPstream::msgType(oldTag); // Restore tag
180 
182 }
183 
184 
186 (
187  Ostream& os
188 ) const
189 {
191  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
192  os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
193  os.writeEntry("nbrPhi", nbrPhiName_);
195 }
196 
197 
198 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
199 
200 namespace Foam
201 {
203  (
205  mappedFlowRateFvPatchVectorField
206  );
207 }
208 
209 
210 // ************************************************************************* //
Foam::surfaceFields.
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
faceListList boundary
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
fvPatchField< vector > fvPatchVectorField
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
Describes a volumetric/mass flow normal vector boundary condition by its magnitude as an integral ove...
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition: UPstream.H:1251
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
A FieldMapper for finite-volume patch fields.
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:336
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
U
Definition: pEqn.H:72
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
mappedFlowRateFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.