Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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"
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 : 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
173 cellToWalls_.setSize(mesh_.nCells());
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;
202 (
203 new mapDistribute
204 (
205 globalWalls,
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 }
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_()
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
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_);
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// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const polyBoundaryMesh & pbm
Base cloud calls templated on particle type.
Definition Cloud.H:64
SubField< vector > subField
Definition Field.H:183
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void setSize(label n)
Alias for resize().
Definition List.H:536
An OFstream that keeps track of vertices and provides convenience output methods for OBJ files.
Definition OBJstream.H:58
Ostream & writeLine(const point &p0, const point &p1)
Write line joining two points.
Definition OBJstream.C:234
virtual const fileName & name() const override
Read/write access to the name of the stream.
Definition OSstream.H:134
label index() const noexcept
Return the hit index.
bool hit() const noexcept
Is there a hit?
const point_type & point() const noexcept
Return point, no checks.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A cloud is a registry collection of lagrangian particles.
Definition cloud.H:56
static word defaultName
The default cloud name: defaultCloud.
Definition cloud.H:84
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
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.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
static int debug
Flag to execute debug content.
virtual bool end()
Called when Time::run() determines that the time-loop exits.
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
Samples near-patch volume fields within an input distance range.
scalar distance_
Distance away from wall.
HashTable< word > reverseFieldMap_
From resulting back to original field.
autoPtr< mapDistribute > getPatchDataMapPtr_
Map from cell based data back to patch based data.
List< Tuple2< word, word > > fieldSet_
Fields to process (input-name output-name).
List< List< point > > cellToSamples_
From cell to tracked end point.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
labelListList cellToWalls_
From cell to seed patch faces.
void calcAddressing()
Calculate addressing from cells back to patch faces.
void createFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &) const
Create the local fields.
void sampleFields(PtrList< GeometricField< Type, fvPatchField, volMesh > > &) const
Sample the fields.
PtrList< volSymmTensorField > vSymmtf_
nearWallFields(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
virtual bool execute()
Execute the function-object operations.
PtrList< volSphericalTensorField > vSpheretf_
virtual bool write()
Write the function-object results.
labelList patchIDs_
Patches to sample.
HashTable< word > fieldMap_
From original field to sampled result.
const Time & time_
Reference to the time database.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
label totalSize() const noexcept
The total addressed size, which corresponds to the end offset and also the sum of all localSizes.
label toGlobal(const label proci, const label i) const
From local to global on proci.
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:
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition tetIndices.H:79
tetPointRef tet(const polyMesh &mesh) const
The tet geometry for this tet, where point0 is the cell centre.
Point centre() const
Return centre (centroid).
A List of wordRe with additional matching capabilities.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
engineTime & runTime
auto & name
wallPoints::trackData td(isBlockedFace, regionToBlockSize)
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
List< label > labelList
A List of labels.
Definition List.H:62
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 POSIX.C:801
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299