mappedFlowRateFvPatchVectorField.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 "volFields.H"
29 #include "fvPatchFieldMapper.H"
30 #include "mappedPatchBase.H"
31 #include "surfaceFields.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
36 (
37  const fvPatch& p,
39 )
40 :
42  nbrPhiName_("none"),
43  phiName_("phi"),
44  rhoName_("rho")
45 {}
46 
47 
49 (
50  const fvPatch& p,
52  const dictionary& dict
53 )
54 :
56  nbrPhiName_(dict.lookupOrDefault<word>("nbrPhi", "phi")),
57  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
58  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
59 {}
60 
61 
63 (
65  const fvPatch& p,
67  const fvPatchFieldMapper& mapper
68 )
69 :
70  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
71  nbrPhiName_(ptf.nbrPhiName_),
72  phiName_(ptf.phiName_),
73  rhoName_(ptf.rhoName_)
74 {}
75 
76 
78 (
81 )
82 :
84  nbrPhiName_(ptf.nbrPhiName_),
85  phiName_(ptf.phiName_),
86  rhoName_(ptf.rhoName_)
87 {}
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  if (updated())
95  {
96  return;
97  }
98 
99  // Since we're inside initEvaluate/evaluate there might be processor
100  // comms underway. Change the tag we use.
101  int oldTag = UPstream::msgType();
102  UPstream::msgType() = oldTag+1;
103 
104  // Get the coupling information from the mappedPatchBase
105  const mappedPatchBase& mpp = refCast<const mappedPatchBase>
106  (
107  patch().patch()
108  );
109  const polyMesh& nbrMesh = mpp.sampleMesh();
110  const fvPatch& nbrPatch = refCast<const fvMesh>
111  (
112  nbrMesh
113  ).boundary()[mpp.samplePolyPatch().index()];
114 
115  scalarList phi =
116  nbrPatch.lookupPatchField<surfaceScalarField, scalar>(nbrPhiName_);
117 
118  mpp.distribute(phi);
119 
120 
121  const surfaceScalarField& phiName =
122  db().lookupObject<surfaceScalarField>(phiName_);
123 
124  scalarField U(-phi/patch().magSf());
125 
126  vectorField n(patch().nf());
127 
128  if (phiName.dimensions() == dimFlux)
129  {
130  // volumetric flow-rate
131  operator==(n*U);
132  }
133  else if (phiName.dimensions() == dimMassFlux)
134  {
135  const fvPatchField<scalar>& rhop =
136  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
137 
138  // mass flow-rate
139  operator==(n*U/rhop);
140 
141  if (debug)
142  {
143  scalar phi = gSum(rhop*(*this) & patch().Sf());
144  Info<< patch().boundaryMesh().mesh().name() << ':'
145  << patch().name() << ':'
146  << this->internalField().name() << " <- "
147  << nbrMesh.name() << ':'
148  << nbrPatch.name() << ':'
149  << this->internalField().name() << " :"
150  << " mass flux[Kg/s]:" << -phi
151  << endl;
152  }
153  }
154  else
155  {
157  << "dimensions of " << phiName_ << " are incorrect" << nl
158  << " on patch " << this->patch().name()
159  << " of field " << this->internalField().name()
160  << " in file " << this->internalField().objectPath()
161  << nl << exit(FatalError);
162  }
163 
164  // Restore tag
165  UPstream::msgType() = oldTag;
166 
168 }
169 
170 
172 (
173  Ostream& os
174 ) const
175 {
177  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
178  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
179  writeEntry(os, "nbrPhi", nbrPhiName_);
180  writeEntry(os, "value", *this);
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
186 namespace Foam
187 {
189  (
192  );
193 }
194 
195 
196 // ************************************************************************* //
Foam::surfaceFields.
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Describes a volumetric/mass flow normal vector boundary condition by its magnitude as an integral ove...
U
Definition: pEqn.H:72
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
const polyMesh & sampleMesh() const
Get the region mesh.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
const dimensionSet & dimensions() const
Return dimensions.
Type gSum(const FieldField< Field, Type > &f)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
const dimensionSet dimFlux
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
phi
Definition: correctPhi.H:3
faceListList boundary(nPatches)
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
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
const dimensionSet dimMassFlux
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
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:216
label index() const
Return the index of this patch in the boundaryMesh.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
label n
mappedFlowRateFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
virtual const word & name() const
Return name.
Definition: fvPatch.H:145
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const polyPatch & samplePolyPatch() const
Get the patch on the region.
Namespace for OpenFOAM.