Loading...
Searching...
No Matches
mapLagrangian.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-2017 OpenFOAM Foundation
9 Copyright (C) 2018-2022 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 "MapLagrangianFields.H"
31#include "meshSearch.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37
38static const scalar perturbFactor = 1e-6;
39
40
41// Special version of findCell that generates a cell guaranteed to be
42// compatible with tracking.
43static label findCell(const Cloud<passiveParticle>& cloud, const point& pt)
44{
45 label celli = -1;
46 label tetFacei = -1;
47 label tetPtI = -1;
48
49 const polyMesh& mesh = cloud.pMesh();
50
51 mesh.findCellFacePt(pt, celli, tetFacei, tetPtI);
52
53 if (celli >= 0)
54 {
55 return celli;
56 }
57 else
58 {
59 // See if particle on face by finding nearest face and shifting
60 // particle.
61
62 meshSearch meshSearcher
63 (
64 mesh,
65 polyMesh::FACE_PLANES // no decomposition needed
66 );
67
68 label facei = meshSearcher.findNearestBoundaryFace(pt);
69
70 if (facei >= 0)
71 {
72 const point& cc = mesh.cellCentres()[mesh.faceOwner()[facei]];
73
74 const point perturbPt = (1-perturbFactor)*pt+perturbFactor*cc;
75
76 mesh.findCellFacePt(perturbPt, celli, tetFacei, tetPtI);
77
78 return celli;
79 }
80 }
81
82 return -1;
83}
84
85
86void mapLagrangian(const meshToMesh& interp)
87{
88 // Determine which particles are in meshTarget
89 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
90
91 const polyMesh& meshSource = interp.srcRegion();
92 const polyMesh& meshTarget = interp.tgtRegion();
93 const labelListList& sourceToTarget = interp.srcToTgtCellAddr();
94
95 const fileNameList cloudDirs
96 (
98 (
99 meshSource.time().timePath()/cloud::prefix,
101 )
102 );
103
104 for (const fileName& cloudDir : cloudDirs)
105 {
106 // Search for list of lagrangian objects for this time
107 IOobjectList objects
108 (
109 meshSource,
110 meshSource.time().timeName(),
111 cloud::prefix/cloudDir
112 );
113
114 if
115 (
117 (
118 objects.found("coordinates") || objects.found("positions")
119 )
120 )
121 {
122 // Has coordinates/positions - so must be a valid cloud
123 Info<< nl << " processing cloud " << cloudDir << endl;
124
125 // Read positions & cell
126 passiveParticleCloud sourceParcels
127 (
128 meshSource,
129 cloudDir,
130 false
131 );
132
133 Info<< " read " << sourceParcels.size()
134 << " parcels from source mesh." << endl;
135
136 // Construct empty target cloud
137 passiveParticleCloud targetParcels
138 (
139 meshTarget,
140 Foam::zero{},
141 cloudDir
142 );
143
144 passiveParticle::trackingData td(targetParcels);
145
146 label sourceParticleI = 0;
147
148 // Indices of source particles that get added to targetParcels
149 DynamicList<label> addParticles(sourceParcels.size());
150
151 // Unmapped particles
152 labelHashSet unmappedSource(sourceParcels.size());
153
154
155 // Initial: track from fine-mesh cell centre to particle position
156 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
157 // This requires there to be no boundary in the way.
158
159
160 for (const passiveParticle& p : sourceParcels)
161 {
162 bool foundCell = false;
163
164 // Assume that cell from read parcel is the correct one...
165 if (p.cell() >= 0)
166 {
167 const labelList& targetCells =
168 sourceToTarget[p.cell()];
169
170 // Particle probably in one of the targetcells. Try
171 // all by tracking from their cell centre to the parcel
172 // position.
173
174 for (const label targetCell : targetCells)
175 {
176 // Track from its cellcentre to position to make sure.
178 (
180 (
181 meshTarget,
182 barycentric(1, 0, 0, 0),
183 targetCell,
184 meshTarget.cells()[targetCell][0],
185 1
186 )
187 );
188 passiveParticle& newP = newPtr();
189
190 newP.track(p.position() - newP.position(), 0);
191
192 if (!newP.onFace())
193 {
194 // Hit position.
195 foundCell = true;
196 addParticles.append(sourceParticleI);
197 targetParcels.addParticle(newPtr.ptr());
198 break;
199 }
200 }
201 }
202
203 if (!foundCell)
204 {
205 // Store for closer analysis
206 unmappedSource.insert(sourceParticleI);
207 }
208
209 sourceParticleI++;
210 }
211
212 Info<< " after meshToMesh addressing found "
213 << targetParcels.size()
214 << " parcels in target mesh." << endl;
215
216
217 // Do closer inspection for unmapped particles
218 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
219
220 if (unmappedSource.size())
221 {
222 sourceParticleI = 0;
223
224 for (passiveParticle& p : sourceParcels)
225 {
226 if (unmappedSource.found(sourceParticleI))
227 {
228 const label targetCell =
229 findCell(targetParcels, p.position());
230
231 if (targetCell >= 0)
232 {
233 unmappedSource.erase(sourceParticleI);
234 addParticles.append(sourceParticleI);
235 targetParcels.addParticle
236 (
238 (
239 meshTarget,
240 p.position(),
241 targetCell
242 )
243 );
244 sourceParcels.remove(&p);
245 }
246 }
247 sourceParticleI++;
248 }
249 }
250 addParticles.shrink();
251
252 Info<< " after additional mesh searching found "
253 << targetParcels.size() << " parcels in target mesh." << endl;
254
255 if (addParticles.size())
256 {
258
259 // addParticles now contains the indices of the sourceMesh
260 // particles that were appended to the target mesh.
261
262 // Map lagrangian fields
263 // ~~~~~~~~~~~~~~~~~~~~~
264
266 (
267 cloudDir,
268 objects,
269 meshTarget,
270 addParticles
271 );
273 (
274 cloudDir,
275 objects,
276 meshTarget,
277 addParticles
278 );
280 (
281 cloudDir,
282 objects,
283 meshTarget,
284 addParticles
285 );
287 (
288 cloudDir,
289 objects,
290 meshTarget,
291 addParticles
292 );
294 (
295 cloudDir,
296 objects,
297 meshTarget,
298 addParticles
299 );
301 (
302 cloudDir,
303 objects,
304 meshTarget,
305 addParticles
306 );
307 }
308 }
309 }
310}
311
312
313// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
314
315} // End namespace Foam
316
317// ************************************************************************* //
Base cloud calls templated on particle type.
Definition Cloud.H:64
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
Helper IO class to read and write particle coordinates (positions).
Definition IOPosition.H:52
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
Definition IOPosition.C:51
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable,...
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
static const word prefix
The prefix to local: lagrangian.
Definition cloud.H:79
A class for handling file names.
Definition fileName.H:75
@ DIRECTORY
A directory.
Definition fileName.H:85
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition meshSearch.H:57
Class to calculate the cell-addressing between two overlapping meshes.
Definition meshToMesh.H:61
A Cloud of passive particles.
Copy of base particle.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition polyMesh.C:1345
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
const vectorField & cellCentres() const
volScalarField & p
dynamicFvMesh & mesh
Gets the indices of (source)particles that have been appended to the target cloud and maps the lagran...
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
Namespace for OpenFOAM.
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition HashSet.H:85
List< fileName > fileNameList
List of fileName.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
vector point
Point is a vector.
Definition point.H:37
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition barycentric.H:45
void mapLagrangian(const meshToMesh0 &meshToMesh0Interp)
Maps lagrangian positions and fields.
void MapLagrangianFields(const string &cloudName, const IOobjectList &objects, const meshToMesh0 &meshToMesh0Interp, const labelList &addParticles, const char *msg)
Gets the indices of (source)particles that have been appended to the.
fileNameList readDir(const fileName &directory, const fileName::Type type=fileName::Type::FILE, const bool filtergz=true, const bool followLink=true)
Read a directory and return the entries as a fileName List.
Definition POSIX.C:965
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & e