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 meshToMesh& interp)
84 {
85  // Determine which particles are in meshTarget
86  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
87 
88  const polyMesh& meshSource = interp.srcRegion();
89  const polyMesh& meshTarget = interp.tgtRegion();
90  const labelListList& sourceToTarget = interp.srcToTgtCellAddr();
91 
92  const pointField& targetCc = meshTarget.cellCentres();
93 
94  fileNameList cloudDirs
95  (
96  readDir
97  (
98  meshSource.time().timePath()/cloud::prefix,
100  )
101  );
102 
103  forAll(cloudDirs, cloudI)
104  {
105  // Search for list of lagrangian objects for this time
106  IOobjectList objects
107  (
108  meshSource,
109  meshSource.time().timeName(),
110  cloud::prefix/cloudDirs[cloudI]
111  );
112 
113  IOobject* positionsPtr = objects.lookup(word("positions"));
114 
115  if (positionsPtr)
116  {
117  Info<< nl << " processing cloud " << cloudDirs[cloudI] << endl;
118 
119  // Read positions & cell
120  passiveParticleCloud sourceParcels
121  (
122  meshSource,
123  cloudDirs[cloudI],
124  false
125  );
126  Info<< " read " << sourceParcels.size()
127  << " parcels from source mesh." << endl;
128 
129  // Construct empty target cloud
130  passiveParticleCloud targetParcels
131  (
132  meshTarget,
133  cloudDirs[cloudI],
134  IDLList<passiveParticle>()
135  );
136 
137  particle::TrackingData<passiveParticleCloud> td(targetParcels);
138 
139  label sourceParticleI = 0;
140 
141  // Indices of source particles that get added to targetParcels
142  DynamicList<label> addParticles(sourceParcels.size());
143 
144  // Unmapped particles
145  labelHashSet unmappedSource(sourceParcels.size());
146 
147 
148  // Initial: track from fine-mesh cell centre to particle position
149  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
150  // This requires there to be no boundary in the way.
151 
152 
153  forAllConstIter(Cloud<passiveParticle>, sourceParcels, iter)
154  {
155  bool foundCell = false;
156 
157  // Assume that cell from read parcel is the correct one...
158  if (iter().cell() >= 0)
159  {
160  const labelList& targetCells =
161  sourceToTarget[iter().cell()];
162 
163  // Particle probably in one of the targetcells. Try
164  // all by tracking from their cell centre to the parcel
165  // position.
166 
167  forAll(targetCells, i)
168  {
169  // Track from its cellcentre to position to make sure.
170  autoPtr<passiveParticle> newPtr
171  (
172  new passiveParticle
173  (
174  meshTarget,
175  targetCc[targetCells[i]],
176  targetCells[i]
177  )
178  );
179  passiveParticle& newP = newPtr();
180 
181  label facei = newP.track(iter().position(), td);
182 
183  if (facei < 0 && newP.cell() >= 0)
184  {
185  // Hit position.
186  foundCell = true;
187  addParticles.append(sourceParticleI);
188  targetParcels.addParticle(newPtr.ptr());
189  break;
190  }
191  }
192  }
193 
194  if (!foundCell)
195  {
196  // Store for closer analysis
197  unmappedSource.insert(sourceParticleI);
198  }
199 
200  sourceParticleI++;
201  }
202 
203  Info<< " after meshToMesh addressing found "
204  << targetParcels.size()
205  << " parcels in target mesh." << endl;
206 
207 
208  // Do closer inspection for unmapped particles
209  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
210 
211  if (unmappedSource.size())
212  {
213  sourceParticleI = 0;
214 
215  forAllIter(Cloud<passiveParticle>, sourceParcels, iter)
216  {
217  if (unmappedSource.found(sourceParticleI))
218  {
219  label targetCell =
220  findCell(targetParcels, iter().position());
221 
222  if (targetCell >= 0)
223  {
224  unmappedSource.erase(sourceParticleI);
225  addParticles.append(sourceParticleI);
226  iter().cell() = targetCell;
227  targetParcels.addParticle
228  (
229  sourceParcels.remove(&iter())
230  );
231  }
232  }
233  sourceParticleI++;
234  }
235  }
236  addParticles.shrink();
237 
238  Info<< " after additional mesh searching found "
239  << targetParcels.size() << " parcels in target mesh." << endl;
240 
241  if (addParticles.size())
242  {
243  IOPosition<passiveParticleCloud>(targetParcels).write();
244 
245  // addParticles now contains the indices of the sourceMesh
246  // particles that were appended to the target mesh.
247 
248  // Map lagrangian fields
249  // ~~~~~~~~~~~~~~~~~~~~~
250 
251  MapLagrangianFields<label>
252  (
253  cloudDirs[cloudI],
254  objects,
255  meshTarget,
256  addParticles
257  );
258  MapLagrangianFields<scalar>
259  (
260  cloudDirs[cloudI],
261  objects,
262  meshTarget,
263  addParticles
264  );
265  MapLagrangianFields<vector>
266  (
267  cloudDirs[cloudI],
268  objects,
269  meshTarget,
270  addParticles
271  );
272  MapLagrangianFields<sphericalTensor>
273  (
274  cloudDirs[cloudI],
275  objects,
276  meshTarget,
277  addParticles
278  );
279  MapLagrangianFields<symmTensor>
280  (
281  cloudDirs[cloudI],
282  objects,
283  meshTarget,
284  addParticles
285  );
286  MapLagrangianFields<tensor>
287  (
288  cloudDirs[cloudI],
289  objects,
290  meshTarget,
291  addParticles
292  );
293  }
294  }
295  }
296 }
297 
298 
299 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
300 
301 } // End namespace Foam
302 
303 // ************************************************************************* //
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
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.