Loading...
Searching...
No Matches
exactPatchDistMethod.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) 2018-2024 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
31#include "volFields.H"
32#include "DynamicField.H"
33#include "OBJstream.H"
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace patchDistMethods
40{
43}
44}
45
47Foam::patchDistMethods::exact::patchSurface() const
48{
49 if (!patchSurfPtr_)
50 {
52
53 Random rndGen(0);
54
55 treeBoundBox localBb(mesh_.points());
56
57 // Determine mesh bounding boxes:
59 (
60 1,
61 localBb.extend(rndGen, 1E-3)
62 );
63
64 // Add any missing properties (but not override existing ones)
65 dict_.add("bounds", meshBb);
66 dict_.add
67 (
68 "distributionType",
70 [
71 //distributedTriSurfaceMesh::FOLLOW // use mesh bb
72 //distributedTriSurfaceMesh::INDEPENDENT // master-only
74 ]
75 );
76 dict_.add("mergeDistance", 1e-6*localBb.mag());
77
78 Info<< "Triangulating local patch faces" << nl << endl;
79
80 labelList mapTriToGlobal;
81
82 patchSurfPtr_.reset
83 (
85 (
87 (
88 "wallSurface.stl",
89 mesh_.time().constant(),// directory
90 "triSurface", // instance
91 mesh_.time(), // registry
94 ),
96 (
97 pbm,
99 mapTriToGlobal
100 ),
101 dict_
102 )
103 );
104
105 // Do redistribution
106 Info<< "Redistributing surface" << nl << endl;
108 autoPtr<mapDistribute> pointMap;
109 patchSurfPtr_().distribute
110 (
111 meshBb,
112 false, //keepNonMapped,
113 faceMap,
114 pointMap
115 );
116 faceMap.clear();
117 pointMap.clear();
119 return patchSurfPtr_();
120}
121
122
123// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
124
125Foam::patchDistMethods::exact::exact
126(
127 const dictionary& dict,
128 const fvMesh& mesh,
131:
133 dict_(dict.subOrEmptyDict(typeName + "Coeffs"))
134{}
135
136
137Foam::patchDistMethods::exact::exact
138(
139 const fvMesh& mesh,
141)
144{}
145
146
147// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
150{
151 return correct(y, volVectorField::null().constCast());
152}
153
154
156(
159)
160{
161 const distributedTriSurfaceMesh& surf = patchSurface();
162
164 surf.findNearest
165 (
166 mesh_.cellCentres(),
167 scalarField(mesh_.nCells(), Foam::sqr(GREAT)),
168 info
169 );
170
171 // Take over hits
172 // label nHits = 0;
173 forAll(info, celli)
174 {
175 if (info[celli].hit())
176 {
177 const point& cc = mesh_.cellCentres()[celli];
178 y[celli] = info[celli].point().dist(cc);
179 // ++nHits;
180 }
181 //else
182 //{
183 // // Miss. Do what? Not possible with GREAT hopefully ...
184 //}
185 }
186 y.correctBoundaryConditions();
187
188 if (debug)
189 {
190 OBJstream str(mesh_.time().timePath()/"wallPoint.obj");
191 Info<< type() << ": dumping nearest wall point to "
192 << str.name() << endl;
193 forAll(mesh_.cellCentres(), celli)
194 {
195 const point& cc = mesh_.cellCentres()[celli];
196 str.writeLine(cc, info[celli].point());
197 }
198 }
199
200
201 // Only calculate n if the field is defined
202 if (notNull(n))
203 {
204 surf.getNormal(info, n.primitiveFieldRef());
205 n.correctBoundaryConditions();
206 }
207
208 return true;
209}
210
211
212// ************************************************************************* //
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
const polyBoundaryMesh & pbm
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
Ostream & writeLine(const point &p0, const point &p1)
Write line joining two points.
Definition OBJstream.C:234
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
Random number generator.
Definition Random.H:56
const word & constant() const noexcept
Return constant name.
Definition TimePathsI.H:131
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
void clear() noexcept
Same as reset(nullptr).
Definition autoPtr.H:255
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
entry * add(entry *entryPtr, bool mergeEntry=false)
Add a new entry.
Definition dictionary.C:625
IOoject and searching on distributed triSurface. All processor hold (possibly overlapping) part of th...
static const Enum< distributionType > distributionTypeNames_
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
virtual void findNearest(const pointField &sample, const scalarField &nearestDistSqr, List< pointIndexHit > &) const
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
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.
Calculation of exact distance to nearest patch for all cells and boundary by constructing a search tr...
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
Standard boundBox with extra functionality for use in octree.
static triSurface triangulate(const polyBoundaryMesh &mBesh, const labelHashSet &includePatches, labelList &faceMap, const bool verbose=false)
Simple triangulation of (selected patches of) boundaryMesh. Needs.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
thermo correct()
dynamicFvMesh & mesh
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
GeometricField< vector, fvPatchField, volMesh > volVectorField
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
vector point
Point is a vector.
Definition point.H:37
bool notNull(const T *ptr) noexcept
True if ptr is not a pointer (of type T) to the nullObject.
Definition nullObject.H:267
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Random rndGen