nearWallFields.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 "nearWallFields.H"
27 #include "wordReList.H"
28 #include "findCellParticle.H"
29 #include "mappedPatchBase.H"
30 #include "OBJstream.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(nearWallFields, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
43 {
44  const fvMesh& mesh = refCast<const fvMesh>(obr_);
45 
46  // Count number of faces
47  label nPatchFaces = 0;
48  forAllConstIter(labelHashSet, patchSet_, iter)
49  {
50  label patchI = iter.key();
51  nPatchFaces += mesh.boundary()[patchI].size();
52  }
53 
54  // Global indexing
55  globalIndex globalWalls(nPatchFaces);
56 
57  if (debug)
58  {
59  Info<< "nearWallFields::calcAddressing() :"
60  << " nPatchFaces:" << globalWalls.size() << endl;
61  }
62 
63  // Construct cloud
65 
66  // Add particles to track to sample locations
67  nPatchFaces = 0;
68 
69  forAllConstIter(labelHashSet, patchSet_, iter)
70  {
71  label patchI = iter.key();
72  const fvPatch& patch = mesh.boundary()[patchI];
73 
74  vectorField nf(patch.nf());
75  vectorField faceCellCentres(patch.patch().faceCellCentres());
76 
77  forAll(patch, patchFaceI)
78  {
79  label meshFaceI = patch.start()+patchFaceI;
80 
81  // Find starting point on face (since faceCentre might not
82  // be on face-diagonal decomposition)
83  pointIndexHit startInfo
84  (
86  (
87  mesh,
88  meshFaceI,
90  )
91  );
92 
93 
94  point start;
95  if (startInfo.hit())
96  {
97  start = startInfo.hitPoint();
98  }
99  else
100  {
101  // Fallback: start tracking from neighbouring cell centre
102  start = faceCellCentres[patchFaceI];
103  }
104 
105  const point end = start-distance_*nf[patchFaceI];
106 
107  // Find tet for starting location
108  label cellI = -1;
109  label tetFaceI = -1;
110  label tetPtI = -1;
111  mesh.findCellFacePt(start, cellI, tetFaceI, tetPtI);
112 
113  // Add to cloud. Add originating face as passive data
114  cloud.addParticle
115  (
116  new findCellParticle
117  (
118  mesh,
119  start,
120  cellI,
121  tetFaceI,
122  tetPtI,
123  end,
124  globalWalls.toGlobal(nPatchFaces) // passive data
125  )
126  );
127 
128  nPatchFaces++;
129  }
130  }
131 
132 
133 
134  if (debug)
135  {
136  // Dump particles
137  OBJstream str
138  (
139  mesh.time().path()
140  /"wantedTracks_" + mesh.time().timeName() + ".obj"
141  );
142  Info<< "nearWallFields::calcAddressing() :"
143  << "Dumping tracks to " << str.name() << endl;
144 
146  {
147  const findCellParticle& tp = iter();
148  str.write(linePointRef(tp.position(), tp.end()));
149  }
150  }
151 
152 
153 
154  // Per cell: empty or global wall index and end location
155  cellToWalls_.setSize(mesh.nCells());
156  cellToSamples_.setSize(mesh.nCells());
157 
158  // Database to pass into findCellParticle::move
159  findCellParticle::trackingData td(cloud, cellToWalls_, cellToSamples_);
160 
161  // Track all particles to their end position.
162  scalar maxTrackLen = 2.0*mesh.bounds().mag();
163 
164 
165  //Debug: collect start points
166  pointField start;
167  if (debug)
168  {
169  start.setSize(nPatchFaces);
170  nPatchFaces = 0;
172  {
173  const findCellParticle& tp = iter();
174  start[nPatchFaces++] = tp.position();
175  }
176  }
177 
178 
179  cloud.move(td, maxTrackLen);
180 
181 
182  // Rework cell-to-globalpatchface into a map
183  List<Map<label> > compactMap;
184  getPatchDataMapPtr_.reset
185  (
186  new mapDistribute
187  (
188  globalWalls,
189  cellToWalls_,
190  compactMap
191  )
192  );
193 
194 
195  // Debug: dump resulting tracks
196  if (debug)
197  {
198  getPatchDataMapPtr_().distribute(start);
199  {
200  OBJstream str
201  (
202  mesh.time().path()
203  /"obtainedTracks_" + mesh.time().timeName() + ".obj"
204  );
205  Info<< "nearWallFields::calcAddressing() :"
206  << "Dumping obtained to " << str.name() << endl;
207 
208  forAll(cellToWalls_, cellI)
209  {
210  const List<point>& ends = cellToSamples_[cellI];
211  const labelList& cData = cellToWalls_[cellI];
212  forAll(cData, i)
213  {
214  str.write(linePointRef(ends[i], start[cData[i]]));
215  }
216  }
217  }
218  }
219 }
220 
221 
222 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
223 
224 Foam::nearWallFields::nearWallFields
225 (
226  const word& name,
227  const objectRegistry& obr,
228  const dictionary& dict,
229  const bool loadFromFiles
230 )
231 :
232  name_(name),
233  obr_(obr),
234  active_(true),
235  fieldSet_()
236 {
237  // Check if the available mesh is an fvMesh otherise deactivate
238  if (isA<fvMesh>(obr_))
239  {
240  read(dict);
241  }
242  else
243  {
244  active_ = false;
245  WarningIn
246  (
247  "nearWallFields::nearWallFields"
248  "("
249  "const word&, "
250  "const objectRegistry&, "
251  "const dictionary&, "
252  "const bool"
253  ")"
254  ) << "No fvMesh available, deactivating " << name_
255  << endl;
256  }
257 
258 }
259 
260 
261 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
262 
264 {
265  if (debug)
266  {
267  Info<< "nearWallFields::~nearWallFields()" << endl;
268  }
269 }
270 
271 
272 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
273 
275 {
276  if (debug)
277  {
278  Info<< "nearWallFields::read(const dictionary&)" << endl;
279  }
280 
281  if (active_)
282  {
283  const fvMesh& mesh = refCast<const fvMesh>(obr_);
284 
285  dict.lookup("fields") >> fieldSet_;
286  patchSet_ =
287  mesh.boundaryMesh().patchSet(wordReList(dict.lookup("patches")));
288  distance_ = readScalar(dict.lookup("distance"));
289 
290 
291  // Clear out any previously loaded fields
292  vsf_.clear();
293  vvf_.clear();
294  vSpheretf_.clear();
295  vSymmtf_.clear();
296  vtf_.clear();
297  fieldMap_.clear();
298  reverseFieldMap_.clear();
299 
300 
301  // Generate fields with mappedField boundary condition
302 
303  // Convert field to map
304  fieldMap_.resize(2*fieldSet_.size());
305  reverseFieldMap_.resize(2*fieldSet_.size());
306  forAll(fieldSet_, setI)
307  {
308  const word& fldName = fieldSet_[setI].first();
309  const word& sampleFldName = fieldSet_[setI].second();
310 
311  fieldMap_.insert(fldName, sampleFldName);
312  reverseFieldMap_.insert(sampleFldName, fldName);
313  }
314 
315  Info<< type() << " " << name_ << ": Sampling " << fieldMap_.size()
316  << " fields" << endl;
317 
318  // Do analysis
319  calcAddressing();
320  }
321 }
322 
323 
325 {
326  if (debug)
327  {
328  Info<< "nearWallFields:execute()" << endl;
329  }
330 
331 
332  if (active_)
333  {
334  if
335  (
336  fieldMap_.size()
337  && vsf_.empty()
338  && vvf_.empty()
339  && vSpheretf_.empty()
340  && vSymmtf_.empty()
341  && vtf_.empty()
342  )
343  {
344  Info<< type() << " " << name_ << ": Creating " << fieldMap_.size()
345  << " fields" << endl;
346 
347  createFields(vsf_);
348  createFields(vvf_);
349  createFields(vSpheretf_);
350  createFields(vSymmtf_);
351  createFields(vtf_);
352 
353  Info<< endl;
354  }
355 
356  Info<< type() << " " << name_ << " output:" << nl;
357 
358  Info<< " Sampling fields to " << obr_.time().timeName()
359  << endl;
360 
361  sampleFields(vsf_);
362  sampleFields(vvf_);
363  sampleFields(vSpheretf_);
364  sampleFields(vSymmtf_);
365  sampleFields(vtf_);
366  }
367 }
368 
369 
371 {
372  if (debug)
373  {
374  Info<< "nearWallFields:end()" << endl;
375  }
376 
377  if (active_)
378  {
379  execute();
380  }
381 }
382 
383 
385 {
386  // Do nothing
387 }
388 
389 
391 {
392  if (debug)
393  {
394  Info<< "nearWallFields:write()" << endl;
395  }
396 
397  if (active_)
398  {
399  Info<< " Writing sampled fields to " << obr_.time().timeName()
400  << endl;
401 
402  forAll(vsf_, i)
403  {
404  vsf_[i].write();
405  }
406  forAll(vvf_, i)
407  {
408  vvf_[i].write();
409  }
410  forAll(vSpheretf_, i)
411  {
412  vSpheretf_[i].write();
413  }
414  forAll(vSymmtf_, i)
415  {
416  vSymmtf_[i].write();
417  }
418  forAll(vtf_, i)
419  {
420  vtf_[i].write();
421  }
422 
423  Info<< endl;
424  }
425 }
426 
427 
428 // ************************************************************************* //
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
label size() const
Global sum of localSizes.
Definition: globalIndexI.H:66
virtual void read(const dictionary &)
Read the field min/max data.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
fileName path() const
Return path.
Definition: Time.H:281
virtual void end()
Execute at the final time-loop, currently does nothing.
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:90
Class containing processor-to-processor mapping information.
const point & end() const
Point to track to.
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
A class for handling words, derived from string.
Definition: word.H:59
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
label nCells() const
messageStream Info
bool hit() const
Is there a hit.
void move(TrackData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:187
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
const Point & hitPoint() const
Return hit point.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
static const char nl
Definition: Ostream.H:260
void setSize(const label)
Reset size of List.
Definition: List.C:318
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
#define WarningIn(functionName)
Report a warning using Foam::Warning.
Particle class that finds cells by tracking.
const vector & position() const
Return current particle position.
Definition: particleI.H:603
tmp< vectorField > faceCellCentres() const
Return face cell centres.
Definition: polyPatch.C:297
#define forAll(list, i)
Definition: UList.H:421
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
virtual ~nearWallFields()
Destructor.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:124
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
Class used to pass tracking data to the trackToFace function.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:427
Registry of regIOobjects.
bool read(const char *, int32_t &)
Definition: int32IO.C:87
virtual void write()
Write.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
void calcAddressing()
Calculate addressing from cells back to patch faces.
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:162
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
virtual void execute()
Execute, currently does nothing.
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1261
static pointIndexHit facePoint(const polyMesh &, const label faceI, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
defineTypeNameAndDebug(combustionModel, 0)
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:552