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-------------------------------------------------------------------------------
11License
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
38namespace Foam
39{
40namespace functionObjects
41{
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.totalSize() << endl;
63
64 // Construct cloud
66 (
67 mesh_,
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 (
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
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;
207 (
208 new mapDistribute
209 (
210 globalWalls,
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:
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// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Base cloud calls templated on particle type.
Definition: Cloud.H:68
Template class for intrusive linked lists.
Definition: ILList.H:69
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void setSize(const label n)
Alias for resize()
Definition: List.H:218
OFstream that keeps track of vertices.
Definition: OBJstream.H:61
virtual Ostream & write(const char c)
Write character.
Definition: OBJstream.C:78
virtual const fileName & name() const
Read/write access to the name of the stream.
Definition: OSstream.H:107
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:66
label index() const noexcept
Return the hit index.
bool hit() const noexcept
Is there a hit?
const point_type & hitPoint() const
Return hit point. Fatal if not hit.
virtual bool read()
Re-read model coefficients if they have changed.
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: SubField.H:62
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
fileName path() const
Return path.
Definition: Time.H:358
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:133
A cloud is a registry collection of lagrangian particles.
Definition: cloud.H:60
static word defaultName
The default cloud name: defaultCloud.
Definition: cloud.H:90
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Class used to pass tracking data to the trackToFace function.
Particle class that finds cells by tracking.
Abstract base-class for Time/database function objects.
static int debug
Flag to execute debug content.
virtual bool end()
Called when Time::run() determines that the time-loop exits.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
Samples near-patch volume fields within an input distance range.
scalar distance_
Distance away from wall.
autoPtr< mapDistribute > getPatchDataMapPtr_
Map from cell based data back to patch based data.
List< List< point > > cellToSamples_
From cell to tracked end point.
virtual bool read(const dictionary &dict)
Read the controls.
labelListList cellToWalls_
From cell to seed patch faces.
void calcAddressing()
Calculate addressing from cells back to patch faces.
labelHashSet patchSet_
Patches to sample.
virtual bool execute()
Calculate the near-wall fields.
virtual bool write()
Write the near-wall fields.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:712
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:260
label totalSize() const
Global sum of localSizes.
Definition: globalIndexI.H:125
Class containing processor-to-processor mapping information.
static pointIndexHit facePoint(const polyMesh &, const label facei, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:462
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:1396
label nCells() const noexcept
Number of mesh cells.
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:84
tetPointRef tet(const polyMesh &mesh) const
Return the geometry corresponding to this tet.
Definition: tetIndicesI.H:134
Point centre() const
Return centre (centroid)
Definition: tetrahedronI.H:165
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
engineTime & runTime
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
Namespace for OpenFOAM.
line< point, const point & > linePointRef
A line using referred points.
Definition: linePointRef.H:47
messageStream Info
Information stream (stdout output on master, null elsewhere)
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333