UList< T > Class Template Reference

A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscript bounds checking, etc. More...

Inheritance diagram for UList< T >:
[legend]
Collaboration diagram for UList< T >:
[legend]

Classes

struct  greater
 A list compare binary predicate for reverse sort. More...
 
struct  Hash
 Hashing function class for UList. More...
 
struct  less
 A list compare binary predicate for normal sort. More...
 

Public Types

typedef T value_type
 The value type the list contains. More...
 
typedef Tpointer
 The pointer type for non-const access to value_type items. More...
 
typedef Treference
 The type used for storing into value_type objects. More...
 
typedef Titerator
 Random access iterator for traversing a UList. More...
 
typedef label size_type
 The type to represent the size of a UList. More...
 
typedef label difference_type
 The difference between iterator objects. More...
 
typedef std::reverse_iterator< iteratorreverse_iterator
 Reverse iterator (non-const access) More...
 
typedef std::reverse_iterator< const_iteratorconst_reverse_iterator
 Reverse iterator (const access) More...
 

Public Member Functions

 UList (const UList< T > &)=default
 Copy construct. More...
 
constexpr UList () noexcept
 Default construct, zero-sized and nullptr. More...
 
 UList (T *__restrict__ v, label size) noexcept
 Construct from components. More...
 
label fcIndex (const label i) const
 
const TfcValue (const label i) const
 Return forward circular value (ie, next value in the list) More...
 
TfcValue (const label i)
 Return forward circular value (ie, next value in the list) More...
 
label rcIndex (const label i) const
 
const TrcValue (const label i) const
 Return reverse circular value (ie, previous value in the list) More...
 
TrcValue (const label i)
 Return reverse circular value (ie, previous value in the list) More...
 
std::streamsize byteSize () const
 
const Tcdata () const
 Return a const pointer to the first data element. More...
 
Tdata ()
 Return a pointer to the first data element. More...
 
Tfirst ()
 Return the first element of the list. More...
 
const Tfirst () const
 Return first element of the list. More...
 
Tlast ()
 Return the last element of the list. More...
 
const Tlast () const
 Return the last element of the list. More...
 
void checkStart (const label start) const
 Check start is within valid range [0,size) More...
 
void checkSize (const label size) const
 Check size is within valid range [0,size]. More...
 
void checkIndex (const label i) const
 Check index is within valid range [0,size) More...
 
bool uniform () const
 True if all entries have identical values, and list is non-empty. More...
 
label find (const T &val, label pos=0) const
 Find index of the first occurrence of the value. More...
 
label rfind (const T &val, label pos=-1) const
 Find index of the last occurrence of the value. More...
 
bool found (const T &val, label pos=0) const
 True if the value if found in the list. More...
 
void moveFirst (const label i)
 Move element to the first position. More...
 
void moveLast (const label i)
 Move element to the last position. More...
 
void swapFirst (const label i)
 Swap element with the first element. Fatal on an empty list. More...
 
void swapLast (const label i)
 Swap element with the last element. Fatal on an empty list. More...
 
void shallowCopy (const UList< T > &list)
 Copy the pointer held by the given UList. More...
 
void deepCopy (const UList< T > &list)
 Copy elements of the given UList. More...
 
Toperator[] (const label i)
 Return element of UList. More...
 
const Toperator[] (const label i) const
 Return element of constant UList. More...
 
UList< Toperator[] (const labelRange &range)
 Return (start,size) subset from UList with non-const access. More...
 
const UList< Toperator[] (const labelRange &range) const
 Return (start,size) subset from UList with const access. More...
 
 operator const Foam::List< T > & () const
 Allow cast to a const List<T>&. More...
 
void operator= (const T &val)
 Assignment of all entries to the given value. More...
 
void operator= (const Foam::zero)
 Assignment of all entries to zero. More...
 
iterator begin ()
 Return an iterator to begin traversing the UList. More...
 
iterator end ()
 Return an iterator to end traversing the UList. More...
 
const_iterator cbegin () const
 Return const_iterator to begin traversing the constant UList. More...
 
const_iterator cend () const
 Return const_iterator to end traversing the constant UList. More...
 
const_iterator begin () const
 Return const_iterator to begin traversing the constant UList. More...
 
const_iterator end () const
 Return const_iterator to end traversing the constant UList. More...
 
reverse_iterator rbegin ()
 Return reverse_iterator to begin reverse traversing the UList. More...
 
reverse_iterator rend ()
 Return reverse_iterator to end reverse traversing the UList. More...
 
const_reverse_iterator crbegin () const
 Return const_reverse_iterator to begin reverse traversing the UList. More...
 
const_reverse_iterator crend () const
 Return const_reverse_iterator to end reverse traversing the UList. More...
 
const_reverse_iterator rbegin () const
 Return const_reverse_iterator to begin reverse traversing the UList. More...
 
const_reverse_iterator rend () const
 Return const_reverse_iterator to end reverse traversing the UList. More...
 
label size () const noexcept
 The number of elements in the UList. More...
 
bool empty () const noexcept
 True if the UList is empty (ie, size() is zero) More...
 
void swap (UList< T > &list)
 Swap content with another UList of the same type in constant time. More...
 
bool operator== (const UList< T > &a) const
 Equality operation on ULists of the same type. More...
 
bool operator!= (const UList< T > &a) const
 The opposite of the equality operation. Takes linear time. More...
 
bool operator< (const UList< T > &list) const
 Compare two ULists lexicographically. Takes linear time. More...
 
bool operator> (const UList< T > &a) const
 Compare two ULists lexicographically. Takes linear time. More...
 
bool operator<= (const UList< T > &a) const
 Return true if !(a > b). Takes linear time. More...
 
bool operator>= (const UList< T > &a) const
 Return true if !(a < b). Takes linear time. More...
 
void writeEntry (const word &keyword, Ostream &os) const
 Write the List as a dictionary entry with keyword. More...
 
OstreamwriteList (Ostream &os, const label shortLen=0) const
 Write List, with line-breaks in ASCII when length exceeds shortLen. More...
 
template<class TypeT = T>
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type test (const label i) const
 A bitSet::test() method for a list of bool. More...
 
template<class TypeT = T>
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type get (const label i) const
 A bitSet::get() method for a list of bool. More...
 
template<class TypeT = T>
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type unset (const label i)
 A bitSet::unset() method for a list of bool. More...
 
template<>
const booloperator[] (const label i) const
 
template<>
Foam::UPstream::commsStructoperator[] (const label procID)
 
template<>
const Foam::UPstream::commsStructoperator[] (const label procID) const
 
template<>
UPstream::commsStructoperator[] (const label)
 
template<>
const UPstream::commsStructoperator[] (const label) const
 

Static Public Member Functions

static const UList< T > & null ()
 Return a UList reference to a nullObject. More...
 
static constexpr label max_size () noexcept
 The size of the largest possible UList. More...
 

Public Attributes

const typedef Tconst_pointer
 The pointer type for const access to value_type items. More...
 
const typedef Tconst_reference
 The type used for reading from constant value_type objects. More...
 
const typedef Tconst_iterator
 Random access iterator for traversing a UList. More...
 

Protected Member Functions

void size (const label n) noexcept
 Override size to be inconsistent with allocated storage. More...
 
void writeEntry (Ostream &os) const
 Write the UList with its compound type. More...
 
labelRange validateRange (const labelRange &range) const
 
UList< T > & operator= (const UList< T > &)=delete
 No copy assignment (default: shallow copy) More...
 

Friends

class List< T >
 Declare friendship with the List class. More...
 
class SubList< T >
 Declare friendship with the SubList class. More...
 
Istreamoperator>> (Istream &os, UList< T > &list)
 Read List contents from Istream. More...
 

Detailed Description

template<class T>
class Foam::UList< T >

A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscript bounds checking, etc.

Storage is not allocated during construction or use but is supplied to the constructor as an argument. This type of list is particularly useful for lists that refer to parts of existing lists such as SubList.

Source files

Definition at line 103 of file HashTable.H.

Member Typedef Documentation

◆ value_type

typedef T value_type

The value type the list contains.

Definition at line 128 of file UList.H.

◆ pointer

typedef T* pointer

The pointer type for non-const access to value_type items.

Definition at line 131 of file UList.H.

◆ reference

typedef T& reference

The type used for storing into value_type objects.

Definition at line 137 of file UList.H.

◆ iterator

typedef T* iterator

Random access iterator for traversing a UList.

Definition at line 143 of file UList.H.

◆ size_type

typedef label size_type

The type to represent the size of a UList.

Definition at line 149 of file UList.H.

◆ difference_type

typedef label difference_type

The difference between iterator objects.

Definition at line 152 of file UList.H.

◆ reverse_iterator

typedef std::reverse_iterator<iterator> reverse_iterator

Reverse iterator (non-const access)

Definition at line 155 of file UList.H.

◆ const_reverse_iterator

typedef std::reverse_iterator<const_iterator> const_reverse_iterator

Reverse iterator (const access)

Definition at line 158 of file UList.H.

Constructor & Destructor Documentation

◆ UList() [1/3]

UList ( const UList< T > &  )
default

Copy construct.

◆ UList() [2/3]

constexpr UList ( )
inlineconstexprnoexcept

Default construct, zero-sized and nullptr.

Definition at line 36 of file UListI.H.

◆ UList() [3/3]

UList ( T *__restrict__  v,
label  size 
)
inlinenoexcept

Construct from components.

Definition at line 44 of file UListI.H.

Member Function Documentation

◆ size() [1/2]

void size ( const label  n)
inlineprotectednoexcept

Override size to be inconsistent with allocated storage.

Use with care

Definition at line 360 of file UListI.H.

Referenced by AABBTree< Type >::AABBTree(), dynamicOversetFvMesh::active(), fvMeshPrimitiveLduAddressing::addAddressing(), dynamicOversetFvMesh::addInterpolation(), faMatrix< Type >::addToInternalField(), fvMatrix< Type >::addToInternalField(), MeshedSurface< Foam::face >::addZones(), pairGAMGAgglomeration::agglomerate(), GAMGAgglomeration::agglomerateLduAddressing(), UPstream::allToAll(), List< substance >::append(), DynamicField< Foam::Vector >::append(), DynamicList< Foam::triPoints >::append(), ITstream::append(), Foam::appendCsvLabels(), Field< Foam::Vector2D >::autoMap(), fvPatchField< vector >::autoMap(), bitSet::bitSet(), isoSurfaceBase::blockCells(), coupledPolyPatch::calcFaceTol(), Foam::checkFaceSizeMatch(), GAMGAgglomeration::checkRestriction(), lduPrimitiveMesh::checkUpperTriangular(), ensightFaces::classify(), ORourkeCollision< CloudType >::collide(), GAMGAgglomeration::compactLevels(), vtuAdaptor::convertField(), Keyed< T >::createList(), Matrix< RectangularMatrix< Type >, Type >::diag(), edgeInterpolationScheme< scalar >::euclidianInterpolate(), Pstream::exchange(), functionObjectList::execute(), regionSizeDistribution::extractData(), FaceCellWave< Type, TrackingData >::FaceCellWave(), NASCore::faceDecomposition(), ABAQUSCore::faceDecomposition(), lduMatrix::faceH(), LduMatrix< Type, DType, LUType >::faceH(), Foam::vtk::Tools::Faces(), fft::fftRenumberRecurse(), FixedList< Type, 3 >::FixedList(), mapDistributeBase::flipAndCombine(), globalIndex::get(), coupledPolyPatch::getAnchorPoints(), Foam::ensightOutput::Detail::getFaceSizes(), weightedPosition::getPoints(), Foam::ensightOutput::Detail::getPolysNFaces(), Foam::Im(), Foam::ImComplexField(), vtuAdaptor::internal(), cyclicACMIGAMGInterface::internalFieldTransfer(), cyclicAMIGAMGInterface::internalFieldTransfer(), edgeInterpolationScheme< scalar >::interpolate(), surfaceInterpolationScheme< GType >::interpolate(), AMIInterpolation::interpolateToSource(), AMIInterpolation::interpolateToTarget(), interpolationTable< scalar >::interpolateValues(), meshToMeshMethod::intersect(), meshToMeshMethod::interVol(), meshToMeshMethod::interVolAndCentroid(), Foam::inv(), cell::labels(), lduPrimitiveMesh::lduPrimitiveMesh(), List< substance >::List(), Field< Foam::Vector2D >::map(), mapDistributeBase::mapDistributeBase(), Foam::matchPoints(), MeshedSurface< Foam::face >::MeshedSurface(), LduMatrix< Type, DType, LUType >::negSumDiag(), lduMatrix::negSumDiag(), noiseFFT::octaves(), noiseModel::octaves(), dlLibraryTable::open(), Hash< List< T > >::operator()(), UList< T >::Hash< HashT >::operator()(), Hash< UList< T > >::operator()(), Foam::operator<<(), BiIndirectList< T >::operator=(), List< substance >::operator=(), FixedList< Type, 3 >::operator=(), HashSet< Foam::string >::operator=(), PackedList< 2 >::PackedList(), vtuAdaptor::points(), Foam::vtk::Tools::Points(), Polynomial< 8 >::Polynomial(), seriesWriter::print(), FIRECore::putFireLabels(), Foam::Re(), hashedWordList::rehash(), Foam::ReImSum(), UnsortedMeshedSurface< Face >::remapFaces(), MeshedSurface< Foam::face >::remapFaces(), UPtrList< Foam::diameterModels::velocityGroup >::reorder(), Foam::reverse(), face::reverseFace(), SortList< T >::reverseSort(), sampledSurface::sampleOnFaces(), sampledSurface::sampleOnPoints(), triSurfaceLoader::select(), weightedPosition::setPoints(), removeCells::setRefinement(), CompactListList< T, Container >::setSize(), Foam::ListOps::setValue(), UnsortedMeshedSurface< Face >::setZones(), directFieldMapper::size(), directPointPatchFieldMapper::size(), directFvPatchFieldMapper::size(), distributedUnallocatedDirectFieldMapper::size(), distributedUnallocatedDirectFvPatchFieldMapper::size(), LLTMatrix< Type >::solve(), SortList< T >::sort(), UPtrList< Foam::diameterModels::velocityGroup >::sortOrder(), extendedEdgeMesh::sortPointsAndEdges(), Foam::stabilise(), faMatrix< Type >::subtractFromInternalField(), fvMatrix< Type >::subtractFromInternalField(), Foam::sumCmptProd(), LduMatrix< Type, DType, LUType >::sumDiag(), lduMatrix::sumDiag(), LduMatrix< Type, DType, LUType >::sumMagOffDiag(), lduMatrix::sumMagOffDiag(), Foam::sumProd(), syncTools::swapBoundaryCellList(), syncTools::swapBoundaryCellPositions(), syncTools::syncBoundaryFaceList(), syncTools::syncPointList(), weightedPosition::syncPoints(), Foam::system(), meshRefinement::testSyncBoundaryFaceList(), fft::transform(), PackedList< 2 >::unpack(), dynamicOversetFvMesh::updateAddressing(), topoBoolSet::updateLabels(), topoBitSet::updateLabels(), topoSet::updateLabels(), Foam::HashSetOps::used(), Foam::vtk::Tools::Vertices(), FLMAsurfaceFormat< Face >::write(), OBJsurfaceFormat< Face >::write(), SMESHsurfaceFormat< Face >::write(), OFFsurfaceFormat< Face >::write(), GTSsurfaceFormat< Face >::write(), ABAQUSsurfaceFormat< Face >::write(), dynamicOversetFvMesh::write(), STARCDsurfaceFormat< Face >::write(), NASsurfaceFormat< Face >::write(), VTPsurfaceFormat< Face >::write(), OBJstream::write(), STARCDsurfaceFormatCore::writeCase(), VTKsurfaceFormatCore::writeCellData(), VTPsurfaceFormatCore::writeCellData(), lumpedPointMovement::writeData(), VTKedgeFormat::writeEdges(), lumpedPointMovement::writeForcesAndMomentsVTP(), VTKsurfaceFormatCore::writeHeader(), VTPsurfaceFormatCore::writeHeader(), AC3DsurfaceFormatCore::writeHeader(), UList< Foam::wordRe >::writeList(), coupledPolyPatch::writeOBJ(), Foam::meshTools::writeOBJ(), Foam::ensightOutput::writePolysPoints(), nastranWriter::writeTemplate(), starcdWriter::writeTemplate(), and Foam::zip().

◆ writeEntry() [1/2]

void writeEntry ( Ostream os) const
protected

Write the UList with its compound type.

Definition at line 38 of file UListIO.C.

Referenced by faceZone::writeDict().

Here is the caller graph for this function:

◆ validateRange()

Foam::labelRange validateRange ( const labelRange range) const
protected

Return a validated (start,size) subset range, which means that it always addresses a valid section of the list.

Definition at line 39 of file UList.C.

◆ operator=() [1/3]

UList<T>& operator= ( const UList< T > &  )
protecteddelete

No copy assignment (default: shallow copy)

Assignment may need to be shallow (copy pointer) or deep (copy elements) depending on context or type of list. Disallow default assignment and provide separate 'shallowCopy' and 'deepCopy' member functions.

Referenced by SortableList< T >::operator=().

Here is the caller graph for this function:

◆ null()

const Foam::UList< T > & null ( )
inlinestatic

Return a UList reference to a nullObject.

Definition at line 54 of file UListI.H.

Referenced by FieldMapper::directAddressing(), pointFieldReconstructor::pointPatchFieldReconstructor::directAddressing(), fvFieldReconstructor::fvPatchFieldReconstructor::directAddressing(), faFieldReconstructor::faPatchFieldReconstructor::directAddressing(), and starcdWriter::write().

Here is the caller graph for this function:

◆ fcIndex()

Foam::label fcIndex ( const label  i) const
inline

Return the forward circular index, i.e. next index which returns to the first at the end of the list

Definition at line 61 of file UListI.H.

◆ fcValue() [1/2]

const T & fcValue ( const label  i) const
inline

Return forward circular value (ie, next value in the list)

Definition at line 68 of file UListI.H.

◆ fcValue() [2/2]

T & fcValue ( const label  i)
inline

Return forward circular value (ie, next value in the list)

Definition at line 75 of file UListI.H.

◆ rcIndex()

Foam::label rcIndex ( const label  i) const
inline

Return the reverse circular index, i.e. previous index which returns to the last at the beginning of the list

Definition at line 82 of file UListI.H.

◆ rcValue() [1/2]

const T & rcValue ( const label  i) const
inline

Return reverse circular value (ie, previous value in the list)

Definition at line 89 of file UListI.H.

◆ rcValue() [2/2]

T & rcValue ( const label  i)
inline

Return reverse circular value (ie, previous value in the list)

Definition at line 96 of file UListI.H.

◆ byteSize()

std::streamsize byteSize ( ) const

Return the binary size in number of characters of the UList if the element is a primitive type

i.e. is_contiguous<T>::value == true. Note that is of type streamsize since used in stream ops

Definition at line 191 of file UList.C.

Referenced by OFstreamCollator::write(), and UList< Foam::wordRe >::writeList().

Here is the caller graph for this function:

◆ cdata()

const T * cdata ( ) const
inline

Return a const pointer to the first data element.

Similar to the STL front() method and the string::data() method This can be used (with caution) when interfacing with C code

Definition at line 198 of file UListI.H.

Referenced by Hash< List< T > >::operator()(), UList< T >::Hash< HashT >::operator()(), Hash< UList< T > >::operator()(), OFstreamCollator::write(), and UList< Foam::wordRe >::writeList().

Here is the caller graph for this function:

◆ data()

T * data ( )
inline

Return a pointer to the first data element.

Similar to the STL front() method and the string::data() method This can be used (with caution) when interfacing with C code

Definition at line 205 of file UListI.H.

◆ first() [1/2]

T & first ( )
inline

Return the first element of the list.

Definition at line 170 of file UListI.H.

Referenced by weightedPosition::getPoints(), and Foam::printTable().

Here is the caller graph for this function:

◆ first() [2/2]

const T & first ( ) const
inline

Return first element of the list.

Definition at line 177 of file UListI.H.

◆ last() [1/2]

T & last ( )
inline

Return the last element of the list.

Definition at line 184 of file UListI.H.

Referenced by globalIndex::gather().

Here is the caller graph for this function:

◆ last() [2/2]

const T & last ( ) const
inline

Return the last element of the list.

Definition at line 191 of file UListI.H.

◆ checkStart()

void checkStart ( const label  start) const
inline

Check start is within valid range [0,size)

Definition at line 103 of file UListI.H.

Referenced by SubList< T >::SubList().

Here is the caller graph for this function:

◆ checkSize()

void checkSize ( const label  size) const
inline

Check size is within valid range [0,size].

Definition at line 116 of file UListI.H.

Referenced by SubList< T >::SubList().

Here is the caller graph for this function:

◆ checkIndex()

void checkIndex ( const label  i) const
inline

Check index is within valid range [0,size)

Definition at line 128 of file UListI.H.

◆ uniform()

bool uniform ( ) const
inline

True if all entries have identical values, and list is non-empty.

Definition at line 146 of file UListI.H.

Referenced by UList< Foam::wordRe >::writeList().

Here is the caller graph for this function:

◆ find()

Foam::label find ( const T val,
label  pos = 0 
) const

Find index of the first occurrence of the value.

Any occurrences before the start pos are ignored. Linear search.

Returns
position in list or -1 if not found.

Definition at line 205 of file UList.C.

◆ rfind()

Foam::label rfind ( const T val,
label  pos = -1 
) const

Find index of the last occurrence of the value.

Any occurrences after the end pos are ignored. Linear search.

Returns
position in list or -1 if not found.

Definition at line 229 of file UList.C.

◆ found()

bool found ( const T val,
label  pos = 0 
) const
inline

True if the value if found in the list.

Any occurrences before the start pos are ignored. Linear search.

Returns
true if found.

Definition at line 212 of file UListI.H.

Referenced by findNearestMaskedOp< PatchType >::operator()(), and foundOp< StringType >::operator()().

Here is the caller graph for this function:

◆ moveFirst()

void moveFirst ( const label  i)

Move element to the first position.

Definition at line 55 of file UList.C.

◆ moveLast()

void moveLast ( const label  i)

Move element to the last position.

Definition at line 67 of file UList.C.

◆ swapFirst()

void swapFirst ( const label  i)

Swap element with the first element. Fatal on an empty list.

Definition at line 79 of file UList.C.

◆ swapLast()

void swapLast ( const label  i)

Swap element with the last element. Fatal on an empty list.

Definition at line 91 of file UList.C.

◆ shallowCopy()

void shallowCopy ( const UList< T > &  list)
inline

Copy the pointer held by the given UList.

Definition at line 219 of file UListI.H.

Referenced by surfMesh::updateRefs().

Here is the caller graph for this function:

◆ deepCopy()

void deepCopy ( const UList< T > &  list)

Copy elements of the given UList.

Definition at line 105 of file UList.C.

Referenced by UPstream::allToAll().

Here is the caller graph for this function:

◆ operator[]() [1/9]

T & operator[] ( const label  i)
inline

Return element of UList.

Definition at line 246 of file UListI.H.

◆ operator[]() [2/9]

const T & operator[] ( const label  i) const
inline

Return element of constant UList.

Note that the bool specialization adds lazy evaluation so reading an out-of-range element returns false without any ill-effects

Definition at line 256 of file UListI.H.

◆ operator[]() [3/9]

Foam::UList< T > operator[] ( const labelRange range)

Return (start,size) subset from UList with non-const access.

The range is subsetted with the list size itself to ensure that the result always addresses a valid section of the list.

Definition at line 143 of file UList.C.

◆ operator[]() [4/9]

const Foam::UList< T > operator[] ( const labelRange range) const

Return (start,size) subset from UList with const access.

The range is subsetted with the list size itself to ensure that the result always addresses a valid section of the list.

Definition at line 152 of file UList.C.

◆ operator const Foam::List< T > &()

operator const Foam::List< T > & ( ) const
inline

Allow cast to a const List<T>&.

Definition at line 266 of file UListI.H.

◆ operator=() [2/3]

void operator= ( const T val)

Assignment of all entries to the given value.

Definition at line 161 of file UList.C.

◆ operator=() [3/3]

void operator= ( const Foam::zero  )

Assignment of all entries to zero.

Definition at line 175 of file UList.C.

◆ begin() [1/2]

◆ end() [1/2]

◆ cbegin()

Foam::UList< T >::const_iterator cbegin ( ) const
inline

Return const_iterator to begin traversing the constant UList.

Definition at line 290 of file UListI.H.

Referenced by IndirectListBase< T, Addr >::begin(), IndirectListBase< T, Addr >::cbegin(), HashTable< Foam::phase * >::erase(), and Foam::isFaceCut().

Here is the caller graph for this function:

◆ cend()

Foam::UList< T >::const_iterator cend ( ) const
inline

Return const_iterator to end traversing the constant UList.

Definition at line 311 of file UListI.H.

Referenced by IndirectListBase< T, Addr >::cend(), IndirectListBase< T, Addr >::end(), HashTable< Foam::phase * >::erase(), and Foam::isFaceCut().

Here is the caller graph for this function:

◆ begin() [2/2]

Foam::UList< T >::const_iterator begin ( ) const
inline

Return const_iterator to begin traversing the constant UList.

Definition at line 283 of file UListI.H.

◆ end() [2/2]

Foam::UList< T >::const_iterator end ( ) const
inline

Return const_iterator to end traversing the constant UList.

Definition at line 304 of file UListI.H.

◆ rbegin() [1/2]

Foam::UList< T >::reverse_iterator rbegin ( )
inline

Return reverse_iterator to begin reverse traversing the UList.

Definition at line 318 of file UListI.H.

◆ rend() [1/2]

Foam::UList< T >::reverse_iterator rend ( )
inline

Return reverse_iterator to end reverse traversing the UList.

Definition at line 339 of file UListI.H.

◆ crbegin()

Foam::UList< T >::const_reverse_iterator crbegin ( ) const
inline

Return const_reverse_iterator to begin reverse traversing the UList.

Definition at line 332 of file UListI.H.

◆ crend()

Foam::UList< T >::const_reverse_iterator crend ( ) const
inline

Return const_reverse_iterator to end reverse traversing the UList.

Definition at line 353 of file UListI.H.

◆ rbegin() [2/2]

Foam::UList< T >::const_reverse_iterator rbegin ( ) const
inline

Return const_reverse_iterator to begin reverse traversing the UList.

Definition at line 325 of file UListI.H.

◆ rend() [2/2]

Foam::UList< T >::const_reverse_iterator rend ( ) const
inline

Return const_reverse_iterator to end reverse traversing the UList.

Definition at line 346 of file UListI.H.

◆ size() [2/2]

Foam::label size ( ) const
inlinenoexcept

The number of elements in the UList.

Definition at line 367 of file UListI.H.

Referenced by UList< Foam::wordRe >::test().

Here is the caller graph for this function:

◆ empty()

bool empty ( ) const
inlinenoexcept

True if the UList is empty (ie, size() is zero)

Definition at line 374 of file UListI.H.

Referenced by AABBTree< Type >::AABBTree(), Foam::printTable(), UnsortedMeshedSurface< Face >::remapFaces(), MeshedSurface< Foam::face >::remapFaces(), Foam::system(), and treeBoundBox::treeBoundBox().

Here is the caller graph for this function:

◆ max_size()

static constexpr label max_size ( )
inlinestaticconstexprnoexcept

The size of the largest possible UList.

Definition at line 427 of file UList.H.

◆ swap()

void swap ( UList< T > &  list)
inline

Swap content with another UList of the same type in constant time.

Definition at line 381 of file UListI.H.

Referenced by Foam::Swap().

Here is the caller graph for this function:

◆ operator==()

bool operator== ( const UList< T > &  a) const

Equality operation on ULists of the same type.

Returns true when the ULists are element-wise equal (using UList::value_type::operator==). Takes linear time

Definition at line 291 of file UList.C.

◆ operator!=()

bool operator!= ( const UList< T > &  a) const

The opposite of the equality operation. Takes linear time.

Definition at line 315 of file UList.C.

◆ operator<()

bool operator< ( const UList< T > &  list) const

Compare two ULists lexicographically. Takes linear time.

Definition at line 322 of file UList.C.

◆ operator>()

bool operator> ( const UList< T > &  a) const

Compare two ULists lexicographically. Takes linear time.

Definition at line 347 of file UList.C.

◆ operator<=()

bool operator<= ( const UList< T > &  a) const

Return true if !(a > b). Takes linear time.

Definition at line 354 of file UList.C.

◆ operator>=()

bool operator>= ( const UList< T > &  a) const

Return true if !(a < b). Takes linear time.

Definition at line 361 of file UList.C.

◆ writeEntry() [2/2]

void writeEntry ( const word keyword,
Ostream os 
) const

Write the List as a dictionary entry with keyword.

Definition at line 66 of file UListIO.C.

◆ writeList()

Foam::Ostream & writeList ( Ostream os,
const label  shortLen = 0 
) const

Write List, with line-breaks in ASCII when length exceeds shortLen.

Using '0' suppresses line-breaks entirely.

Definition at line 79 of file UListIO.C.

Referenced by Foam::operator<<(), and ensightCells::writeDict().

Here is the caller graph for this function:

◆ test()

std::enable_if<std::is_same<bool, TypeT>::value, bool>::type test ( const label  i) const
inline

A bitSet::test() method for a list of bool.

Returns
The element value, or false for out-of-range access

Definition at line 484 of file UList.H.

◆ get()

std::enable_if<std::is_same<bool, TypeT>::value, bool>::type get ( const label  i) const
inline

A bitSet::get() method for a list of bool.

Returns
The element value, or false for out-of-range access

Definition at line 494 of file UList.H.

◆ unset()

std::enable_if<std::is_same<bool, TypeT>::value, bool>::type unset ( const label  i)
inline

A bitSet::unset() method for a list of bool.

Returns
True if value changed and was not out-of-range

Definition at line 504 of file UList.H.

◆ operator[]() [5/9]

const bool & operator[] ( const label  i) const
inline

Definition at line 232 of file UListI.H.

◆ operator[]() [6/9]

Foam::UPstream::commsStruct & operator[] ( const label  procID)

Definition at line 252 of file UPstream.C.

◆ operator[]() [7/9]

const Foam::UPstream::commsStruct & operator[] ( const label  procID) const

Definition at line 353 of file UPstream.C.

◆ operator[]() [8/9]

UPstream::commsStruct & operator[] ( const  label)

◆ operator[]() [9/9]

const UPstream::commsStruct & operator[] ( const  label) const

Friends And Related Function Documentation

◆ List< T >

friend class List< T >
friend

Declare friendship with the List class.

Definition at line 164 of file UList.H.

◆ SubList< T >

friend class SubList< T >
friend

Declare friendship with the SubList class.

Definition at line 167 of file UList.H.

◆ operator>>

Istream& operator>> ( Istream os,
UList< T > &  list 
)
friend

Read List contents from Istream.

Requires size to have been set before

Member Data Documentation

◆ const_pointer

const typedef T* const_pointer

The pointer type for const access to value_type items.

Definition at line 134 of file UList.H.

◆ const_reference

const typedef T& const_reference

The type used for reading from constant value_type objects.

Definition at line 140 of file UList.H.

◆ const_iterator

const typedef T* const_iterator

Random access iterator for traversing a UList.

Definition at line 146 of file UList.H.


The documentation for this class was generated from the following files: