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-2020 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 "phasePair.H"
28 #include "surfaceInterpolate.H"
29 #include "wallFvPatch.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(wallDampingModel, 0);
36  defineRunTimeSelectionTable(wallDampingModel, dictionary);
37 }
38 
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  const phasePair& pair
48 )
49 :
50  wallDependentModel(pair.phase1().mesh()),
51  pair_(pair),
52  Cd_("Cd", dimless, dict),
53  zeroWallDist_
54  (
55  dimensionedScalar::lookupOrDefault
56  (
57  "zeroWallDist",
58  dict,
59  dimLength,
60  0
61  )
62  ),
63  zeroInNearWallCells_
64  (
65  dict.lookupOrDefault<Switch>("zeroInNearWallCells", false)
66  )
67 {}
68 
69 
70 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
71 
73 {}
74 
75 
76 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77 
80 {
81  tmp<volScalarField> tlimiter(limiter().ptr());
82 
84  {
85  volScalarField& limiter = tlimiter.ref();
86 
87  const fvBoundaryMesh& bMesh = limiter.mesh().boundary();
88 
89  forAll(bMesh, patchi)
90  {
91  if (isA<wallFvPatch>(bMesh[patchi]))
92  {
93  const labelUList& faceCells = bMesh[patchi].faceCells();
94 
95  forAll(faceCells, facei)
96  {
97  limiter[faceCells[facei]] = 0;
98  }
99  }
100  }
101 
102  return tlimiter;
103  }
104  else
105  {
106  return tlimiter;
107  }
108 }
109 
110 
113 {
114  return fvc::interpolate(damping());
115 }
116 
117 
118 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< surfaceScalarField > dampingf() const
Return damped face coefficient.
UList< label > labelUList
Definition: UList.H:65
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
static const dimensionSet dimF
Coefficient dimensions.
const Switch zeroInNearWallCells_
Set the value to zero in wall-adjacent cells.
Dimension set for the base types.
Definition: dimensionSet.H:120
dynamicFvMesh & mesh
virtual ~wallDampingModel()
Destructor.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
const Mesh & mesh() const
Return mesh.
defineTypeNameAndDebug(combustionModel, 0)
virtual tmp< volScalarField > limiter() const
Return the force limiter field.
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
Internal & ref()
Return a reference to the dimensioned internal field.
label patchi
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
wallDampingModel(const dictionary &dict, const phasePair &pair)
Construct from components.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > damping() const
Return damped coefficient.
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.