SortList.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) 2019-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26Class
27 Foam::SortList
28
29Description
30 An indirect list with addressing based on sorting.
31 The list is sorted upon construction or when explicitly requested.
32
33 Uses the std::stable_sort() algorithm.
34
35SourceFiles
36 SortListI.H
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef Foam_SortList_H
41#define Foam_SortList_H
42
43#include "IndirectList.H"
44
45// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46
47namespace Foam
48{
49
50/*---------------------------------------------------------------------------*\
51 Class SortList Declaration
52\*---------------------------------------------------------------------------*/
53
54template<class T>
55class SortList
56:
57 public IndirectList<T>
58{
59public:
60
61 // Constructors
62
63 //- Shallow copy values list reference, sort immediately
64 inline explicit SortList(const UList<T>& values);
65
66 //- Shallow copy values list reference,
67 //- sort with given \b value comparator.
68 // \note The comparator is not stored in the class.
69 template<class Compare>
70 inline SortList(const UList<T>& values, const Compare& comp);
71
72
73 // Member Functions
74
75 //- Return the list of sorted indices (updated every sort).
76 // Same as addressing()
77 inline const labelUList& indices() const noexcept;
78
79 //- Return the list of indices (updated every sort).
80 // Same as addressing()
81 inline labelList& indices() noexcept;
82
83 //- Reverse the indices
84 inline void reverse();
85
86 //- Reset list indices to identity
87 inline void reset();
88
89 //- Sort the list using the given \b value comparator
90 template<class Compare>
91 inline void sort(const Compare& comp);
92
93 //- Forward (stable) sort the list.
94 //- Functionally identical to sort with std::less<T>()
95 inline void sort();
96
97 //- Reverse (stable) sort the list.
98 //- Functionally identical to sort with std::greater<T>()
99 inline void reverseSort();
100
101 //- Sort the list, only retaining unique entries
102 inline void uniqueSort();
103
104
105 // Member Operators
106
107 //- Assignment operators
108 using IndirectList<T>::operator=;
109};
110
111
112// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
113
114} // End namespace Foam
115
116// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
117
118#include "SortListI.H"
119
120// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121
122#endif
123
124// ************************************************************************* //
const UList< T > & values() const noexcept
The list of values (without addressing)
A List with indirect addressing.
Definition: IndirectList.H:119
An indirect list with addressing based on sorting. The list is sorted upon construction or when expli...
Definition: SortList.H:57
void reverse()
Reverse the indices.
Definition: SortListI.H:68
const labelUList & indices() const noexcept
Return the list of sorted indices (updated every sort).
Definition: SortListI.H:54
void reverseSort()
Definition: SortListI.H:118
void uniqueSort()
Sort the list, only retaining unique entries.
Definition: SortListI.H:111
void reset()
Reset list indices to identity.
Definition: SortListI.H:75
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: UList.H:94
const volScalarField & T
Namespace for OpenFOAM.