All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
Classes | Public Types | Public Member Functions | Static Public Member Functions | Friends | List of all members
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...

Inherited by List< T >, pTraits< UList< T > >, SubList< T >, List< bool >, List< cell >, List< cellCutType >, List< cellShape >, List< CGAL::indexedVertex >, List< char >, List< const Foam::cellModel *>, List< const Foam::UList *>, List< const label *>, List< const lduInterface *>, List< const lduInterfaceField *>, List< const LduInterfaceField< Type > *>, List< edge >, List< Face >, List< faceList >, List< Field< scalar > >, List< Field< Type > >, List< fileName >, List< Foam::boundBox >, List< Foam::cellModel *>, List< Foam::dictionary >, List< Foam::DynamicList< char > >, List< Foam::DynamicList< DSMCParcel< ParcelType > *> >, List< Foam::DynamicList< Foam::List > >, List< Foam::DynamicList< Foam::molecule *> >, List< Foam::DynamicList< Foam::string > >, List< Foam::DynamicList< Foam::Vector > >, List< Foam::DynamicList< label > >, List< Foam::DynamicList< ParcelType *> >, List< Foam::DynamicList< parcelType *> >, List< Foam::DynamicList< scalar > >, List< Foam::face >, List< Foam::faceSets >, List< Foam::Field >, List< Foam::Field< Foam::Vector > >, List< Foam::Field< Type > *>, List< Foam::FixedList< label, 3 > >, List< Foam::FixedList< label, 4 > >, List< Foam::FixedList< scalar, 2 > >, List< Foam::FixedList< scalar, 3 > >, List< Foam::functionObjects::fieldAverageItem >, List< Foam::fvFieldDecomposer::patchFieldDecomposer *>, List< Foam::fvFieldDecomposer::processorSurfacePatchFieldDecomposer *>, List< Foam::fvFieldDecomposer::processorVolPatchFieldDecomposer *>, List< Foam::GeometricField *>, List< Foam::gradingDescriptors >, List< Foam::ILList< Foam::molecule > >, List< Foam::ILList< ParticleType > >, List< Foam::ILList< typename CloudType::parcelType > >, List< Foam::indexedOctree::node >, List< Foam::List >, List< Foam::List< Foam::face > >, List< Foam::List< Foam::Field< scalar > > >, List< Foam::List< Foam::List< scalar > > >, List< Foam::List< Foam::meshReader::cellFaceIdentifier > >, List< Foam::List< Foam::specieElement > >, List< Foam::List< Foam::treeBoundBox > >, List< Foam::List< Foam::Vector > >, List< Foam::List< Key > >, List< Foam::List< label > >, List< Foam::List< scalar > >, List< Foam::List< T > >, List< Foam::molecule::constantProperties >, List< Foam::objectHit >, List< Foam::objectMap >, List< Foam::PackedBoolList >, List< Foam::Pair >, List< Foam::phaseProperties >, List< Foam::pointConstraint >, List< Foam::pointFieldDecomposer::patchFieldDecomposer *>, List< Foam::Reaction::specieCoeffs >, List< Foam::referredWallFace >, List< Foam::searchableSurface *>, List< Foam::SHA1Digest >, List< Foam::solutionControl::fieldData >, List< Foam::spatialTransform >, List< Foam::SpatialVector >, List< Foam::surfZone >, List< Foam::surfZoneIdentifier >, List< Foam::Tensor >, List< Foam::token >, List< Foam::topoDistanceData >, List< Foam::Tuple2< Foam::Vector, Foam::Vector > >, List< Foam::Tuple2< Foam::word, Foam::Tuple2< bool, scalar > > >, List< Foam::Tuple2< Foam::word, Foam::word > >, List< Foam::Tuple2< Foam::word, label > >, List< Foam::Tuple2< Foam::word, scalar > >, List< Foam::Tuple2< label, Foam::List< Foam::List > > >, List< Foam::Tuple2< label, scalar > >, List< Foam::Tuple2< mapType, Foam::List > >, List< Foam::Tuple2< scalar, Foam::fileName > >, List< Foam::Tuple2< scalar, scalar > >, List< Foam::Tuple2< scalar, Type > >, List< Foam::Tuple2< Type, scalar > >, List< Foam::Tuple2< Type, Type > >, List< Foam::Vector >, List< Foam::VectorSpace >, List< Foam::vectorTensorTransform >, List< Foam::volumeType >, List< Foam::word >, List< geometricSurfacePatch >, List< gradingDescriptor >, List< instant >, List< int >, List< kinematicParcelInjectionData >, List< label >, List< labelledTri >, List< labelList >, List< labelListList >, List< lduScheduleEntry >, List< List< pointIndexHit > >, List< List< scalar > >, List< mergeInfo >, List< patchInteractionData >, List< polyDecomp >, List< reactingMultiphaseParcelInjectionData >, List< reactingParcelInjectionData >, List< refineCell >, List< refineMode >, List< scalar >, List< scalarField >, List< scalarList >, List< scalarRange >, List< sideVolumeType >, List< simpleRegIOobject *>, List< string >, List< substance >, List< surfAndLabel >, List< T *>, List< treeBoundBox >, List< Tuple2< scalar, List< Tuple2< scalar, scalar > > > >, List< Tuple2< scalar, List< Tuple2< scalar, Type > > > >, List< Tuple2< scalar, scalar > >, List< Tuple2< scalar, Type > >, List< typename DSMCParcel< ParcelType > ::constantProperties >, List< typename ParcelType::constantProperties >, List< unsigned int >, List< vector >, List< vectorField >, List< volatileData >, List< vtkTextActor *>, List< word >, List< wordRe >, and SubList< point >.

Classes

class  greater
 Greater function class that can be used for sorting. More...
 
class  less
 Less function class that can be used for sorting. More...
 

Public Types

typedef T value_type
 Type of values the UList contains. More...
 
typedef Treference
 Type that can be used for storing into. More...
 
typedef const Tconst_reference
 Type that can be used for storing into. More...
 
typedef label difference_type
 The type that can represent the difference between any two. More...
 
typedef label size_type
 The type that can represent the size of a UList. More...
 
typedef Titerator
 Random access iterator for traversing UList. More...
 
typedef const Tconst_iterator
 Random access iterator for traversing UList. More...
 
typedef Treverse_iterator
 Reverse iterator for reverse traversal of UList. More...
 
typedef const Tconst_reverse_iterator
 Reverse iterator for reverse traversal of constant UList. More...
 

Public Member Functions

 UList ()
 Null constructor. More...
 
 UList (T *__restrict__ v, label size)
 Construct from components. More...
 
label fcIndex (const label i) const
 Return the forward circular index, i.e. the next index. More...
 
label rcIndex (const label i) const
 Return the reverse circular index, i.e. the previous index. More...
 
std::streamsize byteSize () const
 Return the binary size in number of characters of the UList. More...
 
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-1) More...
 
void checkSize (const label size) const
 Check size is within valid range (0 ... size) More...
 
void checkIndex (const label i) const
 Check index i is within valid range (0 ... size-1) More...
 
void shallowCopy (const UList< T > &)
 Copy the pointer held by the given UList. More...
 
void deepCopy (const UList< T > &)
 Copy elements of the given UList. More...
 
void writeEntry (Ostream &) const
 Write the UList as a dictionary entry. More...
 
void writeEntry (const word &keyword, Ostream &) const
 Write the UList as a dictionary entry with keyword. More...
 
Toperator[] (const label)
 Return element of UList. More...
 
const Toperator[] (const label) const
 Return element of constant UList. More...
 
 operator const Foam::List< T > & () const
 Allow cast to a const List<T>&. More...
 
void operator= (const T &)
 Assignment of all entries to the given value. More...
 
void operator= (const 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
 Return the number of elements in the UList. More...
 
label max_size () const
 Return size of the largest possible UList. More...
 
bool empty () const
 Return true if the UList is empty (ie, size() is zero) More...
 
void swap (UList< T > &)
 Swap two ULists of the same type in constant time. More...
 
bool operator== (const UList< T > &) const
 Equality operation on ULists of the same type. More...
 
bool operator!= (const UList< T > &) const
 The opposite of the equality operation. Takes linear time. More...
 
bool operator< (const UList< T > &) const
 Compare two ULists lexicographically. Takes linear time. More...
 
bool operator> (const UList< T > &) const
 Compare two ULists lexicographically. Takes linear time. More...
 
bool operator<= (const UList< T > &) const
 Return true if !(a > b). Takes linear time. More...
 
bool operator>= (const UList< T > &) const
 Return true if !(a < b). Takes linear time. More...
 
template<>
const bool & operator[] (const label i) const
 

Static Public Member Functions

static const UList< T > & null ()
 Return a null UList. More...
 

Friends

class List< T >
 Declare friendship with the List class. More...
 
class SubList< T >
 Declare friendship with the SubList class. More...
 
Ostreamoperator (Ostream &, const UList< T > &)
 
Istreamoperator>> (Istream &, UList< T > &)
 Read UList contents from Istream. Requires size to have been set. 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 61 of file HashTable.H.

Member Typedef Documentation

◆ value_type

typedef T value_type

Type of values the UList contains.

Definition at line 251 of file UList.H.

◆ reference

typedef T& reference

Type that can be used for storing into.

UList::value_type objects

Definition at line 255 of file UList.H.

◆ const_reference

typedef const T& const_reference

Type that can be used for storing into.

constant UList::value_type objects

Definition at line 259 of file UList.H.

◆ difference_type

The type that can represent the difference between any two.

UList iterator objects

Definition at line 263 of file UList.H.

◆ size_type

typedef label size_type

The type that can represent the size of a UList.

Definition at line 266 of file UList.H.

◆ iterator

typedef T* iterator

Random access iterator for traversing UList.

Definition at line 272 of file UList.H.

◆ const_iterator

typedef const T* const_iterator

Random access iterator for traversing UList.

Definition at line 284 of file UList.H.

◆ reverse_iterator

typedef T* reverse_iterator

Reverse iterator for reverse traversal of UList.

Definition at line 302 of file UList.H.

◆ const_reverse_iterator

typedef const T* const_reverse_iterator

Reverse iterator for reverse traversal of constant UList.

Definition at line 314 of file UList.H.

Constructor & Destructor Documentation

◆ UList() [1/2]

UList ( )
inline

Null constructor.

Definition at line 33 of file UListI.H.

Referenced by UList< T >::greater::operator()().

Here is the caller graph for this function:

◆ UList() [2/2]

UList ( T *__restrict__  v,
label  size 
)
inline

Construct from components.

Definition at line 41 of file UListI.H.

Member Function Documentation

◆ null()

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

◆ fcIndex()

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

Return the forward circular index, i.e. the next index.

which returns to the first at the end of the list

Definition at line 58 of file UListI.H.

Referenced by pointMVCWeight::calcWeights(), tetOverlapVolume::cellCellOverlapMinDecomp(), tetOverlapVolume::cellCellOverlapVolumeMinDecomp(), primitiveMesh::checkCommonOrder(), primitiveMesh::checkEdgeLength(), primitiveMeshGeometry::checkFaceAngles(), polyMeshGeometry::checkFaceAngles(), edgeCollapser::checkMeshQuality(), polyMeshGeometry::checkTriangleTwist(), face::edgeDirection(), primitiveMeshTools::faceConcavity(), face::faceEdge(), primitiveMesh::faceEdges(), hexRef8::faceLevel(), mappedPatchBase::facePoint(), tetIndices::faceTriIs(), searchableSurfacesQueries::findNearest(), isoSurface::generateFaceTriPoints(), edgeFaceCirculator::getMinIndex(), combineFaces::getOutsideFace(), face::inertia(), face::nextLabel(), minEqOpFace::operator()(), listPlusEqOp< T >::operator()(), ifEqEqOp< value >::operator()(), UList< T >::greater::operator()(), triSurfaceTools::otherEdges(), treeDataFace::overlaps(), OBJsurfaceFormat< Face >::read(), VTKsurfaceFormat< Face >::read(), OFFsurfaceFormat< Face >::read(), AC3DsurfaceFormat< Face >::read(), snappySnapDriver::repatchToSurface(), faceCollapser::setRefinement(), createShellMesh::setRefinement(), tetDecomposer::setRefinement(), addPatchCellLayer::setRefinement(), meshRefinement::splitFaces(), syncTools::syncEdgeMap(), MeshedSurface< Foam::face >::triangulate(), triSurfaceTools::triangulateFaceCentre(), patchEdgeFaceRegions::updateEdge(), and patchEdgeFaceRegions::updateFace().

Here is the caller graph for this function:

◆ rcIndex()

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

◆ 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. contiguous<T>() == true. Note that is of type streamsize since used in stream ops

Definition at line 100 of file UList.C.

Referenced by mapDistributeBase::distribute(), globalIndex::gather(), Pstream::gatherList(), List< Field< scalar > >::List(), Pstream::listCombineGather(), Pstream::listCombineScatter(), UList< T >::greater::operator()(), processorLduInterface::receive(), globalIndex::scatter(), processorLduInterface::send(), decomposedBlockData::writeBlocks(), and decomposedBlockData::writeData().

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 142 of file UListI.H.

Referenced by UList< T >::greater::operator()(), and PackedList< 2 >::write().

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 149 of file UListI.H.

Referenced by patchDataWave< TransferType >::correct(), Foam::cwd(), decompositionMethod::decompose(), List< Field< scalar > >::List(), UList< T >::greater::operator()(), PackedList< 2 >::read(), and UList< Foam::wordRe >::writeEntry().

Here is the caller graph for this function:

◆ first() [1/2]

T & first ( )
inline

Return the first element of the list.

Definition at line 114 of file UListI.H.

Referenced by preserveBafflesConstraint::add(), bladeModel::bladeModel(), TableBase< Type >::check(), primitiveMeshGeometry::checkFaceAngles(), polyMeshGeometry::checkFaceAngles(), polyMeshGeometry::checkFaceDotProduct(), polyMeshGeometry::checkFacePyramids(), polyMeshGeometry::checkFaceSkewness(), polyMeshGeometry::checkFaceTets(), polyMeshGeometry::checkFaceWeights(), TableBase< Type >::checkMaxBounds(), edgeCollapser::checkMeshQuality(), TableBase< Type >::checkMinBounds(), interpolation2DTable< scalar >::checkOrder(), polyMeshGeometry::checkVolRatio(), TableBase< Type >::convertTimeBase(), meshRefinement::createZoneBaffles(), Distribution< Type >::cumulativeNormalised(), Distribution< Type >::cumulativeRaw(), displacementInterpolationMotionSolver::displacementInterpolationMotionSolver(), primitiveMeshTools::faceConcavity(), patchProbes::findElements(), mappedPatchBase::findSamples(), TableBase< Type >::interpolator(), cell::labels(), solutionControl::maxTypeResidual(), distribution::median(), Distribution< Type >::median(), meshRefinement::mergeBaffles(), distribution::normalised(), Distribution< Type >::normalised(), distribution::normalisedShifted(), procLess::operator()(), listPlusEqOp< T >::operator()(), UList< T >::greater::operator()(), interpolation2DTable< scalar >::operator()(), Distribution< Type >::operator=(), ParticleCollector< CloudType >::ParticleCollector(), SprayCloud< Foam::DSMCCloud >::penetration(), BSpline::position(), CatmullRomSpline::position(), polyLine::position(), distribution::raw(), Distribution< Type >::raw(), tabulated6DoFAcceleration::read(), tabulated6DoFMotion::read(), ITstream::rewind(), primitiveEntry::startLineNumber(), Distribution< Type >::write(), pairPotential::writeEnergyAndForceTables(), residuals::writeResidual(), TableBase< Type >::x(), and RaviPetersen::~RaviPetersen().

◆ first() [2/2]

const T & first ( ) const
inline

Return first element of the list.

Definition at line 121 of file UListI.H.

◆ last() [1/2]

T & last ( )
inline

Return the last element of the list.

Definition at line 128 of file UListI.H.

Referenced by masterCoarsestGAMGProcAgglomeration::agglomerate(), eagerGAMGProcAgglomeration::agglomerate(), procFacesGAMGProcAgglomeration::agglomerate(), manualGAMGProcAgglomeration::agglomerate(), polyLine::calcParam(), triSurface::checkEdges(), primitiveMeshGeometry::checkFaceAngles(), polyMeshGeometry::checkFaceAngles(), TableBase< Type >::checkMaxBounds(), edgeCollapser::checkMeshQuality(), TableBase< Type >::checkMinBounds(), MeshedSurface< Foam::face >::checkZones(), surfMesh::checkZones(), processorLduInterface::compressedReceive(), processorLduInterface::compressedSend(), decompositionMethod::decompositionMethod(), face::edges(), edgeSurface::edgeSurface(), primitiveEntry::endLineNumber(), primitiveMeshTools::faceConcavity(), Time::findClosestTime(), searchableBox::findLineAll(), searchableExtrudedCircle::findParametricNearest(), globalIndex::gather(), FreeStream< CloudType >::inflow(), distribution::insertMissingKeys(), linearInterpolationWeights::integrationWeights(), lduPrimitiveMesh::lduPrimitiveMesh(), solutionControl::maxTypeResidual(), procLess::operator()(), UList< T >::greater::operator()(), Foam::operator<<(), SprayCloud< Foam::DSMCCloud >::penetration(), BSpline::position(), CatmullRomSpline::position(), projectEdge::position(), projectCurveEdge::position(), polyLine::position(), GAMGAgglomeration::procAgglomerateRestrictAddressing(), ITstream::read(), regIOobject::readHeaderOk(), regIOobject::readIfModified(), Time::setControls(), surfaceFeatures::setFromStatus(), hexRef8::setInstance(), Reaction< ReactionThermo >::setLRhs(), patchInjectionBase::setPositionAndCell(), createShellMesh::setRefinement(), addPatchCellLayer::setRefinement(), globalIndex::size(), TableBase< Type >::value(), and RaviPetersen::~RaviPetersen().

◆ last() [2/2]

const T & last ( ) const
inline

Return the last element of the list.

Definition at line 135 of file UListI.H.

◆ checkStart()

void checkStart ( const label  start) const
inline

Check start is within valid range (0 ... size-1)

Definition at line 72 of file UListI.H.

Referenced by UList< T >::greater::operator()(), and SubList< point >::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 84 of file UListI.H.

Referenced by UList< T >::greater::operator()(), and SubList< point >::SubList().

Here is the caller graph for this function:

◆ checkIndex()

void checkIndex ( const label  i) const
inline

Check index i is within valid range (0 ... size-1)

Definition at line 96 of file UListI.H.

Referenced by UList< T >::greater::operator()().

Here is the caller graph for this function:

◆ shallowCopy()

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

Copy the pointer held by the given UList.

Definition at line 156 of file UListI.H.

Referenced by UList< T >::greater::operator()(), and ODESolver::resizeField().

Here is the caller graph for this function:

◆ deepCopy()

void deepCopy ( const UList< T > &  a)

Copy elements of the given UList.

Definition at line 35 of file UList.C.

Referenced by UPstream::allToAll(), UList< T >::greater::operator()(), and globalIndex::scatter().

Here is the caller graph for this function:

◆ writeEntry() [1/2]

void writeEntry ( Ostream os) const

◆ writeEntry() [2/2]

void writeEntry ( const word keyword,
Ostream os 
) const

Write the UList as a dictionary entry with keyword.

Definition at line 54 of file UListIO.C.

◆ operator[]() [1/3]

T & operator[] ( const label  i)
inline

◆ operator[]() [2/3]

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 196 of file UListI.H.

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

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

Allow cast to a const List<T>&.

Definition at line 206 of file UListI.H.

◆ operator=() [1/2]

void operator= ( const T t)

Assignment of all entries to the given value.

Definition at line 68 of file UList.C.

◆ operator=() [2/2]

void operator= ( const zero  )

Assignment of all entries to zero.

Definition at line 78 of file UList.C.

◆ begin() [1/2]

Foam::UList< T >::iterator begin ( )
inline

Return an iterator to begin traversing the UList.

Definition at line 216 of file UListI.H.

Referenced by lduMatrix::Amul(), TDILUPreconditioner< Type, DType, LUType >::calcInvD(), DICPreconditioner::calcReciprocalD(), DILUPreconditioner::calcReciprocalD(), mapDistributeBase::compact(), processorLduInterface::compressedReceive(), processorLduInterface::compressedSend(), triSurfaceTools::delaunay2D(), mapDistributeBase::distribute(), globalIndex::gather(), Pstream::gatherList(), lduMatrix::H(), lduMatrix::H1(), Pstream::listCombineGather(), Pstream::listCombineScatter(), LUscalarMatrix::LUscalarMatrix(), UList< Foam::wordRe >::operator<(), SortableListDRGEP< Type >::partialSort(), SortableListEFA< Type >::partialSort(), noPreconditioner::precondition(), DICPreconditioner::precondition(), DILUPreconditioner::precondition(), diagonalPreconditioner::precondition(), FDICPreconditioner::precondition(), DILUPreconditioner::preconditionT(), ITstream::print(), masterUncollatedFileOperation::readAndSend(), decomposedBlockData::readMasterHeader(), processorLduInterface::receive(), lduMatrix::residual(), ODESolver::resizeField(), globalIndex::scatter(), processorLduInterface::send(), Foam::shuffle(), GaussSeidelSmoother::smooth(), symGaussSeidelSmoother::smooth(), TGaussSeidelSmoother< Type, DType, LUType >::smooth(), DILUSmoother::smooth(), DICSmoother::smooth(), FDICSmoother::smooth(), nonBlockingGaussSeidelSmoother::smooth(), PBiCICG< Type, DType, LUType >::solve(), PBiCCCG< Type, DType, LUType >::solve(), PCICG< Type, DType, LUType >::solve(), PBiCG::solve(), PCG::solve(), PBiCGStab::solve(), Foam::sort(), Foam::stableSort(), lduMatrix::sumA(), TGaussSeidelSmoother< Type, DType, LUType >::TGaussSeidelSmoother(), lduMatrix::Tmul(), fft::transform(), ensightBinaryStream::write(), decomposedBlockData::writeBlocks(), masterOFstream::~masterOFstream(), and UOPstream::~UOPstream().

◆ end() [1/2]

Foam::UList< T >::iterator end ( )
inline

◆ cbegin()

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

Return const_iterator to begin traversing the constant UList.

Definition at line 230 of file UListI.H.

Referenced by decomposedBlockData::writeData().

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 251 of file UListI.H.

◆ begin() [2/2]

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

Return const_iterator to begin traversing the constant UList.

Definition at line 223 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 244 of file UListI.H.

◆ rbegin() [1/2]

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

Return reverse_iterator to begin reverse traversing the UList.

Definition at line 258 of file UListI.H.

Referenced by ITstream::print().

Here is the caller graph for this function:

◆ rend() [1/2]

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

Return reverse_iterator to end reverse traversing the UList.

Definition at line 279 of file UListI.H.

◆ crbegin()

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

Return const_reverse_iterator to begin reverse traversing the UList.

Definition at line 272 of file UListI.H.

◆ crend()

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

Return const_reverse_iterator to end reverse traversing the UList.

Definition at line 293 of file UListI.H.

◆ rbegin() [2/2]

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

Return const_reverse_iterator to begin reverse traversing the UList.

Definition at line 265 of file UListI.H.

◆ rend() [2/2]

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

Return const_reverse_iterator to end reverse traversing the UList.

Definition at line 286 of file UListI.H.

◆ size()

Foam::label size ( ) const
inline

Return the number of elements in the UList.

Definition at line 299 of file UListI.H.

Referenced by mapDistributeBase::accessAndFlip(), isoSurface::adaptPatchFields(), fvMatrix< Type >::addToInternalField(), MeshedSurface< Foam::face >::addZones(), pairGAMGAgglomeration::agglomerate(), GAMGAgglomeration::agglomerateLduAddressing(), DynamicList< Foam::Tensor >::append(), DynamicField< T, SizeInc, SizeMult, SizeDiv >::append(), List< Field< scalar > >::append(), Field< Foam::Vector2D >::autoMap(), fvPatchField< Type >::autoMap(), Foam::average(), coupledPolyPatch::calcFaceTol(), edgeCollapser::checkMeshQuality(), GAMGAgglomeration::checkRestriction(), Foam::cmptAv(), Foam::cmptMag(), Foam::cmptMax(), Foam::cmptMin(), meshRefinement::collectAndPrint(), ORourkeCollision< CloudType >::collide(), GAMGAgglomeration::compactLevels(), Foam::ComplexField(), processorLduInterface::compressedReceive(), processorLduInterface::compressedSend(), Keyed< T >::createList(), directFvPatchFieldMapper::directFvPatchFieldMapper(), directPointPatchFieldMapper::directPointPatchFieldMapper(), Foam::duplicateOrder(), HashTable< Foam::phase *>::erase(), Pstream::exchange(), regionSizeDistribution::extractData(), FaceCellWave< Type, TrackingData >::FaceCellWave(), lduMatrix::faceH(), LduMatrix< Type, DType, LUType >::faceH(), Foam::facePointN(), FixedList< Type, 3 >::FixedList(), mapDistributeBase::flipAndCombine(), globalIndex::gather(), Foam::gAverage(), meshRefinement::gAverage(), coupledPolyPatch::getAnchorPoints(), Foam::Im(), Foam::ImComplexField(), Foam::inplaceSubset(), cyclicACMIGAMGInterface::internalFieldTransfer(), cyclicAMIGAMGInterface::internalFieldTransfer(), regionCoupledFvPatch::internalFieldTransfer(), regionCoupledWallFvPatch::internalFieldTransfer(), surfaceInterpolationScheme< GType >::interpolate(), AMIInterpolation< SourcePatch, TargetPatch >::interpolateToSource(), AMIInterpolation< SourcePatch, TargetPatch >::interpolateToTarget(), Foam::inv(), cell::labels(), lduPrimitiveMesh::lduPrimitiveMesh(), LduMatrix< Type, DType, LUType >::lower(), Foam::mag(), Foam::magSqr(), Field< Foam::Vector2D >::map(), Foam::matchPoints(), Foam::max(), Foam::maxMagSqr(), Foam::mergePoints(), MeshedSurface< Foam::face >::MeshedSurface(), Foam::min(), Foam::minMagSqr(), LduMatrix< Type, DType, LUType >::negSumDiag(), lduMatrix::negSumDiag(), procLess::operator()(), lessProcPatches::operator()(), UList< T >::greater::operator()(), Foam::operator<<(), UIndirectList< T >::operator=(), BiIndirectList< T >::operator=(), PackedBoolList::operator=(), DynamicList< Foam::Tensor >::operator=(), FixedList< Type, 3 >::operator=(), PackedList< 2 >::operator=(), directFieldMapper::patchFieldSubset(), Polynomial< 8 >::Polynomial(), Foam::pow(), Foam::Re(), OBJsurfaceFormat< Face >::read(), OFFsurfaceFormat< Face >::read(), Foam::ReComplexField(), fvFieldReconstructor::reconstructFvSurfaceField(), fvFieldReconstructor::reconstructFvVolumeField(), Foam::ReImSum(), UnsortedMeshedSurface< Face >::remapFaces(), cuttingPlane::remapFaces(), MeshedSurface< Foam::face >::remapFaces(), ensightPart::renumber(), PtrList< transferModel >::reorder(), UPtrList< Foam::Field< Type > >::reorder(), Foam::reverse(), CompactListList< T, Container >::setSize(), UnsortedMeshedSurface< Face >::setZones(), directFieldMapper::size(), directFvPatchFieldMapper::size(), directPointPatchFieldMapper::size(), Foam::sortedOrder(), Foam::sqr(), Foam::stabilise(), Foam::subset(), fvMatrix< Type >::subtractFromInternalField(), Foam::sum(), Foam::sumCmptMag(), Foam::sumCmptProd(), LduMatrix< Type, DType, LUType >::sumDiag(), lduMatrix::sumDiag(), Foam::sumMag(), LduMatrix< Type, DType, LUType >::sumMagOffDiag(), lduMatrix::sumMagOffDiag(), Foam::sumProd(), Foam::sumSqr(), syncTools::swapBoundaryCellList(), syncTools::swapBoundaryCellPositions(), syncTools::syncBoundaryFaceList(), meshRefinement::testSyncBoundaryFaceList(), Foam::transformList(), UnsortedMeshedSurface< Face >::UnsortedMeshedSurface(), sampledPatch::update(), LduMatrix< Type, DType, LUType >::upper(), decomposedBlockData::writeBlocks(), ensightPartFaces::writeConnectivity(), VTKedgeFormat::writeEdges(), UList< Foam::wordRe >::writeEntry(), ensightPart::writeField(), WRLsurfaceFormatCore::writeHeader(), OFSsurfaceFormatCore::writeHeader(), AC3DsurfaceFormatCore::writeHeader(), coupledPolyPatch::writeOBJ(), Foam::meshTools::writeOBJ(), ensightPart::writeScalarField(), VTKsurfaceFormatCore::writeTail(), ensightPart::writeVectorField(), and ensightPartCells::~ensightPartCells().

◆ max_size()

Foam::label max_size ( ) const
inline

Return size of the largest possible UList.

Definition at line 306 of file UListI.H.

◆ empty()

bool empty ( ) const
inline

Return true if the UList is empty (ie, size() is zero)

Definition at line 313 of file UListI.H.

Referenced by DynamicID< pointZoneMesh >::active(), labelRanges::add(), dynamicIndexedOctree< Foam::dynamicTreeDataPoint >::bb(), indexedOctree< Foam::treeDataFace >::bb(), boundBox::boundBox(), cellVolumeWeightMethod::calculateAddressing(), cellMapper::cellMapper(), blockMesh::cells(), triSurface::checkEdges(), edgeCollapser::checkMeshQuality(), primitiveMesh::checkPoints(), labelRanges::const_iterator::const_iterator(), boundBox::contains(), boundBox::containsAny(), dynamicCode::copyOrCreateFiles(), hexCellLooper::cut(), fvSchemes::d2dt2Scheme(), fvSchemes::ddtScheme(), fvSchemes::divScheme(), edgeSurface::edgeSurface(), primitiveEntry::endLineNumber(), faceMapper::faceMapper(), searchableSurfacesQueries::findAllIntersections(), Foam::findEtcFiles(), coordinateSystems::findIndex(), ZoneMesh< cellZone, polyMesh >::findIndex(), polyBoundaryMesh::findIndex(), triSurfaceRegionSearch::findNearest(), searchableSurfacesQueries::findNearest(), fvSchemes::gradScheme(), Distribution< Type >::index(), DynamicID< pointZoneMesh >::index(), fvSchemes::interpolationScheme(), intersectedSurface::intersectedSurface(), fvSchemes::laplacianScheme(), meshTriangulation::meshTriangulation(), cellCuts::nonAnchorPoints(), Distribution< Type >::normalised(), dlLibraryTable::open(), listPlusEqOp< T >::operator()(), cellToCellStencil::unionEqOp::operator()(), cellToFaceStencil::unionEqOp::operator()(), Foam::operator<<(), triSurfaceMesh::overlaps(), coupledPolyPatch::parallel(), ParticleErosion< CloudType >::ParticleErosion(), blockMesh::patches(), patchInteractionDataList::patchInteractionDataList(), boundaryMesh::patchNames(), PatchPostProcessing< CloudType >::PatchPostProcessing(), polyBoundaryMesh::patchSet(), pointMapper::pointMapper(), blockMesh::points(), polynomialFunction::polynomialFunction(), Distribution< Type >::raw(), boundaryMesh::readTriSurface(), triSurfaceMesh::regions(), labelRanges::remove(), timeSelector::select0(), Reaction< ReactionThermo >::setLRhs(), combineFaces::setUnrefinement(), fvSchemes::snGradScheme(), primitiveEntry::startLineNumber(), surfaceIntersection::surfaceIntersection(), treeBoundBox::treeBoundBox(), sampledPlane::update(), walkPatch::walkPatch(), and writeObjectsBase::write().

◆ swap()

void swap ( UList< T > &  a)

Swap two ULists of the same type in constant time.

Definition at line 90 of file UList.C.

◆ 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 152 of file UList.C.

◆ operator!=()

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

The opposite of the equality operation. Takes linear time.

Definition at line 173 of file UList.C.

◆ operator<()

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

Compare two ULists lexicographically. Takes linear time.

Definition at line 180 of file UList.C.

◆ operator>()

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

Compare two ULists lexicographically. Takes linear time.

Definition at line 211 of file UList.C.

◆ operator<=()

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

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

Definition at line 218 of file UList.C.

◆ operator>=()

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

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

Definition at line 225 of file UList.C.

◆ operator[]() [3/3]

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

Definition at line 180 of file UListI.H.

Friends And Related Function Documentation

◆ List< T >

friend class List< T >
friend

Declare friendship with the List class.

Definition at line 100 of file UList.H.

◆ SubList< T >

friend class SubList< T >
friend

Declare friendship with the SubList class.

Definition at line 103 of file UList.H.

◆ operator

Ostream& operator ( Ostream ,
const UList< T > &   
)
friend

◆ operator>>

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

Read UList contents from Istream. Requires size to have been set.

before


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