Loading...
Searching...
No Matches
acousticDampingSource.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) 2016-2023 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 "fvMesh.H"
30#include "fvMatrices.H"
31#include "fvmSup.H"
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace fv
41{
44}
45}
46
47
48// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49
51{
53
54 const vectorField& Cf = mesh_.C();
55
56 const scalar pi = constant::mathematical::pi;
57
58 forAll(cells_, i)
59 {
60 label celli = cells_[i];
61 scalar d = mag(Cf[celli] - x0_);
62
63 if (d < r1_)
64 {
65 blendFactor_[celli] = 0.0;
66 }
67 else if ((d >= r1_) && (d <= r2_))
68 {
69 blendFactor_[celli] = (1.0 - cos(pi*mag(d - r1_)/(r2_ - r1_)))/2.0;
70 }
71 }
73 blendFactor_.correctBoundaryConditions();
74}
75
76
77// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78
80(
81 const word& name,
82 const word& modelType,
83 const dictionary& dict,
84 const fvMesh& mesh
85)
86:
87 fv::cellSetOption(name, modelType, dict, mesh),
88 blendFactor_
89 (
90 mesh_.newIOobject(IOobject::scopedName(name_, "blend")),
91 mesh_,
92 scalar(1),
93 dimless,
94 fvPatchFieldBase::zeroGradientType()
95 ),
96 frequency_("frequency", dimless/dimTime, Zero),
97 x0_(Zero),
98 r1_(0),
99 r2_(0),
100 URefName_("unknown-URef"),
101 w_(20)
103 read(dict);
104}
105
106
107// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108
110(
111 fvMatrix<vector>& eqn,
112 const label fieldI
113)
114{
115 auto tcoeff = volScalarField::New
116 (
117 IOobject::scopedName(name_, "coeff"),
119 w_*frequency_*blendFactor_
120 );
121 const auto& coeff = tcoeff();
122
123 const volVectorField& U = eqn.psi();
124 const auto& URef = mesh().lookupObject<volVectorField>(URefName_);
125
126 fvMatrix<vector> dampingEqn
128 fvm::Sp(coeff, U) - coeff*URef
129 );
130 eqn -= dampingEqn;
131}
132
133
135(
136 const volScalarField& rho,
137 fvMatrix<vector>& eqn,
138 const label fieldI
139)
140{
141 auto tcoeff = volScalarField::New
142 (
143 IOobject::scopedName(name_, "coeff"),
145 w_*frequency_*blendFactor_
146 );
147 const auto& coeff = tcoeff();
148
149 const volVectorField& U = eqn.psi();
150 const auto& URef = mesh().lookupObject<volVectorField>(URefName_);
151
152 fvMatrix<vector> dampingEqn
154 fvm::Sp(rho*coeff, U) - rho*coeff*URef
155 );
156 eqn -= dampingEqn;
157}
158
159
161(
162 const volScalarField& alpha,
163 const volScalarField& rho,
164 fvMatrix<vector>& eqn,
165 const label fieldI
166)
167{
168 auto tcoeff = volScalarField::New
169 (
170 IOobject::scopedName(name_, "coeff"),
172 w_*frequency_*blendFactor_
173 );
174 const auto& coeff = tcoeff();
175
176 const volVectorField& U = eqn.psi();
177 const auto& URef = mesh().lookupObject<volVectorField>(URefName_);
178
179 fvMatrix<vector> dampingEqn
181 fvm::Sp(alpha*rho*coeff, U) - alpha*rho*coeff*URef
182 );
183 eqn -= dampingEqn;
184}
185
186
188{
190 {
191 if (!coeffs_.readIfPresent("UNames", fieldNames_))
192 {
193 fieldNames_.resize(1);
194 fieldNames_.first() = coeffs_.getOrDefault<word>("U", "U");
195 }
196
198
199 coeffs_.readEntry("frequency", frequency_.value());
200 coeffs_.readEntry("URef", URefName_);
201 coeffs_.readCompat<vector>("origin", {{"centre", -1806}}, x0_);
202 coeffs_.readEntry("radius1", r1_);
203 coeffs_.readEntry("radius2", r2_);
204
205 if (coeffs_.readIfPresent("w", w_))
206 {
207 Info<< name_ << ": Setting stencil width to " << w_ << endl;
208 }
209
210 setBlendingFactor();
211
212 return true;
213 }
214
215 return false;
216}
217
218
219// ************************************************************************* //
constexpr scalar pi(M_PI)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const volVectorField & C() const
Return cell centres as volVectorField.
Template invariant parts for fvPatchField.
Applies sources on velocity (i.e. U) within a specified region to enable acoustic damping.
volScalarField blendFactor_
Blending factor [-].
acousticDampingSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
virtual bool read(const dictionary &dict)
Read dictionary.
word URefName_
Name of reference velocity field.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldI)
Add implicit contribution to momentum equation.
void setBlendingFactor()
Helper function to set the blending factor.
dimensionedScalar frequency_
Frequency [Hz].
Intermediate abstract class for handling cell-set options for the derived fvOptions.
labelList cells_
Set of cells to apply source to.
virtual bool read(const dictionary &dict)
Read source dictionary.
cellSetOption(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition fvOption.H:124
const word & name() const noexcept
Return const access to the source name.
Definition fvOptionI.H:24
const fvMesh & mesh_
Reference to the mesh database.
Definition fvOption.H:142
wordList fieldNames_
Field names to apply source to - populated by derived models.
Definition fvOption.H:157
dictionary coeffs_
Dictionary containing source coefficients.
Definition fvOption.H:152
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition fvOption.C:41
const fvMesh & mesh() const noexcept
Return const access to the mesh database.
Definition fvOptionI.H:30
const word name_
Source name.
Definition fvOption.H:132
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
dynamicFvMesh & mesh
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Calculate the finiteVolume matrix for implicit and explicit sources.
constexpr scalar pi(M_PI)
Namespace for finite-volume.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57
dimensionedScalar cos(const dimensionedScalar &ds)
volScalarField & alpha
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299