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-2022 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 
29 #include "dictionaryContent.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
40 (
41  const fvPatch& p,
43 )
44 :
46  baffle_(nullptr),
47  dict_(),
48  curTimeIndex_(-1),
49  zeroWallVelocity_(true)
50 {
51  refValue() = Zero;
52  refGrad() = Zero;
53  valueFraction() = 1;
54 }
55 
56 
58 (
60  const fvPatch& p,
62  const fvPatchFieldMapper& mapper
63 )
64 :
66  (
67  ptf,
68  p,
69  iF,
70  mapper
71  ),
72  baffle_(nullptr),
73  dict_(ptf.dict_),
74  curTimeIndex_(-1),
75  zeroWallVelocity_(true)
76 {}
77 
78 
80 (
81  const fvPatch& p,
83  const dictionary& dict
84 )
85 :
87  baffle_(nullptr),
88  dict_
89  (
90  // Copy dictionary, but without "heavy" data chunks
91  dictionaryContent::copyDict
92  (
93  dict,
94  wordList(), // allow
95  wordList // deny
96  ({
97  "type", // redundant
98  "value", "refValue", "refGradient", "valueFraction"
99  })
100  )
101  ),
102  curTimeIndex_(-1),
103  zeroWallVelocity_(dict.getOrDefault<bool>("zeroWallVelocity", true))
104 {
105  this->readValueEntry(dict, IOobjectOption::MUST_READ);
106 
107  if (this->readMixedEntries(dict))
108  {
109  // Full restart
110  }
111  else
112  {
113  // Start from user entered data. Assume fixedValue.
114  refValue() = *this;
115  refGrad() = Zero;
116  valueFraction() = 1;
117  }
118 
119  if (!baffle_)
120  {
121  baffle_.reset(baffleType::New(p.boundaryMesh().mesh(), dict_));
122  }
123 }
124 
125 
127 (
130 )
131 :
132  mixedFvPatchField<vector>(ptf, iF),
133  baffle_(nullptr),
134  dict_(ptf.dict_),
135  curTimeIndex_(-1),
136  zeroWallVelocity_(true)
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 {
144  if (this->updated())
145  {
146  return;
147  }
148 
149  // Execute the change only once per time-step
150  if (curTimeIndex_ != this->db().time().timeIndex())
151  {
152  baffle_->evolve();
153 
154  vectorField& pfld = *this;
155 
156  baffle_->vsm().mapToVolumePatch(baffle_->Us(), pfld, patch().index());
157 
158  refGrad() = Zero;
159  valueFraction() = 1;
160 
161  if (zeroWallVelocity_)
162  {
163  refValue() = Zero;
164  }
165  else
166  {
167  refValue() = pfld;
168  }
169 
170  curTimeIndex_ = this->db().time().timeIndex();
171  }
172 
174 }
175 
176 
178 {
180  dict_.write(os, false);
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
185 
187 (
190 );
191 
192 
193 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
194 
195 } // End namespace Foam
196 
197 
198 // ************************************************************************* //
dictionary dict
const objectRegistry & db() const
The associated objectRegistry.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const fvPatch & patch() const noexcept
Return the patch.
Definition: fvPatchField.H:269
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Macros for easy insertion into run-time selection tables.
velocityFilmShellFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual Field< vector > & refValue()
bool updated() const noexcept
True if the boundary condition has already been updated.
Definition: fvPatchField.H:296
const Time & time() const noexcept
Return time registry.
virtual scalarField & valueFraction()
A FieldMapper for finite-volume patch fields.
label timeIndex() const noexcept
Return the current time index.
Definition: TimeStateI.H:43
A wrapper for dictionary content, without operators that could affect inheritance patterns...
const Mesh & mesh() const noexcept
Return mesh.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
OBJstream os(runTime.globalPath()/outputName)
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:204
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:309
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void write(Ostream &) const
Write.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
volScalarField & p
label index() const noexcept
The index of this patch in the boundary mesh.
Definition: fvPatch.H:218
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual Field< vector > & refGrad()
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127