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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2015-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "nearWallFields.H"
30 #include "wordRes.H"
31 #include "findCellParticle.H"
32 #include "mappedPatchBase.H"
33 #include "OBJstream.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace functionObjects
41 {
42  defineTypeNameAndDebug(nearWallFields, 0);
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
51 {
52  // Count number of faces
53  label nPatchFaces = 0;
54  for (const label patchi : patchIDs_)
55  {
56  nPatchFaces += mesh_.boundary()[patchi].size();
57  }
58 
59  // Global indexing
60  const globalIndex globalWalls(nPatchFaces);
61 
62  DebugInFunction << "nPatchFaces: " << globalWalls.totalSize() << endl;
63 
64  // Start with empty cloud
65  Cloud<findCellParticle> cloud(mesh_, Foam::zero{}, cloud::defaultName);
66 
67  // Add particles to track to sample locations
68  nPatchFaces = 0;
69 
70  for (const label patchi : patchIDs_)
71  {
72  const fvPatch& patch = mesh_.boundary()[patchi];
73 
74  const vectorField nf(patch.nf());
75  const vectorField faceCellCentres(patch.patch().faceCellCentres());
76  const labelUList& faceCells = patch.patch().faceCells();
77  const vectorField::subField faceCentres = patch.patch().faceCentres();
78 
79  forAll(patch, patchFacei)
80  {
81  label meshFacei = patch.start()+patchFacei;
82 
83  // Find starting point on face (since faceCentre might not
84  // be on face-diagonal decomposition)
85  pointIndexHit startInfo
86  (
88  (
89  mesh_,
90  meshFacei,
92  )
93  );
94 
95 
96  // Starting point and tet
97  point start;
98  const label celli = faceCells[patchFacei];
99 
100  if (startInfo.hit())
101  {
102  // Move start point slightly in so it is inside the tet
103  const face& f = mesh_.faces()[meshFacei];
104 
105  label tetFacei = meshFacei;
106  label tetPti = (startInfo.index()+1) % f.size();
107 
108  start = startInfo.point();
109 
110  // Uncomment below to shift slightly in:
111  tetIndices tet(celli, tetFacei, tetPti);
112  start =
113  (1.0 - 1e-6)*startInfo.point()
114  + 1e-6*tet.tet(mesh_).centre();
115 
116  // Re-check that we have a valid location
117  mesh_.findTetFacePt(celli, start, tetFacei, tetPti);
118  if (tetFacei == -1)
119  {
120  start = faceCellCentres[patchFacei];
121  }
122  }
123  else
124  {
125  // Fallback: start tracking from neighbouring cell centre
126  start = faceCellCentres[patchFacei];
127  }
128 
129 
130  // TBD: why use start? and not faceCentres[patchFacei]
131  //const point end = start-distance_*nf[patchFacei];
132  const point end = faceCentres[patchFacei]-distance_*nf[patchFacei];
133 
134 
135  // Add a particle to the cloud with originating face as passive data
136  cloud.addParticle
137  (
138  new findCellParticle
139  (
140  mesh_,
141  start,
142  celli,
143  end,
144  globalWalls.toGlobal(nPatchFaces) // passive data
145  )
146  );
147 
148  nPatchFaces++;
149  }
150  }
151 
152 
153 
154  if (debug)
155  {
156  // Dump particles
157  OBJstream str
158  (
159  mesh_.time().path()
160  /"wantedTracks_" + mesh_.time().timeName() + ".obj"
161  );
162  InfoInFunction << "Dumping tracks to " << str.name() << endl;
163 
164  for (const findCellParticle& tp : cloud)
165  {
166  str.writeLine(tp.position(), tp.end());
167  }
168  }
169 
170 
171 
172  // Per cell: empty or global wall index and end location
174  cellToSamples_.setSize(mesh_.nCells());
175 
176  // Database to pass into findCellParticle::move
177  findCellParticle::trackingData td(cloud, cellToWalls_, cellToSamples_);
178 
179  // Track all particles to their end position.
180  scalar maxTrackLen = 2.0*mesh_.bounds().mag();
181 
182 
183  //Debug: collect start points
184  pointField start;
185  if (debug)
186  {
187  start.setSize(nPatchFaces);
188  nPatchFaces = 0;
189  for (const findCellParticle& tp : cloud)
190  {
191  start[nPatchFaces++] = tp.position();
192  }
193  }
194 
195 
196  cloud.move(cloud, td, maxTrackLen);
197 
198 
199  // Rework cell-to-globalpatchface into a map
200  List<Map<label>> compactMap;
201  getPatchDataMapPtr_.reset
202  (
203  new mapDistribute
204  (
205  globalWalls,
206  cellToWalls_,
207  compactMap
208  )
209  );
210 
211 
212  // Debug: dump resulting tracks
213  if (debug)
214  {
215  getPatchDataMapPtr_().distribute(start);
216  {
217  OBJstream str
218  (
219  mesh_.time().path()
220  /"obtainedTracks_" + mesh_.time().timeName() + ".obj"
221  );
222  InfoInFunction << "Dumping obtained to " << str.name() << endl;
223 
224  forAll(cellToWalls_, celli)
225  {
226  const List<point>& ends = cellToSamples_[celli];
227  const labelList& cData = cellToWalls_[celli];
228  forAll(cData, i)
229  {
230  str.writeLine(ends[i], start[cData[i]]);
231  }
232  }
233  }
234  }
235 }
236 
237 
238 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
239 
241 (
242  const word& name,
243  const Time& runTime,
244  const dictionary& dict
245 )
246 :
247  fvMeshFunctionObject(name, runTime, dict),
248  fieldSet_()
249 {
250  read(dict);
251 }
252 
253 
254 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
255 
257 {
258  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
259 
261 
262  dict.readEntry("fields", fieldSet_);
263  dict.readEntry("distance", distance_);
264 
265  // Can also use pbm.indices(), but no warnings...
266  patchIDs_ = pbm.patchSet(dict.get<wordRes>("patches")).sortedToc();
267 
268 
269  // Clear out any previously loaded fields
270  vsf_.clear();
271  vvf_.clear();
272  vSpheretf_.clear();
273  vSymmtf_.clear();
274  vtf_.clear();
275  fieldMap_.clear();
276  reverseFieldMap_.clear();
277 
278 
279  // Generate fields with mappedField boundary condition
280 
281  // Convert field to map
282  fieldMap_.reserve(fieldSet_.size());
283  reverseFieldMap_.reserve(fieldSet_.size());
284  forAll(fieldSet_, seti)
285  {
286  const word& fldName = fieldSet_[seti].first();
287  const word& sampleFldName = fieldSet_[seti].second();
288 
289  fieldMap_.insert(fldName, sampleFldName);
290  reverseFieldMap_.insert(sampleFldName, fldName);
291  }
292 
293  Info<< type() << " " << name()
294  << ": Sampling " << fieldMap_.size() << " fields" << endl;
295 
296  // Do analysis
297  calcAddressing();
298 
299  return true;
300 }
301 
302 
304 {
306 
307  if
308  (
309  fieldMap_.size()
310  && vsf_.empty()
311  && vvf_.empty()
312  && vSpheretf_.empty()
313  && vSymmtf_.empty()
314  && vtf_.empty()
315  )
316  {
317  Log << type() << " " << name()
318  << ": Creating " << fieldMap_.size() << " fields" << endl;
319 
320  createFields(vsf_);
321  createFields(vvf_);
322  createFields(vSpheretf_);
323  createFields(vSymmtf_);
324  createFields(vtf_);
325 
326  Log << endl;
327  }
328 
329  Log << type() << " " << name()
330  << " write:" << nl
331  << " Sampling fields to " << time_.timeName()
332  << endl;
333 
334  sampleFields(vsf_);
335  sampleFields(vvf_);
336  sampleFields(vSpheretf_);
337  sampleFields(vSymmtf_);
338  sampleFields(vtf_);
339 
340  return true;
341 }
342 
343 
345 {
347 
348  Log << " Writing sampled fields to " << time_.timeName()
349  << endl;
350 
351  forAll(vsf_, i)
352  {
353  vsf_[i].write();
354  }
355  forAll(vvf_, i)
356  {
357  vvf_[i].write();
358  }
359  forAll(vSpheretf_, i)
360  {
361  vSpheretf_[i].write();
362  }
363  forAll(vSymmtf_, i)
364  {
365  vSymmtf_[i].write();
366  }
367  forAll(vtf_, i)
368  {
369  vtf_[i].write();
370  }
371 
372  return true;
373 }
374 
375 
376 // ************************************************************************* //
static word defaultName
The default cloud name: defaultCloud.
Definition: cloud.H:84
const polyBoundaryMesh & pbm
dictionary dict
defineTypeNameAndDebug(ObukhovLength, 0)
virtual bool write()
Write the near-wall fields.
fileName path() const
Return path = rootPath/caseName. Same as TimePaths::path()
Definition: Time.H:503
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
labelList patchIDs_
Patches to sample.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
void calcAddressing()
Calculate addressing from cells back to patch faces.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Abstract base-class for Time/database function objects.
List< List< point > > cellToSamples_
From cell to tracked end point.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual bool end()
Called when Time::run() determines that the time-loop exits.
virtual bool execute()
Calculate the near-wall fields.
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual bool read(const dictionary &dict)
Read the controls.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:316
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
autoPtr< mapDistribute > getPatchDataMapPtr_
Map from cell based data back to patch based data.
A class for handling words, derived from Foam::string.
Definition: word.H:63
#define DebugInFunction
Report an information message using Foam::Info.
scalar mag() const
The magnitude/length of the bounding box diagonal.
Definition: boundBoxI.H:198
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
labelListList cellToWalls_
From cell to seed patch faces.
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
nearWallFields(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
labelList f(nPoints)
Samples near-patch volume fields within an input distance range.
vector point
Point is a vector.
Definition: point.H:37
label nCells() const noexcept
Number of mesh cells.
#define Log
Definition: PDRblock.C:28
const std::string patch
OpenFOAM patch number as a std::string.
void findTetFacePt(const label celli, const point &p, label &tetFacei, label &tetPti) const
Find the tetFacei and tetPti for point p in celli.
Definition: polyMesh.C:1385
messageStream Info
Information stream (stdout output on master, null elsewhere)
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
virtual bool read(const dictionary &dict)
Read optional controls.
SubField< vector > subField
Declare type of subField.
Definition: Field.H:128
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
const boundBox & bounds() const noexcept
Return mesh bounding box.
Definition: polyMesh.H:616
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
scalar distance_
Distance away from wall.
const fvMesh & mesh_
Reference to the fvMesh.
Namespace for OpenFOAM.
static int debug
Flag to execute debug content.
#define InfoInFunction
Report an information message using Foam::Info.