inverseFaceDistanceDiffusivity.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-2012 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 "wallPoint.H"
30 #include "MeshWave.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(inverseFaceDistanceDiffusivity, 0);
37 
39  (
40  motionDiffusivity,
41  inverseFaceDistanceDiffusivity,
42  Istream
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
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.size());
75 
76  label nPatchFaces = 0;
77 
78  forAll(patchNames_, i)
79  {
80  const label pID = bdry.findPatchID(patchNames_[i]);
81 
82  if (pID > -1)
83  {
84  patchSet.insert(pID);
85  nPatchFaces += bdry[pID].size();
86  }
87  }
88 
89  List<wallPoint> faceDist(nPatchFaces);
90  labelList changedFaces(nPatchFaces);
91 
92  nPatchFaces = 0;
93 
94  forAllConstIter(labelHashSet, patchSet, iter)
95  {
96  const polyPatch& patch = bdry[iter.key()];
97 
98  const vectorField::subField fc(patch.faceCentres());
99 
100  forAll(fc, patchFaceI)
101  {
102  changedFaces[nPatchFaces] = patch.start() + patchFaceI;
103 
104  faceDist[nPatchFaces] = wallPoint(fc[patchFaceI], 0);
105 
106  nPatchFaces++;
107  }
108  }
109  faceDist.setSize(nPatchFaces);
110  changedFaces.setSize(nPatchFaces);
111 
112  MeshWave<wallPoint> waveInfo
113  (
114  mesh(),
115  changedFaces,
116  faceDist,
117  mesh().globalData().nTotalCells()+1 // max iterations
118  );
119 
120  const List<wallPoint>& faceInfo = waveInfo.allFaceInfo();
121  const List<wallPoint>& cellInfo = waveInfo.allCellInfo();
122 
123  for (label faceI=0; faceI<mesh().nInternalFaces(); faceI++)
124  {
125  scalar dist = faceInfo[faceI].distSqr();
126 
127  faceDiffusivity_[faceI] = 1.0/sqrt(dist);
128  }
129 
130  forAll(faceDiffusivity_.boundaryField(), patchI)
131  {
132  fvsPatchScalarField& bfld = faceDiffusivity_.boundaryField()[patchI];
133 
134  const labelUList& faceCells = bfld.patch().faceCells();
135 
136  if (patchSet.found(patchI))
137  {
138  forAll(bfld, i)
139  {
140  scalar dist = cellInfo[faceCells[i]].distSqr();
141  bfld[i] = 1.0/sqrt(dist);
142  }
143  }
144  else
145  {
146  const label start = bfld.patch().start();
147 
148  forAll(bfld, i)
149  {
150  scalar dist = faceInfo[start+i].distSqr();
151  bfld[i] = 1.0/sqrt(dist);
152  }
153  }
154  }
155 }
156 
157 
158 // ************************************************************************* //
Foam::wallPoint
Holds information regarding nearest wall point. Used in wall distance calculation.
Definition: wallPoint.H:63
Foam::MeshWave::allFaceInfo
const List< Type > & allFaceInfo() const
Get allFaceInfo.
Definition: MeshWave.H:119
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName) const
Find patch index given a name.
Definition: polyBoundaryMesh.C:678
wallPoint.H
Foam::inverseFaceDistanceDiffusivity::inverseFaceDistanceDiffusivity
inverseFaceDistanceDiffusivity(const inverseFaceDistanceDiffusivity &)
Disallow default bitwise copy construct.
Foam::inverseFaceDistanceDiffusivity::correct
virtual void correct()
Correct the motion diffusivity.
Definition: inverseFaceDistanceDiffusivity.C:70
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::uniformDiffusivity
Uniform uniform finite volume mesh motion diffusivity.
Definition: uniformDiffusivity.H:49
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Foam::inverseFaceDistanceDiffusivity::~inverseFaceDistanceDiffusivity
virtual ~inverseFaceDistanceDiffusivity()
Destructor.
Definition: inverseFaceDistanceDiffusivity.C:64
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::HashSet< label, Hash< label > >
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::SubField
Pre-declare related SubField type.
Definition: Field.H:61
Foam::label
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
Foam::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
MeshWave.H
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::MeshWave
FaceCellWave plus data.
Definition: MeshWave.H:56
correct
fvOptions correct(rho)
HashSet.H
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
inverseFaceDistanceDiffusivity.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::cellInfo
Holds information regarding type of cell. Used in inside/outside determination in cellClassification.
Definition: cellInfo.H:54
Foam::fvsPatchField::patch
const fvPatch & patch() const
Return patch.
Definition: fvsPatchField.H:278
Foam::polyPatch::faceCentres
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:308
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::UList
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
Foam::fvPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::MeshWave::allCellInfo
const List< Type > & allCellInfo() const
Get allCellInfo.
Definition: MeshWave.H:125