wallDampingModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "wallDampingModel.H"
27 #include "surfaceInterpolate.H"
28 #include "wallFvPatch.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 (
43  const dictionary& dict,
44  const phaseInterface& interface
45 )
46 :
47  wallDependentModel(interface.mesh()),
48  interface_
49  (
50  interface.modelCast<wallDampingModel, dispersedPhaseInterface>()
51  ),
52  Cd_("Cd", dimless, dict),
53  zeroWallDist_("zeroWallDist", dimLength, dict, 0),
54  zeroInNearWallCells_
55  (
56  dict.lookupOrDefault<Switch>("zeroInNearWallCells", false)
57  )
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
71 {
72  tmp<volScalarField> tlimiter(limiter().ptr());
73 
74  if (zeroInNearWallCells_)
75  {
76  volScalarField& limiter = tlimiter.ref();
77 
78  const fvBoundaryMesh& bMesh = limiter.mesh().boundary();
79 
81  {
82  if (isA<wallFvPatch>(bMesh[patchi]))
83  {
84  const labelUList& faceCells = bMesh[patchi].faceCells();
85 
86  forAll(faceCells, facei)
87  {
88  limiter[faceCells[facei]] = 0;
89  }
90  }
91  }
92 
93  return tlimiter;
94  }
95  else
96  {
97  return tlimiter;
98  }
99 }
100 
101 
104 {
105  return fvc::interpolate(damping());
106 }
107 
108 
109 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic GeometricField class.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Class to represent a interface between phases where one phase is considered dispersed within the othe...
Foam::fvBoundaryMesh.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
Wall damping models can be used to filter interfacial models near the walls. This is particularly use...
virtual tmp< surfaceScalarField > dampingf() const
Return damped face coefficient.
virtual ~wallDampingModel()
Destructor.
virtual tmp< volScalarField > damping() const
Return damped coefficient.
wallDampingModel(const dictionary &dict, const phaseInterface &interface)
Construct from a dictionary and an interface.
A class which provides on-demand creation and caching of wall distance and wall normal fields for use...
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
const polyBoundaryMesh & bMesh
label patchi
void limiter(const control &controls, surfaceScalarField &lambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const scalarField &SuCorr, const surfaceScalarField &phi, const surfaceScalarField &phiCorr, const SpType &Sp, const PsiMaxType &psiMax, const PsiMinType &psiMin)
Definition: MULESlimiter.C:42
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Namespace for OpenFOAM.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
const dimensionSet dimless
const dimensionSet dimLength
defineTypeNameAndDebug(combustionModel, 0)
dictionary dict