flip_zoneGenerator.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) 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 "flip_zoneGenerator.H"
27 #include "polyMesh.H"
28 #include "syncTools.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace zoneGenerators
36  {
39  (
41  flip,
43  );
44  }
45 }
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
50 (
51  const word& name,
52  const polyMesh& mesh,
53  const dictionary& dict
54 )
55 :
57  zoneGenerator_(zoneGenerator::New(mesh, dict))
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
70 {
71  zoneSet zs(zoneGenerator_->generate());
72  const faceZone& fZone = zs.fZone();
73  boolList flipMap(fZone.flipMap());
74 
75  forAll(flipMap, fi)
76  {
77  flipMap[fi] = !flipMap[fi];
78  }
79 
80  return zoneSet
81  (
82  new faceZone
83  (
84  fZone,
85  zoneName_,
86  fZone,
87  flipMap,
88  mesh_.faceZones()
89  )
90  );
91 }
92 
93 
94 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Named list of face indices representing a sub-set of the mesh faces.
Definition: faceZone.H:66
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:262
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
A class for handling words, derived from string.
Definition: word.H:62
Abstract base class for all zoneGenerators, providing runtime selection.
Definition: zoneGenerator.H:57
A zoneGenerator which flips the flipMap of the faceZone returned by the given zoneGenerator.
virtual zoneSet generate() const
Generate and return the zoneSet.
virtual ~flip()
Destructor.
flip(const word &name, const polyMesh &mesh, const dictionary &dict)
Construct from dictionary.
Zone container returned by zoneGenerator::generate.
Definition: zoneSet.H:94
const faceZone & fZone() const
Return a reference to the faceZone if allocated.
Definition: zoneSetI.H:230
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
defineTypeNameAndDebug(cylinderHeadPoints, 0)
addToRunTimeSelectionTable(zoneGenerator, cylinderHeadPoints, dictionary)
Namespace for OpenFOAM.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
dictionary dict