velocityFilmShellFvPatchVectorField.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) 2020 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvPatch& p,
42 )
43 :
45  baffle_(),
46  dict_(dictionary::null),
47  curTimeIndex_(-1),
48  zeroWallVelocity_(true)
49 {
50  refValue() = 0;
51  refGrad() = 0;
52  valueFraction() = 1;
53 }
54 
55 
57 (
59  const fvPatch& p,
61  const fvPatchFieldMapper& mapper
62 )
63 :
65  (
66  ptf,
67  p,
68  iF,
69  mapper
70  ),
71  baffle_(),
72  dict_(ptf.dict_),
73  curTimeIndex_(-1),
74  zeroWallVelocity_(true)
75 {}
76 
77 
79 (
80  const fvPatch& p,
82  const dictionary& dict
83 )
84 :
86  baffle_(nullptr),
87  dict_(dict),
88  curTimeIndex_(-1),
89  zeroWallVelocity_(dict.getOrDefault<bool>("zeroWallVelocity", true))
90 {
92 
94 
95  if (dict.found("refValue"))
96  {
97  // Full restart
98  refValue() = vectorField("refValue", dict, p.size());
99  refGrad() = vectorField("refGradient", dict, p.size());
100  valueFraction() = scalarField("valueFraction", dict, p.size());
101  }
102  else
103  {
104  // Start from user entered data. Assume fixedValue.
105  refValue() = *this;
106  refGrad() = vector::zero;
107  valueFraction() = 1;
108  }
109 
110  if (!baffle_)
111  {
112  baffle_.reset(baffle::New(p, dict).ptr());
113  }
114 }
115 
116 
118 (
121 )
122 :
123  mixedFvPatchField<vector>(ptf, iF),
124  baffle_(),
125  dict_(ptf.dict_),
126  curTimeIndex_(-1),
127  zeroWallVelocity_(true)
128 {}
129 
130 
131 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
132 
133 
135 {
136  if (this->updated())
137  {
138  return;
139  }
140 
141  // Execute the change to the openFraction only once per time-step
142  if (curTimeIndex_ != this->db().time().timeIndex())
143  {
144  baffle_->evolve();
145 
146  volVectorField::Boundary& vfb =
148  (
149  this->internalField().name()
150  ).boundaryFieldRef();
151 
152  baffle_->vsm().mapToVolume(baffle_->Us(), vfb);
153 
154  refGrad() = Zero;
155  valueFraction() = 1;
156 
157  if (zeroWallVelocity_)
158  {
159  refValue() = Zero;
160  }
161  else
162  {
163  refValue() = vfb[patch().index()];
164  }
165  curTimeIndex_ = this->db().time().timeIndex();
166  }
167 
169 }
170 
171 
173 {
175 
176  // Remove value and type already written by mixedFvPatchField
177  dict_.remove("value");
178  dict_.remove("type");
179  dict_.remove("refValue");
180  dict_.remove("refGradient");
181  dict_.remove("valueFraction");
182  dict_.write(os, false);
183 }
184 
185 
186 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
187 
189 (
192 );
193 
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 } // End namespace Foam
198 
199 
200 // ************************************************************************* //
Foam::fvPatchField< vector >
Foam::fvPatchField< vector >::internalField
const DimensionedField< vector, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:363
Foam::velocityFilmShellFvPatchVectorField
Definition: velocityFilmShellFvPatchVectorField.H:211
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::regionModels::areaSurfaceFilmModels::liquidFilmBase
Definition: liquidFilmBase.H:60
Foam::mixedFvPatchField< vector >::valueFraction
virtual scalarField & valueFraction()
Definition: mixedFvPatchField.H:250
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::mixedFvPatchField< vector >::refValue
virtual Field< vector > & refValue()
Definition: mixedFvPatchField.H:230
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Foam::mixedFvPatchField< vector >::refGrad
virtual Field< vector > & refGrad()
Definition: mixedFvPatchField.H:240
Foam::fvPatchField< vector >::updated
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:387
velocityFilmShellFvPatchVectorField.H
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:392
Foam::TimeState::timeIndex
label timeIndex() const noexcept
Return current time index.
Definition: TimeStateI.H:37
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::velocityFilmShellFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: velocityFilmShellFvPatchVectorField.C:134
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::mixedFvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: mixedFvPatchField.C:235
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::fvPatch::index
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:197
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::velocityFilmShellFvPatchVectorField::velocityFilmShellFvPatchVectorField
velocityFilmShellFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: velocityFilmShellFvPatchVectorField.C:39
Foam::mixedFvPatchField< vector >
Foam::fvPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:321
Foam::objectRegistry::lookupObjectRef
Type & lookupObjectRef(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:478
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::fvPatchField< vector >::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:206
Foam::dictionary::remove
bool remove(const word &keyword)
Remove an entry specified by keyword.
Definition: dictionarySearch.C:582
Foam::dictionary::write
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:206
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
timeIndex
label timeIndex
Definition: getTimeIndex.H:30
Foam::fvPatchField< vector >::patch
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:357
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::velocityFilmShellFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: velocityFilmShellFvPatchVectorField.C:172
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::fvPatchField< vector >::operator=
virtual void operator=(const UList< vector > &)
Definition: fvPatchField.C:404
Foam::objectRegistry::time
const Time & time() const noexcept
Return time registry.
Definition: objectRegistry.H:178
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54