Loading...
Searching...
No Matches
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-2025 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
29#include "dictionaryContent.H"
30#include "regionProperties.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37
38// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39
40void velocityFilmShellFvPatchVectorField::create_baffle()
41{
42 if (!baffle_)
43 {
44 baffle_.reset
45 (
46 baffleType::New(this->patch().boundaryMesh().mesh(), dict_)
47 );
48 }
49}
50
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
56 const fvPatch& p,
58)
59:
60 parent_bctype(p, iF),
61 dict_(),
62 curTimeIndex_(-1),
63 zeroWallVelocity_(true)
77)
78:
79 parent_bctype(ptf, p, iF, mapper),
80 dict_(ptf.dict_),
81 curTimeIndex_(-1),
82 zeroWallVelocity_(true)
83{}
84
85
87(
88 const fvPatch& p,
90 const dictionary& dict
91)
92:
93 parent_bctype(p, iF),
94 dict_
95 (
96 // Copy dictionary, but without "heavy" data chunks
97 dictionaryContent::copyDict
98 (
99 dict,
100 wordList(), // allow
101 wordList // deny
102 ({
103 "type", // redundant
104 "value", "refValue", "refGradient", "valueFraction"
105 })
106 )
107 ),
108 curTimeIndex_(-1),
109 zeroWallVelocity_(dict.getOrDefault<bool>("zeroWallVelocity", true))
110{
111 this->readValueEntry(dict, IOobjectOption::MUST_READ);
112
113 if (this->readMixedEntries(dict))
114 {
115 // Full restart
116 }
117 else
118 {
119 // Start from user entered data. Assume fixedValue.
120 refValue() = *this;
121 refGrad() = Zero;
122 valueFraction() = 1;
123 }
124
125 // Create baffle now.
126 // Lazy evaluation has issues with loading the finite-area fields
138)
139:
140 parent_bctype(ptf, iF),
141 dict_(ptf.dict_),
142 curTimeIndex_(-1),
143 zeroWallVelocity_(true)
144{}
145
146
147// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
148
150{
151 if (this->updated())
152 {
153 return;
154 }
155
156 // Create baffle if needed, but ignore if regionFaModels are disabled
157 if (!baffle_)
158 {
160 {
161 create_baffle();
162 }
163 else
164 {
165 static bool warned = false;
166 if (!warned)
167 {
168 warned = true;
170 << "Ignoring, regionFaModels are disabled" << endl;
171 }
172 return;
173 }
174 }
175
176 // Execute the change only once per time-step
177 if (curTimeIndex_ != this->db().time().timeIndex())
178 {
179 vectorField& pfld = *this;
180
181 auto& baffle = baffle_();
182
183 baffle.evolve();
184
185 baffle.vsm().mapToVolumePatch(baffle.Us(), pfld, patch().index());
186
187 refGrad() = Zero;
188 valueFraction() = 1;
189
190 if (zeroWallVelocity_)
191 {
192 refValue() = Zero;
193 }
194 else
195 {
196 refValue() = pfld;
197 }
198
199 curTimeIndex_ = this->db().time().timeIndex();
200 }
201
203}
204
205
207{
209 dict_.write(os, false);
210}
211
212
213// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214
216(
219);
220
221
222// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223
224} // End namespace Foam
225
226
227// ************************************************************************* //
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...
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
void write(Ostream &os, const bool subDict=true) const
Write dictionary, normally with sub-dictionary formatting.
const objectRegistry & db() const
The associated objectRegistry.
const fvPatch & patch() const noexcept
Return the patch.
bool updated() const noexcept
True if the boundary condition has already been updated.
A FieldMapper for finite-volume patch fields.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual void write(Ostream &) const
virtual Field< vector > & refGrad()
virtual Field< vector > & refValue()
virtual scalarField & valueFraction()
const Time & time() const noexcept
Return time registry.
static autoPtr< liquidFilmBase > New(const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected model using dictionary.
velocityFilmShellFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
volScalarField & p
dynamicFvMesh & mesh
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
#define InfoInFunction
Report an information message using Foam::Info.
bool allowFaModels() noexcept
The enable/disable state for regionFaModel (default: true).
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
fvPatchField< vector > fvPatchVectorField
label timeIndex
dictionary dict