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  Info<< " read " << sourceParcels.size()
137  << " parcels from source mesh." << endl;
138 
139  // Construct empty target cloud
140  passiveParticleCloud targetParcels
141  (
142  meshTarget,
143  cloudDir,
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  for (const passiveParticle& p : sourceParcels)
164  {
165  bool foundCell = false;
166 
167  // Assume that cell from read parcel is the correct one...
168  if (p.cell() >= 0)
169  {
170  const labelList& targetCells =
171  sourceToTargets[p.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  for (const label targetCell : targetCells)
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  targetCell,
187  meshTarget.cells()[targetCell][0],
188  1
189  )
190  );
191  passiveParticle& newP = newPtr();
192 
193  newP.track(p.position() - newP.position(), 0);
194 
195  if (!newP.onFace())
196  {
197  // Hit position.
198  foundCell = true;
199  addParticles.append(sourceParticleI);
200  targetParcels.addParticle(newPtr.ptr());
201  break;
202  }
203  }
204  }
205 
206  if (!foundCell)
207  {
208  // Store for closer analysis
209  unmappedSource.insert(sourceParticleI);
210  }
211 
212  sourceParticleI++;
213  }
214 
215  Info<< " after meshToMesh0 addressing found "
216  << targetParcels.size()
217  << " parcels in target mesh." << endl;
218 
219 
220  // Do closer inspection for unmapped particles
221  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
222 
223  if (unmappedSource.size())
224  {
225  sourceParticleI = 0;
226 
227  for (passiveParticle& p : sourceParcels)
228  {
229  if (unmappedSource.found(sourceParticleI))
230  {
231  const label targetCell =
232  findCell(targetParcels, p.position());
233 
234  if (targetCell >= 0)
235  {
236  unmappedSource.erase(sourceParticleI);
237  addParticles.append(sourceParticleI);
238  targetParcels.addParticle
239  (
240  new passiveParticle
241  (
242  meshTarget,
243  p.position(),
244  targetCell
245  )
246  );
247  sourceParcels.remove(&p);
248  }
249  }
250  sourceParticleI++;
251  }
252  }
253  addParticles.shrink();
254 
255  Info<< " after additional mesh searching found "
256  << targetParcels.size() << " parcels in target mesh." << endl;
257 
258  if (addParticles.size())
259  {
260  IOPosition<passiveParticleCloud>(targetParcels).write();
261 
262  // addParticles now contains the indices of the sourceMesh
263  // particles that were appended to the target mesh.
264 
265  // Map lagrangian fields
266  // ~~~~~~~~~~~~~~~~~~~~~
267 
268  MapLagrangianFields<label>
269  (cloudDir, objects, meshToMesh0Interp, addParticles);
270  MapLagrangianFields<scalar>
271  (cloudDir, objects, meshToMesh0Interp, addParticles);
272  MapLagrangianFields<vector>
273  (cloudDir, objects, meshToMesh0Interp, addParticles);
274  MapLagrangianFields<sphericalTensor>
275  (cloudDir, objects, meshToMesh0Interp, addParticles);
276  MapLagrangianFields<symmTensor>
277  (cloudDir, objects, meshToMesh0Interp, addParticles);
278  MapLagrangianFields<tensor>
279  (cloudDir, objects, meshToMesh0Interp, addParticles);
280  }
281  }
282  }
283 }
284 
285 
286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
287 
288 } // End namespace Foam
289 
290 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::cloud::prefix
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:87
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
Foam::polyMesh::findCellFacePt
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1356
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::fileNameList
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:58
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
meshSearch.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
passiveParticleCloud.H
Foam::mapLagrangian
void mapLagrangian(const meshToMesh0 &meshToMesh0Interp)
Maps lagrangian positions and fields.
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::fileName::DIRECTORY
A directory.
Definition: fileName.H:84
Foam::polyMesh::FACE_PLANES
Definition: polyMesh.H:102
Foam::barycentric
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition: barycentric.H:47
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::invertOneToMany
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:114
Foam::readDir
fileNameList readDir(const fileName &directory, const fileName::Type type=fileName::FILE, const bool filtergz=true, const bool followLink=true)
Read a directory and return the entries as a fileName List.
Definition: MSwindows.C:707