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-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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"
30 #include "passiveParticleCloud.H"
31 #include "meshSearch.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 
38 static const scalar perturbFactor = 1e-6;
39 
40 
41 // Special version of findCell that generates a cell guaranteed to be
42 // compatible with tracking.
43 static 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 
86 void mapLagrangian(const meshToMesh0& meshToMesh0Interp)
87 {
88  // Determine which particles are in meshTarget
89  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
90 
91  // Target to source cell map
92  const labelList& cellAddressing = meshToMesh0Interp.cellAddressing();
93 
94  // Invert celladdressing to get source to target(s).
95  // Note: could use sparse addressing but that is too storage inefficient
96  // (Map<labelList>)
97  const labelListList sourceToTargets
98  (
99  invertOneToMany(meshToMesh0Interp.fromMesh().nCells(), cellAddressing)
100  );
101 
102  const fvMesh& meshSource = meshToMesh0Interp.fromMesh();
103  const fvMesh& meshTarget = meshToMesh0Interp.toMesh();
104 
105  const fileNameList cloudDirs
106  (
107  readDir
108  (
109  meshSource.time().timePath()/cloud::prefix,
111  )
112  );
113 
114  for (const fileName& cloudDir : cloudDirs)
115  {
116  // Search for list of lagrangian objects for this time
117  IOobjectList objects
118  (
119  meshSource,
120  meshSource.time().timeName(),
121  cloud::prefix/cloudDir
122  );
123 
124  if (objects.found("coordinates") || objects.found("positions"))
125  {
126  // Has coordinates/positions - so must be a valid cloud
127  Info<< nl << " processing cloud " << cloudDir << endl;
128 
129  // Read positions & cell
130  passiveParticleCloud sourceParcels
131  (
132  meshSource,
133  cloudDir,
134  false
135  );
136 
137  Info<< " read " << sourceParcels.size()
138  << " parcels from source mesh." << endl;
139 
140  // Construct empty target cloud
141  passiveParticleCloud targetParcels
142  (
143  meshTarget,
144  Foam::zero{},
145  cloudDir
146  );
147 
148  passiveParticle::trackingData td(targetParcels);
149 
150  label sourceParticleI = 0;
151 
152  // Indices of source particles that get added to targetParcels
153  DynamicList<label> addParticles(sourceParcels.size());
154 
155  // Unmapped particles
156  labelHashSet unmappedSource(sourceParcels.size());
157 
158 
159  // Initial: track from fine-mesh cell centre to particle position
160  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
161  // This requires there to be no boundary in the way.
162 
163 
164  for (const passiveParticle& p : sourceParcels)
165  {
166  bool foundCell = false;
167 
168  // Assume that cell from read parcel is the correct one...
169  if (p.cell() >= 0)
170  {
171  const labelList& targetCells =
172  sourceToTargets[p.cell()];
173 
174  // Particle probably in one of the targetcells. Try
175  // all by tracking from their cell centre to the parcel
176  // position.
177 
178  for (const label targetCell : targetCells)
179  {
180  // Track from its cellcentre to position to make sure.
181  autoPtr<passiveParticle> newPtr
182  (
183  new passiveParticle
184  (
185  meshTarget,
186  barycentric(1, 0, 0, 0),
187  targetCell,
188  meshTarget.cells()[targetCell][0],
189  1
190  )
191  );
192  passiveParticle& newP = newPtr();
193 
194  newP.track(p.position() - newP.position(), 0);
195 
196  if (!newP.onFace())
197  {
198  // Hit position.
199  foundCell = true;
200  addParticles.append(sourceParticleI);
201  targetParcels.addParticle(newPtr.ptr());
202  break;
203  }
204  }
205  }
206 
207  if (!foundCell)
208  {
209  // Store for closer analysis
210  unmappedSource.insert(sourceParticleI);
211  }
212 
213  sourceParticleI++;
214  }
215 
216  Info<< " after meshToMesh0 addressing found "
217  << targetParcels.size()
218  << " parcels in target mesh." << endl;
219 
220 
221  // Do closer inspection for unmapped particles
222  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
223 
224  if (unmappedSource.size())
225  {
226  sourceParticleI = 0;
227 
228  for (passiveParticle& p : sourceParcels)
229  {
230  if (unmappedSource.found(sourceParticleI))
231  {
232  const label targetCell =
233  findCell(targetParcels, p.position());
234 
235  if (targetCell >= 0)
236  {
237  unmappedSource.erase(sourceParticleI);
238  addParticles.append(sourceParticleI);
239  targetParcels.addParticle
240  (
241  new passiveParticle
242  (
243  meshTarget,
244  p.position(),
245  targetCell
246  )
247  );
248  sourceParcels.remove(&p);
249  }
250  }
251  sourceParticleI++;
252  }
253  }
254  addParticles.shrink();
255 
256  Info<< " after additional mesh searching found "
257  << targetParcels.size() << " parcels in target mesh." << endl;
258 
259  if (addParticles.size())
260  {
261  IOPosition<passiveParticleCloud>(targetParcels).write();
262 
263  // addParticles now contains the indices of the sourceMesh
264  // particles that were appended to the target mesh.
265 
266  // Map lagrangian fields
267  // ~~~~~~~~~~~~~~~~~~~~~
268 
269  MapLagrangianFields<label>
270  (cloudDir, objects, meshToMesh0Interp, addParticles);
271  MapLagrangianFields<scalar>
272  (cloudDir, objects, meshToMesh0Interp, addParticles);
273  MapLagrangianFields<vector>
274  (cloudDir, objects, meshToMesh0Interp, addParticles);
275  MapLagrangianFields<sphericalTensor>
276  (cloudDir, objects, meshToMesh0Interp, addParticles);
277  MapLagrangianFields<symmTensor>
278  (cloudDir, objects, meshToMesh0Interp, addParticles);
279  MapLagrangianFields<tensor>
280  (cloudDir, objects, meshToMesh0Interp, addParticles);
281  }
282  }
283  }
284 }
285 
286 
287 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
288 
289 } // End namespace Foam
290 
291 // ************************************************************************* //
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition: barycentric.H:41
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
void write(vtk::formatter &fmt, const Type &val, const label n=1)
Component-wise write of a value (N times)
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:107
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1360
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
const vectorField & cellCentres() const
void mapLagrangian(const meshToMesh0 &meshToMesh0Interp)
Maps lagrangian positions and fields.
vector point
Point is a vector.
Definition: point.H:37
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
List< fileName > fileNameList
List of fileName.
Definition: fileNameList.H:32
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:963
Namespace for OpenFOAM.
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:79