nearWallFields.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 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 "OBJstream.H"
30 #include "globalIndex.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
49  // Count number of faces
50  label nPatchFaces = 0;
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 
74  {
75  label patchi = iter.key();
76  const fvPatch& patch = mesh_.boundary()[patchi];
77 
78  vectorField nf(patch.nf());
79 
80  forAll(patch, patchFacei)
81  {
82  cloud.addParticle
83  (
85  (
86  mesh_,
87  patch.Cf()[patchFacei],
88  patch.faceCells()[patchFacei],
89  - distance_*nf[patchFacei],
90  globalWalls.toGlobal(nPatchFaces) // passive data
91  )
92  );
93 
94  nPatchFaces++;
95  }
96  }
97 
98 
99 
100  if (debug)
101  {
102  // Dump particles
103  OBJstream str
104  (
105  time_.path()
106  /"wantedTracks_" + time_.name() + ".obj"
107  );
108  InfoInFunction << "Dumping tracks to " << str.name() << endl;
109 
111  {
112  const vector p = iter().position(mesh_);
113  str.write(linePointRef(p, p + iter().displacement()));
114  }
115  }
116 
117 
118 
119  // Per cell: empty or global wall index and end location
121  cellToSamples_.setSize(mesh_.nCells());
122 
123  // Database to pass into findCellParticle::move
125 
126  // Debug: collect start points
127  pointField start;
128  if (debug)
129  {
130  start.setSize(nPatchFaces);
131  nPatchFaces = 0;
133  {
134  const findCellParticle& tp = iter();
135  start[nPatchFaces++] = tp.position(mesh_);
136  }
137  }
138 
139 
140  // Track
141  cloud.move(cloud, td);
142 
143 
144  // Rework cell-to-globalpatchface into a map
145  List<Map<label>> compactMap;
146  getPatchDataMapPtr_.reset
147  (
148  new distributionMap
149  (
150  globalWalls,
151  cellToWalls_,
152  compactMap
153  )
154  );
155 
156 
157  // Debug: dump resulting tracks
158  if (debug)
159  {
160  getPatchDataMapPtr_().distribute(start);
161  {
162  OBJstream str
163  (
164  time_.path()
165  /"obtainedTracks_" + time_.name() + ".obj"
166  );
167  InfoInFunction << "Dumping obtained to " << str.name() << endl;
168 
169  forAll(cellToWalls_, celli)
170  {
171  const List<point>& ends = cellToSamples_[celli];
172  const labelList& cData = cellToWalls_[celli];
173  forAll(cData, i)
174  {
175  str.write(linePointRef(ends[i], start[cData[i]]));
176  }
177  }
178  }
179  }
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
184 
186 (
187  const word& name,
188  const Time& runTime,
189  const dictionary& dict
190 )
191 :
192  fvMeshFunctionObject(name, runTime, dict),
193  fieldSet_()
194 {
195  read(dict);
196 }
197 
198 
199 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
200 
202 {
204 }
205 
206 
207 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
208 
210 {
212 
213  dict.lookup("fields") >> fieldSet_;
214  patchSet_ =
215  mesh_.boundaryMesh().patchSet(wordReList(dict.lookup("patches")));
216  distance_ = dict.lookup<scalar>("distance");
217 
218 
219  // Clear out any previously loaded fields
220  vsf_.clear();
221  vvf_.clear();
222  vSpheretf_.clear();
223  vSymmtf_.clear();
224  vtf_.clear();
225  fieldMap_.clear();
226  reverseFieldMap_.clear();
227 
228 
229  // Generate fields with mappedField boundary condition
230 
231  // Convert field to map
232  fieldMap_.resize(2*fieldSet_.size());
233  reverseFieldMap_.resize(2*fieldSet_.size());
234  forAll(fieldSet_, seti)
235  {
236  const word& fldName = fieldSet_[seti].first();
237  const word& sampleFldName = fieldSet_[seti].second();
238 
239  fieldMap_.insert(fldName, sampleFldName);
240  reverseFieldMap_.insert(sampleFldName, fldName);
241  }
242 
243  Log << type() << " " << name()
244  << ": Sampling " << fieldMap_.size() << " fields" << endl;
245 
246  // Do analysis
247  calcAddressing();
248 
249  return true;
250 }
251 
252 
254 {
255  wordList fields(fieldSet_.size());
256 
257  forAll(fieldSet_, fieldi)
258  {
259  fields[fieldi] = fieldSet_[fieldi].first();
260  }
261 
262  return fields;
263 }
264 
265 
267 {
269 
270  if
271  (
272  fieldMap_.size()
273  && vsf_.empty()
274  && vvf_.empty()
275  && vSpheretf_.empty()
276  && vSymmtf_.empty()
277  && vtf_.empty()
278  )
279  {
280  Log << type() << " " << name()
281  << ": Creating " << fieldMap_.size() << " fields" << endl;
282 
283  createFields(vsf_);
284  createFields(vvf_);
285  createFields(vSpheretf_);
286  createFields(vSymmtf_);
287  createFields(vtf_);
288 
289  Log << endl;
290  }
291 
292  Log << type() << " " << name()
293  << " write:" << nl
294  << " Sampling fields to " << time_.name()
295  << endl;
296 
297  sampleFields(vsf_);
298  sampleFields(vvf_);
299  sampleFields(vSpheretf_);
300  sampleFields(vSymmtf_);
301  sampleFields(vtf_);
302 
303  return true;
304 }
305 
306 
308 {
310 
311  Log << " Writing sampled fields to " << time_.name()
312  << endl;
313 
314  forAll(vsf_, i)
315  {
316  vsf_[i].write();
317  }
318  forAll(vvf_, i)
319  {
320  vvf_[i].write();
321  }
322  forAll(vSpheretf_, i)
323  {
324  vSpheretf_[i].write();
325  }
326  forAll(vSymmtf_, i)
327  {
328  vSymmtf_[i].write();
329  }
330  forAll(vtf_, i)
331  {
332  vtf_[i].write();
333  }
334 
335  return true;
336 }
337 
338 
339 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
Base cloud calls templated on particle type.
Definition: Cloud.H:74
Template class for intrusive linked lists.
Definition: ILList.H:67
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:91
void setSize(const label)
Reset size of List.
Definition: List.C:281
OFstream which keeps track of vertices.
Definition: OBJstream.H:56
virtual Ostream & write(const char)
Write character.
Definition: OBJstream.C:81
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:120
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
fileName path() const
Explicitly inherit path from TimePaths to disambiguate from.
Definition: TimePaths.H:138
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
A cloud is a collection of lagrangian particles.
Definition: cloud.H:55
static word defaultName
The default cloud name: defaultCloud.
Definition: cloud.H:66
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const word & name() const
Return const reference to name.
Class containing processor-to-processor mapping information.
Class used to pass tracking data to the trackToFace function.
Particle class that finds cells by tracking.
Abstract base-class for Time/database functionObjects.
const Time & time_
Reference to time.
Specialisation 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.
scalar distance_
Distance away from wall.
virtual wordList fields() const
Return the list of fields required.
List< List< point > > cellToSamples_
From cell to tracked end point.
labelListList cellToWalls_
From cell to seed patch faces.
void calcAddressing()
Calculate addressing from cells back to patch faces.
labelHashSet patchSet_
Patches to sample.
nearWallFields(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
autoPtr< distributionMap > getPatchDataMapPtr_
Map from cell based data back to patch based data.
virtual bool execute()
Calculate the near-wall fields.
virtual bool write()
Write the near-wall fields.
virtual bool read(const dictionary &)
Read the controls.
virtual bool read(const dictionary &)
Read optional controls.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:130
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:105
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:99
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
label size() const
Global sum of localSizes.
Definition: globalIndexI.H:66
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
vector position(const polyMesh &mesh) const
Return current particle position.
Definition: particleI.H:255
label nCells() const
A class for handling words, derived from string.
Definition: word.H:62
label patchi
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
#define Log
Report write to Foam::Info if the local log switch is true.
#define DebugInFunction
Report an information message using Foam::Info.
#define InfoInFunction
Report an information message using Foam::Info.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
volScalarField & p