inversePointDistanceDiffusivity.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-2023 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 "HashSet.H"
29 #include "pointEdgeDist.H"
30 #include "PointEdgeWave.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
38 
40  (
43  Istream
44  );
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const fvMesh& mesh,
53  Istream& mdData
54 )
55 :
56  motionDiffusivity(mesh),
57  patchNames_(mdData)
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
71 {
72  tmp<surfaceScalarField> tfaceDiffusivity
73  (
75  (
76  IOobject
77  (
78  "faceDiffusivity",
79  mesh().time().name(),
80  mesh(),
83  ),
84  mesh(),
86  )
87  );
88 
89  surfaceScalarField& faceDiffusivity = tfaceDiffusivity.ref();
90 
91  const polyBoundaryMesh& bdry = mesh().boundaryMesh();
92 
93  const labelHashSet patchSet(bdry.patchSet(patchNames_));
94 
95  label nPatchEdges = 0;
96 
97  forAllConstIter(labelHashSet, patchSet, iter)
98  {
99  nPatchEdges += bdry[iter.key()].nEdges();
100  }
101 
102  // Distance to wall on points and edges.
103  List<pointEdgeDist> pointWallDist(mesh().nPoints());
104  List<pointEdgeDist> edgeWallDist(mesh().nEdges());
105 
106  pointEdgeDist::data pointEdgeData(mesh().points());
107 
108  {
109  // Seeds
110  List<pointEdgeDist> seedInfo(nPatchEdges);
111  labelList seedPoints(nPatchEdges);
112 
113  nPatchEdges = 0;
114 
115  forAllConstIter(labelHashSet, patchSet, iter)
116  {
117  const polyPatch& patch = bdry[iter.key()];
118 
119  const labelList& meshPoints = patch.meshPoints();
120 
121  forAll(meshPoints, i)
122  {
123  const label pointi = meshPoints[i];
124 
125  if (!pointWallDist[pointi].valid(pointEdgeData))
126  {
127  // Not yet seeded
128  seedInfo[nPatchEdges] = pointEdgeDist
129  (
130  mesh().points()[pointi],
131  0
132  );
133  seedPoints[nPatchEdges] = pointi;
134  pointWallDist[pointi] = seedInfo[nPatchEdges];
135 
136  nPatchEdges++;
137  }
138  }
139  }
140  seedInfo.setSize(nPatchEdges);
141  seedPoints.setSize(nPatchEdges);
142 
143  // Do calculations
145  <
148  > waveInfo
149  (
150  mesh(),
151  seedPoints,
152  seedInfo,
153 
154  pointWallDist,
155  edgeWallDist,
156  mesh().globalData().nTotalPoints(),// max iterations
157  pointEdgeData
158  );
159  }
160 
161 
162  for (label facei=0; facei<mesh().nInternalFaces(); facei++)
163  {
164  const face& f = mesh().faces()[facei];
165 
166  scalar dist = 0;
167 
168  forAll(f, fp)
169  {
170  dist += sqrt(pointWallDist[f[fp]].distSqr());
171  }
172  dist /= f.size();
173 
174  faceDiffusivity[facei] = 1.0/dist;
175  }
176 
177 
178  surfaceScalarField::Boundary& faceDiffusivityBf =
179  faceDiffusivity.boundaryFieldRef();
180 
181  forAll(faceDiffusivityBf, patchi)
182  {
183  fvsPatchScalarField& bfld = faceDiffusivityBf[patchi];
184 
185  if (patchSet.found(patchi))
186  {
187  const labelUList& faceCells = bfld.patch().faceCells();
188 
189  forAll(bfld, i)
190  {
191  const cell& ownFaces = mesh().cells()[faceCells[i]];
192 
193  labelHashSet cPoints(4*ownFaces.size());
194 
195  scalar dist = 0;
196 
197  forAll(ownFaces, ownFacei)
198  {
199  const face& f = mesh().faces()[ownFaces[ownFacei]];
200 
201  forAll(f, fp)
202  {
203  if (cPoints.insert(f[fp]))
204  {
205  dist += sqrt(pointWallDist[f[fp]].distSqr());
206  }
207  }
208  }
209  dist /= cPoints.size();
210 
211  bfld[i] = 1.0/dist;
212  }
213  }
214  else
215  {
216  const label start = bfld.patch().start();
217 
218  forAll(bfld, i)
219  {
220  const face& f = mesh().faces()[start+i];
221 
222  scalar dist = 0;
223 
224  forAll(f, fp)
225  {
226  dist += sqrt(pointWallDist[f[fp]].distSqr());
227  }
228  dist /= f.size();
229 
230  bfld[i] = 1.0/dist;
231  }
232  }
233  }
234 
235  return tfaceDiffusivity;
236 }
237 
238 
239 // ************************************************************************* //
#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
Macros for easy insertion into run-time selection tables.
Generic GeometricBoundaryField class.
Generic GeometricField class.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:89
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:132
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
const fvPatch & patch() const
Return patch.
Inverse distance to the given patches motion diffusivity.
virtual tmp< surfaceScalarField > operator()() const
Return diffusivity field.
inversePointDistanceDiffusivity(const fvMesh &mesh, Istream &mdData)
Construct for the given fvMesh and data Istream.
Abstract base class for cell-centre mesh motion diffusivity.
Class used to pass data into container.
Definition: pointEdgeDist.H:80
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
Definition: pointEdgeDist.H:66
Foam::polyBoundaryMesh.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch indices corresponding to the given names.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
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:181
label patchi
const pointField & points
label nPoints
bool valid(const PtrList< ModelType > &l)
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimless
dimensionedScalar sqrt(const dimensionedScalar &ds)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
labelList f(nPoints)
Foam::surfaceFields.