Loading...
Searching...
No Matches
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 Copyright (C) 2020,2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "directions.H"
30#include "polyMesh.H"
31#include "twoDPointCorrector.H"
32#include "directionInfo.H"
33#include "MeshWave.H"
34#include "OFstream.H"
35#include "meshTools.H"
36#include "hexMatcher.H"
37#include "globalMeshData.H"
38#include "coordinateSystem.H"
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42const Foam::Enum
43<
45>
46Foam::directions::directionTypeNames_
47({
48 { directionType::TAN1, "tan1" },
49 { directionType::TAN2, "tan2" },
50 { directionType::NORMAL, "normal" },
51});
52
53
54// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55
56void Foam::directions::writeOBJ(Ostream& os, const point& pt)
57{
58 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
59}
60
61
62void Foam::directions::writeOBJ
63(
64 Ostream& os,
65 const point& pt0,
66 const point& pt1,
67 label& vertI
68)
69{
70 writeOBJ(os, pt0);
71 writeOBJ(os, pt1);
72
73 os << "l " << vertI + 1 << ' ' << vertI + 2 << endl;
74
75 vertI += 2;
76}
77
78
79void Foam::directions::writeOBJ
80(
81 const fileName& fName,
82 const primitiveMesh& mesh,
83 const vectorField& dirs
84)
85{
86 Pout<< "Writing cell info to " << fName << " as vectors at the cellCentres"
87 << endl << endl;
88
89 OFstream xDirStream(fName);
90
91 label vertI = 0;
92
93 forAll(dirs, celli)
94 {
95 const point& ctr = mesh.cellCentres()[celli];
96
97 // Calculate local length scale
98 scalar minDist = GREAT;
99
100 const labelList& nbrs = mesh.cellCells()[celli];
101
102 forAll(nbrs, nbrI)
103 {
104 minDist = min(minDist, mag(mesh.cellCentres()[nbrs[nbrI]] - ctr));
105 }
106
107 scalar scale = 0.5*minDist;
108
109 writeOBJ(xDirStream, ctr, ctr + scale*dirs[celli], vertI);
110 }
111}
112
113
114void Foam::directions::check2D
115(
116 const twoDPointCorrector* correct2DPtr,
117 const vector& vec
118)
119{
120 if (correct2DPtr)
121 {
122 if (mag(correct2DPtr->planeNormal() & vec) > 1e-6)
123 {
125 << "is not normal to plane defined in dynamicMeshDict."
126 << endl
127 << "Either make case 3D or adjust vector."
128 << exit(FatalError);
129 }
130 }
131}
132
133
134Foam::vectorField Foam::directions::propagateDirection
135(
136 const polyMesh& mesh,
137 const bool useTopo,
138 const polyPatch& pp,
139 const vectorField& ppField,
140 const vector& defaultDir
141)
142{
143 // Seed all faces on patch
144 labelList changedFaces(pp.size());
145 List<directionInfo> changedFacesInfo(pp.size());
146
147 if (useTopo)
148 {
149 forAll(pp, patchFacei)
150 {
151 label meshFacei = pp.start() + patchFacei;
152
153 label celli = mesh.faceOwner()[meshFacei];
154
155 if (!hexMatcher::test(mesh, celli))
156 {
158 << "useHexTopology specified but cell " << celli
159 << " on face " << patchFacei << " of patch " << pp.name()
160 << " is not a hex" << exit(FatalError);
161 }
162
163 const vector& cutDir = ppField[patchFacei];
164
165 // Get edge(bundle) on cell most in direction of cutdir
166 label edgeI = meshTools::cutDirToEdge(mesh, celli, cutDir);
167
168 // Convert edge into index on face
169 label faceIndex =
171 (
172 mesh,
173 celli,
174 meshFacei,
175 edgeI
176 );
177
178 // Set initial face and direction
179 changedFaces[patchFacei] = meshFacei;
180 changedFacesInfo[patchFacei] =
182 (
183 faceIndex,
184 cutDir
185 );
186 }
187 }
188 else
189 {
190 forAll(pp, patchFacei)
191 {
192 changedFaces[patchFacei] = pp.start() + patchFacei;
193 changedFacesInfo[patchFacei] =
195 (
196 -2, // Geometric information only
197 ppField[patchFacei]
198 );
199 }
200 }
201
202 MeshWave<directionInfo> directionCalc
203 (
204 mesh,
205 changedFaces,
206 changedFacesInfo,
208 );
209
210 const List<directionInfo>& cellInfo = directionCalc.allCellInfo();
211
212 vectorField dirField(cellInfo.size());
213
214 label nUnset = 0;
215 label nGeom = 0;
216 label nTopo = 0;
217
218 forAll(cellInfo, celli)
219 {
220 label index = cellInfo[celli].index();
221
222 if (index == -3)
223 {
224 // Never visited
226 << "Cell " << celli << " never visited to determine "
227 << "local coordinate system" << endl
228 << "Using direction " << defaultDir << " instead" << endl;
229
230 dirField[celli] = defaultDir;
231
232 nUnset++;
233 }
234 else if (index == -2)
235 {
236 // Geometric direction
237 dirField[celli] = cellInfo[celli].n();
238
239 nGeom++;
240 }
241 else if (index == -1)
242 {
244 << "Illegal index " << index << endl
245 << "Value is only allowed on faces" << abort(FatalError);
246 }
247 else
248 {
249 // Topological edge cut. Convert into average cut direction.
250 dirField[celli] = meshTools::edgeToCutDir(mesh, celli, index);
251
252 nTopo++;
253 }
254 }
255
256 reduce(nGeom, sumOp<label>());
257 reduce(nTopo, sumOp<label>());
258 reduce(nUnset, sumOp<label>());
259
260 Info<< "Calculated local coords for " << defaultDir
261 << endl
262 << " Geometric cut cells : " << nGeom << endl
263 << " Topological cut cells : " << nTopo << endl
264 << " Unset cells : " << nUnset << endl
265 << endl;
267 return dirField;
268}
269
270
271// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
272
273Foam::directions::directions
274(
275 const polyMesh& mesh,
276 const dictionary& dict,
277 const twoDPointCorrector* correct2DPtr
278)
279:
281{
282 const wordList wantedDirs(dict.get<wordList>("directions"));
283 const word coordSystem(dict.get<word>("coordinateSystem"));
284
285 List<vectorField>::resize(wantedDirs.size());
286
287 bool wantNormal = false;
288 bool wantTan1 = false;
289 bool wantTan2 = false;
290 label nDirs = 0;
291
292 if (coordSystem != "fieldBased")
293 {
294 for (const word& wantedName : wantedDirs)
295 {
296 directionType wantedDir = directionTypeNames_[wantedName];
297
298 if (wantedDir == NORMAL)
299 {
300 wantNormal = true;
301 }
302 else if (wantedDir == TAN1)
303 {
304 wantTan1 = true;
305 }
306 else if (wantedDir == TAN2)
307 {
308 wantTan2 = true;
309 }
310 }
311 }
312
313
314 if (coordSystem == "global")
315 {
316 const dictionary& globalDict = dict.subDict("globalCoeffs");
317
318 vector tan1(globalDict.get<vector>("tan1"));
319 check2D(correct2DPtr, tan1);
320
321 vector tan2(globalDict.get<vector>("tan2"));
322 check2D(correct2DPtr, tan2);
323
324 const vector normal = normalised(tan1 ^ tan2);
325
326 Info<< "Global Coordinate system:" << endl
327 << " normal : " << normal << endl
328 << " tan1 : " << tan1 << endl
329 << " tan2 : " << tan2
330 << endl << endl;
331
332 if (wantNormal)
333 {
334 operator[](nDirs++) = vectorField(1, normal);
335 }
336 if (wantTan1)
337 {
338 operator[](nDirs++) = vectorField(1, tan1);
339 }
340 if (wantTan2)
341 {
342 operator[](nDirs++) = vectorField(1, tan2);
343 }
344 }
345 else if (coordSystem == "user")
346 {
347 const dictionary& globalDict = dict.subDict("userCoeffs");
348
349 auto csysPtr(coordinateSystem::New(mesh.thisDb(), globalDict));
350 const auto& cs = csysPtr();
351
352 if (wantNormal)
353 {
354 // 'Normal' is usually e3.
355 vectorField& result = operator[](nDirs++);
356 result = cs.transform(mesh.cellCentres(), vector(0, 0, 1));
357 }
358 if (wantTan1)
359 {
360 vectorField& result = operator[](nDirs++);
361 result = cs.transform(mesh.cellCentres(), vector(1, 0, 0));
362 }
363 if (wantTan2)
364 {
365 vectorField& result = operator[](nDirs++);
366 result = cs.transform(mesh.cellCentres(), vector(0, 1, 0));
367 }
368 }
369 else if (coordSystem == "patchLocal")
370 {
371 const dictionary& patchDict = dict.subDict("patchLocalCoeffs");
372
373 const word patchName(patchDict.get<word>("patch"));
374
375 const label patchi = mesh.boundaryMesh().findPatchID(patchName);
376
377 if (patchi == -1)
378 {
380 << "Cannot find patch "
381 << patchName
382 << exit(FatalError);
383 }
384
385 // Take zeroth face on patch
386 const polyPatch& pp = mesh.boundaryMesh()[patchi];
387
388 vector tan1(patchDict.get<vector>("tan1"));
389
390 const vector& n0 = pp.faceNormals()[0];
391
392 if (correct2DPtr)
393 {
394 tan1 = correct2DPtr->planeNormal() ^ n0;
395
397 << "Discarding user specified tan1 since 2D case." << endl
398 << "Recalculated tan1 from face normal and planeNormal as "
399 << tan1 << endl << endl;
400 }
401
402 const bool useTopo(dict.get<bool>("useHexTopology"));
403
404 vectorField normalDirs;
405 vectorField tan1Dirs;
406
407 if (wantNormal || wantTan2)
408 {
409 normalDirs =
410 propagateDirection
411 (
412 mesh,
413 useTopo,
414 pp,
415 pp.faceNormals(),
416 n0
417 );
418
419 if (wantNormal)
420 {
421 this->operator[](nDirs++) = normalDirs;
422 }
423 }
424
425 if (wantTan1 || wantTan2)
426 {
427 tan1Dirs =
428 propagateDirection
429 (
430 mesh,
431 useTopo,
432 pp,
433 vectorField(pp.size(), tan1),
434 tan1
435 );
436
437
438 if (wantTan1)
439 {
440 this->operator[](nDirs++) = tan1Dirs;
441 }
442 }
443 if (wantTan2)
444 {
445 tmp<vectorField> tan2Dirs = normalDirs ^ tan1Dirs;
446
447 this->operator[](nDirs++) = tan2Dirs;
448 }
449 }
450 else if (coordSystem == "fieldBased")
451 {
452 forAll(wantedDirs, i)
453 {
454 operator[](nDirs++) =
456 (
457 IOobject
458 (
459 mesh.instance()/wantedDirs[i],
460 mesh,
463 )
464 );
465 }
466 }
467 else
468 {
470 << "Unknown coordinate system "
471 << coordSystem << endl
472 << "Known types are global, user, patchLocal and fieldBased"
473 << exit(FatalError);
474 }
475}
476
477
478// ************************************************************************* //
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
label size() const noexcept
The number of elements in the list.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
FaceCellWave plus data.
Definition MeshWave.H:58
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
vectorField & operator[](const label i)
Definition UListI.H:363
Holds information regarding type of cell. Used in inside/outside determination in cellClassification.
Definition cellInfo.H:61
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Holds direction in which to split cell (in fact a local coordinate axes). Information is a label and ...
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.
directionType
Enumeration listing the possible coordinate directions.
Definition directions.H:76
A class for handling file names.
Definition fileName.H:75
label nTotalCells() const noexcept
Total global number of mesh cells.
static bool test(const UList< face > &faces)
Test if given list of faces satisfies criteria for HEX. (6 quad).
Definition hexMatcher.C:80
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual const labelList & faceOwner() const
Return face owner.
Definition polyMesh.C:1101
const globalMeshData & globalData() const
Return parallel info (demand-driven).
Definition polyMesh.C:1296
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Cell-face mesh analysis engine.
const vectorField & cellCentres() const
const labelListList & cellCells() const
A class for managing temporary objects.
Definition tmp.H:75
Class applies a two-dimensional correction to mesh motion point field.
const vector & planeNormal() const
Return plane normal.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
Definition word.H:66
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for coordinate systems.
Definition cartesianCS.C:30
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:756
label cutDirToEdge(const primitiveMesh &mesh, const label celli, const vector &cutDir)
Reverse of edgeToCutDir: given direction find edge bundle and.
Definition meshTools.C:803
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition meshTools.C:196
List< word > wordList
List of word.
Definition fileName.H:60
List< label > labelList
A List of labels.
Definition List.H:62
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
IOField< vector > vectorIOField
IO for a Field of vector.
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299