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-2021 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  (
91  (
92  IOobject
93  (
94  name_ + ":blend",
95  mesh_.time().timeName(),
96  mesh_,
97  IOobject::NO_READ,
98  IOobject::NO_WRITE
99  ),
100  mesh_,
101  scalar(1),
102  dimless,
104  )
105  ),
106  frequency_("frequency", dimless/dimTime, 0),
107  x0_(Zero),
108  r1_(0),
109  r2_(0),
110  URefName_("unknown-URef"),
111  w_(20)
112 {
114 }
115 
116 
117 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118 
120 (
121  fvMatrix<vector>& eqn,
122  const label fieldI
123 )
124 {
125  const volVectorField& U = eqn.psi();
126  const volScalarField coeff(name_ + ":coeff", w_*frequency_*blendFactor_);
127  const auto& URef = mesh().lookupObject<volVectorField>(URefName_);
128 
129  fvMatrix<vector> dampingEqn
130  (
131  fvm::Sp(coeff, U) - coeff*URef
132  );
133  eqn -= dampingEqn;
134 }
135 
136 
138 (
139  const volScalarField& rho,
140  fvMatrix<vector>& eqn,
141  const label fieldI
142 )
143 {
144  const volVectorField& U = eqn.psi();
145  const volScalarField coeff(name_ + ":coeff", w_*frequency_*blendFactor_);
146  const auto& URef = mesh().lookupObject<volVectorField>(URefName_);
147 
148  fvMatrix<vector> dampingEqn
149  (
150  fvm::Sp(rho*coeff, U) - rho*coeff*URef
151  );
152  eqn -= dampingEqn;
153 }
154 
155 
157 (
158  const volScalarField& alpha,
159  const volScalarField& rho,
160  fvMatrix<vector>& eqn,
161  const label fieldI
162 )
163 {
164  const volVectorField& U = eqn.psi();
165  const volScalarField coeff(name_ + ":coeff", w_*frequency_*blendFactor_);
166  const auto& URef = mesh().lookupObject<volVectorField>(URefName_);
167 
168  fvMatrix<vector> dampingEqn
169  (
170  fvm::Sp(alpha*rho*coeff, U) - alpha*rho*coeff*URef
171  );
172  eqn -= dampingEqn;
173 }
174 
175 
177 {
179  {
180  if (!coeffs_.readIfPresent("UNames", fieldNames_))
181  {
182  fieldNames_.resize(1);
183  fieldNames_.first() = coeffs_.getOrDefault<word>("U", "U");
184  }
185 
187 
188  coeffs_.readEntry("frequency", frequency_.value());
189  coeffs_.readEntry("URef", URefName_);
190  coeffs_.readCompat<vector>("origin", {{"centre", -1806}}, x0_);
191  coeffs_.readEntry("radius1", r1_);
192  coeffs_.readEntry("radius2", r2_);
193 
194  if (coeffs_.readIfPresent("w", w_))
195  {
196  Info<< name_ << ": Setting stencil width to " << w_ << endl;
197  }
198 
199  setBlendingFactor();
200 
201  return true;
202  }
203 
204  return false;
205 }
206 
207 
208 // ************************************************************************* //
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
word timeName
Definition: getTimeIndex.H:3
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.
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:172
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