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 -------------------------------------------------------------------------------
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 
28 #include "acousticDampingSource.H"
29 #include "fvMesh.H"
30 #include "fvMatrices.H"
31 #include "fvmSup.H"
33 #include "mathematicalConstants.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace fv
41 {
42  defineTypeNameAndDebug(acousticDampingSource, 0);
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  }
72 
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,
95  ),
96  frequency_("frequency", dimless/dimTime, Zero),
97  x0_(Zero),
98  r1_(0),
99  r2_(0),
100  URefName_("unknown-URef"),
101  w_(20)
102 {
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
127  (
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
153  (
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
180  (
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 // ************************************************************************* //
dictionary dict
const word zeroGradientType
A zeroGradient patch field type.
virtual bool read(const dictionary &dict)
Read dictionary.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:142
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual bool read(const dictionary &dict)
Read source dictionary.
Template invariant parts for fvPatchField.
Definition: fvPatchField.H:77
const dimensionSet dimless
Dimensionless.
void setBlendingFactor()
Helper function to set the blending factor.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
labelList fv(nPoints)
constexpr scalar pi(M_PI)
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
acousticDampingSource(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
Applies sources on velocity (i.e. U) within a specified region to enable acoustic damping...
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
labelList cells_
Set of cells to apply source to.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
U
Definition: pEqn.H:72
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void resetApplied()
Resize/reset applied flag list for all fieldNames_ entries.
Definition: fvOption.C:41
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
A special matrix type and solver, designed for finite volume solutions of scalar equations.
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
Information stream (stdout output on master, null elsewhere)
Intermediate abstract class for handling cell-set options for the derived fvOptions.
virtual void addSup(fvMatrix< vector > &eqn, const label fieldI)
Add implicit contribution to momentum equation.
const volVectorField & C() const
Return cell centres as volVectorField.
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
Do not request registration (bool: false)
volScalarField blendFactor_
Blending factor [-].
Calculate the finiteVolume matrix for implicit and explicit sources.
Namespace for OpenFOAM.
Base abstract class for handling finite volume options (i.e. fvOption).
Definition: fvOption.H:123
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127