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-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 "nearWallFields.H"
27 #include "wordReList.H"
28 #include "findCellParticle.H"
29 #include "mappedPatchBase.H"
30 #include "OBJstream.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
39  defineTypeNameAndDebug(nearWallFields, 0);
40  addToRunTimeSelectionTable(functionObject, nearWallFields, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
46 
48 {
49  // Count number of faces
50  label nPatchFaces = 0;
51  forAllConstIter(labelHashSet, patchSet_, iter)
52  {
53  label patchi = iter.key();
54  nPatchFaces += mesh_.boundary()[patchi].size();
55  }
56 
57  // Global indexing
58  globalIndex globalWalls(nPatchFaces);
59 
60  DebugInFunction << "nPatchFaces: " << globalWalls.size() << endl;
61 
62  // Construct cloud
64 
65  // Add particles to track to sample locations
66  nPatchFaces = 0;
67 
68  forAllConstIter(labelHashSet, patchSet_, iter)
69  {
70  label patchi = iter.key();
71  const fvPatch& patch = mesh_.boundary()[patchi];
72 
73  vectorField nf(patch.nf());
74  vectorField faceCellCentres(patch.patch().faceCellCentres());
75 
76  forAll(patch, patchFacei)
77  {
78  label meshFacei = patch.start()+patchFacei;
79 
80  // Find starting point on face (since faceCentre might not
81  // be on face-diagonal decomposition)
82  pointIndexHit startInfo
83  (
85  (
86  mesh_,
87  meshFacei,
89  )
90  );
91 
92 
93  point start;
94  if (startInfo.hit())
95  {
96  start = startInfo.hitPoint();
97  }
98  else
99  {
100  // Fallback: start tracking from neighbouring cell centre
101  start = faceCellCentres[patchFacei];
102  }
103 
104  const point end = start-distance_*nf[patchFacei];
105 
106  // Find tet for starting location
107  label celli = -1;
108  label tetFacei = -1;
109  label tetPtI = -1;
110  mesh_.findCellFacePt(start, celli, tetFacei, tetPtI);
111 
112  // Add to cloud. Add originating face as passive data
113  cloud.addParticle
114  (
115  new findCellParticle
116  (
117  mesh_,
118  start,
119  celli,
120  tetFacei,
121  tetPtI,
122  end,
123  globalWalls.toGlobal(nPatchFaces) // passive data
124  )
125  );
126 
127  nPatchFaces++;
128  }
129  }
130 
131 
132 
133  if (debug)
134  {
135  // Dump particles
136  OBJstream str
137  (
138  mesh_.time().path()
139  /"wantedTracks_" + mesh_.time().timeName() + ".obj"
140  );
141  InfoInFunction << "Dumping tracks to " << str.name() << endl;
142 
144  {
145  const findCellParticle& tp = iter();
146  str.write(linePointRef(tp.position(), tp.end()));
147  }
148  }
149 
150 
151 
152  // Per cell: empty or global wall index and end location
153  cellToWalls_.setSize(mesh_.nCells());
154  cellToSamples_.setSize(mesh_.nCells());
155 
156  // Database to pass into findCellParticle::move
157  findCellParticle::trackingData td(cloud, cellToWalls_, cellToSamples_);
158 
159  // Track all particles to their end position.
160  scalar maxTrackLen = 2.0*mesh_.bounds().mag();
161 
162 
163  //Debug: collect start points
164  pointField start;
165  if (debug)
166  {
167  start.setSize(nPatchFaces);
168  nPatchFaces = 0;
170  {
171  const findCellParticle& tp = iter();
172  start[nPatchFaces++] = tp.position();
173  }
174  }
175 
176 
177  cloud.move(td, maxTrackLen);
178 
179 
180  // Rework cell-to-globalpatchface into a map
181  List<Map<label>> compactMap;
182  getPatchDataMapPtr_.reset
183  (
184  new mapDistribute
185  (
186  globalWalls,
187  cellToWalls_,
188  compactMap
189  )
190  );
191 
192 
193  // Debug: dump resulting tracks
194  if (debug)
195  {
196  getPatchDataMapPtr_().distribute(start);
197  {
198  OBJstream str
199  (
200  mesh_.time().path()
201  /"obtainedTracks_" + mesh_.time().timeName() + ".obj"
202  );
203  InfoInFunction << "Dumping obtained to " << str.name() << endl;
204 
205  forAll(cellToWalls_, celli)
206  {
207  const List<point>& ends = cellToSamples_[celli];
208  const labelList& cData = cellToWalls_[celli];
209  forAll(cData, i)
210  {
211  str.write(linePointRef(ends[i], start[cData[i]]));
212  }
213  }
214  }
215  }
216 }
217 
218 
219 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
220 
221 Foam::functionObjects::nearWallFields::nearWallFields
222 (
223  const word& name,
224  const Time& runTime,
225  const dictionary& dict
226 )
227 :
228  fvMeshFunctionObject(name, runTime, dict),
229  fieldSet_()
230 {
231  read(dict);
232 }
233 
234 
235 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
236 
238 {
240 }
241 
242 
243 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
244 
246 {
248 
249  dict.lookup("fields") >> fieldSet_;
250  patchSet_ =
251  mesh_.boundaryMesh().patchSet(wordReList(dict.lookup("patches")));
252  distance_ = readScalar(dict.lookup("distance"));
253 
254 
255  // Clear out any previously loaded fields
256  vsf_.clear();
257  vvf_.clear();
258  vSpheretf_.clear();
259  vSymmtf_.clear();
260  vtf_.clear();
261  fieldMap_.clear();
262  reverseFieldMap_.clear();
263 
264 
265  // Generate fields with mappedField boundary condition
266 
267  // Convert field to map
268  fieldMap_.resize(2*fieldSet_.size());
269  reverseFieldMap_.resize(2*fieldSet_.size());
270  forAll(fieldSet_, setI)
271  {
272  const word& fldName = fieldSet_[setI].first();
273  const word& sampleFldName = fieldSet_[setI].second();
274 
275  fieldMap_.insert(fldName, sampleFldName);
276  reverseFieldMap_.insert(sampleFldName, fldName);
277  }
278 
279  Log << type() << " " << name()
280  << ": Sampling " << fieldMap_.size() << " fields" << endl;
281 
282  // Do analysis
283  calcAddressing();
284 
285  return true;
286 }
287 
288 
290 {
292 
293  if
294  (
295  fieldMap_.size()
296  && vsf_.empty()
297  && vvf_.empty()
298  && vSpheretf_.empty()
299  && vSymmtf_.empty()
300  && vtf_.empty()
301  )
302  {
303  Log << type() << " " << name()
304  << ": Creating " << fieldMap_.size() << " fields" << endl;
305 
306  createFields(vsf_);
307  createFields(vvf_);
308  createFields(vSpheretf_);
309  createFields(vSymmtf_);
310  createFields(vtf_);
311 
312  Log << endl;
313  }
314 
315  Log << type() << " " << name()
316  << " write:" << nl
317  << " Sampling fields to " << time_.timeName()
318  << endl;
319 
320  sampleFields(vsf_);
321  sampleFields(vvf_);
322  sampleFields(vSpheretf_);
323  sampleFields(vSymmtf_);
324  sampleFields(vtf_);
325 
326  return true;
327 }
328 
329 
331 {
333 
334  Log << " Writing sampled fields to " << time_.timeName()
335  << endl;
336 
337  forAll(vsf_, i)
338  {
339  vsf_[i].write();
340  }
341  forAll(vvf_, i)
342  {
343  vvf_[i].write();
344  }
345  forAll(vSpheretf_, i)
346  {
347  vSpheretf_[i].write();
348  }
349  forAll(vSymmtf_, i)
350  {
351  vSymmtf_[i].write();
352  }
353  forAll(vtf_, i)
354  {
355  vtf_[i].write();
356  }
357 
358  Log << endl;
359 
360  return true;
361 }
362 
363 
364 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual bool write()
Write the near-wall fields.
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
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
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
bool hit() const
Is there a hit.
void calcAddressing()
Calculate addressing from cells back to patch faces.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
Macros for easy insertion into run-time selection tables.
const point & end() const
Point to track to.
label size() const
Global sum of localSizes.
Definition: globalIndexI.H:66
virtual bool execute()
Calculate the near-wall fields.
const vector & position() const
Return current particle position.
Definition: particleI.H:586
virtual bool read(const dictionary &)
Read optional controls.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:162
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
A class for handling words, derived from string.
Definition: word.H:59
#define DebugInFunction
Report an information message using Foam::Info.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
void move(TrackData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:187
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
static const char nl
Definition: Ostream.H:262
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:124
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void setSize(const label)
Reset size of List.
Definition: List.C:295
label patchi
Class containing processor-to-processor mapping information.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
addToRunTimeSelectionTable(functionObject, blendingFactor, dictionary)
Particle class that finds cells by tracking.
const Point & hitPoint() const
Return hit point.
#define Log
Report write to Foam::Info if the local log switch is true.
virtual bool read(const dictionary &)
Read the controls.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
tmp< vectorField > faceCellCentres() const
Return face cell centres.
Definition: polyPatch.C:296
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
Class used to pass tracking data to the trackToFace function.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
#define InfoInFunction
Report an information message using Foam::Info.