patchEdgeSet.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) 2019-2020 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
28#include "patchEdgeSet.H"
29#include "polyMesh.H"
31#include "searchableSurface.H"
32#include "Time.H"
33#include "mergePoints.H"
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
41}
42
43// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44
45void Foam::patchEdgeSet::genSamples()
46{
47 // Storage for sample points
48 DynamicList<point> samplingPts;
49 DynamicList<label> samplingCells;
50 DynamicList<label> samplingFaces;
51 DynamicList<label> samplingSegments;
52 DynamicList<scalar> samplingCurveDist;
53
54 const labelList patchIDs(patchSet_.sortedToc());
55
56 for (const label patchi : patchIDs)
57 {
58 const polyPatch& pp = mesh().boundaryMesh()[patchi];
59 const edgeList& edges = pp.edges();
60 const labelList& mp = pp.meshPoints();
61 const pointField& pts = pp.points();
62
63 pointField start(edges.size());
64 pointField end(edges.size());
65 forAll(edges, edgei)
66 {
67 const edge& e = edges[edgei];
68 start[edgei] = pts[mp[e[0]]];
69 end[edgei] = pts[mp[e[1]]];
70 }
71
72
73 List<List<pointIndexHit>> hits;
74 surfPtr_->findLineAll(start, end, hits);
75
76 forAll(hits, edgei)
77 {
78 const List<pointIndexHit>& eHits = hits[edgei];
79
80 if (eHits.size())
81 {
82 const label meshFacei = pp.start()+pp.edgeFaces()[edgei][0];
83 const label celli = mesh().faceOwner()[meshFacei];
84 for (const auto& hit : eHits)
85 {
86 const point& pt = hit.hitPoint();
87 samplingPts.append(pt);
88 samplingCells.append(celli);
89 samplingFaces.append(meshFacei);
90 samplingSegments.append(0);
91 samplingCurveDist.append(mag(pt-origin_));
92 }
93 }
94 }
95
97 //forAll(edges, edgei)
98 //{
99 // const edge& e = edges[edgei];
100 // const point& p0 = pts[mp[e[0]]];
101 // const point& p1 = pts[mp[e[1]]];
102 // const vector dir(p1-p0);
103 // const scalar s = pl_.normalIntersect(p0, dir);
104 //
105 // if (s >= 0.0 && s < 1.0)
106 // {
107 // const point pt(p0+s*dir);
108 //
109 // // Calculate distance on curve
110 // //const scalar dist = (pt-start_)&curve;
111 // //if (dist >= 0 && dist < magCurve)
112 // const scalar dist = mag(pt-pl_.origin());
113 //
114 // // Take any face using this edge
115 // const label meshFacei = pp.start()+pp.edgeFaces()[edgei][0];
116 //
117 // samplingPts.append(pt);
118 // samplingCells.append(mesh().faceOwner()[meshFacei]);
119 // samplingFaces.append(meshFacei);
120 // samplingSegments.append(0);
121 // samplingCurveDist.append(dist);
122 // }
123 //}
124 }
125
126 samplingPts.shrink();
127 samplingCells.shrink();
128 samplingFaces.shrink();
129 samplingSegments.shrink();
130 samplingCurveDist.shrink();
131
132
133 labelList pointMap;
134 const label nMerged = Foam::mergePoints
135 (
136 samplingPts,
137 SMALL, // mergeTol
138 false, // verbose = false
139 pointMap
140 );
141
142 if (nMerged == samplingPts.size())
143 {
144 // Nothing merged
146 (
147 samplingPts,
148 samplingCells,
149 samplingFaces,
150 samplingSegments,
151 samplingCurveDist
152 );
153 }
154 else
155 {
156 // Compress out duplicates
157
158 List<point> newSamplingPts(nMerged);
159 List<label> newSamplingCells(nMerged);
160 List<label> newSamplingFaces(nMerged);
161 List<label> newSamplingSegments(nMerged);
162 List<scalar> newSamplingCurveDist(nMerged);
163
164 forAll(pointMap, i)
165 {
166 const label newi = pointMap[i];
167 newSamplingPts[newi] = samplingPts[i];
168 newSamplingCells[newi] = samplingCells[i];
169 newSamplingFaces[newi] = samplingFaces[i];
170 newSamplingSegments[newi] = samplingSegments[i];
171 newSamplingCurveDist[newi] = samplingCurveDist[i];
172 }
173
175 (
176 newSamplingPts,
177 newSamplingCells,
178 newSamplingFaces,
179 newSamplingSegments,
180 newSamplingCurveDist
181 );
182 }
183
184 if (debug)
185 {
186 write(Info);
187 }
188}
189
190
191// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
192
194(
195 const word& name,
196 const polyMesh& mesh,
197 const meshSearch& searchEngine,
198 const dictionary& dict
199)
200:
201 sampledSet(name, mesh, searchEngine, dict),
202 surfPtr_
203 (
205 (
206 dict.get<word>("surfaceType"),
208 (
209 dict.getOrDefault("surfaceName", name),
210 mesh.time().constant(), // directory
211 "triSurface", // instance
212 mesh.time(), // registry
213 IOobject::MUST_READ,
214 IOobject::NO_WRITE
215 ),
216 dict
217 )
218 ),
219 origin_(dict.get<point>("origin")),
220 //end_(dict.get<point>("end")),
221 //pl_(dict),
222 patchSet_
223 (
224 mesh.boundaryMesh().patchSet(dict.get<wordRes>("patches"))
225 )
226{
227 genSamples();
228}
229
230
231// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
List< Key > sortedToc() const
The table of contents (the keys) in sorted order.
Definition: HashTable.C:137
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition: UListI.H:350
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
Definition: boundaryMesh.H:63
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:61
Like Foam::uniformSet but samples patch edges.
Definition: patchEdgeSet.H:162
constant condensation/saturation model.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:86
void setSamples(const List< point > &samplingPts, const labelList &samplingCells, const labelList &samplingFaces, const labelList &samplingSegments, const scalarList &samplingDistance)
Set sample data. Copy list contents.
Definition: sampledSet.C:379
const polyMesh & mesh() const noexcept
Definition: sampledSet.H:319
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
dynamicFvMesh & mesh
Geometric merging of points. See below.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: List.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
label mergePoints(const PointList &points, labelList &pointToUnique, labelList &uniquePoints, const scalar mergeTol=SMALL, const bool verbose=false)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
runTime write()
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333