directions.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-2025 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 "directions.H"
27 #include "polyMesh.H"
28 #include "twoDPointCorrector.H"
29 #include "directionInfo.H"
30 #include "FaceCellWave.H"
31 #include "OFstream.H"
32 #include "meshTools.H"
33 #include "hexMatcher.H"
34 #include "Switch.H"
35 #include "globalMeshData.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
40 Foam::directions::directionTypeNames_
41 {
42  "e1",
43  "e2",
44  "e3"
45 };
46 
47 
48 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
49 
50 void Foam::directions::writeOBJ(Ostream& os, const point& pt)
51 {
52  os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
53 }
54 
55 
56 void Foam::directions::writeOBJ
57 (
58  Ostream& os,
59  const point& pt0,
60  const point& pt1,
61  label& vertI
62 )
63 {
64  writeOBJ(os, pt0);
65  writeOBJ(os, pt1);
66 
67  os << "l " << vertI + 1 << ' ' << vertI + 2 << endl;
68 
69  vertI += 2;
70 }
71 
72 
73 void Foam::directions::writeOBJ
74 (
75  const fileName& fName,
76  const primitiveMesh& mesh,
77  const vectorField& dirs
78 )
79 {
80  Pout<< "Writing cell info to " << fName << " as vectors at the cellCentres"
81  << endl << endl;
82 
83  OFstream xDirStream(fName);
84 
85  label vertI = 0;
86 
87  forAll(dirs, celli)
88  {
89  const point& ctr = mesh.cellCentres()[celli];
90 
91  // Calculate local length scale
92  scalar minDist = great;
93 
94  const labelList& nbrs = mesh.cellCells()[celli];
95 
96  forAll(nbrs, nbrI)
97  {
98  minDist = min(minDist, mag(mesh.cellCentres()[nbrs[nbrI]] - ctr));
99  }
100 
101  scalar scale = 0.5*minDist;
102 
103  writeOBJ(xDirStream, ctr, ctr + scale*dirs[celli], vertI);
104  }
105 }
106 
107 
108 void Foam::directions::check2D
109 (
110  const twoDPointCorrector* correct2DPtr,
111  const vector& vec
112 )
113 {
114  if (correct2DPtr)
115  {
116  if (mag(correct2DPtr->planeNormal() & vec) > 1e-6)
117  {
119  << "is not normal to plane defined in dynamicMeshDict."
120  << endl
121  << "Either make case 3D or adjust vector."
122  << exit(FatalError);
123  }
124  }
125 }
126 
127 
128 Foam::vectorField Foam::directions::propagateDirection
129 (
130  const polyMesh& mesh,
131  const bool useTopo,
132  const polyPatch& pp,
133  const vectorField& ppField,
134  const vector& defaultDir
135 )
136 {
137  // Seed all faces on patch
138  labelList changedFaces(pp.size());
139  List<directionInfo> changedFacesInfo(pp.size());
140 
141  if (useTopo)
142  {
143  forAll(pp, patchFacei)
144  {
145  label meshFacei = pp.start() + patchFacei;
146 
147  label celli = mesh.faceOwner()[meshFacei];
148 
149  if (!hexMatcher().isA(mesh, celli))
150  {
152  << "useHexTopology specified but cell " << celli
153  << " on face " << patchFacei << " of patch " << pp.name()
154  << " is not a hex" << exit(FatalError);
155  }
156 
157  const vector& cutDir = ppField[patchFacei];
158 
159  // Get edge(bundle) on cell most in direction of cutdir
160  label edgeI = meshTools::cutDirToEdge(mesh, celli, cutDir);
161 
162  // Convert edge into index on face
163  label faceIndex =
165  (
166  mesh,
167  celli,
168  meshFacei,
169  edgeI
170  );
171 
172  // Set initial face and direction
173  changedFaces[patchFacei] = meshFacei;
174  changedFacesInfo[patchFacei] =
175  directionInfo
176  (
177  faceIndex,
178  cutDir
179  );
180  }
181  }
182  else
183  {
184  forAll(pp, patchFacei)
185  {
186  changedFaces[patchFacei] = pp.start() + patchFacei;
187  changedFacesInfo[patchFacei] =
188  directionInfo
189  (
190  -2, // Geometric information only
191  ppField[patchFacei]
192  );
193  }
194  }
195 
196  List<directionInfo> faceInfo(mesh.nFaces()), cellInfo(mesh.nCells());
197  FaceCellWave<directionInfo> directionCalc
198  (
199  mesh,
200  changedFaces,
201  changedFacesInfo,
202  faceInfo,
203  cellInfo,
204  mesh.globalData().nTotalCells() + 1
205  );
206 
207  vectorField dirField(cellInfo.size());
208 
209  label nUnset = 0;
210  label nGeom = 0;
211  label nTopo = 0;
212 
213  forAll(cellInfo, celli)
214  {
215  label index = cellInfo[celli].index();
216 
217  if (index == -3)
218  {
219  // Never visited
221  << "Cell " << celli << " never visited to determine "
222  << "local coordinate system" << endl
223  << "Using direction " << defaultDir << " instead" << endl;
224 
225  dirField[celli] = defaultDir;
226 
227  nUnset++;
228  }
229  else if (index == -2)
230  {
231  // Geometric direction
232  dirField[celli] = cellInfo[celli].n();
233 
234  nGeom++;
235  }
236  else if (index == -1)
237  {
239  << "Illegal index " << index << endl
240  << "Value is only allowed on faces" << abort(FatalError);
241  }
242  else
243  {
244  // Topological edge cut. Convert into average cut direction.
245  dirField[celli] = meshTools::edgeToCutDir(mesh, celli, index);
246 
247  nTopo++;
248  }
249  }
250 
251  reduce(nGeom, sumOp<label>());
252  reduce(nTopo, sumOp<label>());
253  reduce(nUnset, sumOp<label>());
254 
255  Info<< "Calculated local coords for " << defaultDir
256  << endl
257  << " Geometric cut cells : " << nGeom << endl
258  << " Topological cut cells : " << nTopo << endl
259  << " Unset cells : " << nUnset << endl
260  << endl;
261 
262  return dirField;
263 }
264 
265 
266 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
267 
269 (
270  const polyMesh& mesh,
271  const dictionary& dict,
272  const twoDPointCorrector* correct2DPtr
273 )
274 {
275  word coordSystem;
276  wordList wantedDirs;
277 
278  if (dict.found("type"))
279  {
280  coordSystem = dict.lookup<word>("type");
281  wantedDirs = dict.lookup<wordList>("directions");
282  }
283  else
284  {
285  coordSystem = dict.lookup<word>("coordinateSystem");
286  wantedDirs = dict.lookup<wordList>("directions");
287  }
288 
289  setSize(wantedDirs.size());
290 
291  bool wantE3 = false;
292  bool wantE1 = false;
293  bool wantE2 = false;
294  label nDirs = 0;
295 
296  if (coordSystem != "fieldBased")
297  {
298  forAll(wantedDirs, i)
299  {
300  directionType wantedDir = directionTypeNames_[wantedDirs[i]];
301 
302  if (wantedDir == e3)
303  {
304  wantE3 = true;
305  }
306  else if (wantedDir == e1)
307  {
308  wantE1 = true;
309  }
310  else if (wantedDir == e2)
311  {
312  wantE2 = true;
313  }
314  }
315  }
316 
317 
318  if (coordSystem == "global")
319  {
320  const dictionary& globalDict = dict.optionalSubDict("globalCoeffs");
321 
322  vector e1(globalDict.lookup("e1"));
323  check2D(correct2DPtr, e1);
324 
325  vector e2(globalDict.lookup("e2"));
326  check2D(correct2DPtr, e2);
327 
328  vector e3 = e1 ^ e2;
329  e3 /= mag(e3);
330 
331  Info<< "Global Coordinate system:" << endl
332  << " e1 : " << e1 << nl
333  << " e2 : " << e2 << nl
334  << " e3 : " << e3 << nl
335  << endl;
336 
337  if (wantE3)
338  {
339  operator[](nDirs++) = vectorField(1, e3);
340  }
341  if (wantE1)
342  {
343  operator[](nDirs++) = vectorField(1, e1);
344  }
345  if (wantE2)
346  {
347  operator[](nDirs++) = vectorField(1, e2);
348  }
349  }
350  else if (coordSystem == "patchLocal")
351  {
352  const dictionary& patchDict = dict.optionalSubDict("patchLocalCoeffs");
353 
354  const word patchName(patchDict.lookup("patch"));
355 
356  const label patchi = mesh.boundaryMesh().findIndex(patchName);
357 
358  if (patchi == -1)
359  {
361  << "Cannot find patch "
362  << patchName
363  << exit(FatalError);
364  }
365 
366  // Take zeroth face on patch
367  const polyPatch& pp = mesh.boundaryMesh()[patchi];
368 
369  vector e1(patchDict.lookup("e1"));
370 
371  const vector& n0 = pp.faceNormals()[0];
372 
373  if (correct2DPtr)
374  {
375  e1 = correct2DPtr->planeNormal() ^ n0;
376 
378  << "Discarding user specified e1 since 2D case." << endl
379  << "Recalculated e1 from face normal and planeNormal as "
380  << e1 << endl << endl;
381  }
382 
383  Switch useTopo(dict.lookupOrDefault("useHexTopology", false));
384 
385  vectorField e3Dirs;
386  vectorField e1Dirs;
387 
388  if (wantE3 || wantE2)
389  {
390  e3Dirs =
391  propagateDirection
392  (
393  mesh,
394  useTopo,
395  pp,
396  pp.faceNormals(),
397  n0
398  );
399 
400  if (wantE3)
401  {
402  this->operator[](nDirs++) = e3Dirs;
403  }
404  }
405 
406  if (wantE1 || wantE2)
407  {
408  e1Dirs =
409  propagateDirection
410  (
411  mesh,
412  useTopo,
413  pp,
414  vectorField(pp.size(), e1),
415  e1
416  );
417 
418 
419  if (wantE1)
420  {
421  this->operator[](nDirs++) = e1Dirs;
422  }
423  }
424  if (wantE2)
425  {
426  tmp<vectorField> e2Dirs = e3Dirs ^ e1Dirs;
427 
428  this->operator[](nDirs++) = e2Dirs;
429  }
430  }
431  else if (coordSystem == "fieldBased")
432  {
433  forAll(wantedDirs, i)
434  {
435  operator[](nDirs++) =
437  (
438  IOobject
439  (
440  mesh.instance()/wantedDirs[i],
441  mesh,
444  )
445  );
446  }
447  }
448  else
449  {
451  << "Unknown coordinate system "
452  << coordSystem << endl
453  << "Known types are global, patchLocal and fieldBased"
454  << exit(FatalError);
455  }
456 }
457 
458 
459 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:352
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:55
const Field< PointType > & faceNormals() const
Return face normals for patch.
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
const dictionary & optionalSubDict(const word &) const
Find and return a sub-dictionary if found.
Definition: dictionary.C:927
T lookupOrDefault(const word &, const T &, const bool writeDefault=writeOptionalEntries > 0) const
Find and return a T, if not found return the given default.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
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:93
directions(const polyMesh &mesh, const dictionary &dict, const twoDPointCorrector *correct2DPtr=nullptr)
Construct from mesh and dictionary and optional 2D corrector.
Definition: directions.C:269
directionType
Enumeration listing the possible coordinate directions.
Definition: directions.H:75
label nTotalCells() const
Return total number of cells in decomposed mesh.
label findIndex(const word &patchName) const
Find patch index given a name.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1357
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1521
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:401
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const vectorField & cellCentres() const
const labelListList & cellCells() const
label nCells() const
label nFaces() const
A class for managing temporary objects.
Definition: tmp.H:55
Class applies a two-dimensional correction to mesh motion point field.
const vector & planeNormal() const
Return plane normal.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
label cutDirToEdge(const primitiveMesh &, const label celli, const vector &cutDir)
Reverse of edgeToCutDir: given direction find edge bundle and.
Definition: meshTools.C:795
vector edgeToCutDir(const primitiveMesh &, const label celli, const label edgeI)
Given edge on hex find all 'parallel' (i.e. non-connected)
Definition: meshTools.C:747
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:106
scalar minDist(const List< pointIndexHit > &hitList)
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
messageStream Info
vector point
Point is a vector.
Definition: point.H:41
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
IOField< vector > vectorIOField
vectorField with IO.
Definition: vectorIOField.H:43
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:171
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Field< vector > vectorField
Specialisation of Field<T> for vector.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
static const char nl
Definition: Ostream.H:267
points setSize(newPointi)
dictionary dict