1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by
13 the Free Software Foundation, either version 3 of the License, or
14 (at your option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "uniformSet.H"
27 #include "meshSearch.H"
28 #include "DynamicList.H"
31 #include "addToRunTimeSelectionTable.H"
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 defineTypeNameAndDebug(uniformSet, 0);
38 addToRunTimeSelectionTable(sampledSet, uniformSet, word);
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 bool Foam::uniformSet::nextSample
46 const point& currentPt,
48 const scalar smallDist,
53 bool pointFound = false;
55 const vector normOffset = offset/mag(offset);
60 for (; sampleI < nPoints_; sampleI++)
62 scalar s = (samplePt - currentPt) & normOffset;
66 // samplePt is close to or beyond currentPt -> use it
78 bool Foam::uniformSet::trackToBoundary
80 passiveParticle& singleParticle,
83 DynamicList<point>& samplingPts,
84 DynamicList<label>& samplingCells,
85 DynamicList<label>& samplingFaces,
86 DynamicList<scalar>& samplingCurveDist
89 // distance vector between sampling points
90 const vector offset = (end_ - start_)/(nPoints_ - 1);
91 const vector smallVec = tol*offset;
92 const scalar smallDist = mag(smallVec);
95 const point& trackPt = singleParticle.position();
97 passiveParticleCloud particleCloud(mesh());
98 particle::TrackingData<passiveParticleCloud> trackData(particleCloud);
102 // Find next samplePt on/after trackPt. Update samplePt, sampleI
103 if (!nextSample(trackPt, offset, smallDist, samplePt, sampleI))
108 Info<< "trackToBoundary : Reached end : samplePt now:"
109 << samplePt << " sampleI now:" << sampleI << endl;
114 if (mag(samplePt - trackPt) < smallDist)
116 // trackPt corresponds with samplePt. Store and use next samplePt
119 Info<< "trackToBoundary : samplePt corresponds to trackPt : "
120 << " trackPt:" << trackPt << " samplePt:" << samplePt
124 samplingPts.append(trackPt);
125 samplingCells.append(singleParticle.cell());
126 samplingFaces.append(-1);
127 samplingCurveDist.append(mag(trackPt - start_));
129 // go to next samplePt
130 if (!nextSample(trackPt, offset, smallDist, samplePt, sampleI))
135 Info<< "trackToBoundary : Reached end : "
136 << " samplePt now:" << samplePt
137 << " sampleI now:" << sampleI
148 Info<< "Searching along trajectory from "
149 << " trackPt:" << trackPt
150 << " trackCellI:" << singleParticle.cell()
151 << " to:" << samplePt << endl;
154 point oldPos = trackPt;
158 singleParticle.stepFraction() = 0;
159 singleParticle.track(samplePt, trackData);
163 Info<< "Result of tracking "
164 << " trackPt:" << trackPt
165 << " trackCellI:" << singleParticle.cell()
166 << " trackFaceI:" << singleParticle.face()
167 << " onBoundary:" << singleParticle.onBoundary()
168 << " samplePt:" << samplePt
169 << " smallDist:" << smallDist
175 !singleParticle.onBoundary()
176 && (mag(trackPt - oldPos) < smallDist)
179 if (singleParticle.onBoundary())
181 //Info<< "trackToBoundary : reached boundary" << endl;
182 if (mag(trackPt - samplePt) < smallDist)
184 //Info<< "trackToBoundary : boundary is also sampling point"
186 // Reached samplePt on boundary
187 samplingPts.append(trackPt);
188 samplingCells.append(singleParticle.cell());
189 samplingFaces.append(facei);
190 samplingCurveDist.append(mag(trackPt - start_));
196 //Info<< "trackToBoundary : reached internal sampling point" << endl;
197 // Reached samplePt in cell or on internal face
198 samplingPts.append(trackPt);
199 samplingCells.append(singleParticle.cell());
200 samplingFaces.append(-1);
201 samplingCurveDist.append(mag(trackPt - start_));
203 // go to next samplePt
208 void Foam::uniformSet::calcSamples
210 DynamicList<point>& samplingPts,
211 DynamicList<label>& samplingCells,
212 DynamicList<label>& samplingFaces,
213 DynamicList<label>& samplingSegments,
214 DynamicList<scalar>& samplingCurveDist
217 // distance vector between sampling points
218 if ((nPoints_ < 2) || (mag(end_ - start_) < SMALL))
220 FatalErrorIn("uniformSet::calcSamples()")
221 << "Incorrect sample specification. Either too few points or"
222 << " start equals end point." << endl
223 << "nPoints:" << nPoints_
224 << " start:" << start_
229 const vector offset = (end_ - start_)/(nPoints_ - 1);
230 const vector normOffset = offset/mag(offset);
231 const vector smallVec = tol*offset;
232 const scalar smallDist = mag(smallVec);
234 // Force calculation of minimum-tet decomposition.
235 (void) mesh().tetBasePtIs();
237 // Get all boundary intersections
238 List<pointIndexHit> bHits = searchEngine().intersections
244 point bPoint(GREAT, GREAT, GREAT);
249 bPoint = bHits[0].hitPoint();
250 bFaceI = bHits[0].index();
253 // Get first tracking point. Use bPoint, bFaceI if provided.
256 label trackCellI = -1;
257 label trackFaceI = -1;
272 if (trackCellI == -1)
274 // Line start_ - end_ does not intersect domain at all.
275 // (or is along edge)
276 // Set points and cell/face labels to empty lists
283 samplingPts.append(start_);
284 samplingCells.append(trackCellI);
285 samplingFaces.append(trackFaceI);
286 samplingCurveDist.append(0.0);
290 // Track until hit end of all boundary intersections
293 // current segment number
296 // starting index of current segment in samplePts
297 label startSegmentI = 0;
300 point samplePt = start_;
302 // index in bHits; current boundary intersection
307 // Initialize tracking starting from trackPt
308 passiveParticle singleParticle(mesh(), trackPt, trackCellI);
310 bool reachedBoundary = trackToBoundary
321 // fill sampleSegments
322 for (label i = samplingPts.size() - 1; i >= startSegmentI; --i)
324 samplingSegments.append(segmentI);
328 if (!reachedBoundary)
332 Info<< "calcSamples : Reached end of samples: "
333 << " samplePt now:" << samplePt
334 << " sampleI now:" << sampleI
341 bool foundValidB = false;
343 while (bHitI < bHits.size())
346 (bHits[bHitI].hitPoint() - singleParticle.position())
351 Info<< "Finding next boundary : "
352 << "bPoint:" << bHits[bHitI].hitPoint()
353 << " tracking:" << singleParticle.position()
358 if (dist > smallDist)
360 // hitpoint is past tracking position
372 // No valid boundary intersection found beyond tracking position
376 // Update starting point for tracking
378 trackPt = pushIn(bPoint, trackFaceI);
379 trackCellI = getBoundaryCell(trackFaceI);
383 startSegmentI = samplingPts.size();
388 void Foam::uniformSet::genSamples()
390 // Storage for sample points
391 DynamicList<point> samplingPts;
392 DynamicList<label> samplingCells;
393 DynamicList<label> samplingFaces;
394 DynamicList<label> samplingSegments;
395 DynamicList<scalar> samplingCurveDist;
406 samplingPts.shrink();
407 samplingCells.shrink();
408 samplingFaces.shrink();
409 samplingSegments.shrink();
410 samplingCurveDist.shrink();
422 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
424 Foam::uniformSet::uniformSet
427 const polyMesh& mesh,
428 meshSearch& searchEngine,
435 sampledSet(name, mesh, searchEngine, axis),
449 Foam::uniformSet::uniformSet
452 const polyMesh& mesh,
453 meshSearch& searchEngine,
454 const dictionary& dict
457 sampledSet(name, mesh, searchEngine, dict),
458 start_(dict.lookup("start")),
459 end_(dict.lookup("end")),
460 nPoints_(readLabel(dict.lookup("nPoints")))
471 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
473 Foam::uniformSet::~uniformSet()
477 // ************************************************************************* //