pressureInletOutletVelocityFvPatchVectorField.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-2018 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 "volFields.H"
29 #include "surfaceFields.H"
30 
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39 )
40 :
41  directionMixedFvPatchVectorField(p, iF),
42  phiName_("phi")
43 {
44  refValue() = Zero;
45  refGrad() = Zero;
46  valueFraction() = Zero;
47 }
48 
49 
52 (
54  const fvPatch& p,
56  const fvPatchFieldMapper& mapper
57 )
58 :
59  directionMixedFvPatchVectorField(ptf, p, iF, mapper),
60  phiName_(ptf.phiName_)
61 {
62  if (ptf.tangentialVelocity_.size())
63  {
64  tangentialVelocity_ = mapper(ptf.tangentialVelocity_);
65  }
66 }
67 
68 
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
77  directionMixedFvPatchVectorField(p, iF),
78  phiName_(dict.lookupOrDefault<word>("phi", "phi"))
79 {
80  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
81 
82  if (dict.found("tangentialVelocity"))
83  {
84  setTangentialVelocity
85  (
86  vectorField("tangentialVelocity", dict, p.size())
87  );
88  }
89  else
90  {
91  refValue() = Zero;
92  }
93 
94  refGrad() = Zero;
95  valueFraction() = Zero;
96 }
97 
98 
101 (
103 )
104 :
105  directionMixedFvPatchVectorField(pivpvf),
106  phiName_(pivpvf.phiName_),
107  tangentialVelocity_(pivpvf.tangentialVelocity_)
108 {}
109 
110 
113 (
116 )
117 :
118  directionMixedFvPatchVectorField(pivpvf, iF),
119  phiName_(pivpvf.phiName_),
120  tangentialVelocity_(pivpvf.tangentialVelocity_)
121 {}
122 
123 
124 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
125 
127 setTangentialVelocity(const vectorField& tangentialVelocity)
128 {
129  tangentialVelocity_ = tangentialVelocity;
130  const vectorField n(patch().nf());
131  refValue() = tangentialVelocity_ - n*(n & tangentialVelocity_);
132 }
133 
134 
136 (
137  const fvPatchFieldMapper& m
138 )
139 {
140  directionMixedFvPatchVectorField::autoMap(m);
141  if (tangentialVelocity_.size())
142  {
143  tangentialVelocity_.autoMap(m);
144  }
145 }
146 
147 
149 (
150  const fvPatchVectorField& ptf,
151  const labelList& addr
152 )
153 {
154  directionMixedFvPatchVectorField::rmap(ptf, addr);
155 
156  if (tangentialVelocity_.size())
157  {
159  refCast<const pressureInletOutletVelocityFvPatchVectorField>(ptf);
160 
161  tangentialVelocity_.rmap(tiptf.tangentialVelocity_, addr);
162  }
163 }
164 
165 
167 {
168  if (updated())
169  {
170  return;
171  }
172 
173  const fvsPatchField<scalar>& phip =
174  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
175 
176  valueFraction() = neg(phip)*(I - sqr(patch().nf()));
177 
178  directionMixedFvPatchVectorField::updateCoeffs();
179  directionMixedFvPatchVectorField::evaluate();
180 }
181 
182 
184 (
185  Ostream& os
186 )
187 const
188 {
190  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
191  if (tangentialVelocity_.size())
192  {
193  tangentialVelocity_.writeEntry("tangentialVelocity", os);
194  }
195  writeEntry("value", os);
196 }
197 
198 
199 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
200 
201 void Foam::pressureInletOutletVelocityFvPatchVectorField::operator=
202 (
203  const fvPatchField<vector>& pvf
204 )
205 {
206  tmp<vectorField> normalValue = transform(valueFraction(), refValue());
207  tmp<vectorField> transformGradValue = transform(I - valueFraction(), pvf);
208  fvPatchField<vector>::operator=(normalValue + transformGradValue);
209 }
210 
211 
212 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
213 
214 namespace Foam
215 {
217  (
220  );
221 }
222 
223 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensionedSymmTensor sqr(const dimensionedVector &dv)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:726
volVectorField vectorField(fieldObject, mesh)
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
pressureInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
const vectorField & tangentialVelocity() const
Return the tangential velocity.
static const Identity< scalar > I
Definition: Identity.H:93
This velocity inlet/outlet boundary condition is applied to pressure boundaries where the pressure is...
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:97
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:521
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:584
void setTangentialVelocity(const vectorField &tangentialVelocity)
Reset the tangential velocity.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:395
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
label n
A class for managing temporary objects.
Definition: PtrList.H:53
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.
dimensionSet transform(const dimensionSet &)
Definition: dimensionSet.C:477