Loading...
Searching...
No Matches
vibrationShellFvPatchScalarField.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) 2019-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 vibrationShellFvPatchScalarField::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,
57 const DimensionedField<scalar, volMesh>& iF
58)
59:
60 parent_bctype(p, iF),
61 dict_()
74 const fvPatchFieldMapper& mapper
76:
77 parent_bctype(ptf, p, iF, mapper),
78 dict_(ptf.dict_)
79{}
80
81
83(
84 const fvPatch& p,
86 const dictionary& dict
87)
88:
89 parent_bctype(p, iF),
90 dict_
91 (
92 // Copy dictionary, but without "heavy" data chunks
93 dictionaryContent::copyDict
94 (
95 dict,
96 wordList(), // allow
97 wordList // deny
98 ({
99 "type", // redundant
100 "value", "refValue", "refGradient", "valueFraction"
101 })
102 )
103 )
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 // Create baffle now.
120 // Lazy evaluation has issues with loading the finite-area fields
131 const DimensionedField<scalar, volMesh>& iF
132)
133:
134 parent_bctype(ptf, iF),
135 dict_(ptf.dict_)
136{}
137
138
139// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
140
142{
143 if (this->updated())
144 {
145 return;
146 }
147
148 // Create baffle if needed, but ignore if regionFaModels are disabled
149 if (!baffle_)
150 {
152 {
153 create_baffle();
154 }
155 else
156 {
157 static bool warned = false;
158 if (!warned)
159 {
160 warned = true;
162 << "Ignoring, regionFaModels are disabled" << endl;
163 }
164 return;
165 }
166
167 }
168 auto& baffle = baffle_();
169
170 const auto& transportProperties =
171 db().lookupObject<IOdictionary>("transportProperties");
172
174
175 baffle.evolve();
176
177 // rho * acceleration
178
179 refGrad() = Zero; // safety (for any unmapped values)
180
181 baffle.vsm().mapToVolumePatch(baffle.a(), refGrad(), patch().index());
182
183 refGrad() *= rho.value();
184
190
191
193{
195 dict_.write(os, false);
196}
197
198
199// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200
202(
205);
206
207// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
208
209} // End namespace Foam
210
211
212// ************************************************************************* //
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...
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ MUST_READ
Reading required.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
This boundary condition provides a coupled acoustic pressure condition between a primary region (3D m...
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< scalar > & refGrad()
virtual Field< scalar > & refValue()
virtual scalarField & valueFraction()
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
static autoPtr< vibrationShellModel > New(const fvMesh &mesh, const dictionary &dict)
Return a reference to the selected model using dictionary.
vibrationShellFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, 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
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
fvPatchField< scalar > fvPatchScalarField
dictionary dict
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))