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-2017 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 // * * * * * * * * * * * * Protected 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  mesh_,
68  );
69 
70  // Add particles to track to sample locations
71  nPatchFaces = 0;
72 
73  forAllConstIter(labelHashSet, patchSet_, iter)
74  {
75  label patchi = iter.key();
76  const fvPatch& patch = mesh_.boundary()[patchi];
77 
78  vectorField nf(patch.nf());
79  vectorField faceCellCentres(patch.patch().faceCellCentres());
80 
81  forAll(patch, patchFacei)
82  {
83  label meshFacei = patch.start()+patchFacei;
84 
85  // Find starting point on face (since faceCentre might not
86  // be on face-diagonal decomposition)
87  pointIndexHit startInfo
88  (
90  (
91  mesh_,
92  meshFacei,
94  )
95  );
96 
97 
98  point start;
99  if (startInfo.hit())
100  {
101  start = startInfo.hitPoint();
102  }
103  else
104  {
105  // Fallback: start tracking from neighbouring cell centre
106  start = faceCellCentres[patchFacei];
107  }
108 
109  const point end = start-distance_*nf[patchFacei];
110 
111  // Add a particle to the cloud with originating face as passive data
112  cloud.addParticle
113  (
114  new findCellParticle
115  (
116  mesh_,
117  start,
118  -1,
119  end,
120  globalWalls.toGlobal(nPatchFaces) // passive data
121  )
122  );
123 
124  nPatchFaces++;
125  }
126  }
127 
128 
129 
130  if (debug)
131  {
132  // Dump particles
133  OBJstream str
134  (
135  mesh_.time().path()
136  /"wantedTracks_" + mesh_.time().timeName() + ".obj"
137  );
138  InfoInFunction << "Dumping tracks to " << str.name() << endl;
139 
141  {
142  const findCellParticle& tp = iter();
143  str.write(linePointRef(tp.position(), tp.end()));
144  }
145  }
146 
147 
148 
149  // Per cell: empty or global wall index and end location
150  cellToWalls_.setSize(mesh_.nCells());
151  cellToSamples_.setSize(mesh_.nCells());
152 
153  // Database to pass into findCellParticle::move
154  findCellParticle::trackingData td(cloud, cellToWalls_, cellToSamples_);
155 
156  // Track all particles to their end position.
157  scalar maxTrackLen = 2.0*mesh_.bounds().mag();
158 
159 
160  //Debug: collect start points
161  pointField start;
162  if (debug)
163  {
164  start.setSize(nPatchFaces);
165  nPatchFaces = 0;
167  {
168  const findCellParticle& tp = iter();
169  start[nPatchFaces++] = tp.position();
170  }
171  }
172 
173 
174  cloud.move(td, maxTrackLen);
175 
176 
177  // Rework cell-to-globalpatchface into a map
178  List<Map<label>> compactMap;
179  getPatchDataMapPtr_.reset
180  (
181  new mapDistribute
182  (
183  globalWalls,
184  cellToWalls_,
185  compactMap
186  )
187  );
188 
189 
190  // Debug: dump resulting tracks
191  if (debug)
192  {
193  getPatchDataMapPtr_().distribute(start);
194  {
195  OBJstream str
196  (
197  mesh_.time().path()
198  /"obtainedTracks_" + mesh_.time().timeName() + ".obj"
199  );
200  InfoInFunction << "Dumping obtained to " << str.name() << endl;
201 
202  forAll(cellToWalls_, celli)
203  {
204  const List<point>& ends = cellToSamples_[celli];
205  const labelList& cData = cellToWalls_[celli];
206  forAll(cData, i)
207  {
208  str.write(linePointRef(ends[i], start[cData[i]]));
209  }
210  }
211  }
212  }
213 }
214 
215 
216 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
217 
218 Foam::functionObjects::nearWallFields::nearWallFields
219 (
220  const word& name,
221  const Time& runTime,
222  const dictionary& dict
223 )
224 :
225  fvMeshFunctionObject(name, runTime, dict),
226  fieldSet_()
227 {
228  read(dict);
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
233 
235 {
237 }
238 
239 
240 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
241 
243 {
245 
246  dict.lookup("fields") >> fieldSet_;
247  patchSet_ =
248  mesh_.boundaryMesh().patchSet(wordReList(dict.lookup("patches")));
249  distance_ = readScalar(dict.lookup("distance"));
250 
251 
252  // Clear out any previously loaded fields
253  vsf_.clear();
254  vvf_.clear();
255  vSpheretf_.clear();
256  vSymmtf_.clear();
257  vtf_.clear();
258  fieldMap_.clear();
259  reverseFieldMap_.clear();
260 
261 
262  // Generate fields with mappedField boundary condition
263 
264  // Convert field to map
265  fieldMap_.resize(2*fieldSet_.size());
266  reverseFieldMap_.resize(2*fieldSet_.size());
267  forAll(fieldSet_, setI)
268  {
269  const word& fldName = fieldSet_[setI].first();
270  const word& sampleFldName = fieldSet_[setI].second();
271 
272  fieldMap_.insert(fldName, sampleFldName);
273  reverseFieldMap_.insert(sampleFldName, fldName);
274  }
275 
276  Log << type() << " " << name()
277  << ": Sampling " << fieldMap_.size() << " fields" << endl;
278 
279  // Do analysis
280  calcAddressing();
281 
282  return true;
283 }
284 
285 
287 {
289 
290  if
291  (
292  fieldMap_.size()
293  && vsf_.empty()
294  && vvf_.empty()
295  && vSpheretf_.empty()
296  && vSymmtf_.empty()
297  && vtf_.empty()
298  )
299  {
300  Log << type() << " " << name()
301  << ": Creating " << fieldMap_.size() << " fields" << endl;
302 
303  createFields(vsf_);
304  createFields(vvf_);
305  createFields(vSpheretf_);
306  createFields(vSymmtf_);
307  createFields(vtf_);
308 
309  Log << endl;
310  }
311 
312  Log << type() << " " << name()
313  << " write:" << nl
314  << " Sampling fields to " << time_.timeName()
315  << endl;
316 
317  sampleFields(vsf_);
318  sampleFields(vvf_);
319  sampleFields(vSpheretf_);
320  sampleFields(vSymmtf_);
321  sampleFields(vtf_);
322 
323  return true;
324 }
325 
326 
328 {
330 
331  Log << " Writing sampled fields to " << time_.timeName()
332  << endl;
333 
334  forAll(vsf_, i)
335  {
336  vsf_[i].write();
337  }
338  forAll(vvf_, i)
339  {
340  vvf_[i].write();
341  }
342  forAll(vSpheretf_, i)
343  {
344  vSpheretf_[i].write();
345  }
346  forAll(vSymmtf_, i)
347  {
348  vSymmtf_[i].write();
349  }
350  forAll(vtf_, i)
351  {
352  vtf_[i].write();
353  }
354 
355  return true;
356 }
357 
358 
359 // ************************************************************************* //
Template class for intrusive linked lists.
Definition: ILList.H:50
#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
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:60
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:124
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
Macros for easy insertion into run-time selection tables.
virtual bool execute()
Calculate the near-wall fields.
const Point & hitPoint() const
Return hit point.
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:141
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
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
#define DebugInFunction
Report an information message using Foam::Info.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
label size() const
Global sum of localSizes.
Definition: globalIndexI.H:66
bool hit() const
Is there a hit.
tmp< vectorField > faceCellCentres() const
Return face cell centres.
Definition: polyPatch.C:296
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:185
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
static const char nl
Definition: Ostream.H:262
static word defaultName
The default cloud name: defaultCloud.
Definition: cloud.H:74
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
const point & end() const
Point to track to.
void setSize(const label)
Reset size of List.
Definition: List.C:281
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
label patchi
Class containing processor-to-processor mapping information.
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
Particle class that finds cells by tracking.
#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
Class used to pass tracking data to the trackToFace function.
vector position() const
Return current particle position.
Definition: particleI.H:204
addToRunTimeSelectionTable(functionObject, add, dictionary)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
#define InfoInFunction
Report an information message using Foam::Info.