inverseFaceDistanceDiffusivity.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) 2011-2022 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 
27 #include "surfaceFields.H"
28 #include "fvPatchDistWave.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(inverseFaceDistanceDiffusivity, 0);
36 
38  (
39  motionDiffusivity,
40  inverseFaceDistanceDiffusivity,
41  Istream
42  );
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const fvMesh& mesh,
51  Istream& mdData
52 )
53 :
54  motionDiffusivity(mesh),
55  patchNames_(mdData)
56 {}
57 
58 
59 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
60 
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
66 
69 {
71  (
72  IOobject
73  (
74  "y",
75  mesh().time().timeName(),
76  mesh()
77  ),
78  mesh(),
80  );
81 
82  const labelHashSet patchIDs(mesh().boundaryMesh().patchSet(patchNames_));
83 
84  if (patchNames_.size())
85  {
87  (
88  mesh(),
89  patchIDs,
90  -vGreat,
91  y
92  );
93 
94  // Use cell distance on faces that are part of the patch set. This
95  // avoids divide-by-zero issues.
96  forAllConstIter(labelHashSet, patchIDs, iter)
97  {
98  const label patchi = iter.key();
99 
100  const labelUList& patchCells =
101  mesh().boundary()[patchi].faceCells();
102 
103  forAll(patchCells, patchFacei)
104  {
105  y.boundaryFieldRef()[patchi][patchFacei] =
106  mag
107  (
108  mesh().Cf().boundaryField()[patchi][patchFacei]
109  - mesh().C()[patchCells[patchFacei]]
110  );
111  }
112  }
113  }
114 
115  return surfaceScalarField::New("faceDiffusivity", 1/y);
116 }
117 
118 
119 // ************************************************************************* //
Foam::surfaceFields.
inverseFaceDistanceDiffusivity(const fvMesh &mesh, Istream &mdData)
Construct for the given fvMesh and data Istream.
Graphite solid properties.
Definition: C.H:48
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
static tmp< GeometricField< scalar, fvsPatchField, surfaceMesh > > New(const word &name, const Internal &, const PtrList< fvsPatchField< scalar >> &)
Return a temporary field constructed from name,.
virtual tmp< surfaceScalarField > operator()() const
Return diffusivity field.
const dimensionSet dimless
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
scalar y
label calculate(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, GeometricField< scalar, PatchField, GeoMesh > &distance)
Calculate distance data from patches.
word timeName
Definition: getTimeIndex.H:3
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Abstract base class for cell-centre mesh motion diffusivity.
defineTypeNameAndDebug(combustionModel, 0)
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Namespace for OpenFOAM.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:800