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-2024 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 
73  label nLocateBoundaryHits = 0;
74 
76  {
77  label patchi = iter.key();
78  const fvPatch& patch = mesh_.boundary()[patchi];
79 
80  vectorField nf(patch.nf());
81 
82  forAll(patch, patchFacei)
83  {
84  cloud.addParticle
85  (
87  (
88  mesh_,
89  patch.Cf()[patchFacei],
90  patch.faceCells()[patchFacei],
91  nLocateBoundaryHits,
92  - distance_*nf[patchFacei],
93  globalWalls.toGlobal(nPatchFaces) // passive data
94  )
95  );
96 
97  nPatchFaces++;
98  }
99  }
100 
101 
102 
103  if (debug)
104  {
105  // Dump particles
106  OBJstream str
107  (
108  time_.path()
109  /"wantedTracks_" + time_.name() + ".obj"
110  );
111  InfoInFunction << "Dumping tracks to " << str.name() << endl;
112 
114  {
115  const vector p = iter().position(mesh_);
116  str.write(linePointRef(p, p + iter().displacement()));
117  }
118  }
119 
120 
121 
122  // Per cell: empty or global wall index and end location
124  cellToSamples_.setSize(mesh_.nCells());
125 
126  // Database to pass into findCellParticle::move
128 
129  // Debug: collect start points
130  pointField start;
131  if (debug)
132  {
133  start.setSize(nPatchFaces);
134  nPatchFaces = 0;
136  {
137  const findCellParticle& tp = iter();
138  start[nPatchFaces++] = tp.position(mesh_);
139  }
140  }
141 
142 
143  // Track
144  cloud.move(cloud, td);
145 
146 
147  // Rework cell-to-globalpatchface into a map
148  List<Map<label>> compactMap;
149  getPatchDataMapPtr_.reset
150  (
151  new distributionMap
152  (
153  globalWalls,
154  cellToWalls_,
155  compactMap
156  )
157  );
158 
159 
160  // Debug: dump resulting tracks
161  if (debug)
162  {
163  getPatchDataMapPtr_().distribute(start);
164  {
165  OBJstream str
166  (
167  time_.path()
168  /"obtainedTracks_" + time_.name() + ".obj"
169  );
170  InfoInFunction << "Dumping obtained to " << str.name() << endl;
171 
172  forAll(cellToWalls_, celli)
173  {
174  const List<point>& ends = cellToSamples_[celli];
175  const labelList& cData = cellToWalls_[celli];
176  forAll(cData, i)
177  {
178  str.write(linePointRef(ends[i], start[cData[i]]));
179  }
180  }
181  }
182  }
183 }
184 
185 
186 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
187 
189 (
190  const word& name,
191  const Time& runTime,
192  const dictionary& dict
193 )
194 :
195  fvMeshFunctionObject(name, runTime, dict),
196  fieldSet_()
197 {
198  read(dict);
199 }
200 
201 
202 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
203 
205 {
207 }
208 
209 
210 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
211 
213 {
215 
216  dict.lookup("fields") >> fieldSet_;
217  patchSet_ = patchSet(dict);
218  distance_ = dict.lookup<scalar>("distance");
219 
220 
221  // Clear out any previously loaded fields
222  vsf_.clear();
223  vvf_.clear();
224  vSpheretf_.clear();
225  vSymmtf_.clear();
226  vtf_.clear();
227  fieldMap_.clear();
228  reverseFieldMap_.clear();
229 
230 
231  // Generate fields with mappedField boundary condition
232 
233  // Convert field to map
234  fieldMap_.resize(2*fieldSet_.size());
235  reverseFieldMap_.resize(2*fieldSet_.size());
236  forAll(fieldSet_, seti)
237  {
238  const word& fieldName = fieldSet_[seti].first();
239  const word& sampleFieldName = fieldSet_[seti].second();
240 
241  fieldMap_.insert(fieldName, sampleFieldName);
242  reverseFieldMap_.insert(sampleFieldName, fieldName);
243  }
244 
245  Log << type() << " " << name()
246  << ": Sampling " << fieldMap_.size() << " fields" << endl;
247 
248  // Do analysis
249  calcAddressing();
250 
251  return true;
252 }
253 
254 
256 {
257  wordList fields(fieldSet_.size());
258 
259  forAll(fieldSet_, fieldi)
260  {
261  fields[fieldi] = fieldSet_[fieldi].first();
262  }
263 
264  return fields;
265 }
266 
267 
269 {
271 
272  if
273  (
274  fieldMap_.size()
275  && vsf_.empty()
276  && vvf_.empty()
277  && vSpheretf_.empty()
278  && vSymmtf_.empty()
279  && vtf_.empty()
280  )
281  {
282  Log << type() << " " << name()
283  << ": Creating " << fieldMap_.size() << " fields" << endl;
284 
285  createFields(vsf_);
286  createFields(vvf_);
287  createFields(vSpheretf_);
288  createFields(vSymmtf_);
289  createFields(vtf_);
290 
291  Log << endl;
292  }
293 
294  Log << type() << " " << name()
295  << " write:" << nl
296  << " Sampling fields to " << time_.name()
297  << endl;
298 
299  sampleFields(vsf_);
300  sampleFields(vvf_);
301  sampleFields(vSpheretf_);
302  sampleFields(vSymmtf_);
303  sampleFields(vtf_);
304 
305  return true;
306 }
307 
308 
310 {
312 
313  Log << " Writing sampled fields to " << time_.name()
314  << endl;
315 
316  forAll(vsf_, i)
317  {
318  vsf_[i].write();
319  }
320  forAll(vvf_, i)
321  {
322  vvf_[i].write();
323  }
324  forAll(vSpheretf_, i)
325  {
326  vSpheretf_[i].write();
327  }
328  forAll(vSymmtf_, i)
329  {
330  vSymmtf_[i].write();
331  }
332  forAll(vtf_, i)
333  {
334  vtf_[i].write();
335  }
336 
337  return true;
338 }
339 
340 
341 // ************************************************************************* //
#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:162
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:857
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:228
#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:257
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
static const char nl
Definition: Ostream.H:266
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