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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(nearWallFields, 0);
37 }
38 
39 
40 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
41 
43 {
44  const fvMesh& mesh = refCast<const fvMesh>(obr_);
45 
46  // Count number of faces
47  label nPatchFaces = 0;
49  {
50  label patchI = iter.key();
51  nPatchFaces += mesh.boundary()[patchI].size();
52  }
53 
54  // Global indexing
55  globalIndex globalWalls(nPatchFaces);
56 
57  if (debug)
58  {
59  Info<< "nearWallFields::calcAddressing() :"
60  << " nPatchFaces:" << globalWalls.size() << endl;
61  }
62 
63  // Construct cloud
65 
66  // Add particles to track to sample locations
67  nPatchFaces = 0;
68 
70  {
71  label patchI = iter.key();
72  const fvPatch& patch = mesh.boundary()[patchI];
73 
74  vectorField nf(patch.nf());
75  vectorField faceCellCentres(patch.patch().faceCellCentres());
76 
77  forAll(patch, patchFaceI)
78  {
79  label meshFaceI = patch.start()+patchFaceI;
80 
81  // Find starting point on face (since faceCentre might not
82  // be on face-diagonal decomposition)
83  pointIndexHit startInfo
84  (
86  (
87  mesh,
88  meshFaceI,
90  )
91  );
92 
93 
94  point start;
95  if (startInfo.hit())
96  {
97  start = startInfo.hitPoint();
98  }
99  else
100  {
101  // Fallback: start tracking from neighbouring cell centre
102  start = faceCellCentres[patchFaceI];
103  }
104 
105  const point end = start-distance_*nf[patchFaceI];
106 
107  // Find tet for starting location
108  label cellI = -1;
109  label tetFaceI = -1;
110  label tetPtI = -1;
111  mesh.findCellFacePt(start, cellI, tetFaceI, tetPtI);
112 
113  // Add to cloud. Add originating face as passive data
114  cloud.addParticle
115  (
116  new findCellParticle
117  (
118  mesh,
119  start,
120  cellI,
121  tetFaceI,
122  tetPtI,
123  end,
124  globalWalls.toGlobal(nPatchFaces) // passive data
125  )
126  );
127 
128  nPatchFaces++;
129  }
130  }
131 
132 
133 
134  if (debug)
135  {
136  // Dump particles
137  OBJstream str
138  (
139  mesh.time().path()
140  /"wantedTracks_" + mesh.time().timeName() + ".obj"
141  );
142  Info<< "nearWallFields::calcAddressing() :"
143  << "Dumping tracks to " << str.name() << endl;
144 
146  {
147  const findCellParticle& tp = iter();
148  str.write(linePointRef(tp.position(), tp.end()));
149  }
150  }
151 
152 
153 
154  // Per cell: empty or global wall index and end location
155  cellToWalls_.setSize(mesh.nCells());
156  cellToSamples_.setSize(mesh.nCells());
157 
158  // Database to pass into findCellParticle::move
160 
161  // Track all particles to their end position.
162  scalar maxTrackLen = 2.0*mesh.bounds().mag();
163 
164 
165  //Debug: collect start points
166  pointField start;
167  if (debug)
168  {
169  start.setSize(nPatchFaces);
170  nPatchFaces = 0;
172  {
173  const findCellParticle& tp = iter();
174  start[nPatchFaces++] = tp.position();
175  }
176  }
177 
178 
179  cloud.move(td, maxTrackLen);
180 
181 
182  // Rework cell-to-globalpatchface into a map
183  List<Map<label> > compactMap;
184  getPatchDataMapPtr_.reset
185  (
186  new mapDistribute
187  (
188  globalWalls,
189  cellToWalls_,
190  compactMap
191  )
192  );
193 
194 
195  // Debug: dump resulting tracks
196  if (debug)
197  {
198  getPatchDataMapPtr_().distribute(start);
199  {
200  OBJstream str
201  (
202  mesh.time().path()
203  /"obtainedTracks_" + mesh.time().timeName() + ".obj"
204  );
205  Info<< "nearWallFields::calcAddressing() :"
206  << "Dumping obtained to " << str.name() << endl;
207 
208  forAll(cellToWalls_, cellI)
209  {
210  const List<point>& ends = cellToSamples_[cellI];
211  const labelList& cData = cellToWalls_[cellI];
212  forAll(cData, i)
213  {
214  str.write(linePointRef(ends[i], start[cData[i]]));
215  }
216  }
217  }
218  }
219 }
220 
221 
222 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
223 
225 (
226  const word& name,
227  const objectRegistry& obr,
228  const dictionary& dict,
229  const bool loadFromFiles
230 )
231 :
232  name_(name),
233  obr_(obr),
234  active_(true),
235  fieldSet_()
236 {
237  // Check if the available mesh is an fvMesh otherise deactivate
238  if (isA<fvMesh>(obr_))
239  {
240  read(dict);
241  }
242  else
243  {
244  active_ = false;
246  << "No fvMesh available, deactivating " << name_
247  << endl;
248  }
249 
250 }
251 
252 
253 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
254 
256 {
257  if (debug)
258  {
259  Info<< "nearWallFields::~nearWallFields()" << endl;
260  }
261 }
262 
263 
264 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
265 
267 {
268  if (debug)
269  {
270  Info<< "nearWallFields::read(const dictionary&)" << endl;
271  }
272 
273  if (active_)
274  {
275  const fvMesh& mesh = refCast<const fvMesh>(obr_);
276 
277  dict.lookup("fields") >> fieldSet_;
278  patchSet_ =
280  distance_ = readScalar(dict.lookup("distance"));
281 
282 
283  // Clear out any previously loaded fields
284  vsf_.clear();
285  vvf_.clear();
286  vSpheretf_.clear();
287  vSymmtf_.clear();
288  vtf_.clear();
289  fieldMap_.clear();
290  reverseFieldMap_.clear();
291 
292 
293  // Generate fields with mappedField boundary condition
294 
295  // Convert field to map
296  fieldMap_.resize(2*fieldSet_.size());
297  reverseFieldMap_.resize(2*fieldSet_.size());
298  forAll(fieldSet_, setI)
299  {
300  const word& fldName = fieldSet_[setI].first();
301  const word& sampleFldName = fieldSet_[setI].second();
302 
303  fieldMap_.insert(fldName, sampleFldName);
304  reverseFieldMap_.insert(sampleFldName, fldName);
305  }
306 
307  Info<< type() << " " << name_ << ": Sampling " << fieldMap_.size()
308  << " fields" << endl;
309 
310  // Do analysis
311  calcAddressing();
312  }
313 }
314 
315 
317 {
318  if (debug)
319  {
320  Info<< "nearWallFields:execute()" << endl;
321  }
322 
323 
324  if (active_)
325  {
326  if
327  (
328  fieldMap_.size()
329  && vsf_.empty()
330  && vvf_.empty()
331  && vSpheretf_.empty()
332  && vSymmtf_.empty()
333  && vtf_.empty()
334  )
335  {
336  Info<< type() << " " << name_ << ": Creating " << fieldMap_.size()
337  << " fields" << endl;
338 
339  createFields(vsf_);
340  createFields(vvf_);
341  createFields(vSpheretf_);
342  createFields(vSymmtf_);
343  createFields(vtf_);
344 
345  Info<< endl;
346  }
347 
348  Info<< type() << " " << name_ << " output:" << nl;
349 
350  Info<< " Sampling fields to " << obr_.time().timeName()
351  << endl;
352 
353  sampleFields(vsf_);
354  sampleFields(vvf_);
355  sampleFields(vSpheretf_);
356  sampleFields(vSymmtf_);
357  sampleFields(vtf_);
358  }
359 }
360 
361 
363 {
364  // Do nothing
365 }
366 
367 
369 {
370  // Do nothing
371 }
372 
373 
375 {
376  if (debug)
377  {
378  Info<< "nearWallFields:write()" << endl;
379  }
380 
381  if (active_)
382  {
383  Info<< " Writing sampled fields to " << obr_.time().timeName()
384  << endl;
385 
386  forAll(vsf_, i)
387  {
388  vsf_[i].write();
389  }
390  forAll(vvf_, i)
391  {
392  vvf_[i].write();
393  }
394  forAll(vSpheretf_, i)
395  {
396  vSpheretf_[i].write();
397  }
398  forAll(vSymmtf_, i)
399  {
400  vSymmtf_[i].write();
401  }
402  forAll(vtf_, i)
403  {
404  vtf_[i].write();
405  }
406 
407  Info<< endl;
408  }
409 }
410 
411 
412 // ************************************************************************* //
Foam::polyMesh::FACE_DIAG_TRIS
@ FACE_DIAG_TRIS
Definition: polyMesh.H:105
Foam::nearWallFields::nearWallFields
nearWallFields(const nearWallFields &)
Disallow default bitwise copy construct.
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::OBJstream
OFstream which keeps track of vertices.
Definition: OBJstream.H:53
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::OBJstream::write
virtual Ostream & write(const char)
Write character.
Definition: OBJstream.C:82
Foam::nearWallFields::cellToWalls_
labelListList cellToWalls_
From cell to seed patch faces.
Definition: nearWallFields.H:163
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::wordReList
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::PointIndexHit::hit
bool hit() const
Is there a hit.
Definition: PointIndexHit.H:105
tp
volScalarField tp(IOobject("tp", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar("tp", dimensionSet(0, 0, 1, 0, 0, 0, 0), 0))
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::nearWallFields::cellToSamples_
List< List< point > > cellToSamples_
From cell to tracked end point.
Definition: nearWallFields.H:166
Foam::HashSet< label, Hash< label > >
Foam::nearWallFields::getPatchDataMapPtr_
autoPtr< mapDistribute > getPatchDataMapPtr_
Map from cell based data back to patch based data.
Definition: nearWallFields.H:169
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
nearWallFields.H
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::fvPatch::nf
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:124
Foam::PointIndexHit
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
Definition: PointIndexHit.H:53
Foam::nearWallFields::obr_
const objectRegistry & obr_
Definition: nearWallFields.H:134
Foam::label
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
Foam::nearWallFields::patchSet_
labelHashSet patchSet_
Patches to sample.
Definition: nearWallFields.H:148
Foam::IDLList
Intrusive doubly-linked list.
Definition: IDLList.H:47
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::PointIndexHit::hitPoint
const Point & hitPoint() const
Return hit point.
Definition: PointIndexHit.H:117
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:152
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::nearWallFields::end
virtual void end()
Execute at the final time-loop, currently does nothing.
Definition: nearWallFields.C:362
findCellParticle.H
Foam::findCellParticle
Particle class that finds cells by tracking.
Definition: findCellParticle.H:51
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::nearWallFields::~nearWallFields
virtual ~nearWallFields()
Destructor.
Definition: nearWallFields.C:255
Foam::polyPatch::faceCellCentres
tmp< vectorField > faceCellCentres() const
Return face cell centres.
Definition: polyPatch.C:321
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::cloud
A cloud is a collection of lagrangian particles.
Definition: cloud.H:51
Foam::linePointRef
line< point, const point & > linePointRef
Line using referred points.
Definition: linePointRef.H:45
Foam::readScalar
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
Foam::mappedPatchBase::facePoint
static pointIndexHit facePoint(const polyMesh &, const label faceI, const polyMesh::cellDecomposition)
Get a point on the face given a face decomposition method:
Definition: mappedPatchBase.C:1307
Foam::globalIndex::size
label size() const
Global sum of localSizes.
Definition: globalIndexI.H:66
Foam::Vector< scalar >
wordReList.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::nearWallFields::distance_
scalar distance_
Distance away from wall.
Definition: nearWallFields.H:151
Foam::Cloud< findCellParticle >
Foam::fvPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:155
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
Foam::nearWallFields::write
virtual void write()
Write.
Definition: nearWallFields.C:374
Foam::nearWallFields::timeSet
virtual void timeSet()
Called when time was set at the end of the Time::operator++.
Definition: nearWallFields.C:368
Foam::OFstream::name
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:118
Foam::polyBoundaryMesh::patchSet
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
Definition: polyBoundaryMesh.C:750
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::nearWallFields::read
virtual void read(const dictionary &)
Read the field min/max data.
Definition: nearWallFields.C:266
Foam::nearWallFields::calcAddressing
void calcAddressing()
Calculate addressing from cells back to patch faces.
Definition: nearWallFields.C:42
OBJstream.H
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::nearWallFields::execute
virtual void execute()
Execute, currently does nothing.
Definition: nearWallFields.C:316
Foam::findCellParticle::trackingData
Class used to pass tracking data to the trackToFace function.
Definition: findCellParticle.H:69
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
mappedPatchBase.H