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-------------------------------------------------------------------------------
11License
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"
32#include "fvPatchFieldMapper.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 int oldTag = UPstream::msgType();
117 UPstream::msgType() = oldTag+1;
118
119 // Get the coupling information from the mappedPatchBase
120 const mappedPatchBase& mpp = refCast<const mappedPatchBase>
121 (
122 patch().patch()
123 );
124 const polyMesh& nbrMesh = mpp.sampleMesh();
125 const fvPatch& nbrPatch = refCast<const fvMesh>
126 (
127 nbrMesh
128 ).boundary()[mpp.samplePolyPatch().index()];
129
131 nbrPatch.lookupPatchField<surfaceScalarField, scalar>(nbrPhiName_);
132
133 mpp.distribute(phi);
134
135
136 const surfaceScalarField& phiName =
137 db().lookupObject<surfaceScalarField>(phiName_);
138
139 scalarField U(-phi/patch().magSf());
140
141 vectorField n(patch().nf());
142
143 if (phiName.dimensions() == dimVelocity*dimArea)
144 {
145 // volumetric flow-rate
146 operator==(n*U);
147 }
148 else if (phiName.dimensions() == dimDensity*dimVelocity*dimArea)
149 {
150 const fvPatchField<scalar>& rhop =
151 patch().lookupPatchField<volScalarField, scalar>(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 // Restore tag
180 UPstream::msgType() = oldTag;
181
183}
184
185
187(
188 Ostream& os
189) const
190{
192 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
193 os.writeEntryIfDifferent<word>("rho", "rho", rhoName_);
194 os.writeEntry("nbrPhi", nbrPhiName_);
195 writeEntry("value", os);
196}
197
198
199// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200
201namespace Foam
202{
204 (
207 );
208}
209
210
211// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
surfaceScalarField & phi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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:251
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
virtual const word & name() const
Return name.
Definition: fvPatch.H:173
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Describes a volumetric/mass flow normal vector boundary condition by its magnitude as an integral ove...
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
const polyMesh & sampleMesh() const
Get the region mesh.
const polyPatch & samplePolyPatch() const
Get the patch on the region.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
label index() const noexcept
The index of this patch in the boundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
A class for handling words, derived from Foam::string.
Definition: word.H:68
U
Definition: pEqn.H:72
volScalarField & p
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
Type gSum(const FieldField< Field, Type > &f)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
error FatalError
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
Foam::surfaceFields.