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 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
39// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40
41const Foam::Enum
42<
44>
45Foam::directions::directionTypeNames_
46({
47 { directionType::TAN1, "tan1" },
48 { directionType::TAN2, "tan2" },
49 { directionType::NORMAL, "normal" },
50});
51
52
53// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54
55void Foam::directions::writeOBJ(Ostream& os, const point& pt)
56{
57 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
58}
59
60
61void Foam::directions::writeOBJ
62(
63 Ostream& os,
64 const point& pt0,
65 const point& pt1,
66 label& vertI
67)
68{
69 writeOBJ(os, pt0);
70 writeOBJ(os, pt1);
71
72 os << "l " << vertI + 1 << ' ' << vertI + 2 << endl;
73
74 vertI += 2;
75}
76
77
78void Foam::directions::writeOBJ
79(
80 const fileName& fName,
81 const primitiveMesh& mesh,
82 const vectorField& dirs
83)
84{
85 Pout<< "Writing cell info to " << fName << " as vectors at the cellCentres"
86 << endl << endl;
87
88 OFstream xDirStream(fName);
89
90 label vertI = 0;
91
92 forAll(dirs, celli)
93 {
94 const point& ctr = mesh.cellCentres()[celli];
95
96 // Calculate local length scale
97 scalar minDist = GREAT;
98
99 const labelList& nbrs = mesh.cellCells()[celli];
100
101 forAll(nbrs, nbrI)
102 {
103 minDist = min(minDist, mag(mesh.cellCentres()[nbrs[nbrI]] - ctr));
104 }
105
106 scalar scale = 0.5*minDist;
107
108 writeOBJ(xDirStream, ctr, ctr + scale*dirs[celli], vertI);
109 }
110}
111
112
113void Foam::directions::check2D
114(
115 const twoDPointCorrector* correct2DPtr,
116 const vector& vec
117)
118{
119 if (correct2DPtr)
120 {
121 if (mag(correct2DPtr->planeNormal() & vec) > 1e-6)
122 {
124 << "is not normal to plane defined in dynamicMeshDict."
125 << endl
126 << "Either make case 3D or adjust vector."
127 << exit(FatalError);
128 }
129 }
130}
131
132
133Foam::vectorField Foam::directions::propagateDirection
134(
135 const polyMesh& mesh,
136 const bool useTopo,
137 const polyPatch& pp,
138 const vectorField& ppField,
139 const vector& defaultDir
140)
141{
142 // Seed all faces on patch
143 labelList changedFaces(pp.size());
144 List<directionInfo> changedFacesInfo(pp.size());
145
146 if (useTopo)
147 {
148 forAll(pp, patchFacei)
149 {
150 label meshFacei = pp.start() + patchFacei;
151
152 label celli = mesh.faceOwner()[meshFacei];
153
154 if (!hexMatcher::test(mesh, celli))
155 {
157 << "useHexTopology specified but cell " << celli
158 << " on face " << patchFacei << " of patch " << pp.name()
159 << " is not a hex" << exit(FatalError);
160 }
161
162 const vector& cutDir = ppField[patchFacei];
163
164 // Get edge(bundle) on cell most in direction of cutdir
165 label edgeI = meshTools::cutDirToEdge(mesh, celli, cutDir);
166
167 // Convert edge into index on face
168 label faceIndex =
170 (
171 mesh,
172 celli,
173 meshFacei,
174 edgeI
175 );
176
177 // Set initial face and direction
178 changedFaces[patchFacei] = meshFacei;
179 changedFacesInfo[patchFacei] =
180 directionInfo
181 (
182 faceIndex,
183 cutDir
184 );
185 }
186 }
187 else
188 {
189 forAll(pp, patchFacei)
190 {
191 changedFaces[patchFacei] = pp.start() + patchFacei;
192 changedFacesInfo[patchFacei] =
193 directionInfo
194 (
195 -2, // Geometric information only
196 ppField[patchFacei]
197 );
198 }
199 }
200
201 MeshWave<directionInfo> directionCalc
202 (
203 mesh,
204 changedFaces,
205 changedFacesInfo,
207 );
208
209 const List<directionInfo>& cellInfo = directionCalc.allCellInfo();
210
211 vectorField dirField(cellInfo.size());
212
213 label nUnset = 0;
214 label nGeom = 0;
215 label nTopo = 0;
216
217 forAll(cellInfo, celli)
218 {
219 label index = cellInfo[celli].index();
220
221 if (index == -3)
222 {
223 // Never visited
225 << "Cell " << celli << " never visited to determine "
226 << "local coordinate system" << endl
227 << "Using direction " << defaultDir << " instead" << endl;
228
229 dirField[celli] = defaultDir;
230
231 nUnset++;
232 }
233 else if (index == -2)
234 {
235 // Geometric direction
236 dirField[celli] = cellInfo[celli].n();
237
238 nGeom++;
239 }
240 else if (index == -1)
241 {
243 << "Illegal index " << index << endl
244 << "Value is only allowed on faces" << abort(FatalError);
245 }
246 else
247 {
248 // Topological edge cut. Convert into average cut direction.
249 dirField[celli] = meshTools::edgeToCutDir(mesh, celli, index);
250
251 nTopo++;
252 }
253 }
254
255 reduce(nGeom, sumOp<label>());
256 reduce(nTopo, sumOp<label>());
257 reduce(nUnset, sumOp<label>());
258
259 Info<< "Calculated local coords for " << defaultDir
260 << endl
261 << " Geometric cut cells : " << nGeom << endl
262 << " Topological cut cells : " << nTopo << endl
263 << " Unset cells : " << nUnset << endl
264 << endl;
265
266 return dirField;
267}
268
269
270// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
271
273(
274 const polyMesh& mesh,
275 const dictionary& dict,
276 const twoDPointCorrector* correct2DPtr
277)
278:
280{
281 const wordList wantedDirs(dict.get<wordList>("directions"));
282 const word coordSystem(dict.get<word>("coordinateSystem"));
283
284 List<vectorField>::resize(wantedDirs.size());
285
286 bool wantNormal = false;
287 bool wantTan1 = false;
288 bool wantTan2 = false;
289 label nDirs = 0;
290
291 if (coordSystem != "fieldBased")
292 {
293 for (const word& wantedName : wantedDirs)
294 {
295 directionType wantedDir = directionTypeNames_[wantedName];
296
297 if (wantedDir == NORMAL)
298 {
299 wantNormal = true;
300 }
301 else if (wantedDir == TAN1)
302 {
303 wantTan1 = true;
304 }
305 else if (wantedDir == TAN2)
306 {
307 wantTan2 = true;
308 }
309 }
310 }
311
312
313 if (coordSystem == "global")
314 {
315 const dictionary& globalDict = dict.subDict("globalCoeffs");
316
317 vector tan1(globalDict.get<vector>("tan1"));
318 check2D(correct2DPtr, tan1);
319
320 vector tan2(globalDict.get<vector>("tan2"));
321 check2D(correct2DPtr, tan2);
322
323 const vector normal = normalised(tan1 ^ tan2);
324
325 Info<< "Global Coordinate system:" << endl
326 << " normal : " << normal << endl
327 << " tan1 : " << tan1 << endl
328 << " tan2 : " << tan2
329 << endl << endl;
330
331 if (wantNormal)
332 {
333 operator[](nDirs++) = vectorField(1, normal);
334 }
335 if (wantTan1)
336 {
337 operator[](nDirs++) = vectorField(1, tan1);
338 }
339 if (wantTan2)
340 {
341 operator[](nDirs++) = vectorField(1, tan2);
342 }
343 }
344 else if (coordSystem == "patchLocal")
345 {
346 const dictionary& patchDict = dict.subDict("patchLocalCoeffs");
347
348 const word patchName(patchDict.get<word>("patch"));
349
350 const label patchi = mesh.boundaryMesh().findPatchID(patchName);
351
352 if (patchi == -1)
353 {
355 << "Cannot find patch "
356 << patchName
357 << exit(FatalError);
358 }
359
360 // Take zeroth face on patch
361 const polyPatch& pp = mesh.boundaryMesh()[patchi];
362
363 vector tan1(patchDict.get<vector>("tan1"));
364
365 const vector& n0 = pp.faceNormals()[0];
366
367 if (correct2DPtr)
368 {
369 tan1 = correct2DPtr->planeNormal() ^ n0;
370
372 << "Discarding user specified tan1 since 2D case." << endl
373 << "Recalculated tan1 from face normal and planeNormal as "
374 << tan1 << endl << endl;
375 }
376
377 const bool useTopo(dict.get<bool>("useHexTopology"));
378
379 vectorField normalDirs;
380 vectorField tan1Dirs;
381
382 if (wantNormal || wantTan2)
383 {
384 normalDirs =
385 propagateDirection
386 (
387 mesh,
388 useTopo,
389 pp,
390 pp.faceNormals(),
391 n0
392 );
393
394 if (wantNormal)
395 {
396 this->operator[](nDirs++) = normalDirs;
397 }
398 }
399
400 if (wantTan1 || wantTan2)
401 {
402 tan1Dirs =
403 propagateDirection
404 (
405 mesh,
406 useTopo,
407 pp,
408 vectorField(pp.size(), tan1),
409 tan1
410 );
411
412
413 if (wantTan1)
414 {
415 this->operator[](nDirs++) = tan1Dirs;
416 }
417 }
418 if (wantTan2)
419 {
420 tmp<vectorField> tan2Dirs = normalDirs ^ tan1Dirs;
421
422 this->operator[](nDirs++) = tan2Dirs;
423 }
424 }
425 else if (coordSystem == "fieldBased")
426 {
427 forAll(wantedDirs, i)
428 {
429 operator[](nDirs++) =
431 (
433 (
434 mesh.instance()/wantedDirs[i],
435 mesh,
438 )
439 );
440 }
441 }
442 else
443 {
445 << "Unknown coordinate system "
446 << coordSystem << endl
447 << "Known types are global, patchLocal and fieldBased"
448 << exit(FatalError);
449 }
450}
451
452
453// ************************************************************************* //
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
virtual bool resize()
Resize the ODE solver.
Definition: Euler.C:53
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const fileName & instance() const noexcept
Read access to instance path component.
Definition: IOobjectI.H:196
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
T & operator[](const label i)
Return element of UList.
Definition: UListI.H:299
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
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:96
Set of directions for each cell in the mesh. Either uniform and size=1 or one set of directions per c...
Definition: directions.H:70
directionType
Enumeration listing the possible coordinate directions.
Definition: directions.H:77
label nTotalCells() const noexcept
Return total number of cells in decomposed mesh.
static bool test(const UList< face > &faces)
Definition: hexMatcher.C:87
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1121
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1310
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
const vectorField & cellCentres() const
const labelListList & cellCells() const
A class for managing temporary objects.
Definition: tmp.H:65
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:68
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
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
label cutDirToEdge(const primitiveMesh &mesh, const label celli, const vector &cutDir)
Reverse of edgeToCutDir: given direction find edge bundle and.
Definition: meshTools.C:810
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
List< label > labelList
A List of labels.
Definition: List.H:66
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:680
messageStream Info
Information stream (stdout output on master, null elsewhere)
vector point
Point is a vector.
Definition: point.H:43
IOField< vector > vectorIOField
vectorField with IO.
Definition: vectorIOField.H:44
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Field< vector > vectorField
Specialisation of Field<T> for vector.
error FatalError
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333