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-2020 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 : patchSet_)
55  {
56  nPatchFaces += mesh_.boundary()[patchi].size();
57  }
58 
59  // Global indexing
60  globalIndex globalWalls(nPatchFaces);
61 
62  DebugInFunction << "nPatchFaces: " << globalWalls.size() << endl;
63 
64  // Construct cloud
65  Cloud<findCellParticle> cloud
66  (
67  mesh_,
69  IDLList<findCellParticle>()
70  );
71 
72  // Add particles to track to sample locations
73  nPatchFaces = 0;
74 
75  for (const label patchi : patchSet_)
76  {
77  const fvPatch& patch = mesh_.boundary()[patchi];
78 
79  const vectorField nf(patch.nf());
80  const vectorField faceCellCentres(patch.patch().faceCellCentres());
81  const labelUList& faceCells = patch.patch().faceCells();
82  const vectorField::subField& faceCentres = patch.patch().faceCentres();
83 
84  forAll(patch, patchFacei)
85  {
86  label meshFacei = patch.start()+patchFacei;
87 
88  // Find starting point on face (since faceCentre might not
89  // be on face-diagonal decomposition)
90  pointIndexHit startInfo
91  (
93  (
94  mesh_,
95  meshFacei,
97  )
98  );
99 
100 
101  // Starting point and tet
102  point start;
103  const label celli = faceCells[patchFacei];
104 
105  if (startInfo.hit())
106  {
107  // Move start point slightly in so it is inside the tet
108  const face& f = mesh_.faces()[meshFacei];
109 
110  label tetFacei = meshFacei;
111  label tetPti = (startInfo.index()+1) % f.size();
112 
113  start = startInfo.hitPoint();
114 
115  // Uncomment below to shift slightly in:
116  tetIndices tet(celli, tetFacei, tetPti);
117  start =
118  (1.0 - 1e-6)*startInfo.hitPoint()
119  + 1e-6*tet.tet(mesh_).centre();
120 
121  // Re-check that we have a valid location
122  mesh_.findTetFacePt(celli, start, tetFacei, tetPti);
123  if (tetFacei == -1)
124  {
125  start = faceCellCentres[patchFacei];
126  }
127  }
128  else
129  {
130  // Fallback: start tracking from neighbouring cell centre
131  start = faceCellCentres[patchFacei];
132  }
133 
134 
135  // TBD: why use start? and not faceCentres[patchFacei]
136  //const point end = start-distance_*nf[patchFacei];
137  const point end = faceCentres[patchFacei]-distance_*nf[patchFacei];
138 
139 
140  // Add a particle to the cloud with originating face as passive data
141  cloud.addParticle
142  (
143  new findCellParticle
144  (
145  mesh_,
146  start,
147  celli,
148  end,
149  globalWalls.toGlobal(nPatchFaces) // passive data
150  )
151  );
152 
153  nPatchFaces++;
154  }
155  }
156 
157 
158 
159  if (debug)
160  {
161  // Dump particles
162  OBJstream str
163  (
164  mesh_.time().path()
165  /"wantedTracks_" + mesh_.time().timeName() + ".obj"
166  );
167  InfoInFunction << "Dumping tracks to " << str.name() << endl;
168 
169  for (const findCellParticle& tp : cloud)
170  {
171  str.write(linePointRef(tp.position(), tp.end()));
172  }
173  }
174 
175 
176 
177  // Per cell: empty or global wall index and end location
179  cellToSamples_.setSize(mesh_.nCells());
180 
181  // Database to pass into findCellParticle::move
182  findCellParticle::trackingData td(cloud, cellToWalls_, cellToSamples_);
183 
184  // Track all particles to their end position.
185  scalar maxTrackLen = 2.0*mesh_.bounds().mag();
186 
187 
188  //Debug: collect start points
189  pointField start;
190  if (debug)
191  {
192  start.setSize(nPatchFaces);
193  nPatchFaces = 0;
194  for (const findCellParticle& tp : cloud)
195  {
196  start[nPatchFaces++] = tp.position();
197  }
198  }
199 
200 
201  cloud.move(cloud, td, maxTrackLen);
202 
203 
204  // Rework cell-to-globalpatchface into a map
205  List<Map<label>> compactMap;
206  getPatchDataMapPtr_.reset
207  (
208  new mapDistribute
209  (
210  globalWalls,
211  cellToWalls_,
212  compactMap
213  )
214  );
215 
216 
217  // Debug: dump resulting tracks
218  if (debug)
219  {
220  getPatchDataMapPtr_().distribute(start);
221  {
222  OBJstream str
223  (
224  mesh_.time().path()
225  /"obtainedTracks_" + mesh_.time().timeName() + ".obj"
226  );
227  InfoInFunction << "Dumping obtained to " << str.name() << endl;
228 
229  forAll(cellToWalls_, celli)
230  {
231  const List<point>& ends = cellToSamples_[celli];
232  const labelList& cData = cellToWalls_[celli];
233  forAll(cData, i)
234  {
235  str.write(linePointRef(ends[i], start[cData[i]]));
236  }
237  }
238  }
239  }
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
244 
246 (
247  const word& name,
248  const Time& runTime,
249  const dictionary& dict
250 )
251 :
252  fvMeshFunctionObject(name, runTime, dict),
253  fieldSet_()
254 {
255  read(dict);
256 }
257 
258 
259 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
260 
262 {
264 
265  dict.readEntry("fields", fieldSet_);
266  dict.readEntry("distance", distance_);
267 
268  patchSet_ =
269  mesh_.boundaryMesh().patchSet
270  (
271  dict.get<wordRes>("patches")
272  );
273 
274 
275  // Clear out any previously loaded fields
276  vsf_.clear();
277  vvf_.clear();
278  vSpheretf_.clear();
279  vSymmtf_.clear();
280  vtf_.clear();
281  fieldMap_.clear();
282  reverseFieldMap_.clear();
283 
284 
285  // Generate fields with mappedField boundary condition
286 
287  // Convert field to map
288  fieldMap_.resize(2*fieldSet_.size());
289  reverseFieldMap_.resize(2*fieldSet_.size());
290  forAll(fieldSet_, seti)
291  {
292  const word& fldName = fieldSet_[seti].first();
293  const word& sampleFldName = fieldSet_[seti].second();
294 
295  fieldMap_.insert(fldName, sampleFldName);
296  reverseFieldMap_.insert(sampleFldName, fldName);
297  }
298 
299  Info<< type() << " " << name()
300  << ": Sampling " << fieldMap_.size() << " fields" << endl;
301 
302  // Do analysis
303  calcAddressing();
304 
305  return true;
306 }
307 
308 
310 {
312 
313  if
314  (
315  fieldMap_.size()
316  && vsf_.empty()
317  && vvf_.empty()
318  && vSpheretf_.empty()
319  && vSymmtf_.empty()
320  && vtf_.empty()
321  )
322  {
323  Log << type() << " " << name()
324  << ": Creating " << fieldMap_.size() << " fields" << endl;
325 
326  createFields(vsf_);
327  createFields(vvf_);
328  createFields(vSpheretf_);
329  createFields(vSymmtf_);
330  createFields(vtf_);
331 
332  Log << endl;
333  }
334 
335  Log << type() << " " << name()
336  << " write:" << nl
337  << " Sampling fields to " << time_.timeName()
338  << endl;
339 
340  sampleFields(vsf_);
341  sampleFields(vvf_);
342  sampleFields(vSpheretf_);
343  sampleFields(vSymmtf_);
344  sampleFields(vtf_);
345 
346  return true;
347 }
348 
349 
351 {
353 
354  Log << " Writing sampled fields to " << time_.timeName()
355  << endl;
356 
357  forAll(vsf_, i)
358  {
359  vsf_[i].write();
360  }
361  forAll(vvf_, i)
362  {
363  vvf_[i].write();
364  }
365  forAll(vSpheretf_, i)
366  {
367  vSpheretf_[i].write();
368  }
369  forAll(vSymmtf_, i)
370  {
371  vSymmtf_[i].write();
372  }
373  forAll(vtf_, i)
374  {
375  vtf_[i].write();
376  }
377 
378  return true;
379 }
380 
381 
382 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:63
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
runTime
engineTime & runTime
Definition: createEngineTime.H:13
wordRes.H
InfoInFunction
#define InfoInFunction
Definition: messageStream.H:387
Foam::tetIndices::tet
tetPointRef tet(const polyMesh &mesh) const
Definition: tetIndicesI.H:148
Foam::functionObjects::nearWallFields::cellToSamples_
List< List< point > > cellToSamples_
Definition: nearWallFields.H:191
Foam::polyMesh::FACE_DIAG_TRIS
@ FACE_DIAG_TRIS
Definition: polyMesh.H:103
Log
#define Log
Definition: PDRblock.C:28
Foam::boundBox::mag
scalar mag() const
Definition: boundBoxI.H:126
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::OBJstream
OFstream that keeps track of vertices.
Definition: OBJstream.H:54
Foam::OFstream::name
virtual const fileName & name() const
Definition: OSstream.H:103
Foam::read
bool read(const char *buf, int32_t &val)
Definition: int32.H:125
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::functionObjects::nearWallFields::write
virtual bool write()
Definition: nearWallFields.C:343
Foam::functionObjects::nearWallFields::getPatchDataMapPtr_
autoPtr< mapDistribute > getPatchDataMapPtr_
Definition: nearWallFields.H:194
Foam::functionObjects::nearWallFields
Samples near-patch volume fields within an input distance range.
Definition: nearWallFields.H:159
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::PointIndexHit::hitPoint
const point_type & hitPoint() const
Definition: PointIndexHit.H:150
Foam::OBJstream::write
virtual Ostream & write(const char c)
Definition: OBJstream.C:71
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:60
Foam::functionObject
Abstract base-class for Time/database function objects.
Definition: functionObject.H:329
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:48
nearWallFields.H
Foam::SubField
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: Field.H:60
Foam::primitiveMesh::nCells
label nCells() const noexcept
Definition: primitiveMeshI.H:89
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:48
Foam::PointIndexHit::hit
bool hit() const noexcept
Definition: PointIndexHit.H:126
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::functionObjects::nearWallFields::calcAddressing
void calcAddressing()
Definition: nearWallFields.C:43
Foam::Info
messageStream Info
DebugInFunction
#define DebugInFunction
Definition: messageStream.H:445
Foam::List::setSize
void setSize(const label n)
Definition: List.H:218
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:159
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::functionObjects::nearWallFields::patchSet_
labelHashSet patchSet_
Definition: nearWallFields.H:173
findCellParticle.H
Foam::findCellParticle
Particle class that finds cells by tracking.
Definition: findCellParticle.H:55
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::ILList
Template class for intrusive linked lists.
Definition: ILList.H:48
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Definition: regionFunctionObject.C:166
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::Ostream::write
virtual bool write(const token &tok)=0
Foam
Definition: atmBoundaryLayer.C:26
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
Foam::cloud::defaultName
static word defaultName
Definition: cloud.H:86
Foam::functionObject::end
virtual bool end()
Definition: functionObject.C:185
Foam::functionObjects::nearWallFields::execute
virtual bool execute()
Definition: nearWallFields.C:302
Foam::functionObjects::nearWallFields::distance_
scalar distance_
Definition: nearWallFields.H:176
Foam::functionObjects::nearWallFields::read
virtual bool read(const dictionary &dict)
Definition: nearWallFields.C:254
Foam::polyMesh::bounds
const boundBox & bounds() const
Definition: polyMesh.H:446
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Definition: fvMesh.C:678
Foam::tetIndices
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:79
Foam::cloud
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:53
Foam::linePointRef
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:41
Foam::polyMesh::faces
virtual const faceList & faces() const
Definition: polyMesh.C:1087
Foam::pointIndexHit
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:40
Foam::nl
constexpr char nl
Definition: Ostream.H:424
Foam::Time::path
fileName path() const
Definition: Time.H:354
Foam::functionObject::debug
static int debug
Definition: functionObject.H:364
f
labelList f(nPoints)
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::globalIndex::size
label size() const
Definition: globalIndexI.H:144
Foam::foamVersion::patch
const std::string patch
Foam::Vector< scalar >
Foam::functionObjects::nearWallFields::cellToWalls_
labelListList cellToWalls_
Definition: nearWallFields.H:188
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::PointIndexHit::index
label index() const noexcept
Definition: PointIndexHit.H:132
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Definition: POSIX.C:717
Foam::functionObjects::fvMeshFunctionObject::mesh_
const fvMesh & mesh_
Definition: fvMeshFunctionObject.H:69
Foam::tetrahedron::centre
Point centre() const
Definition: tetrahedronI.H:158
Foam::functionObjects::nearWallFields::nearWallFields
nearWallFields(const word &name, const Time &runTime, const dictionary &dict)
Definition: nearWallFields.C:239
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:99
Foam::Cloud
Base cloud calls templated on particle type.
Definition: Cloud.H:51
Foam::constant::electromagnetic::e
const dimensionedScalar e
Definition: createFields.H:11
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:47
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::name
word name(const expressions::valueTypeCode typeCode)
Definition: exprTraits.C:52
Foam::fvMesh::time
const Time & time() const
Definition: fvMesh.H:276
Foam::Field::subField
SubField< Type > subField
Definition: Field.H:85
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
Foam::point
vector point
Point is a vector.
Definition: point.H:37
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:81
Foam::polyMesh::findTetFacePt
void findTetFacePt(const label celli, const point &p, label &tetFacei, label &tetPti) const
Definition: polyMesh.C:1374
Foam::mappedPatchBase::facePoint
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Definition: mappedPatchBase.C:1716
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
OBJstream.H
Foam::findCellParticle::trackingData
Definition: findCellParticle.H:76
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
Definition: globalIndexI.H:233
mappedPatchBase.H