directions.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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "directions.H"
29 #include "polyMesh.H"
30 #include "twoDPointCorrector.H"
31 #include "directionInfo.H"
32 #include "MeshWave.H"
33 #include "OFstream.H"
34 #include "meshTools.H"
35 #include "hexMatcher.H"
36 #include "globalMeshData.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 const Foam::Enum
41 <
43 >
44 Foam::directions::directionTypeNames_
45 ({
46  { directionType::TAN1, "tan1" },
47  { directionType::TAN2, "tan2" },
48  { directionType::NORMAL, "normal" },
49 });
50 
51 
52 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
53 
54 void Foam::directions::writeOBJ(Ostream& os, const point& pt)
55 {
56  os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
57 }
58 
59 
60 void Foam::directions::writeOBJ
61 (
62  Ostream& os,
63  const point& pt0,
64  const point& pt1,
65  label& vertI
66 )
67 {
68  writeOBJ(os, pt0);
69  writeOBJ(os, pt1);
70 
71  os << "l " << vertI + 1 << ' ' << vertI + 2 << endl;
72 
73  vertI += 2;
74 }
75 
76 
77 void Foam::directions::writeOBJ
78 (
79  const fileName& fName,
80  const primitiveMesh& mesh,
81  const vectorField& dirs
82 )
83 {
84  Pout<< "Writing cell info to " << fName << " as vectors at the cellCentres"
85  << endl << endl;
86 
87  OFstream xDirStream(fName);
88 
89  label vertI = 0;
90 
91  forAll(dirs, celli)
92  {
93  const point& ctr = mesh.cellCentres()[celli];
94 
95  // Calculate local length scale
96  scalar minDist = GREAT;
97 
98  const labelList& nbrs = mesh.cellCells()[celli];
99 
100  forAll(nbrs, nbrI)
101  {
102  minDist = min(minDist, mag(mesh.cellCentres()[nbrs[nbrI]] - ctr));
103  }
104 
105  scalar scale = 0.5*minDist;
106 
107  writeOBJ(xDirStream, ctr, ctr + scale*dirs[celli], vertI);
108  }
109 }
110 
111 
112 void Foam::directions::check2D
113 (
114  const twoDPointCorrector* correct2DPtr,
115  const vector& vec
116 )
117 {
118  if (correct2DPtr)
119  {
120  if (mag(correct2DPtr->planeNormal() & vec) > 1e-6)
121  {
123  << "is not normal to plane defined in dynamicMeshDict."
124  << endl
125  << "Either make case 3D or adjust vector."
126  << exit(FatalError);
127  }
128  }
129 }
130 
131 
132 Foam::vectorField Foam::directions::propagateDirection
133 (
134  const polyMesh& mesh,
135  const bool useTopo,
136  const polyPatch& pp,
137  const vectorField& ppField,
138  const vector& defaultDir
139 )
140 {
141  // Seed all faces on patch
142  labelList changedFaces(pp.size());
143  List<directionInfo> changedFacesInfo(pp.size());
144 
145  if (useTopo)
146  {
147  forAll(pp, patchFacei)
148  {
149  label meshFacei = pp.start() + patchFacei;
150 
151  label celli = mesh.faceOwner()[meshFacei];
152 
153  if (!hexMatcher().isA(mesh, celli))
154  {
156  << "useHexTopology specified but cell " << celli
157  << " on face " << patchFacei << " of patch " << pp.name()
158  << " is not a hex" << exit(FatalError);
159  }
160 
161  const vector& cutDir = ppField[patchFacei];
162 
163  // Get edge(bundle) on cell most in direction of cutdir
164  label edgeI = meshTools::cutDirToEdge(mesh, celli, cutDir);
165 
166  // Convert edge into index on face
167  label faceIndex =
169  (
170  mesh,
171  celli,
172  meshFacei,
173  edgeI
174  );
175 
176  // Set initial face and direction
177  changedFaces[patchFacei] = meshFacei;
178  changedFacesInfo[patchFacei] =
179  directionInfo
180  (
181  faceIndex,
182  cutDir
183  );
184  }
185  }
186  else
187  {
188  forAll(pp, patchFacei)
189  {
190  changedFaces[patchFacei] = pp.start() + patchFacei;
191  changedFacesInfo[patchFacei] =
192  directionInfo
193  (
194  -2, // Geometric information only
195  ppField[patchFacei]
196  );
197  }
198  }
199 
200  MeshWave<directionInfo> directionCalc
201  (
202  mesh,
203  changedFaces,
204  changedFacesInfo,
206  );
207 
208  const List<directionInfo>& cellInfo = directionCalc.allCellInfo();
209 
210  vectorField dirField(cellInfo.size());
211 
212  label nUnset = 0;
213  label nGeom = 0;
214  label nTopo = 0;
215 
216  forAll(cellInfo, celli)
217  {
218  label index = cellInfo[celli].index();
219 
220  if (index == -3)
221  {
222  // Never visited
224  << "Cell " << celli << " never visited to determine "
225  << "local coordinate system" << endl
226  << "Using direction " << defaultDir << " instead" << endl;
227 
228  dirField[celli] = defaultDir;
229 
230  nUnset++;
231  }
232  else if (index == -2)
233  {
234  // Geometric direction
235  dirField[celli] = cellInfo[celli].n();
236 
237  nGeom++;
238  }
239  else if (index == -1)
240  {
242  << "Illegal index " << index << endl
243  << "Value is only allowed on faces" << abort(FatalError);
244  }
245  else
246  {
247  // Topological edge cut. Convert into average cut direction.
248  dirField[celli] = meshTools::edgeToCutDir(mesh, celli, index);
249 
250  nTopo++;
251  }
252  }
253 
254  reduce(nGeom, sumOp<label>());
255  reduce(nTopo, sumOp<label>());
256  reduce(nUnset, sumOp<label>());
257 
258  Info<< "Calculated local coords for " << defaultDir
259  << endl
260  << " Geometric cut cells : " << nGeom << endl
261  << " Topological cut cells : " << nTopo << endl
262  << " Unset cells : " << nUnset << endl
263  << endl;
264 
265  return dirField;
266 }
267 
268 
269 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
270 
271 Foam::directions::directions
272 (
273  const polyMesh& mesh,
274  const dictionary& dict,
275  const twoDPointCorrector* correct2DPtr
276 )
277 :
278  List<vectorField>(wordList(dict.lookup("directions")).size())
279 {
280  const wordList wantedDirs(dict.lookup("directions"));
281  const word coordSystem(dict.get<word>("coordinateSystem"));
282 
283  bool wantNormal = false;
284  bool wantTan1 = false;
285  bool wantTan2 = false;
286  label nDirs = 0;
287 
288  if (coordSystem != "fieldBased")
289  {
290  for (const word& wantedName : wantedDirs)
291  {
292  directionType wantedDir = directionTypeNames_[wantedName];
293 
294  if (wantedDir == NORMAL)
295  {
296  wantNormal = true;
297  }
298  else if (wantedDir == TAN1)
299  {
300  wantTan1 = true;
301  }
302  else if (wantedDir == TAN2)
303  {
304  wantTan2 = true;
305  }
306  }
307  }
308 
309 
310  if (coordSystem == "global")
311  {
312  const dictionary& globalDict = dict.subDict("globalCoeffs");
313 
314  vector tan1(globalDict.lookup("tan1"));
315  check2D(correct2DPtr, tan1);
316 
317  vector tan2(globalDict.lookup("tan2"));
318  check2D(correct2DPtr, tan2);
319 
320  const vector normal = normalised(tan1 ^ tan2);
321 
322  Info<< "Global Coordinate system:" << endl
323  << " normal : " << normal << endl
324  << " tan1 : " << tan1 << endl
325  << " tan2 : " << tan2
326  << endl << endl;
327 
328  if (wantNormal)
329  {
330  operator[](nDirs++) = vectorField(1, normal);
331  }
332  if (wantTan1)
333  {
334  operator[](nDirs++) = vectorField(1, tan1);
335  }
336  if (wantTan2)
337  {
338  operator[](nDirs++) = vectorField(1, tan2);
339  }
340  }
341  else if (coordSystem == "patchLocal")
342  {
343  const dictionary& patchDict = dict.subDict("patchLocalCoeffs");
344 
345  const word patchName(patchDict.get<word>("patch"));
346 
347  const label patchi = mesh.boundaryMesh().findPatchID(patchName);
348 
349  if (patchi == -1)
350  {
352  << "Cannot find patch "
353  << patchName
354  << exit(FatalError);
355  }
356 
357  // Take zeroth face on patch
358  const polyPatch& pp = mesh.boundaryMesh()[patchi];
359 
360  vector tan1(patchDict.lookup("tan1"));
361 
362  const vector& n0 = pp.faceNormals()[0];
363 
364  if (correct2DPtr)
365  {
366  tan1 = correct2DPtr->planeNormal() ^ n0;
367 
369  << "Discarding user specified tan1 since 2D case." << endl
370  << "Recalculated tan1 from face normal and planeNormal as "
371  << tan1 << endl << endl;
372  }
373 
374  const bool useTopo(dict.get<bool>("useHexTopology"));
375 
376  vectorField normalDirs;
377  vectorField tan1Dirs;
378 
379  if (wantNormal || wantTan2)
380  {
381  normalDirs =
382  propagateDirection
383  (
384  mesh,
385  useTopo,
386  pp,
387  pp.faceNormals(),
388  n0
389  );
390 
391  if (wantNormal)
392  {
393  this->operator[](nDirs++) = normalDirs;
394  }
395  }
396 
397  if (wantTan1 || wantTan2)
398  {
399  tan1Dirs =
400  propagateDirection
401  (
402  mesh,
403  useTopo,
404  pp,
405  vectorField(pp.size(), tan1),
406  tan1
407  );
408 
409 
410  if (wantTan1)
411  {
412  this->operator[](nDirs++) = tan1Dirs;
413  }
414  }
415  if (wantTan2)
416  {
417  tmp<vectorField> tan2Dirs = normalDirs ^ tan1Dirs;
418 
419  this->operator[](nDirs++) = tan2Dirs;
420  }
421  }
422  else if (coordSystem == "fieldBased")
423  {
424  forAll(wantedDirs, i)
425  {
426  operator[](nDirs++) =
428  (
429  IOobject
430  (
431  mesh.instance()/wantedDirs[i],
432  mesh,
435  )
436  );
437  }
438  }
439  else
440  {
442  << "Unknown coordinate system "
443  << coordSystem << endl
444  << "Known types are global, patchLocal and fieldBased"
445  << exit(FatalError);
446  }
447 }
448 
449 
450 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::vectorIOField
IOField< vector > vectorIOField
vectorField with IO.
Definition: vectorIOField.H:44
Foam::directionInfo::edgeToFaceIndex
static label edgeToFaceIndex(const primitiveMesh &mesh, const label celli, const label facei, const label edgeI)
Given edge on hex cell find corresponding edge on face. Is either.
Definition: directionInfo.C:95
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
meshTools.H
directionInfo.H
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:51
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::twoDPointCorrector::planeNormal
const vector & planeNormal() const
Return plane normal.
Definition: twoDPointCorrector.C:248
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
globalMeshData.H
Foam::IOobject::instance
const fileName & instance() const
Definition: IOobjectI.H:167
hexMatcher.H
Foam::globalMeshData::nTotalCells
label nTotalCells() const
Return total number of cells in decomposed mesh.
Definition: globalMeshData.H:394
Foam::directions::directionType
directionType
Enumeration listing the possible coordinate directions.
Definition: directions.H:76
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
polyMesh.H
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
OFstream.H
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:59
Foam::twoDPointCorrector
Class applies a two-dimensional correction to mesh motion point field.
Definition: twoDPointCorrector.H:63
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
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:62
Foam::Field< vector >
MeshWave.H
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
directions.H
Foam::meshTools::cutDirToEdge
label cutDirToEdge(const primitiveMesh &mesh, const label celli, const vector &cutDir)
Reverse of edgeToCutDir: given direction find edge bundle and.
Definition: meshTools.C:810
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1076
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:523
Foam::primitiveMesh::cellCells
const labelListList & cellCells() const
Definition: primitiveMeshCellCells.C:102
Foam::polyBoundaryMesh::findPatchID
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Definition: polyBoundaryMesh.C:766
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:419
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:496
Foam::isA
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:176
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::PrimitivePatch::faceNormals
const Field< PointType > & faceNormals() const
Return face unit normals for patch.
Definition: PrimitivePatch.C:577
Foam::Vector< scalar >
Foam::List< vectorField >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::meshTools::edgeToCutDir
vector edgeToCutDir(const primitiveMesh &mesh, const label celli, const label edgeI)
Given edge on hex find all 'parallel' (i.e. non-connected)
Definition: meshTools.C:763
twoDPointCorrector.H
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1241
Foam::point
vector point
Point is a vector.
Definition: point.H:43
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::IOobject::MUST_READ
Definition: IOobject.H:120