pressureDirectedInletVelocityFvPatchVectorField.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 
28 #include "fvPatchFieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
37 (
38  const fvPatch& p,
40 )
41 :
42  fixedValueFvPatchVectorField(p, iF),
43  phiName_("phi"),
44  rhoName_("rho"),
45  inletDir_(p.size())
46 {}
47 
48 
51 (
53  const fvPatch& p,
55  const fvPatchFieldMapper& mapper
56 )
57 :
58  fixedValueFvPatchVectorField(ptf, p, iF, mapper),
59  phiName_(ptf.phiName_),
60  rhoName_(ptf.rhoName_),
61  inletDir_(ptf.inletDir_, mapper)
62 {}
63 
64 
67 (
68  const fvPatch& p,
70  const dictionary& dict
71 )
72 :
73  fixedValueFvPatchVectorField(p, iF),
74  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
75  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
76  inletDir_("inletDirection", dict, p.size())
77 {
78  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
79 }
80 
81 
84 (
86 )
87 :
88  fixedValueFvPatchVectorField(pivpvf),
89  phiName_(pivpvf.phiName_),
90  rhoName_(pivpvf.rhoName_),
91  inletDir_(pivpvf.inletDir_)
92 {}
93 
94 
97 (
100 )
101 :
102  fixedValueFvPatchVectorField(pivpvf, iF),
103  phiName_(pivpvf.phiName_),
104  rhoName_(pivpvf.rhoName_),
105  inletDir_(pivpvf.inletDir_)
106 {}
107 
108 
109 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
110 
112 (
113  const fvPatchFieldMapper& m
114 )
115 {
116  fixedValueFvPatchVectorField::autoMap(m);
117  inletDir_.autoMap(m);
118 }
119 
120 
122 (
123  const fvPatchVectorField& ptf,
124  const labelList& addr
125 )
126 {
127  fixedValueFvPatchVectorField::rmap(ptf, addr);
128 
130  refCast<const pressureDirectedInletVelocityFvPatchVectorField>(ptf);
131 
132  inletDir_.rmap(tiptf.inletDir_, addr);
133 }
134 
135 
137 {
138  if (updated())
139  {
140  return;
141  }
142 
143  const surfaceScalarField& phi =
144  db().lookupObject<surfaceScalarField>(phiName_);
145 
146  const fvsPatchField<scalar>& phip =
147  patch().patchField<surfaceScalarField, scalar>(phi);
148 
149  tmp<vectorField> n = patch().nf();
150  tmp<scalarField> ndmagS = (n & inletDir_)*patch().magSf();
151 
152  if (phi.dimensions() == dimVelocity*dimArea)
153  {
154  operator==(inletDir_*phip/ndmagS);
155  }
156  else if (phi.dimensions() == dimDensity*dimVelocity*dimArea)
157  {
158  const fvPatchField<scalar>& rhop =
159  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
160 
161  operator==(inletDir_*phip/(rhop*ndmagS));
162  }
163  else
164  {
166  << "dimensions of phi are not correct"
167  << "\n on patch " << this->patch().name()
168  << " of field " << this->internalField().name()
169  << " in file " << this->internalField().objectPath()
170  << exit(FatalError);
171  }
172 
173  fixedValueFvPatchVectorField::updateCoeffs();
174 }
175 
176 
178 (
179  Ostream& os
180 ) const
181 {
183  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
184  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
185  inletDir_.writeEntry("inletDirection", os);
186  writeEntry("value", os);
187 }
188 
189 
190 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
191 
192 void Foam::pressureDirectedInletVelocityFvPatchVectorField::operator=
193 (
194  const fvPatchField<vector>& pvf
195 )
196 {
197  fvPatchField<vector>::operator=(inletDir_*(inletDir_ & pvf));
198 }
199 
200 
201 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
202 
203 namespace Foam
204 {
206  (
209  );
210 }
211 
212 // ************************************************************************* //
Foam::surfaceFields.
surfaceScalarField & phi
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:719
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volVectorField vectorField(fieldObject, mesh)
pressureDirectedInletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Macros for easy insertion into run-time selection tables.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:161
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
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
const dimensionSet dimDensity
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:396
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
This velocity inlet boundary condition is applied to patches where the pressure is specified...
label n
A class for managing temporary objects.
Definition: PtrList.H:54
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
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