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-2020 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  option,
46  acousticDampingSource,
47  dictionary
48  );
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
54 
56 {
58 
59  const vectorField& Cf = mesh_.C();
60 
61  const scalar pi = constant::mathematical::pi;
62 
63  forAll(cells_, i)
64  {
65  label celli = cells_[i];
66  scalar d = mag(Cf[celli] - x0_);
67 
68  if (d < r1_)
69  {
70  blendFactor_[celli] = 0.0;
71  }
72  else if ((d >= r1_) && (d <= r2_))
73  {
74  blendFactor_[celli] = (1.0 - cos(pi*mag(d - r1_)/(r2_ - r1_)))/2.0;
75  }
76  }
77 
79 }
80 
81 
82 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
83 
84 Foam::fv::acousticDampingSource::acousticDampingSource
85 (
86  const word& name,
87  const word& modelType,
88  const dictionary& dict,
89  const fvMesh& mesh
90 )
91 :
92  cellSetOption(name, modelType, dict, mesh),
93  frequency_("frequency", dimless/dimTime, 0),
94  blendFactor_
95  (
97  (
98  IOobject
99  (
100  name_ + ":blend",
101  mesh_.time().timeName(),
102  mesh_,
105  ),
106  mesh_,
107  dimensionedScalar("blend0", dimless, 1.0),
109  )
110  ),
111  URefName_("unknown-URef"),
112  x0_(Zero),
113  r1_(0),
114  r2_(0),
115  w_(20)
116 {
117  read(dict);
118 }
119 
120 
121 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
122 
124 (
125  fvMatrix<vector>& eqn,
126  const label fieldI
127 )
128 {
129  const volVectorField& U = eqn.psi();
130  const volScalarField coeff(name_ + ":coeff", w_*frequency_*blendFactor_);
131  const volVectorField& URef(mesh().lookupObject<volVectorField>(URefName_));
132 
133  fvMatrix<vector> dampingEqn
134  (
135  fvm::Sp(coeff, U) - coeff*URef
136  );
137  eqn -= dampingEqn;
138 }
139 
140 
142 (
143  const volScalarField& rho,
144  fvMatrix<vector>& eqn,
145  const label fieldI
146 )
147 {
148  const volVectorField& U = eqn.psi();
149  const volScalarField coeff(name_ + ":coeff", w_*frequency_*blendFactor_);
150  const volVectorField& 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  const volVectorField& U = eqn.psi();
169  const volScalarField coeff(name_ + ":coeff", w_*frequency_*blendFactor_);
170  const volVectorField& URef(mesh().lookupObject<volVectorField>(URefName_));
171 
172  fvMatrix<vector> dampingEqn
173  (
174  fvm::Sp(alpha*rho*coeff, U) - alpha*rho*coeff*URef
175  );
176  eqn -= dampingEqn;
177 }
178 
179 
181 {
183  {
184  if (!coeffs_.readIfPresent("UNames", fieldNames_))
185  {
186  fieldNames_.resize(1);
187  fieldNames_.first() = coeffs_.getOrDefault<word>("U", "U");
188  }
189 
190  applied_.setSize(fieldNames_.size(), false);
191 
192  coeffs_.readEntry("frequency", frequency_.value());
193  coeffs_.readEntry("URef", URefName_);
194  coeffs_.readCompat<vector>("origin", {{"centre", -1806}}, x0_);
195  coeffs_.readEntry("radius1", r1_);
196  coeffs_.readEntry("radius2", r2_);
197 
198  if (coeffs_.readIfPresent("w", w_))
199  {
200  Info<< name_ << ": Setting stencil width to " << w_ << endl;
201  }
202 
203  setBlendingFactor();
204 
205  return true;
206  }
207 
208  return false;
209 }
210 
211 
212 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
mathematicalConstants.H
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fv::cellSetOption
Cell-set options abstract base class. Provides a base set of controls, e.g.:
Definition: cellSetOption.H:72
Foam::fv::acousticDampingSource::addSup
virtual void addSup(fvMatrix< vector > &eqn, const label fieldI)
Add implicit contribution to momentum equation.
Definition: acousticDampingSource.C:124
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::fv::acousticDampingSource::r2_
scalar r2_
Definition: acousticDampingSource.H:102
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::fv::acousticDampingSource::setBlendingFactor
void setBlendingFactor()
Helper function to set the blending factor.
Definition: acousticDampingSource.C:55
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:107
Foam::fv::acousticDampingSource::x0_
point x0_
Definition: acousticDampingSource.H:96
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
Foam::fv::acousticDampingSource::r1_
scalar r1_
Definition: acousticDampingSource.H:99
rho
rho
Definition: readInitialConditions.H:88
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:285
Foam::fv::option::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: fvOption.H:82
Foam::zeroGradientFvPatchField
This boundary condition applies a zero-gradient condition from the patch internal field onto the patc...
Definition: zeroGradientFvPatchField.H:64
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fv::acousticDampingSource::read
virtual bool read(const dictionary &dict)
Read dictionary.
Definition: acousticDampingSource.C:180
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::fv::cellSetOption::read
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: cellSetOption.C:255
Foam::fv::cellSetOption::cells_
labelList cells_
Set of cells to apply source to.
Definition: cellSetOption.H:113
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:341
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
fvmSup.H
Calculate the matrix for implicit and explicit sources.
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:940
acousticDampingSource.H
fv
labelList fv(nPoints)
U
U
Definition: pEqn.H:72
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Vector< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::fv::acousticDampingSource::blendFactor_
volScalarField blendFactor_
Blending factor [].
Definition: acousticDampingSource.H:90
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::fv::defineTypeNameAndDebug
defineTypeNameAndDebug(atmAmbientTurbSource, 0)
Foam::fv::addToRunTimeSelectionTable
addToRunTimeSelectionTable(option, atmAmbientTurbSource, dictionary)
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:123
zeroGradientFvPatchFields.H
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265