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-2022 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(const Cloud<passiveParticle>& cloud, const point& pt)
42 {
43  label celli = -1;
44  label tetFacei = -1;
45  label tetPtI = -1;
46 
47  const polyMesh& mesh = cloud.pMesh();
48 
49  mesh.findCellFacePt(pt, celli, tetFacei, tetPtI);
50 
51  if (celli >= 0)
52  {
53  return celli;
54  }
55  else
56  {
57  // See if particle on face by finding nearest face and shifting
58  // particle.
59 
60  meshSearch meshSearcher
61  (
62  mesh,
63  polyMesh::FACE_PLANES // no decomposition needed
64  );
65 
66  label facei = meshSearcher.findNearestBoundaryFace(pt);
67 
68  if (facei >= 0)
69  {
70  const point& cc = mesh.cellCentres()[mesh.faceOwner()[facei]];
71 
72  const point perturbPt = (1-perturbFactor)*pt+perturbFactor*cc;
73 
74  mesh.findCellFacePt(perturbPt, celli, tetFacei, tetPtI);
75 
76  return celli;
77  }
78  }
79 
80  return -1;
81 }
82 
83 
84 void mapLagrangian(const meshToMesh0& meshToMesh0Interp)
85 {
86  // Determine which particles are in meshTarget
87  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
88 
89  // target to source cell map
90  const labelList& cellAddressing = meshToMesh0Interp.cellAddressing();
91 
92  // Invert celladdressing to get source to target(s).
93  // Note: could use sparse addressing but that is too storage inefficient
94  // (Map<labelList>)
95  labelListList sourceToTargets
96  (
97  invertOneToMany(meshToMesh0Interp.fromMesh().nCells(), cellAddressing)
98  );
99 
100  const fvMesh& meshSource = meshToMesh0Interp.fromMesh();
101  const fvMesh& meshTarget = meshToMesh0Interp.toMesh();
102 
103 
104  fileNameList cloudDirs
105  (
106  readDir
107  (
108  meshSource.time().timePath()/cloud::prefix,
110  )
111  );
112 
113  forAll(cloudDirs, cloudI)
114  {
115  // Search for list of lagrangian objects for this time
116  IOobjectList objects
117  (
118  meshSource,
119  meshSource.time().name(),
120  cloud::prefix/cloudDirs[cloudI]
121  );
122 
123  IOobject* positionsPtr = objects.lookup("positions");
124 
125  if (positionsPtr)
126  {
127  Info<< nl << " processing cloud " << cloudDirs[cloudI] << endl;
128 
129  // Read positions & cell
130  passiveParticleCloud sourceParcels
131  (
132  meshSource,
133  cloudDirs[cloudI],
134  false
135  );
136  Info<< " read " << sourceParcels.size()
137  << " parcels from source mesh." << endl;
138 
139  // Construct empty target cloud
140  passiveParticleCloud targetParcels
141  (
142  meshTarget,
143  cloudDirs[cloudI],
144  IDLList<passiveParticle>()
145  );
146 
147  passiveParticle::trackingData td(targetParcels);
148 
149  label sourceParticleI = 0;
150 
151  // Indices of source particles that get added to targetParcels
152  DynamicList<label> addParticles(sourceParcels.size());
153 
154  // Unmapped particles
155  labelHashSet unmappedSource(sourceParcels.size());
156 
157 
158  // Initial: track from fine-mesh cell centre to particle position
159  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
160  // This requires there to be no boundary in the way.
161 
162 
163  forAllConstIter(Cloud<passiveParticle>, sourceParcels, iter)
164  {
165  bool foundCell = false;
166 
167  // Assume that cell from read parcel is the correct one...
168  if (iter().cell() >= 0)
169  {
170  const labelList& targetCells =
171  sourceToTargets[iter().cell()];
172 
173  // Particle probably in one of the targetcells. Try
174  // all by tracking from their cell centre to the parcel
175  // position.
176 
177  forAll(targetCells, i)
178  {
179  // Track from its cellcentre to position to make sure.
180  autoPtr<passiveParticle> newPtr
181  (
182  new passiveParticle
183  (
184  meshTarget,
185  barycentric(1, 0, 0, 0),
186  targetCells[i],
187  meshTarget.cells()[targetCells[i]][0],
188  1
189  )
190  );
191  passiveParticle& newP = newPtr();
192 
193  newP.track
194  (
195  meshTarget,
196  iter().position(meshSource)
197  - newP.position(meshTarget),
198  0
199  );
200 
201  if (!newP.onFace())
202  {
203  // Hit position.
204  foundCell = true;
205  addParticles.append(sourceParticleI);
206  targetParcels.addParticle(newPtr.ptr());
207  break;
208  }
209  }
210  }
211 
212  if (!foundCell)
213  {
214  // Store for closer analysis
215  unmappedSource.insert(sourceParticleI);
216  }
217 
218  sourceParticleI++;
219  }
220 
221  Info<< " after meshToMesh0 addressing found "
222  << targetParcels.size()
223  << " parcels in target mesh." << endl;
224 
225 
226  // Do closer inspection for unmapped particles
227  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
228 
229  if (unmappedSource.size())
230  {
231  sourceParticleI = 0;
232 
233  forAllIter(Cloud<passiveParticle>, sourceParcels, iter)
234  {
235  if (unmappedSource.found(sourceParticleI))
236  {
237  label targetCell =
238  findCell
239  (
240  targetParcels,
241  iter().position(meshSource)
242  );
243 
244  if (targetCell >= 0)
245  {
246  unmappedSource.erase(sourceParticleI);
247  addParticles.append(sourceParticleI);
248  targetParcels.addParticle
249  (
250  new passiveParticle
251  (
252  meshTarget,
253  iter().position(meshSource),
254  targetCell
255  )
256  );
257  sourceParcels.remove(&iter());
258  }
259  }
260  sourceParticleI++;
261  }
262  }
263  addParticles.shrink();
264 
265  Info<< " after additional mesh searching found "
266  << targetParcels.size() << " parcels in target mesh." << endl;
267 
268  if (addParticles.size())
269  {
270  IOPosition<passiveParticleCloud>(targetParcels).write();
271 
272  // addParticles now contains the indices of the sourceMesh
273  // particles that were appended to the target mesh.
274 
275  // Map lagrangian fields
276  // ~~~~~~~~~~~~~~~~~~~~~
277 
278  MapLagrangianFields<label>
279  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
280  MapLagrangianFields<scalar>
281  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
282  MapLagrangianFields<vector>
283  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
284  MapLagrangianFields<sphericalTensor>
285  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
286  MapLagrangianFields<symmTensor>
287  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
288  MapLagrangianFields<tensor>
289  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
290  }
291  }
292  }
293 }
294 
295 
296 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
297 
298 } // End namespace Foam
299 
300 // ************************************************************************* //
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:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
virtual Ostream & write(const char)=0
Write character.
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:63
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:105
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:251
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:211
static const char nl
Definition: Ostream.H:260
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