structuredRenumber.H
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) 2012-2016 OpenFOAM Foundation
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 Class
28  Foam::structuredRenumber
29 
30 Description
31  Renumbering according to mesh layers.
32  depthFirst = true:
33  first column gets ids 0..nLayer-1,
34  second nLayers..2*nLayers-1 etc.
35  depthFirst = false:
36  first layer gets ids 0,1,2 etc.
37 
38 SourceFiles
39  structuredRenumber.C
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #ifndef Foam_structuredRenumber_H
44 #define Foam_structuredRenumber_H
45 
46 #include "renumberMethod.H"
47 
48 namespace Foam
49 {
50 
51 // Forward Declarations
52 template<class Type> class topoDistanceData;
53 
54 /*---------------------------------------------------------------------------*\
55  Class structuredRenumber Declaration
56 \*---------------------------------------------------------------------------*/
57 
59 :
60  public renumberMethod
61 {
62 public:
63 
64  // Public Classes
65 
66  //- Less function class that can be used for sorting according to
67  // column and layer
68  class layerLess
69  {
70  const bool depthFirst_;
71  const labelList& order_;
72  const List<topoDistanceData<label>>& distance_;
73 
74  public:
75 
76  layerLess
77  (
78  const bool depthFirst,
79  const labelList& order,
81  )
82  :
83  depthFirst_(depthFirst),
84  order_(order),
85  distance_(distance)
86  {}
87 
88  bool operator()(const label a, const label b);
89  };
90 
91 
92  // Private Data
93 
95 
97 
98  const label nLayers_;
99 
100  const bool depthFirst_;
101 
102  const bool reverse_;
103 
105 
106 
107  // Private Member Functions
108 
109  //- No copy construct
110  structuredRenumber(const structuredRenumber&) = delete;
111 
112  //- No copy assignment
113  void operator=(const structuredRenumber&) = delete;
114 
115 
116 public:
117 
118  //- Runtime type information
119  TypeName("structured");
120 
121 
122  // Constructors
123 
124  //- Construct given the renumber dictionary
125  explicit structuredRenumber(const dictionary& dict);
126 
127 
128  //- Destructor
129  virtual ~structuredRenumber() = default;
130 
131 
132  // Member Functions
133 
134  //- Return the order in which cells need to be visited
135  //- (ie. from ordered back to original cell label).
136  // This is only defined for geometric renumberMethods.
137  virtual labelList renumber(const pointField&) const
138  {
140  return labelList();
141  }
142 
143  //- Return the order in which cells need to be visited
144  //- (ie. from ordered back to original cell label).
145  // Use the mesh connectivity (if needed)
146  virtual labelList renumber
147  (
148  const polyMesh& mesh,
149  const pointField& cc
150  ) const;
151 
152  //- Return the order in which cells need to be visited
153  //- (ie. from ordered back to original cell label).
154  virtual labelList renumber
155  (
156  const CompactListList<label>& cellCells,
157  const pointField& cellCentres
158  ) const
159  {
161  return labelList();
162  }
163 
164  //- Return the order in which cells need to be visited
165  //- (ie. from ordered back to original cell label).
166  virtual labelList renumber
167  (
168  const labelListList& cellCells,
169  const pointField& cellCentres
170  ) const
171  {
173  return labelList();
174  }
175 };
176 
177 
178 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179 
180 } // End namespace Foam
181 
182 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
183 
184 #endif
186 // ************************************************************************* //
Renumbering according to mesh layers. depthFirst = true: first column gets ids 0..nLayer-1, second nLayers..2*nLayers-1 etc. depthFirst = false: first layer gets ids 0,1,2 etc.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
bool operator()(const label a, const label b)
Less function class that can be used for sorting according to.
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Abstract base class for renumbering.
TypeName("structured")
Runtime type information.
dynamicFvMesh & mesh
layerLess(const bool depthFirst, const labelList &order, const List< topoDistanceData< label >> &distance)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
structuredRenumber(const structuredRenumber &)=delete
No copy construct.
const autoPtr< renumberMethod > method_
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:53
A packed storage of objects of type <T> using an offset table for access.
virtual ~structuredRenumber()=default
Destructor.
For use with FaceCellWave. Determines topological distance to starting faces. Templated on passive tr...
void operator=(const structuredRenumber &)=delete
No copy assignment.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< label > labelList
A List of labels.
Definition: List.H:62
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited (ie. from ordered back to original cell label)...
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:686
Namespace for OpenFOAM.
const dictionary & coeffsDict_