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 -------------------------------------------------------------------------------
10 License
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 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(patchEdgeSet, 0);
40  addToRunTimeSelectionTable(sampledSet, patchEdgeSet, word);
41 }
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 void 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 = mergePoints
135  (
136  samplingPts,
137  SMALL, //const scalar mergeTol,
138  false, //const bool verbose,
139  pointMap,
140  origin_
141  );
142 
143  if (nMerged == samplingPts.size())
144  {
145  // Nothing merged
146  setSamples
147  (
148  samplingPts,
149  samplingCells,
150  samplingFaces,
151  samplingSegments,
152  samplingCurveDist
153  );
154  }
155  else
156  {
157  // Compress out duplicates
158 
159  List<point> newSamplingPts(nMerged);
160  List<label> newSamplingCells(nMerged);
161  List<label> newSamplingFaces(nMerged);
162  List<label> newSamplingSegments(nMerged);
163  List<scalar> newSamplingCurveDist(nMerged);
164 
165  forAll(pointMap, i)
166  {
167  const label newi = pointMap[i];
168  newSamplingPts[newi] = samplingPts[i];
169  newSamplingCells[newi] = samplingCells[i];
170  newSamplingFaces[newi] = samplingFaces[i];
171  newSamplingSegments[newi] = samplingSegments[i];
172  newSamplingCurveDist[newi] = samplingCurveDist[i];
173  }
174 
175  setSamples
176  (
177  newSamplingPts,
178  newSamplingCells,
179  newSamplingFaces,
180  newSamplingSegments,
181  newSamplingCurveDist
182  );
183  }
184 
185  if (debug)
186  {
187  write(Info);
188  }
189 }
190 
191 
192 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
193 
195 (
196  const word& name,
197  const polyMesh& mesh,
198  const meshSearch& searchEngine,
199  const dictionary& dict
200 )
201 :
202  sampledSet(name, mesh, searchEngine, dict),
203  surfPtr_
204  (
206  (
207  dict.get<word>("surfaceType"),
208  IOobject
209  (
210  dict.getOrDefault("surfaceName", name),
211  mesh.time().constant(), // directory
212  "triSurface", // instance
213  mesh.time(), // registry
216  ),
217  dict
218  )
219  ),
220  origin_(dict.get<point>("origin")),
221  //end_(dict.get<point>("end")),
222  //pl_(dict),
223  patchSet_
224  (
225  mesh.boundaryMesh().patchSet(dict.get<wordRes>("patches"))
226  )
227 {
228  genSamples();
229 }
230 
231 
232 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::sampledSet
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition: sampledSet.H:83
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::constant::atomic::mp
const dimensionedScalar mp
Proton mass.
Foam::sampledSet::write
Ostream & write(Ostream &) const
Output for debugging.
Definition: sampledSet.C:554
Foam::edgeList
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::meshSearch
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:60
searchableSurface.H
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
patchEdgeSet.H
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
stdFoam::end
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:121
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::searchableSurface::New
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
Definition: searchableSurface.C:43
Foam::mergePoints
label mergePoints(const PointList &points, const scalar mergeTol, const bool verbose, labelList &pointMap, typename PointList::const_reference origin=PointList::value_type::zero)
Sorts and merges points. All points closer than/equal mergeTol get merged.
Time.H
Foam::sampledSet::mesh
const polyMesh & mesh() const
Definition: sampledSet.H:281
Foam::sampledSet::setSamples
void setSamples(const List< point > &samplingPts, const labelList &samplingCells, const labelList &samplingFaces, const labelList &samplingSegments, const scalarList &samplingCurveDist)
Set sample data. Copy list contents.
Definition: sampledSet.C:382
Foam::Vector< scalar >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::polyBoundaryMesh::patchSet
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
Definition: polyBoundaryMesh.C:864
mergePoints.H
Merge points. See below.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::patchEdgeSet::patchEdgeSet
patchEdgeSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const dictionary &dict)
Construct from dictionary.
Definition: patchEdgeSet.C:195
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::IOobject::MUST_READ
Definition: IOobject.H:185