Loading...
Searching...
No Matches
meshWaveAddressingPatchDistMethod.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) 2023-2025 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
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
29#include "wallDistAddressing.H"
30#include "volFields.H"
31#include "surfaceFields.H"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace patchDistMethods
39{
42}
43}
44
45
46// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47
49(
50 const dictionary& dict,
51 const fvMesh& mesh,
54:
56 correctWalls_(dict.getOrDefault("correctWalls", true))
57{}
58
59
61(
62 const fvMesh& mesh,
64 const bool correctWalls
65)
66:
68 correctWalls_(correctWalls)
69{}
70
71
72// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
75{
76 return correct(y, volVectorField::null().constCast());
77}
78
79
81(
84)
85{
86 // If not yet constructed, construct with our options.
87 // Note that registration name is still 'wallDistAddressing'
88 // since the supplied options (wall patchIDs and correctWalls)
89 // are not assumed to be different from the defaults of other usage.
90
91 const labelList patchIds(patchIDs_.sortedToc());
92
93 // Calculate distance and addressing
94 const auto& wDist = wallDistAddressing::New
95 (
96 "wall",
97 mesh_,
99 correctWalls_
100 );
101
102 y = wDist.y();
103
104 // Note: copying value only so might not be consistent with supplied
105 // patch types (e.g. zeroGradient when called from wallDist). Assume
106 // only affected ones are the supplied patches ...
107
108 y.boundaryFieldRef().evaluateSelected(patchIds);
109
110
111 // Only calculate n if the field is defined
112 if (notNull(n))
113 {
114 auto& bfld = n.boundaryFieldRef();
115
116 for (const label patchi : patchIds)
117 {
118 bfld[patchi] == bfld[patchi].patch().nf();
119 }
120
121 // No problem with inconsistency as for y (see above) since doing
122 // correctBoundaryConditions on actual n field.
123 wDist.map(n, mapDistribute::transform());
124 }
125
126 return true;
127}
128
129
130// ************************************************************************* //
scalar y
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
labelList patchIDs
static FOAM_NO_DANGLING_REFERENCE const wallDistAddressing & New(const fvMesh &mesh, Args &&... args)
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Default transformation behaviour.
Specialisation of patchDist for wall distance calculation.
const fvMesh & mesh_
Reference to the mesh.
patchDistMethod(const patchDistMethod &)=delete
No copy construct.
const labelHashSet patchIDs_
Set of patch IDs.
const labelHashSet & patchIDs() const
Return the patchIDs.
Variant of patchDistMethods::meshWave that stores nearest-wall addressing instead of directly transpo...
const bool correctWalls_
Do accurate distance calculation for near-wall cells.
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
meshWaveAddressing(const meshWaveAddressing &)=delete
No copy construct.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
thermo correct()
labelList patchIds
dynamicFvMesh & mesh
Namespace for OpenFOAM.
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
bool notNull(const T *ptr) noexcept
True if ptr is not a pointer (of type T) to the nullObject.
Definition nullObject.H:267
dictionary dict
Foam::surfaceFields.