inversePointDistanceDiffusivity.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 | Copyright (C) 2011-2016 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 
28 #include "HashSet.H"
29 #include "pointEdgePoint.H"
30 #include "PointEdgeWave.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(inversePointDistanceDiffusivity, 0);
37 
39  (
40  motionDiffusivity,
41  inversePointDistanceDiffusivity,
42  Istream
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 Foam::inversePointDistanceDiffusivity::inversePointDistanceDiffusivity
50 (
51  const fvMesh& mesh,
52  Istream& mdData
53 )
54 :
55  uniformDiffusivity(mesh, mdData),
56  patchNames_(mdData)
57 {
58  correct();
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
71 {
72  const polyBoundaryMesh& bdry = mesh().boundaryMesh();
73 
74  labelHashSet patchSet(bdry.patchSet(patchNames_));
75 
76  label nPatchEdges = 0;
77 
78  forAllConstIter(labelHashSet, patchSet, iter)
79  {
80  nPatchEdges += bdry[iter.key()].nEdges();
81  }
82 
83  // Distance to wall on points and edges.
84  List<pointEdgePoint> pointWallDist(mesh().nPoints());
85  List<pointEdgePoint> edgeWallDist(mesh().nEdges());
86 
87  int dummyTrackData = 0;
88 
89 
90  {
91  // Seeds
92  List<pointEdgePoint> seedInfo(nPatchEdges);
93  labelList seedPoints(nPatchEdges);
94 
95  nPatchEdges = 0;
96 
97  forAllConstIter(labelHashSet, patchSet, iter)
98  {
99  const polyPatch& patch = bdry[iter.key()];
100 
101  const labelList& meshPoints = patch.meshPoints();
102 
103  forAll(meshPoints, i)
104  {
105  label pointi = meshPoints[i];
106 
107  if (!pointWallDist[pointi].valid(dummyTrackData))
108  {
109  // Not yet seeded
110  seedInfo[nPatchEdges] = pointEdgePoint
111  (
112  mesh().points()[pointi],
113  0.0
114  );
115  seedPoints[nPatchEdges] = pointi;
116  pointWallDist[pointi] = seedInfo[nPatchEdges];
117 
118  nPatchEdges++;
119  }
120  }
121  }
122  seedInfo.setSize(nPatchEdges);
123  seedPoints.setSize(nPatchEdges);
124 
125  // Do calculations
127  (
128  mesh(),
129  seedPoints,
130  seedInfo,
131 
132  pointWallDist,
133  edgeWallDist,
134  mesh().globalData().nTotalPoints(),// max iterations
135  dummyTrackData
136  );
137  }
138 
139 
140  for (label facei=0; facei<mesh().nInternalFaces(); facei++)
141  {
142  const face& f = mesh().faces()[facei];
143 
144  scalar dist = 0;
145 
146  forAll(f, fp)
147  {
148  dist += sqrt(pointWallDist[f[fp]].distSqr());
149  }
150  dist /= f.size();
151 
152  faceDiffusivity_[facei] = 1.0/dist;
153  }
154 
155 
156  surfaceScalarField::Boundary& faceDiffusivityBf =
157  faceDiffusivity_.boundaryFieldRef();
158 
159  forAll(faceDiffusivityBf, patchi)
160  {
161  fvsPatchScalarField& bfld = faceDiffusivityBf[patchi];
162 
163  if (patchSet.found(patchi))
164  {
165  const labelUList& faceCells = bfld.patch().faceCells();
166 
167  forAll(bfld, i)
168  {
169  const cell& ownFaces = mesh().cells()[faceCells[i]];
170 
171  labelHashSet cPoints(4*ownFaces.size());
172 
173  scalar dist = 0;
174 
175  forAll(ownFaces, ownFacei)
176  {
177  const face& f = mesh().faces()[ownFaces[ownFacei]];
178 
179  forAll(f, fp)
180  {
181  if (cPoints.insert(f[fp]))
182  {
183  dist += sqrt(pointWallDist[f[fp]].distSqr());
184  }
185  }
186  }
187  dist /= cPoints.size();
188 
189  bfld[i] = 1.0/dist;
190  }
191  }
192  else
193  {
194  const label start = bfld.patch().start();
195 
196  forAll(bfld, i)
197  {
198  const face& f = mesh().faces()[start+i];
199 
200  scalar dist = 0;
201 
202  forAll(f, fp)
203  {
204  dist += sqrt(pointWallDist[f[fp]].distSqr());
205  }
206  dist /= f.size();
207 
208  bfld[i] = 1.0/dist;
209  }
210  }
211  }
212 }
213 
214 
215 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
Uniform uniform finite volume mesh motion diffusivity.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
label nInternalFaces() const
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
virtual void correct()
Correct the motion diffusivity.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const cellList & cells() const
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
Macros for easy insertion into run-time selection tables.
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
dynamicFvMesh & mesh
const pointField & points
label nPoints
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1036
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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:61
Foam::polyBoundaryMesh.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-5.0/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
const fvPatch & patch() const
Return patch.
void setSize(const label)
Reset size of List.
Definition: List.C:281
label patchi
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:85
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.