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-2019 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 (
53  const fvPatch& p,
55  const dictionary& dict
56 )
57 :
58  directionMixedFvPatchVectorField(p, iF),
59  phiName_(dict.lookupOrDefault<word>("phi", "phi"))
60 {
61  fvPatchVectorField::operator=(vectorField("value", dict, p.size()));
62 
63  if (dict.found("tangentialVelocity"))
64  {
65  tangentialVelocity_ =
66  Function1<vector>::New("tangentialVelocity", dict);
67  }
68 
69  refValue() = Zero;
70  refGrad() = Zero;
71  valueFraction() = Zero;
72 }
73 
74 
77 (
79  const fvPatch& p,
81  const fvPatchFieldMapper& mapper
82 )
83 :
84  directionMixedFvPatchVectorField(ptf, p, iF, mapper),
85  phiName_(ptf.phiName_),
86  tangentialVelocity_(ptf.tangentialVelocity_, false)
87 {}
88 
89 
92 (
94 )
95 :
96  directionMixedFvPatchVectorField(pivpvf),
97  phiName_(pivpvf.phiName_),
98  tangentialVelocity_(pivpvf.tangentialVelocity_, false)
99 {}
100 
101 
104 (
107 )
108 :
109  directionMixedFvPatchVectorField(pivpvf, iF),
110  phiName_(pivpvf.phiName_),
111  tangentialVelocity_(pivpvf.tangentialVelocity_, false)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 {
119  if (updated())
120  {
121  return;
122  }
123 
124  if (tangentialVelocity_.valid())
125  {
126  const scalar t = this->db().time().timeOutputValue();
127  const vector tangentialVelocity = tangentialVelocity_->value(t);
128  const vectorField n(patch().nf());
129  refValue() = tangentialVelocity - n*(n & tangentialVelocity);
130  }
131 
132  const fvsPatchField<scalar>& phip =
133  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
134 
135  valueFraction() = neg(phip)*(I - sqr(patch().nf()));
136 
137  directionMixedFvPatchVectorField::updateCoeffs();
138  directionMixedFvPatchVectorField::evaluate();
139 }
140 
141 
143 (
144  Ostream& os
145 )
146 const
147 {
149  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
150  if (tangentialVelocity_.valid())
151  {
152  writeEntry(os, tangentialVelocity_());
153  }
154  writeEntry(os, "value", *this);
155 }
156 
157 
158 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
159 
160 void Foam::pressureInletOutletVelocityFvPatchVectorField::operator=
161 (
162  const fvPatchField<vector>& pvf
163 )
164 {
165  tmp<vectorField> normalValue = transform(valueFraction(), refValue());
166  tmp<vectorField> transformGradValue = transform(I - valueFraction(), pvf);
167  fvPatchField<vector>::operator=(normalValue + transformGradValue);
168 }
169 
170 
171 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
172 
173 namespace Foam
174 {
176  (
179  );
180 }
181 
182 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensionedSymmTensor sqr(const dimensionedVector &dv)
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
volVectorField vectorField(fieldObject, mesh)
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
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.
virtual Type value(const scalar x) const =0
Return value as a function of (scalar) independent variable.
static const Identity< scalar > I
Definition: Identity.H:93
This velocity inlet/outlet boundary condition is applied to pressure boundaries where the pressure is...
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:97
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
virtual label size() const
Return size.
Definition: fvPatch.H:155
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
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:295
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