mappedFlowRateVelocityFvPatchVectorField.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-2023 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 
37 (
38  const fvPatch& p,
40  const dictionary& dict
41 )
42 :
44  nbrPhiName_(dict.lookupOrDefault<word>("nbrPhi", "phi")),
45  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
46  rhoName_(dict.lookupOrDefault<word>("rho", "rho"))
47 {
49  (
50  *this,
51  iF,
52  dict,
54  );
55 }
56 
57 
60 (
62  const fvPatch& p,
64  const fvPatchFieldMapper& mapper
65 )
66 :
67  fixedValueFvPatchField<vector>(ptf, p, iF, mapper),
68  nbrPhiName_(ptf.nbrPhiName_),
69  phiName_(ptf.phiName_),
70  rhoName_(ptf.rhoName_)
71 {}
72 
73 
76 (
79 )
80 :
82  nbrPhiName_(ptf.nbrPhiName_),
83  phiName_(ptf.phiName_),
84  rhoName_(ptf.rhoName_)
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 {
92  if (updated())
93  {
94  return;
95  }
96 
97  // Since we're inside initEvaluate/evaluate there might be processor
98  // comms underway. Change the tag we use.
99  int oldTag = UPstream::msgType();
100  UPstream::msgType() = oldTag+1;
101 
102  const mappedPatchBase& mapper = mappedPatchBase::getMap(patch().patch());
103  const fvMesh& nbrMesh = refCast<const fvMesh>(mapper.nbrMesh());
104  const label nbrPatchi = mapper.nbrPolyPatch().index();
105  const fvPatch& nbrPatch = nbrMesh.boundary()[nbrPatchi];
106 
107  const surfaceScalarField& phi =
108  db().lookupObject<surfaceScalarField>(phiName_);
109 
110  const scalarField phip
111  (
112  mapper.fromNeighbour
113  (
114  nbrPatch.lookupPatchField<surfaceScalarField, scalar>(nbrPhiName_)
115  )
116  );
117 
118  const scalarField Unp(-phip/patch().magSf());
119 
120  const vectorField np(patch().nf());
121 
122  if (phi.dimensions() == dimFlux)
123  {
124  // volumetric flow-rate
125  operator==(np*Unp);
126  }
127  else if (phi.dimensions() == dimMassFlux)
128  {
129  const fvPatchField<scalar>& rhop =
130  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
131 
132  // mass flow-rate
133  operator==(np*Unp/rhop);
134  }
135  else
136  {
138  << "dimensions of " << phiName_ << " are incorrect" << nl
139  << " on patch " << this->patch().name()
140  << " of field " << this->internalField().name()
141  << " in file " << this->internalField().objectPath()
142  << nl << exit(FatalError);
143  }
144 
145  // Restore tag
146  UPstream::msgType() = oldTag;
147 
149 }
150 
151 
153 (
154  Ostream& os
155 ) const
156 {
158  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
159  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
160  writeEntry(os, "nbrPhi", nbrPhiName_);
161  writeEntry(os, "value", *this);
162 }
163 
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 namespace Foam
168 {
170  (
173  );
174 }
175 
176 
177 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
This boundary condition supplies a fixed value constraint, and is the base class for a number of othe...
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
const GeometricField::Patch & lookupPatchField(const word &name) const
Lookup and return the patchField of the named field from the.
This boundary condition maps the flow rate from a neighbouring patch to this patch,...
mappedFlowRateVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Engine which provides mapping between two patches.
const polyPatch & nbrPolyPatch() const
Get the patch to map from.
static const mappedPatchBase & getMap(const polyPatch &patch)
Cast the given polyPatch to a mappedPatchBase. Handle errors.
const polyMesh & nbrMesh() const
Get the mesh for the region to map from.
static void validateMapForField(const PatchFieldType &field, const FieldType &iF, const dictionary &context, const label froms=from::any)
Validate that the map exists and is appropriate for the given.
tmp< Field< Type > > fromNeighbour(const Field< Type > &nbrFld) const
Map/interpolate the neighbour patch field to this patch.
label index() const
Return the index of this patch in the boundaryMesh.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const dimensionSet dimMassFlux
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:260
const dimensionSet dimFlux
dictionary dict
volScalarField & p
static const label differentPatch
Foam::surfaceFields.