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 | Copyright (C) 2011-2016 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 (
51  const fvPatch& p,
53  const fvPatchFieldMapper& mapper
54 )
55 :
56  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
57  nbrPhiName_(ptf.nbrPhiName_),
58  phiName_(ptf.phiName_),
59  rhoName_(ptf.rhoName_)
60 {}
61 
62 
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
71  nbrPhiName_(dict.lookupOrDefault<word>("nbrPhi", "phi")),
72  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
73  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
74 {}
75 
76 
78 (
80 )
81 :
83  nbrPhiName_(ptf.nbrPhiName_),
84  phiName_(ptf.phiName_),
85  rhoName_(ptf.rhoName_)
86 {}
87 
88 
90 (
93 )
94 :
96  nbrPhiName_(ptf.nbrPhiName_),
97  phiName_(ptf.phiName_),
98  rhoName_(ptf.rhoName_)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 {
106  if (updated())
107  {
108  return;
109  }
110 
111  // Since we're inside initEvaluate/evaluate there might be processor
112  // comms underway. Change the tag we use.
113  int oldTag = UPstream::msgType();
114  UPstream::msgType() = oldTag+1;
115 
116  // Get the coupling information from the mappedPatchBase
117  const mappedPatchBase& mpp = refCast<const mappedPatchBase>
118  (
119  patch().patch()
120  );
121  const polyMesh& nbrMesh = mpp.sampleMesh();
122  const fvPatch& nbrPatch = refCast<const fvMesh>
123  (
124  nbrMesh
125  ).boundary()[mpp.samplePolyPatch().index()];
126 
127  scalarList phi =
128  nbrPatch.lookupPatchField<surfaceScalarField, scalar>(nbrPhiName_);
129 
130  mpp.distribute(phi);
131 
132 
133  const surfaceScalarField& phiName =
134  db().lookupObject<surfaceScalarField>(phiName_);
135 
136  scalarField U(-phi/patch().magSf());
137 
138  vectorField n(patch().nf());
139 
140  if (phiName.dimensions() == dimVelocity*dimArea)
141  {
142  // volumetric flow-rate
143  operator==(n*U);
144  }
145  else if (phiName.dimensions() == dimDensity*dimVelocity*dimArea)
146  {
147  const fvPatchField<scalar>& rhop =
148  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
149 
150  // mass flow-rate
151  operator==(n*U/rhop);
152 
153  if (debug)
154  {
155  scalar phi = gSum(rhop*(*this) & patch().Sf());
156  Info<< patch().boundaryMesh().mesh().name() << ':'
157  << patch().name() << ':'
158  << this->internalField().name() << " <- "
159  << nbrMesh.name() << ':'
160  << nbrPatch.name() << ':'
161  << this->internalField().name() << " :"
162  << " mass flux[Kg/s]:" << -phi
163  << endl;
164  }
165  }
166  else
167  {
169  << "dimensions of " << phiName_ << " are incorrect" << nl
170  << " on patch " << this->patch().name()
171  << " of field " << this->internalField().name()
172  << " in file " << this->internalField().objectPath()
173  << nl << exit(FatalError);
174  }
175 
176  // Restore tag
177  UPstream::msgType() = oldTag;
178 
180 }
181 
182 
184 (
185  Ostream& os
186 ) const
187 {
189  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
190  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
191  os.writeKeyword("nbrPhi") << nbrPhiName_ << token::END_STATEMENT << nl;
192  writeEntry("value", os);
193 }
194 
195 
196 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
197 
198 namespace Foam
199 {
201  (
204  );
205 }
206 
207 
208 // ************************************************************************* //
Foam::surfaceFields.
surfaceScalarField & phi
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
dictionary dict
U
Definition: pEqn.H:83
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...
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const word & name() const
Return name.
Definition: fvPatch.H:149
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
const polyMesh & sampleMesh() const
Get the region mesh.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:464
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
Determines a mapping between patch face centres and mesh cell or face centres and processors they&#39;re ...
faceListList boundary(nPatches)
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=NULL, const Type *=NULL) const
Lookup and return the patchField of the named field from the.
const dimensionSet & dimensions() const
Return dimensions.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
const dimensionSet dimDensity
const polyPatch & samplePolyPatch() const
Get the patch on the region.
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:312
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:74
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
label index() const
Return the index of this patch in the boundaryMesh.
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Namespace for OpenFOAM.
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:363
const dimensionSet dimVelocity