structuredRenumber.H
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) 2012-2019 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 Class
25  Foam::structuredRenumber
26 
27 Description
28  Renumbering according to mesh layers.
29  depthFirst = true:
30  first column gets ids 0..nLayer-1,
31  second nLayers..2*nLayers-1 etc.
32  depthFirst = false:
33  first layer gets ids 0,1,2 etc.
34 
35 SourceFiles
36  structuredRenumber.C
37 
38 \*---------------------------------------------------------------------------*/
39 
40 #ifndef structuredRenumber_H
41 #define structuredRenumber_H
42 
43 #include "renumberMethod.H"
44 #include "topoDistanceData.H"
45 #include "Switch.H"
46 
47 namespace Foam
48 {
49 
50 /*---------------------------------------------------------------------------*\
51  Class structuredRenumber Declaration
52 \*---------------------------------------------------------------------------*/
53 
55 :
56  public renumberMethod
57 {
58 public:
59 
60  // Public classes
61 
62  //- Less function class that can be used for sorting according to
63  // column and layer
64  class layerLess
65  {
66  const Switch depthFirst_;
67  const labelList& order_;
68  const List<topoDistanceData>& distance_;
69 
70  public:
71 
73  (
74  const Switch depthFirst,
75  const labelList& order,
77  )
78  :
79  depthFirst_(depthFirst),
80  order_(order),
81  distance_(distance)
82  {}
83 
84  bool operator()(const label a, const label b);
85  };
86 
87 
88  // Private Data
89 
90  const dictionary methodDict_;
91 
92  const wordReList patches_;
93 
94  const label nLayers_;
95 
96  const Switch depthFirst_;
97 
99 
100  const Switch reverse_;
101 
102 
103 public:
104 
105  //- Runtime type information
106  TypeName("structured");
107 
108 
109  // Constructors
110 
111  //- Construct given the renumber dictionary
112  structuredRenumber(const dictionary& renumberDict);
113 
114  //- Disallow default bitwise copy construction
115  structuredRenumber(const structuredRenumber&) = delete;
116 
117 
118  //- Destructor
119  virtual ~structuredRenumber()
120  {}
121 
122 
123  // Member Functions
124 
125  //- Return the order in which cells need to be visited, i.e.
126  // from ordered back to original cell label.
127  // This is only defined for geometric renumberMethods.
128  virtual labelList renumber(const pointField&) const
129  {
131  return labelList(0);
132  }
133 
134  //- Return the order in which cells need to be visited, i.e.
135  // from ordered back to original cell label.
136  // Use the mesh connectivity (if needed)
137  virtual labelList renumber
138  (
139  const polyMesh& mesh,
140  const pointField& cc
141  ) const;
142 
143  //- Return the order in which cells need to be visited, i.e.
144  // from ordered back to original cell label.
145  // The connectivity is equal to mesh.cellCells() except
146  // - the connections are across coupled patches
148  (
149  const labelListList& cellCells,
150  const pointField& cc
151  ) const
152  {
154  return labelList(0);
155  }
156 
157 
158  // Member Operators
159 
160  //- Disallow default bitwise assignment
161  void operator=(const structuredRenumber&) = delete;
162 };
163 
164 
165 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
166 
167 } // End namespace Foam
168 
169 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
170 
171 #endif
172 
173 // ************************************************************************* //
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.
virtual ~structuredRenumber()
Destructor.
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
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
structuredRenumber(const dictionary &renumberDict)
Construct given the renumber dictionary.
bool operator()(const label a, const label b)
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none/any.
Definition: Switch.H:60
Less function class that can be used for sorting according to.
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
const dictionary methodDict_
Abstract base class for renumbering.
TypeName("structured")
Runtime type information.
dynamicFvMesh & mesh
const dimensionedScalar & b
Wien displacement law constant: default SI units: [m K].
Definition: createFields.H:27
const autoPtr< renumberMethod > method_
List< label > labelList
A List of labels.
Definition: labelList.H:56
void operator=(const structuredRenumber &)=delete
Disallow default bitwise assignment.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
virtual labelList renumber(const pointField &) const
Return the order in which cells need to be visited, i.e.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
layerLess(const Switch depthFirst, const labelList &order, const List< topoDistanceData > &distance)
Namespace for OpenFOAM.