mapLagrangian.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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
19  for more details.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "MapLagrangianFields.H"
27 #include "passiveParticleCloud.H"
28 #include "meshSearch.H"
29 #include "OSspecific.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 static const scalar perturbFactor = 1e-6;
37 
38 
39 // Special version of findCell that generates a cell guaranteed to be
40 // compatible with tracking.
41 static label findCell
42 (
43  const lagrangian::Cloud<passiveParticle>& cloud,
44  const point& pt
45 )
46 {
47  label celli = -1;
48  label tetFacei = -1;
49  label tetPtI = -1;
50 
51  const polyMesh& mesh = cloud.pMesh();
52 
53  mesh.findCellFacePt(pt, celli, tetFacei, tetPtI);
54 
55  if (celli >= 0)
56  {
57  return celli;
58  }
59  else
60  {
61  // See if particle on face by finding nearest face and shifting
62  // particle.
63 
64  meshSearch meshSearcher
65  (
66  mesh,
67  polyMesh::FACE_PLANES // no decomposition needed
68  );
69 
70  label facei = meshSearcher.findNearestBoundaryFace(pt);
71 
72  if (facei >= 0)
73  {
74  const point& cc = mesh.cellCentres()[mesh.faceOwner()[facei]];
75 
76  const point perturbPt = (1-perturbFactor)*pt+perturbFactor*cc;
77 
78  mesh.findCellFacePt(perturbPt, celli, tetFacei, tetPtI);
79 
80  return celli;
81  }
82  }
83 
84  return -1;
85 }
86 
87 
88 void mapLagrangian(const meshToMesh0& meshToMesh0Interp)
89 {
90  // Determine which particles are in meshTarget
91  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
92 
93  // target to source cell map
94  const labelList& cellAddressing = meshToMesh0Interp.cellAddressing();
95 
96  // Invert celladdressing to get source to target(s).
97  // Note: could use sparse addressing but that is too storage inefficient
98  // (Map<labelList>)
99  labelListList sourceToTargets
100  (
101  invertOneToMany(meshToMesh0Interp.fromMesh().nCells(), cellAddressing)
102  );
103 
104  const fvMesh& meshSource = meshToMesh0Interp.fromMesh();
105  const fvMesh& meshTarget = meshToMesh0Interp.toMesh();
106 
107 
108  fileNameList cloudDirs
109  (
110  readDir
111  (
112  meshSource.time().timePath()/lagrangian::cloud::prefix,
114  )
115  );
116 
117  forAll(cloudDirs, cloudI)
118  {
119  // Search for list of lagrangian objects for this time
120  IOobjectList objects
121  (
122  meshSource,
123  meshSource.time().name(),
124  lagrangian::cloud::prefix/cloudDirs[cloudI]
125  );
126 
127  IOobject* positionsPtr = objects.lookup("positions");
128 
129  if (positionsPtr)
130  {
131  Info<< nl << " processing cloud " << cloudDirs[cloudI] << endl;
132 
133  // Read positions & cell
134  passiveParticleCloud sourceParcels
135  (
136  meshSource,
137  cloudDirs[cloudI],
138  false
139  );
140  Info<< " read " << sourceParcels.size()
141  << " parcels from source mesh." << endl;
142 
143  // Construct empty target cloud
144  passiveParticleCloud targetParcels
145  (
146  meshTarget,
147  cloudDirs[cloudI],
148  IDLList<passiveParticle>()
149  );
150 
151  passiveParticle::trackingData td(targetParcels);
152 
153  label sourceParticleI = 0;
154 
155  // Indices of source particles that get added to targetParcels
156  DynamicList<label> addParticles(sourceParcels.size());
157 
158  // Unmapped particles
159  labelHashSet unmappedSource(sourceParcels.size());
160 
161 
162  // Initial: track from fine-mesh cell centre to particle position
163  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
164  // This requires there to be no boundary in the way.
165 
166 
168  (
169  lagrangian::Cloud<passiveParticle>,
170  sourceParcels,
171  iter
172  )
173  {
174  bool foundCell = false;
175 
176  // Assume that cell from read parcel is the correct one...
177  if (iter().cell() >= 0)
178  {
179  const labelList& targetCells =
180  sourceToTargets[iter().cell()];
181 
182  // Particle probably in one of the targetcells. Try
183  // all by tracking from their cell centre to the parcel
184  // position.
185 
186  forAll(targetCells, i)
187  {
188  // Track from its cellcentre to position to make sure.
189  autoPtr<passiveParticle> newPtr
190  (
191  new passiveParticle
192  (
193  meshTarget,
194  barycentric(1, 0, 0, 0),
195  targetCells[i],
196  meshTarget.cells()[targetCells[i]][0],
197  1
198  )
199  );
200  passiveParticle& newP = newPtr();
201 
202  newP.track
203  (
204  meshTarget,
205  iter().position(meshSource)
206  - newP.position(meshTarget),
207  0
208  );
209 
210  if (!newP.onFace())
211  {
212  // Hit position.
213  foundCell = true;
214  addParticles.append(sourceParticleI);
215  targetParcels.addParticle(newPtr.ptr());
216  break;
217  }
218  }
219  }
220 
221  if (!foundCell)
222  {
223  // Store for closer analysis
224  unmappedSource.insert(sourceParticleI);
225  }
226 
227  sourceParticleI++;
228  }
229 
230  Info<< " after meshToMesh0 addressing found "
231  << targetParcels.size()
232  << " parcels in target mesh." << endl;
233 
234 
235  // Do closer inspection for unmapped particles
236  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
237 
238  if (unmappedSource.size())
239  {
240  sourceParticleI = 0;
241 
242  forAllIter
243  (
244  lagrangian::Cloud<passiveParticle>,
245  sourceParcels,
246  iter
247  )
248  {
249  if (unmappedSource.found(sourceParticleI))
250  {
251  label targetCell =
252  findCell
253  (
254  targetParcels,
255  iter().position(meshSource)
256  );
257 
258  if (targetCell >= 0)
259  {
260  label nLocateBoundaryHits = 0;
261  autoPtr<passiveParticle> pPtr
262  (
263  new passiveParticle
264  (
265  meshTarget,
266  iter().position(meshSource),
267  targetCell,
268  nLocateBoundaryHits
269  )
270  );
271 
272  if (nLocateBoundaryHits == 0)
273  {
274  unmappedSource.erase(sourceParticleI);
275  addParticles.append(sourceParticleI);
276  targetParcels.addParticle(pPtr.ptr());
277  sourceParcels.remove(&iter());
278  }
279  }
280  }
281  sourceParticleI++;
282  }
283  }
284  addParticles.shrink();
285 
286  Info<< " after additional mesh searching found "
287  << targetParcels.size() << " parcels in target mesh." << endl;
288 
289  if (addParticles.size())
290  {
291  IOPosition<passiveParticleCloud>(targetParcels).write();
292 
293  // addParticles now contains the indices of the sourceMesh
294  // particles that were appended to the target mesh.
295 
296  // Map lagrangian fields
297  // ~~~~~~~~~~~~~~~~~~~~~
298 
299  MapLagrangianFields<label>
300  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
301  MapLagrangianFields<scalar>
302  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
303  MapLagrangianFields<vector>
304  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
305  MapLagrangianFields<sphericalTensor>
306  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
307  MapLagrangianFields<symmTensor>
308  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
309  MapLagrangianFields<tensor>
310  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
311  }
312  }
313  }
314 }
315 
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 } // End namespace Foam
320 
321 // ************************************************************************* //
Gets the indices of (source)particles that have been appended to the target cloud and maps the lagran...
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:458
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
virtual Ostream & write(const char)=0
Write character.
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:69
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1582
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
const vectorField & cellCentres() const
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
point position(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label facei, const label faceTrii, const scalar stepFraction)
Return the position given the coordinates and tet topology.
Definition: trackingI.H:224
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:106
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition: barycentric.H:45
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
vector point
Point is a vector.
Definition: point.H:41
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
void mapLagrangian(const meshToMesh0 &meshToMesh0Interp)
Maps lagrangian positions and fields.
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:213
static const char nl
Definition: Ostream.H:267
fileNameList readDir(const fileName &, const fileType=fileType::file, const bool filterVariants=true, const bool followLink=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:662
objects