Loading...
Searching...
No Matches
pointPatchDist.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) 2013-2016 OpenFOAM Foundation
9 Copyright (C) 2018 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
31#include "pointMesh.H"
32#include "PointEdgeWave.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36Foam::pointPatchDist::pointPatchDist
37(
38 const pointMesh& pMesh,
40 const pointField& points
41)
42:
44 (
46 (
47 "pointDistance",
48 pMesh.db().time().timeName(),
49 pMesh.db()
50 ),
51 pMesh,
52 dimensionedScalar("y", dimLength, GREAT)
53 ),
54 points_(points),
55 patchIDs_(patchIDs),
56 nUnset_(0)
58 correct();
59}
60
61
62// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63
65{
66 const pointBoundaryMesh& pbm = mesh().boundary();
67
68 label nPoints = 0;
69
70 for (const label patchi : patchIDs_)
71 {
72 nPoints += pbm[patchi].meshPoints().size();
73 }
74
75 externalPointEdgePoint::trackingData td(points_);
76
77 // Set initial changed points to all the patch points(if patch present)
78 List<externalPointEdgePoint> wallInfo(nPoints);
79 labelList wallPoints(nPoints);
80 nPoints = 0;
81
82 for (const label patchi : patchIDs_)
83 {
84 // Retrieve the patch now we have its index in patches.
85
86 const labelList& mp = pbm[patchi].meshPoints();
87
88 forAll(mp, ppI)
89 {
90 label meshPointi = mp[ppI];
91 wallPoints[nPoints] = meshPointi;
92 wallInfo[nPoints] = externalPointEdgePoint
93 (
94 td.points_[meshPointi],
95 0.0
96 );
97 nPoints++;
98 }
99 }
100
101 // Current info on points
102 List<externalPointEdgePoint> allPointInfo(mesh()().nPoints());
103
104 // Current info on edges
105 List<externalPointEdgePoint> allEdgeInfo(mesh()().nEdges());
106
107 PointEdgeWave
108 <
109 externalPointEdgePoint,
110 externalPointEdgePoint::trackingData
111 > wallCalc
112 (
113 mesh()(),
114 wallPoints,
115 wallInfo,
116
117 allPointInfo,
118 allEdgeInfo,
119 mesh().globalData().nTotalPoints(), // max iterations
120 td
121 );
122
123 pointScalarField& psf = *this;
124
125
126 forAll(allPointInfo, pointi)
127 {
128 if (allPointInfo[pointi].valid(td))
129 {
130 psf[pointi] = Foam::sqrt(allPointInfo[pointi].distSqr());
131 }
132 else
133 {
134 nUnset_++;
135 }
136 }
137}
138
139
140// ************************************************************************* //
labelList patchIDs
const polyBoundaryMesh & pbm
const Mesh & mesh() const noexcept
Return const reference to mesh.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const Time & time() const noexcept
Return Time associated with the objectRegistry.
Definition IOobject.C:456
const objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
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
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Class used to pass data into container.
Holds information regarding nearest wall point. Used in PointEdgeWave. (so not standard FaceCellWave)...
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
A pointBoundaryMesh is a pointPatch list with registered IO, a reference to the associated pointMesh,...
Mesh representing a set of points created from polyMesh.
Definition pointMesh.H:49
void correct()
Correct for mesh geom/topo changes.
For use with FaceCellWave. Determines topological distance to starting faces.
Definition wallPoints.H:60
thermo correct()
dynamicFvMesh & mesh
const pointField & points
label nPoints
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
word timeName
Definition getTimeIndex.H:3
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
dimensionedScalar sqrt(const dimensionedScalar &ds)
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299