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-------------------------------------------------------------------------------
10License
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
33namespace 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;
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
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// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
label timeIndex() const noexcept
Return current time index.
Definition: TimeStateI.H:37
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool remove(const word &keyword)
Remove an entry specified by keyword.
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
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
Definition: dictionaryIO.C:206
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:210
virtual void operator=(const UList< vector > &)
Definition: fvPatchField.C:408
const DimensionedField< Type, volMesh > & internalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:368
bool updated() const
Return true if the boundary condition has already been updated.
Definition: fvPatchField.H:392
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:362
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:203
This boundary condition provides a base class for 'mixed' type boundary conditions,...
virtual Field< vector > & refGrad()
virtual Field< vector > & refValue()
virtual scalarField & valueFraction()
const Time & time() const noexcept
Return time registry.
Type & lookupObjectRef(const word &name, const bool recursive=false) const
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
volScalarField & p
bool
Definition: EEqn.H:20
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
label timeIndex
Definition: getTimeIndex.H:30
dictionary dict
A non-counting (dummy) refCount.
Definition: refCount.H:59