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-2015 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  bool foundPositions =
114  returnReduce(objects.found("positions"), orOp<bool>());;
115 
116  if (foundPositions)
117  {
118  Info<< nl << " processing cloud " << cloudDirs[cloudI] << endl;
119 
120  // Read positions & cell
121  passiveParticleCloud sourceParcels
122  (
123  meshSource,
124  cloudDirs[cloudI],
125  false
126  );
127  Info<< " read " << sourceParcels.size()
128  << " parcels from source mesh." << endl;
129 
130  // Construct empty target cloud
131  passiveParticleCloud targetParcels
132  (
133  meshTarget,
134  cloudDirs[cloudI],
135  IDLList<passiveParticle>()
136  );
137 
138  particle::TrackingData<passiveParticleCloud> td(targetParcels);
139 
140  label sourceParticleI = 0;
141 
142  // Indices of source particles that get added to targetParcels
143  DynamicList<label> addParticles(sourceParcels.size());
144 
145  // Unmapped particles
146  labelHashSet unmappedSource(sourceParcels.size());
147 
148 
149  // Initial: track from fine-mesh cell centre to particle position
150  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
151  // This requires there to be no boundary in the way.
152 
153 
154  forAllConstIter(Cloud<passiveParticle>, sourceParcels, iter)
155  {
156  bool foundCell = false;
157 
158  // Assume that cell from read parcel is the correct one...
159  if (iter().cell() >= 0)
160  {
161  const labelList& targetCells =
162  sourceToTarget[iter().cell()];
163 
164  // Particle probably in one of the targetcells. Try
165  // all by tracking from their cell centre to the parcel
166  // position.
167 
168  forAll(targetCells, i)
169  {
170  // Track from its cellcentre to position to make sure.
171  autoPtr<passiveParticle> newPtr
172  (
173  new passiveParticle
174  (
175  meshTarget,
176  targetCc[targetCells[i]],
177  targetCells[i]
178  )
179  );
180  passiveParticle& newP = newPtr();
181 
182  label faceI = newP.track(iter().position(), td);
183 
184  if (faceI < 0 && newP.cell() >= 0)
185  {
186  // Hit position.
187  foundCell = true;
188  addParticles.append(sourceParticleI);
189  targetParcels.addParticle(newPtr.ptr());
190  break;
191  }
192  }
193  }
194 
195  if (!foundCell)
196  {
197  // Store for closer analysis
198  unmappedSource.insert(sourceParticleI);
199  }
200 
201  sourceParticleI++;
202  }
203 
204  Info<< " after meshToMesh addressing found "
205  << targetParcels.size()
206  << " parcels in target mesh." << endl;
207 
208 
209  // Do closer inspection for unmapped particles
210  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
211 
212  if (unmappedSource.size())
213  {
214  sourceParticleI = 0;
215 
216  forAllIter(Cloud<passiveParticle>, sourceParcels, iter)
217  {
218  if (unmappedSource.found(sourceParticleI))
219  {
220  label targetCell =
221  findCell(targetParcels, iter().position());
222 
223  if (targetCell >= 0)
224  {
225  unmappedSource.erase(sourceParticleI);
226  addParticles.append(sourceParticleI);
227  iter().cell() = targetCell;
228  targetParcels.addParticle
229  (
230  sourceParcels.remove(&iter())
231  );
232  }
233  }
234  sourceParticleI++;
235  }
236  }
237  addParticles.shrink();
238 
239  Info<< " after additional mesh searching found "
240  << targetParcels.size() << " parcels in target mesh." << endl;
241 
242  if (addParticles.size())
243  {
244  IOPosition<passiveParticleCloud>(targetParcels).write();
245 
246  // addParticles now contains the indices of the sourceMesh
247  // particles that were appended to the target mesh.
248 
249  // Map lagrangian fields
250  // ~~~~~~~~~~~~~~~~~~~~~
251 
252  MapLagrangianFields<label>
253  (
254  cloudDirs[cloudI],
255  objects,
256  meshTarget,
257  addParticles
258  );
259  MapLagrangianFields<scalar>
260  (
261  cloudDirs[cloudI],
262  objects,
263  meshTarget,
264  addParticles
265  );
266  MapLagrangianFields<vector>
267  (
268  cloudDirs[cloudI],
269  objects,
270  meshTarget,
271  addParticles
272  );
273  MapLagrangianFields<sphericalTensor>
274  (
275  cloudDirs[cloudI],
276  objects,
277  meshTarget,
278  addParticles
279  );
280  MapLagrangianFields<symmTensor>
281  (
282  cloudDirs[cloudI],
283  objects,
284  meshTarget,
285  addParticles
286  );
287  MapLagrangianFields<tensor>
288  (
289  cloudDirs[cloudI],
290  objects,
291  meshTarget,
292  addParticles
293  );
294  }
295  }
296  }
297 }
298 
299 
300 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
301 
302 } // End namespace Foam
303 
304 // ************************************************************************* //
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::label
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
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
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:1204
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::fileNameList
List< fileName > fileNameList
A List of fileNames.
Definition: fileNameList.H:50
Foam::cloud::prefix
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:71
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:211
meshSearch.H
passiveParticleCloud.H
Foam::mapLagrangian
void mapLagrangian(const meshToMesh0 &meshToMesh0Interp)
Maps lagrangian positions and fields.
Foam::fileName::DIRECTORY
@ DIRECTORY
Definition: fileName.H:86
Foam::polyMesh::FACE_PLANES
@ FACE_PLANES
Definition: polyMesh.H:100
Foam::readDir
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:660
Foam::point
vector point
Point is a vector.
Definition: point.H:41
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210