nearWallDist.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 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 "fvPatchDistWave.H"
28 #include "wallPolyPatch.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
35 }
36 
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
41 :
43  <
44  fvMesh,
47  >(mesh),
48  y_
49  (
50  mesh.boundary(),
51  volScalarField::Internal::null(),
52  calculatedFvPatchScalarField::typeName
53  )
54 {
56 
58  (
59  mesh,
61  -vGreat,
62  2,
63  yVf
64  );
65 
66  forAll(y_, patchi)
67  {
68  const labelUList& faceCells = mesh.boundary()[patchi].faceCells();
69  forAll(y_[patchi], patchFacei)
70  {
71  y_[patchi][patchFacei] = yVf[faceCells[patchFacei]];
72  }
73  }
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
78 
80 {}
81 
82 
83 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
MeshObject types:
Definition: MeshObjects.H:67
Templated abstract base-class for demand-driven mesh objects used to automate their allocation to the...
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:857
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
Definition: nearWallDist.H:58
virtual ~nearWallDist()
Destructor.
Definition: nearWallDist.C:79
nearWallDist(const fvMesh &mesh)
Construct from mesh.
Definition: nearWallDist.C:40
labelList findIndices(const wordRe &, const bool usePatchGroups=true) const
Return patch indices for all matches. Optionally matches patchGroups.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
Foam::wallPolyPatch.
Definition: wallPolyPatch.H:51
label patchi
void correct(const fvMesh &mesh, const labelHashSet &patchIDs, const scalar minFaceFraction, const label nCorrections, GeometricField< scalar, PatchField, GeoMesh > &distance)
Correct distance data from patches.
Namespace for OpenFOAM.
const dimensionSet dimLength
defineTypeNameAndDebug(combustionModel, 0)
faceListList boundary(nPatches)