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