pressureInletOutletVelocityFvPatchVectorField.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) 2017-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 
34 
35 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
39 (
40  const fvPatch& p,
42 )
43 :
44  directionMixedFvPatchVectorField(p, iF),
45  phiName_("phi")
46 {
47  refValue() = Zero;
48  refGrad() = Zero;
49  valueFraction() = Zero;
50 }
51 
52 
55 (
57  const fvPatch& p,
59  const fvPatchFieldMapper& mapper
60 )
61 :
62  directionMixedFvPatchVectorField(ptf, p, iF, mapper),
63  phiName_(ptf.phiName_)
64 {
65  if (ptf.tangentialVelocity_.size())
66  {
67  tangentialVelocity_ = mapper(ptf.tangentialVelocity_);
68  }
69 }
70 
71 
74 (
75  const fvPatch& p,
77  const dictionary& dict
78 )
79 :
80  directionMixedFvPatchVectorField(p, iF),
81  phiName_(dict.getOrDefault<word>("phi", "phi"))
82 {
83  patchType() = dict.getOrDefault<word>("patchType", word::null);
84  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
85 
86  if (dict.found("tangentialVelocity"))
87  {
88  setTangentialVelocity
89  (
90  vectorField("tangentialVelocity", dict, p.size())
91  );
92  }
93  else
94  {
95  refValue() = Zero;
96  }
97 
98  refGrad() = Zero;
99  valueFraction() = Zero;
100 }
101 
102 
105 (
107 )
108 :
109  directionMixedFvPatchVectorField(pivpvf),
110  phiName_(pivpvf.phiName_),
111  tangentialVelocity_(pivpvf.tangentialVelocity_)
112 {}
113 
114 
117 (
120 )
121 :
122  directionMixedFvPatchVectorField(pivpvf, iF),
123  phiName_(pivpvf.phiName_),
124  tangentialVelocity_(pivpvf.tangentialVelocity_)
125 {}
126 
127 
128 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
129 
131 setTangentialVelocity(const vectorField& tangentialVelocity)
132 {
133  tangentialVelocity_ = tangentialVelocity;
134  const vectorField n(patch().nf());
135  refValue() = tangentialVelocity_ - n*(n & tangentialVelocity_);
136 }
137 
138 
140 (
141  const fvPatchFieldMapper& m
142 )
143 {
144  directionMixedFvPatchVectorField::autoMap(m);
145  if (tangentialVelocity_.size())
146  {
147  tangentialVelocity_.autoMap(m);
148  }
149 }
150 
151 
153 (
154  const fvPatchVectorField& ptf,
155  const labelList& addr
156 )
157 {
158  directionMixedFvPatchVectorField::rmap(ptf, addr);
159 
160  if (tangentialVelocity_.size())
161  {
163  refCast<const pressureInletOutletVelocityFvPatchVectorField>(ptf);
164 
165  tangentialVelocity_.rmap(tiptf.tangentialVelocity_, addr);
166  }
167 }
168 
169 
171 {
172  if (updated())
173  {
174  return;
175  }
176 
177  const fvsPatchField<scalar>& phip =
178  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
179 
180  valueFraction() = neg(phip)*(I - sqr(patch().nf()));
181 
182  directionMixedFvPatchVectorField::updateCoeffs();
184 }
185 
186 
188 (
189  Ostream& os
190 )
191 const
192 {
194  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
195  if (tangentialVelocity_.size())
196  {
197  tangentialVelocity_.writeEntry("tangentialVelocity", os);
198  }
199  writeEntry("value", os);
200 }
201 
202 
203 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
204 
205 void Foam::pressureInletOutletVelocityFvPatchVectorField::operator=
206 (
207  const fvPatchField<vector>& pvf
208 )
209 {
210  tmp<vectorField> normalValue = transform(valueFraction(), refValue());
211  tmp<vectorField> transformGradValue = transform(I - valueFraction(), pvf);
212  fvPatchField<vector>::operator=(normalValue + transformGradValue);
213 }
214 
215 
216 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
217 
218 namespace Foam
219 {
221  (
224  );
225 }
226 
227 // ************************************************************************* //
Foam::fvPatchField< vector >
volFields.H
Foam::fvPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::Ostream::writeEntryIfDifferent
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:248
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
surfaceFields.H
Foam::surfaceFields.
Foam::pressureInletOutletVelocityFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: pressureInletOutletVelocityFvPatchVectorField.C:188
Foam::pressureInletOutletVelocityFvPatchVectorField::tangentialVelocity
const vectorField & tangentialVelocity() const
Return the tangential velocity.
Definition: pressureInletOutletVelocityFvPatchVectorField.H:205
Foam::pressureInletOutletVelocityFvPatchVectorField::pressureInletOutletVelocityFvPatchVectorField
pressureInletOutletVelocityFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: pressureInletOutletVelocityFvPatchVectorField.C:39
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
Foam::pressureInletOutletVelocityFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: pressureInletOutletVelocityFvPatchVectorField.C:170
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::Field< vector >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::pressureInletOutletVelocityFvPatchVectorField::rmap
virtual void rmap(const fvPatchVectorField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: pressureInletOutletVelocityFvPatchVectorField.C:153
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
pressureInletOutletVelocityFvPatchVectorField.H
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:54
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::pressureInletOutletVelocityFvPatchVectorField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: pressureInletOutletVelocityFvPatchVectorField.C:140
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::pressureInletOutletVelocityFvPatchVectorField::setTangentialVelocity
void setTangentialVelocity(const vectorField &tangentialVelocity)
Reset the tangential velocity.
Definition: pressureInletOutletVelocityFvPatchVectorField.C:131
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::fvPatchField::operator=
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:404
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199
Foam::I
static const Identity< scalar > I
Definition: Identity.H:95
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::pressureInletOutletVelocityFvPatchVectorField
This velocity inlet/outlet boundary condition is applied to velocity boundaries where the pressure is...
Definition: pressureInletOutletVelocityFvPatchVectorField.H:98
Foam::stringOps::evaluate
string evaluate(label fieldWidth, const std::string &s, size_t pos=0, size_t len=std::string::npos)
String evaluation with specified (positive, non-zero) field width.
Definition: stringOpsEvaluate.C:37