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-2018 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 "mappedPatchBase.H"
30 #include "OBJstream.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
39  defineTypeNameAndDebug(nearWallFields, 0);
40  addToRunTimeSelectionTable(functionObject, nearWallFields, dictionary);
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
49  // Count number of faces
50  label nPatchFaces = 0;
51  forAllConstIter(labelHashSet, patchSet_, iter)
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  forAllConstIter(labelHashSet, patchSet_, iter)
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  mesh_.time().path()
106  /"wantedTracks_" + mesh_.time().timeName() + ".obj"
107  );
108  InfoInFunction << "Dumping tracks to " << str.name() << endl;
109 
111  {
112  const vector p = iter().position();
113  str.write(linePointRef(p, p + iter().displacement()));
114  }
115  }
116 
117 
118 
119  // Per cell: empty or global wall index and end location
120  cellToWalls_.setSize(mesh_.nCells());
121  cellToSamples_.setSize(mesh_.nCells());
122 
123  // Database to pass into findCellParticle::move
124  findCellParticle::trackingData td(cloud, cellToWalls_, cellToSamples_);
125 
126  // Track all particles to their end position.
127  scalar maxTrackLen = 2.0*mesh_.bounds().mag();
128 
129 
130  // Debug: collect start points
131  pointField start;
132  if (debug)
133  {
134  start.setSize(nPatchFaces);
135  nPatchFaces = 0;
137  {
138  const findCellParticle& tp = iter();
139  start[nPatchFaces++] = tp.position();
140  }
141  }
142 
143 
144  cloud.move(cloud, td, maxTrackLen);
145 
146 
147  // Rework cell-to-globalpatchface into a map
148  List<Map<label>> compactMap;
149  getPatchDataMapPtr_.reset
150  (
151  new mapDistribute
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  mesh_.time().path()
168  /"obtainedTracks_" + mesh_.time().timeName() + ".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 
188 Foam::functionObjects::nearWallFields::nearWallFields
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_ =
218  mesh_.boundaryMesh().patchSet(wordReList(dict.lookup("patches")));
219  distance_ = readScalar(dict.lookup("distance"));
220 
221 
222  // Clear out any previously loaded fields
223  vsf_.clear();
224  vvf_.clear();
225  vSpheretf_.clear();
226  vSymmtf_.clear();
227  vtf_.clear();
228  fieldMap_.clear();
229  reverseFieldMap_.clear();
230 
231 
232  // Generate fields with mappedField boundary condition
233 
234  // Convert field to map
235  fieldMap_.resize(2*fieldSet_.size());
236  reverseFieldMap_.resize(2*fieldSet_.size());
237  forAll(fieldSet_, setI)
238  {
239  const word& fldName = fieldSet_[setI].first();
240  const word& sampleFldName = fieldSet_[setI].second();
241 
242  fieldMap_.insert(fldName, sampleFldName);
243  reverseFieldMap_.insert(sampleFldName, fldName);
244  }
245 
246  Log << type() << " " << name()
247  << ": Sampling " << fieldMap_.size() << " fields" << endl;
248 
249  // Do analysis
250  calcAddressing();
251 
252  return true;
253 }
254 
255 
257 {
259 
260  if
261  (
262  fieldMap_.size()
263  && vsf_.empty()
264  && vvf_.empty()
265  && vSpheretf_.empty()
266  && vSymmtf_.empty()
267  && vtf_.empty()
268  )
269  {
270  Log << type() << " " << name()
271  << ": Creating " << fieldMap_.size() << " fields" << endl;
272 
273  createFields(vsf_);
274  createFields(vvf_);
275  createFields(vSpheretf_);
276  createFields(vSymmtf_);
277  createFields(vtf_);
278 
279  Log << endl;
280  }
281 
282  Log << type() << " " << name()
283  << " write:" << nl
284  << " Sampling fields to " << time_.timeName()
285  << endl;
286 
287  sampleFields(vsf_);
288  sampleFields(vvf_);
289  sampleFields(vSpheretf_);
290  sampleFields(vSymmtf_);
291  sampleFields(vtf_);
292 
293  return true;
294 }
295 
296 
298 {
300 
301  Log << " Writing sampled fields to " << time_.timeName()
302  << endl;
303 
304  forAll(vsf_, i)
305  {
306  vsf_[i].write();
307  }
308  forAll(vvf_, i)
309  {
310  vvf_[i].write();
311  }
312  forAll(vSpheretf_, i)
313  {
314  vSpheretf_[i].write();
315  }
316  forAll(vSymmtf_, i)
317  {
318  vSymmtf_[i].write();
319  }
320  forAll(vtf_, i)
321  {
322  vtf_[i].write();
323  }
324 
325  return true;
326 }
327 
328 
329 // ************************************************************************* //
Template class for intrusive linked lists.
Definition: ILList.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual bool write()
Write the near-wall fields.
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:124
void calcAddressing()
Calculate addressing from cells back to patch faces.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
virtual bool execute()
Calculate the near-wall fields.
virtual bool read(const dictionary &)
Read optional controls.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void addParticle(ParticleType *pPtr)
Transfer particle to cloud.
Definition: Cloud.C:96
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
A class for handling words, derived from string.
Definition: word.H:59
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
#define DebugInFunction
Report an information message using Foam::Info.
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
label size() const
Global sum of localSizes.
Definition: globalIndexI.H:66
const vectorField & Cf() const
Return face centres.
Definition: fvPatch.C:99
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Base cloud calls templated on particle type.
Definition: Cloud.H:52
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
static const char nl
Definition: Ostream.H:265
static word defaultName
The default cloud name: defaultCloud.
Definition: cloud.H:74
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
void setSize(const label)
Reset size of List.
Definition: List.C:281
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
void move(TrackCloudType &cloud, typename ParticleType::trackingData &td, const scalar trackTime)
Move the particles.
Definition: Cloud.C:141
label patchi
Class containing processor-to-processor mapping information.
Particle class that finds cells by tracking.
#define Log
Report write to Foam::Info if the local log switch is true.
virtual bool read(const dictionary &)
Read the controls.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
volScalarField & p
Class used to pass tracking data to the trackToFace function.
vector position() const
Return current particle position.
Definition: particleI.H:283
addToRunTimeSelectionTable(functionObject, add, dictionary)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
#define InfoInFunction
Report an information message using Foam::Info.