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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "nearWallDist.H"
29 #include "fvMesh.H"
30 #include "cellDistFuncs.H"
31 #include "wallFvPatch.H"
32 #include "surfaceFields.H"
33 
34 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
35 
36 void Foam::nearWallDist::calculate()
37 {
38  cellDistFuncs wallUtils(mesh_);
39 
40  // Get patch ids of walls
41  labelHashSet wallPatchIDs(wallUtils.getPatchIDs<wallPolyPatch>());
42 
43  // Size neighbours array for maximum possible
44 
45  labelList neighbours(wallUtils.maxPatchSize(wallPatchIDs));
46 
47 
48  // Correct all cells with face on wall
49 
50  const volVectorField& cellCentres = mesh_.C();
51 
52  forAll(mesh_.boundary(), patchi)
53  {
54  fvPatchScalarField& ypatch = operator[](patchi);
55 
56  const fvPatch& patch = mesh_.boundary()[patchi];
57 
58  if (isA<wallFvPatch>(patch))
59  {
60  const polyPatch& pPatch = patch.patch();
61 
62  const labelUList& faceCells = patch.faceCells();
63 
64  // Check cells with face on wall
65  forAll(patch, patchFacei)
66  {
67  label nNeighbours = wallUtils.getPointNeighbours
68  (
69  pPatch,
70  patchFacei,
71  neighbours
72  );
73 
74  label minFacei = -1;
75 
76  ypatch[patchFacei] = wallUtils.smallestDist
77  (
78  cellCentres[faceCells[patchFacei]],
79  pPatch,
80  nNeighbours,
81  neighbours,
82  minFacei
83  );
84  }
85  }
86  else
87  {
88  ypatch = 0.0;
89  }
90  }
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
96 Foam::nearWallDist::nearWallDist(const Foam::fvMesh& mesh)
97 :
98  volScalarField::Boundary
99  (
100  mesh.boundary(),
101  mesh.V(), // Dummy internal field,
102  calculatedFvPatchScalarField::typeName
103  ),
104  mesh_(mesh)
105 {
106  calculate();
107 }
108 
109 
110 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
111 
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 {
120  if (mesh_.topoChanging())
121  {
122  const DimensionedField<scalar, volMesh>& V = mesh_.V();
123  const fvBoundaryMesh& bnd = mesh_.boundary();
124 
125  this->setSize(bnd.size());
126  forAll(*this, patchi)
127  {
128  this->set
129  (
130  patchi,
132  (
133  calculatedFvPatchScalarField::typeName,
134  bnd[patchi],
135  V
136  )
137  );
138  }
139  }
140 
141  calculate();
142 }
143 
144 
145 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:50
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:40
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
setSize
points setSize(newPointi)
wallFvPatch.H
surfaceFields.H
Foam::surfaceFields.
cellDistFuncs.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::nearWallDist::correct
virtual void correct()
Correct for mesh geom/topo changes.
Definition: nearWallDist.C:118
Foam::fvBoundaryMesh
Foam::fvBoundaryMesh.
Definition: fvBoundaryMesh.H:57
nearWallDist.H
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:62
Foam::nearWallDist::~nearWallDist
virtual ~nearWallDist()
Destructor.
Definition: nearWallDist.C:112
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:358
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:84
fvMesh.H
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:60
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:589
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:79
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
Definition: HashSet.H:415
Foam::GeometricField< scalar, fvPatchField, volMesh >
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54