Loading...
Searching...
No Matches
SortListI.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-2023 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
26\*---------------------------------------------------------------------------*/
27
28#include "ListOps.H" // For uniqueOrder()
29
30// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31
32template<class T>
34:
37 sort();
38}
39
40
41template<class T>
42template<class Compare>
43inline Foam::SortList<T>::SortList(const UList<T>& values, const Compare& comp)
44:
45 IndirectList<T>(values, labelList())
46{
48}
49
50
51// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
52
53template<class T>
55{
56 return this->addressing();
57}
58
59
60template<class T>
65
66
67template<class T>
69{
70 Foam::reverse(this->indices());
71}
72
73
74template<class T>
75inline void Foam::SortList<T>::reset()
76{
77 auto& addr = this->indices();
78
79 addr.resize_nocopy(this->values().size());
80 Foam::identity(addr, 0);
81}
82
83
84template<class T>
85template<class Compare>
86inline void Foam::SortList<T>::sort(const Compare& comp)
87{
88 auto& vals = this->values();
89 auto& addr = this->indices();
90
91 addr.resize_nocopy(vals.size());
92 Foam::identity(addr, 0);
93
94 std::stable_sort
95 (
96 addr.begin(),
97 addr.end(),
98 [&](label a, label b) -> bool { return comp(vals[a], vals[b]); }
99 );
100}
101
102
103template<class T>
105{
106 Foam::sortedOrder(this->values(), this->indices());
107}
108
109
110template<class T>
112{
113 Foam::uniqueOrder(this->values(), this->indices());
114}
115
116
117template<class T>
119{
120 // Reverse sorted order for indices
122 (
123 this->values(),
124 this->indices(),
125 typename UList<T>::greater(this->values())
126 );
127}
128
129
130// ************************************************************************* //
Various functions to operate on Lists.
const UList< T > & values() const noexcept
A List with indirect addressing.
IndirectList(const UList< T > &values, const labelUList &addr)
Copy construct addressing, shallow copy values reference.
const Addr & addressing() const noexcept
The list addressing.
void reverse()
Reverse the indices.
Definition SortListI.H:61
void sort()
Forward (stable) sort the list. Functionally identical to sort with std::less<T>().
Definition SortListI.H:97
const labelUList & indices() const noexcept
Return the list of sorted indices (updated every sort).
Definition SortListI.H:47
void reverseSort()
Reverse (stable) sort the list. Functionally identical to sort with std::greater<T>().
Definition SortListI.H:111
void uniqueSort()
Sort the list, only retaining unique entries.
Definition SortListI.H:104
void sort(const Compare &comp)
Sort the list using the given value comparator.
Definition SortListI.H:79
void reset()
Reset list indices to identity.
Definition SortListI.H:68
SortList(const UList< T > &values)
Shallow copy values list reference, sort immediately.
Definition SortListI.H:26
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:89
List< label > labelList
A List of labels.
Definition List.H:62
void reverse(UList< T > &list, const label n)
Reverse the first n elements of the list.
Definition UListI.H:539
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
const direction noexcept
Definition scalarImpl.H:265
labelList uniqueOrder(const UList< T > &input)
Return (sorted) indices corresponding to unique list values.
UList< label > labelUList
A UList of labels.
Definition UList.H:75
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
volScalarField & b
A list compare binary predicate for reverse sort.
Definition UList.H:255