nearWallDist.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-2015 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 "nearWallDist.H"
27 #include "fvMesh.H"
28 #include "cellDistFuncs.H"
29 #include "wallFvPatch.H"
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
33 
35 {
36  cellDistFuncs wallUtils(mesh_);
37 
38  // Get patch ids of walls
39  labelHashSet wallPatchIDs(wallUtils.getPatchIDs<wallPolyPatch>());
40 
41  // Size neighbours array for maximum possible
42 
43  labelList neighbours(wallUtils.maxPatchSize(wallPatchIDs));
44 
45 
46  // Correct all cells with face on wall
47 
48  const volVectorField& cellCentres = mesh_.C();
49 
50  forAll(mesh_.boundary(), patchI)
51  {
52  fvPatchScalarField& ypatch = operator[](patchI);
53 
54  const fvPatch& patch = mesh_.boundary()[patchI];
55 
56  if (isA<wallFvPatch>(patch))
57  {
58  const polyPatch& pPatch = patch.patch();
59 
60  const labelUList& faceCells = patch.faceCells();
61 
62  // Check cells with face on wall
63  forAll(patch, patchFaceI)
64  {
65  label nNeighbours = wallUtils.getPointNeighbours
66  (
67  pPatch,
68  patchFaceI,
69  neighbours
70  );
71 
72  label minFaceI = -1;
73 
74  ypatch[patchFaceI] = wallUtils.smallestDist
75  (
76  cellCentres[faceCells[patchFaceI]],
77  pPatch,
78  nNeighbours,
79  neighbours,
80  minFaceI
81  );
82  }
83  }
84  else
85  {
86  ypatch = 0.0;
87  }
88  }
89 }
90 
91 
92 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
93 
95 :
96  volScalarField::GeometricBoundaryField
97  (
98  mesh.boundary(),
99  mesh.V(), // Dummy internal field,
100  calculatedFvPatchScalarField::typeName
101  ),
102  mesh_(mesh)
103 {
104  calculate();
105 }
106 
107 
108 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
109 
111 {}
112 
113 
114 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
115 
117 {
118  if (mesh_.topoChanging())
119  {
120  const DimensionedField<scalar, volMesh>& V = mesh_.V();
121  const fvBoundaryMesh& bnd = mesh_.boundary();
122 
123  this->setSize(bnd.size());
124  forAll(*this, patchI)
125  {
126  this->set
127  (
128  patchI,
130  (
131  calculatedFvPatchScalarField::typeName,
132  bnd[patchI],
133  V
134  )
135  );
136  }
137  }
138 
139  calculate();
140 }
141 
142 
143 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
setSize
points setSize(newPointi)
Foam::nearWallDist::nearWallDist
nearWallDist(const nearWallDist &)
Disallow default bitwise copy construct.
Foam::nearWallDist::mesh_
const fvMesh & mesh_
Reference to mesh.
Definition: nearWallDist.H:58
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::cellDistFuncs::maxPatchSize
label maxPatchSize(const labelHashSet &patchIDs) const
Size of largest patch (out of supplied subset of patches)
Definition: cellDistFuncs.C:229
boundary
faceListList boundary(nPatches)
wallFvPatch.H
surfaceFields.H
Foam::surfaceFields.
cellDistFuncs.H
Foam::HashSet< label, Hash< label > >
Foam::cellDistFuncs::getPatchIDs
labelHashSet getPatchIDs(const wordReList &patchNames) const
Return the set of patch IDs corresponding to the given names.
Definition: cellDistFuncs.C:71
Foam::nearWallDist::correct
virtual void correct()
Correct for mesh geom/topo changes.
Definition: nearWallDist.C:116
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:52
nearWallDist.H
Foam::nearWallDist::calculate
void calculate()
Do all calculations.
Definition: nearWallDist.C:34
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::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::nearWallDist::~nearWallDist
virtual ~nearWallDist()
Destructor.
Definition: nearWallDist.C:110
Foam::cellDistFuncs::getPointNeighbours
label getPointNeighbours(const primitivePatch &, const label patchFaceI, labelList &) const
Get faces sharing point with face on patch.
Definition: cellDistFuncs.C:117
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::wallPolyPatch
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:48
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam::cellDistFuncs::smallestDist
scalar smallestDist(const point &p, const polyPatch &patch, const label nWallFaces, const labelList &wallFaces, label &meshFaceI) const
Calculate smallest true distance (and face index)
Definition: cellDistFuncs.C:83
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:550
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:369
Foam::cellDistFuncs
Collection of functions used in wall distance calculation.
Definition: cellDistFuncs.H:61
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::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::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51