Loading...
Searching...
No Matches
cloudSet.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2025 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
29#include "cloudSet.H"
30#include "sampledSet.H"
31#include "meshSearch.H"
32#include "DynamicList.H"
33#include "polyMesh.H"
35#include "word.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40namespace Foam
41{
44}
45
46
47// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48
49void Foam::cloudSet::calcSamples
50(
51 DynamicList<point>& samplingPts,
52 DynamicList<label>& samplingCells,
53 DynamicList<label>& samplingFaces,
54 DynamicList<label>& samplingSegments,
55 DynamicList<scalar>& samplingCurveDist
56) const
57{
58 const meshSearch& queryMesh = searchEngine();
59
60 labelList foundProc(sampleCoords_.size(), -1);
61 forAll(sampleCoords_, sampleI)
62 {
63 label celli = queryMesh.findCell(sampleCoords_[sampleI]);
64
65 if (celli != -1)
66 {
67 samplingPts.append(sampleCoords_[sampleI]);
68 samplingCells.append(celli);
69 samplingFaces.append(-1);
70 samplingSegments.append(0);
71 samplingCurveDist.append(1.0 * sampleI);
72
73 foundProc[sampleI] = Pstream::myProcNo();
74 }
75 }
76
77 // Check that all have been found
78 labelList maxFoundProc(foundProc);
79 labelList minFoundProc(foundProc.size(), labelMax);
80 forAll(foundProc, i)
81 {
82 if (foundProc[i] != -1)
83 {
84 minFoundProc[i] = foundProc[i];
85 }
86 }
87 Pstream::listReduce(minFoundProc, minOp<label>());
88 Pstream::listReduce(maxFoundProc, maxOp<label>());
89
90
91 DynamicField<point> missingPoints(sampleCoords_.size());
92
93 forAll(sampleCoords_, sampleI)
94 {
95 if (maxFoundProc[sampleI] == -1)
96 {
97 // No processor has found the location.
98 missingPoints.append(sampleCoords_[sampleI]);
99 }
100 else if (minFoundProc[sampleI] != maxFoundProc[sampleI])
101 {
103 << "For sample set " << name()
104 << " location " << sampleCoords_[sampleI]
105 << " seems to be on multiple domains: "
106 << minFoundProc[sampleI] << " and " << maxFoundProc[sampleI]
107 << nl
108 << "This might happen if the location is on"
109 << " a processor patch. Change the location slightly"
110 << " to prevent this." << endl;
111 }
112 }
113
114
115 if (missingPoints.size() > 0)
116 {
117 if (missingPoints.size() < 100 || debug)
118 {
120 << "For sample set " << name()
121 << " did not found " << missingPoints.size()
122 << " points out of " << sampleCoords_.size()
123 << nl
124 << "Missing points:" << missingPoints << endl;
125 }
126 else
127 {
129 << "For sample set " << name()
130 << " did not found " << missingPoints.size()
131 << " points out of " << sampleCoords_.size()
132 << nl
133 << "Print missing points by setting the debug flag"
134 << " for " << cloudSet::typeName << endl;
135 }
136 }
137}
138
139
140void Foam::cloudSet::genSamples()
141{
142 // Storage for sample points
143 DynamicList<point> samplingPts;
144 DynamicList<label> samplingCells;
145 DynamicList<label> samplingFaces;
146 DynamicList<label> samplingSegments;
147 DynamicList<scalar> samplingCurveDist;
148
149 calcSamples
150 (
151 samplingPts,
152 samplingCells,
153 samplingFaces,
154 samplingSegments,
155 samplingCurveDist
156 );
157
158 samplingPts.shrink();
159 samplingCells.shrink();
160 samplingFaces.shrink();
161 samplingSegments.shrink();
162 samplingCurveDist.shrink();
163
164 // Move into *this
165 setSamples
166 (
167 std::move(samplingPts),
168 std::move(samplingCells),
169 std::move(samplingFaces),
170 std::move(samplingSegments),
171 std::move(samplingCurveDist)
172 );
173
174 if (debug)
175 {
177 }
178}
179
180
181// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
182
184(
185 const word& name,
186 const polyMesh& mesh,
187 const meshSearch& searchEngine,
188 const word& axis,
189 const List<point>& sampleCoords
190)
191:
193 sampleCoords_(sampleCoords)
194{
195 genSamples();
196}
197
198
200(
201 const word& name,
202 const polyMesh& mesh,
203 const meshSearch& searchEngine,
204 const dictionary& dict
205)
206:
207 sampledSet(name, mesh, searchEngine, dict),
208 sampleCoords_(dict.get<pointField>("points"))
209{
210 genSamples();
211}
212
213
214// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
static const char *const typeName
Typename for Field.
Definition Field.H:93
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
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
bool get(const label i) const
Definition UList.H:868
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
Samples at arbitrary locations with a volume mesh.
Definition cloudSet.H:82
cloudSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const word &axis, const List< point > &sampleCoords)
Construct from components.
Definition cloudSet.C:177
const word & axis() const
The sort axis name.
Definition coordSet.H:160
const word & name() const noexcept
The coord-set name.
Definition coordSet.H:152
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition meshSearch.H:57
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
Holds list of sampling points which is filled at construction time. Various implementations of this b...
Definition sampledSet.H:82
const meshSearch & searchEngine() const noexcept
Definition sampledSet.H:378
sampledSet(const word &name, const polyMesh &mesh, const meshSearch &searchEngine, const coordSet::coordFormat axisType)
Construct from components.
Definition sampledSet.C:405
const polyMesh & mesh() const noexcept
Definition sampledSet.H:373
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
auto & name
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
constexpr label labelMax
Definition label.H:55
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
runTime write()
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299