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-2022 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 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 void Foam::nearWallDist::resize()
41 {
42  y_.setSize(mesh().boundary().size());
43 
44  forAll(y_, patchi)
45  {
46  y_.set
47  (
48  patchi,
50  (
51  calculatedFvPatchScalarField::typeName,
52  mesh().boundary()[patchi],
54  )
55  );
56  }
57 }
58 
59 
60 void Foam::nearWallDist::correct()
61 {
62  volScalarField yVf(volScalarField::New("y", mesh(), dimLength));
63 
65  (
66  mesh(),
67  mesh().boundaryMesh().findPatchIDs<wallPolyPatch>(),
68  -vGreat,
69  2,
70  yVf
71  );
72 
73  forAll(y_, patchi)
74  {
75  const labelUList& faceCells = mesh().boundary()[patchi].faceCells();
76  forAll(y_[patchi], patchFacei)
77  {
78  y_[patchi][patchFacei] = yVf[faceCells[patchFacei]];
79  }
80  }
81 }
82 
83 
84 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
85 
87 :
89  <
90  fvMesh,
93  >(mesh),
94  y_
95  (
96  mesh.boundary(),
97  volScalarField::Internal::null(),
98  calculatedFvPatchScalarField::typeName
99  )
100 {
101  correct();
102 }
103 
104 
105 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
106 
108 {}
109 
110 
111 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
112 
114 {
115  resize();
116 
117  correct();
118  return true;
119 }
120 
121 
123 {
124  resize();
125  correct();
126 }
127 
128 
130 {
131  resize();
132  correct();
133 }
134 
135 
137 {
138  resize();
139  correct();
140 }
141 
142 
143 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Templated abstract base-class for demand-driven mesh objects used to automate their allocation to the...
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &)
Return a temporary field constructed from name,.
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
static tmp< fvPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
Distance calculation for cells with face on a wall. Searches pointNeighbours to find closest.
Definition: nearWallDist.H:55
virtual bool movePoints()
Update the y-field when the mesh moves.
Definition: nearWallDist.C:113
virtual void topoChange(const polyTopoChangeMap &)
Update the y-field when the mesh changes.
Definition: nearWallDist.C:122
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
Definition: nearWallDist.C:136
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
Definition: nearWallDist.C:129
virtual ~nearWallDist()
Destructor.
Definition: nearWallDist.C:107
nearWallDist(const fvMesh &mesh)
Construct from mesh.
Definition: nearWallDist.C:86
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
label patchi
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
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
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
defineTypeNameAndDebug(combustionModel, 0)
UList< label > labelUList
Definition: UList.H:65
faceListList boundary(nPatches)
triSurfaceToAgglom resize(surfacesMesh.size())