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 | Copyright (C) 2011-2016 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 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 static const scalar perturbFactor = 1e-6;
36 
37 
38 // Special version of findCell that generates a cell guaranteed to be
39 // compatible with tracking.
40 static label findCell(const Cloud<passiveParticle>& cloud, const point& pt)
41 {
42  label celli = -1;
43  label tetFacei = -1;
44  label tetPtI = -1;
45 
46  const polyMesh& mesh = cloud.pMesh();
47 
48  mesh.findCellFacePt(pt, celli, tetFacei, tetPtI);
49 
50  if (celli >= 0)
51  {
52  return celli;
53  }
54  else
55  {
56  // See if particle on face by finding nearest face and shifting
57  // particle.
58 
59  meshSearch meshSearcher
60  (
61  mesh,
62  polyMesh::FACE_PLANES // no decomposition needed
63  );
64 
65  label facei = meshSearcher.findNearestBoundaryFace(pt);
66 
67  if (facei >= 0)
68  {
69  const point& cc = mesh.cellCentres()[mesh.faceOwner()[facei]];
70 
71  const point perturbPt = (1-perturbFactor)*pt+perturbFactor*cc;
72 
73  mesh.findCellFacePt(perturbPt, celli, tetFacei, tetPtI);
74 
75  return celli;
76  }
77  }
78 
79  return -1;
80 }
81 
82 
83 void mapLagrangian(const meshToMesh0& meshToMesh0Interp)
84 {
85  // Determine which particles are in meshTarget
86  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
87 
88  // target to source cell map
89  const labelList& cellAddressing = meshToMesh0Interp.cellAddressing();
90 
91  // Invert celladdressing to get source to target(s).
92  // Note: could use sparse addressing but that is too storage inefficient
93  // (Map<labelList>)
94  labelListList sourceToTargets
95  (
96  invertOneToMany(meshToMesh0Interp.fromMesh().nCells(), cellAddressing)
97  );
98 
99  const fvMesh& meshSource = meshToMesh0Interp.fromMesh();
100  const fvMesh& meshTarget = meshToMesh0Interp.toMesh();
101  const pointField& targetCc = meshTarget.cellCentres();
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().timeName(),
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  particle::TrackingData<passiveParticleCloud> 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  targetCc[targetCells[i]],
186  targetCells[i]
187  )
188  );
189  passiveParticle& newP = newPtr();
190 
191  label facei = newP.track(iter().position(), td);
192 
193  if (facei < 0 && newP.cell() >= 0)
194  {
195  // Hit position.
196  foundCell = true;
197  addParticles.append(sourceParticleI);
198  targetParcels.addParticle(newPtr.ptr());
199  break;
200  }
201  }
202  }
203 
204  if (!foundCell)
205  {
206  // Store for closer analysis
207  unmappedSource.insert(sourceParticleI);
208  }
209 
210  sourceParticleI++;
211  }
212 
213  Info<< " after meshToMesh0 addressing found "
214  << targetParcels.size()
215  << " parcels in target mesh." << endl;
216 
217 
218  // Do closer inspection for unmapped particles
219  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
220 
221  if (unmappedSource.size())
222  {
223  sourceParticleI = 0;
224 
225  forAllIter(Cloud<passiveParticle>, sourceParcels, iter)
226  {
227  if (unmappedSource.found(sourceParticleI))
228  {
229  label targetCell =
230  findCell(targetParcels, iter().position());
231 
232  if (targetCell >= 0)
233  {
234  unmappedSource.erase(sourceParticleI);
235  addParticles.append(sourceParticleI);
236  iter().cell() = targetCell;
237  targetParcels.addParticle
238  (
239  sourceParcels.remove(&iter())
240  );
241  }
242  }
243  sourceParticleI++;
244  }
245  }
246  addParticles.shrink();
247 
248  Info<< " after additional mesh searching found "
249  << targetParcels.size() << " parcels in target mesh." << endl;
250 
251  if (addParticles.size())
252  {
253  IOPosition<passiveParticleCloud>(targetParcels).write();
254 
255  // addParticles now contains the indices of the sourceMesh
256  // particles that were appended to the target mesh.
257 
258  // Map lagrangian fields
259  // ~~~~~~~~~~~~~~~~~~~~~
260 
261  MapLagrangianFields<label>
262  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
263  MapLagrangianFields<scalar>
264  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
265  MapLagrangianFields<vector>
266  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
267  MapLagrangianFields<sphericalTensor>
268  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
269  MapLagrangianFields<symmTensor>
270  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
271  MapLagrangianFields<tensor>
272  (cloudDirs[cloudI], objects, meshToMesh0Interp, addParticles);
273  }
274  }
275  }
276 }
277 
278 
279 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280 
281 } // End namespace Foam
282 
283 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
const double e
Elementary charge.
Definition: doubleFloat.H:78
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:453
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
static const char nl
Definition: Ostream.H:262
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:71
void mapLagrangian(const meshToMesh0 &meshToMesh0Interp)
Maps lagrangian positions and fields.
vector point
Point is a vector.
Definition: point.H:41
messageStream Info
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
fileNameList readDir(const fileName &, const fileName::Type=fileName::FILE, const bool filtergz=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:527
runTime write()
Namespace for OpenFOAM.