MOAB: Mesh Oriented datABase  (version 5.4.1)
moab::TupleList Class Reference

#include <TupleList.hpp>

Classes

class  buffer
struct  SortData

Public Member Functions

 TupleList (uint mi, uint ml, uint mul, uint mr, uint max)
 TupleList ()
 ~TupleList ()
void initialize (uint mi, uint ml, uint mul, uint mr, uint max)
ErrorCode resize (uint max)
ErrorCode sort (uint key, TupleList::buffer *buf)
void reset ()
void reserve ()
int find (unsigned int key_num, sint value)
int find (unsigned int key_num, slong value)
int find (unsigned int key_num, Ulong value)
int find (unsigned int key_num, realType value)
sint get_sint (unsigned int index, unsigned int m)
slong get_int (unsigned int index, unsigned int m)
Ulong get_ulong (unsigned int index, unsigned int m)
realType get_double (unsigned int index, unsigned int m)
ErrorCode get (unsigned int index, const sint *&sp, const slong *&ip, const Ulong *&lp, const realType *&dp)
unsigned int push_back (sint *sp, slong *ip, Ulong *lp, realType *dp)
void enableWriteAccess ()
void disableWriteAccess ()
void getTupleSize (uint &mi_out, uint &ml_out, uint &mul_out, uint &mr_out) const
void set_n (uint n_in)
uint get_n () const
uint get_max () const
bool get_writeEnabled () const
uint inc_n ()
void print (const char *) const
void print_to_file (const char *) const

Public Attributes

sintvi_wr
slongvl_wr
Ulongvul_wr
realTypevr_wr
const sintvi_rd
slongvl_rd
Ulongvul_rd
realTypevr_rd

Private Types

typedef uint Index

Private Member Functions

void permute (uint *perm, void *work)

Static Private Member Functions

template<class Value >
static Value radix_count (const Value *A, const Value *end, Index stride, Index count[DIGITS][DIGIT_VALUES])
static void radix_offsets (Index *c)
template<class Value >
static unsigned radix_zeros (Value bitorkey, Index count[DIGITS][DIGIT_VALUES], unsigned *shift, Index **offsets)
template<class Value >
static void radix_index_pass_b (const Value *A, Index n, Index stride, unsigned sh, Index *off, SortData< Value > *out)
template<class Value >
static void radix_index_pass_m (const SortData< Value > *src, const SortData< Value > *end, unsigned sh, Index *off, SortData< Value > *out)
template<class Value >
static void radix_index_pass_e (const SortData< Value > *src, const SortData< Value > *end, unsigned sh, Index *off, Index *out)
template<class Value >
static void radix_index_pass_be (const Value *A, Index n, Index stride, unsigned sh, Index *off, Index *out)
template<class Value >
static void radix_index_sort (const Value *A, Index n, Index stride, Index *idx, SortData< Value > *work)
template<class Value >
static void merge_index_sort (const Value *A, const Index An, Index stride, Index *idx, SortData< Value > *work)
template<class Value >
static void index_sort (const Value *A, Index n, Index stride, Index *idx, SortData< Value > *work)

Private Attributes

uint mi
uint ml
uint mul
uint mr
uint n
uint max
sintvi
slongvl
Ulongvul
realTypevr
int last_sorted
bool writeEnabled

Detailed Description

Examples:
CrystalRouterExample.cpp.

Definition at line 78 of file TupleList.hpp.


Member Typedef Documentation

typedef uint moab::TupleList::Index [private]

Definition at line 314 of file TupleList.hpp.


Constructor & Destructor Documentation

moab::TupleList::TupleList ( uint  mi,
uint  ml,
uint  mul,
uint  mr,
uint  max 
)

Constructor that takes all parameters and initializes the TupleList

Definition at line 65 of file TupleList.cpp.

References initialize().

    : vi( NULL ), vl( NULL ), vul( NULL ), vr( NULL ), last_sorted( -1 )
{
    initialize( p_mi, p_ml, p_mul, p_mr, p_max );
}

Default constructor (Note: TupleList must be initialized before use!)

Definition at line 71 of file TupleList.cpp.

References disableWriteAccess().

    : vi_rd( NULL ), vl_rd( NULL ), vul_rd( NULL ), vr_rd( NULL ), mi( 0 ), ml( 0 ), mul( 0 ), mr( 0 ), n( 0 ),
      max( 0 ), vi( NULL ), vl( NULL ), vul( NULL ), vr( NULL ), last_sorted( -1 )
{
    disableWriteAccess();
}

Definition at line 135 of file TupleList.hpp.

References reset().

    {
        reset();
    };

Member Function Documentation

int moab::TupleList::find ( unsigned int  key_num,
sint  value 
)

Finds index of the tuple containing 'value' at the key_numth index of said tuple; return -1 if key_num is out of bounds or if 'value' not found Uses binary search if TupleList is sorted by the key_numth field, seqential otherwise (very slow for large TupleLists; please sort before you search)

param key_num index of the tuple where to search for value param value value to search for at the given key_num return the index of the tuple that contains value

Definition at line 236 of file TupleList.cpp.

References last_sorted, mi, n, and vi.

{
    // we are passing an int, no issue, leave it at long
    long uvalue = (long)value;
    if( !( key_num > mi ) )
    {
        // Binary search: only if the tuple_list is sorted
        if( last_sorted == (int)key_num )
        {
            int lb = 0, ub = n, index;  // lb=lower bound, ub=upper bound, index=mid
            for( ; lb <= ub; )
            {
                index = ( lb + ub ) / 2;
                if( vi[index * mi + key_num] == uvalue )
                    return index;
                else if( vi[index * mi + key_num] > uvalue )
                    ub = index - 1;
                else if( vi[index * mi + key_num] < uvalue )
                    lb = index + 1;
            }
        }
        else
        {
            // Sequential search: if tuple_list is not sorted
            for( uint index = 0; index < n; index++ )
            {
                if( vi[index * mi + key_num] == uvalue ) return index;
            }
        }
    }
    return -1;  // If the value wasn't present or an invalid key was given
}
int moab::TupleList::find ( unsigned int  key_num,
slong  value 
)

Definition at line 269 of file TupleList.cpp.

References last_sorted, mi, ml, n, and vl.

{
    long uvalue = (long)value;
    if( !( key_num > ml ) )
    {
        if( last_sorted - mi == key_num )
        {
            int lb = 0, ub = n, index;  // lb=lower bound, ub=upper bound, index=mid
            for( ; lb <= ub; )
            {
                index = ( lb + ub ) / 2;
                if( vl[index * ml + key_num] == uvalue )
                    return index;
                else if( vl[index * ml + key_num] > uvalue )
                    ub = index - 1;
                else if( vl[index * ml + key_num] < uvalue )
                    lb = index + 1;
            }
        }
        else
        {
            // Sequential search: if tuple_list is not sorted
            for( uint index = 0; index < n; index++ )
            {
                if( vl[index * ml + key_num] == uvalue ) return index;
            }
        }
    }
    return -1;  // If the value wasn't present or an invalid key was given
}
int moab::TupleList::find ( unsigned int  key_num,
Ulong  value 
)

Definition at line 300 of file TupleList.cpp.

References last_sorted, mi, ml, mul, n, and vul.

{
    if( !( key_num > mul ) )
    {
        if( last_sorted - mi - ml == key_num )
        {
            int lb = 0, ub = n - 1, index;  // lb=lower bound, ub=upper bound, index=mid
            for( ; lb <= ub; )
            {
                index = ( lb + ub ) / 2;
                if( vul[index * mul + key_num] == value )
                    return index;
                else if( vul[index * mul + key_num] > value )
                    ub = index - 1;
                else if( vul[index * mul + key_num] < value )
                    lb = index + 1;
            }
        }
        else
        {
            // Sequential search: if tuple_list is not sorted
            for( uint index = 0; index < n; index++ )
            {
                if( vul[index * mul + key_num] == value ) return index;
            }
        }
    }
    return -1;  // If the value wasn't present or an invalid key was given
}
int moab::TupleList::find ( unsigned int  key_num,
realType  value 
)

Definition at line 330 of file TupleList.cpp.

References mr, n, and vr.

{
    if( !( key_num > mr ) )
    {
        // Sequential search: TupleList cannot be sorted by reals
        for( uint index = 0; index < n; index++ )
        {
            if( vr[index * mr + key_num] == value ) return index;
        }
    }
    return -1;  // If the value wasn't present or an invalid key was given
}
ErrorCode moab::TupleList::get ( unsigned int  index,
const sint *&  sp,
const slong *&  ip,
const Ulong *&  lp,
const realType *&  dp 
)

get pointers to the data for the index'th tuple; ptr is NULL if that type is not part of this tuple

param index index of the tuple needed param *&sp, *&ip, *&lp, *&dp pointers to each piece of the tuple

Definition at line 367 of file TupleList.cpp.

References MB_SUCCESS, mi, ml, mr, mul, n, vi, vl, vr, and vul.

{
    if( index <= n )
    {
        if( mi )
            *&sp = &vi[index * mi];
        else
            *&sp = NULL;
        if( ml )
            *&ip = &vl[index * ml];
        else
            *&ip = NULL;
        if( mul )
            *&lp = &vul[index * mul];
        else
            *&lp = NULL;
        if( mr )
            *&dp = &vr[index * mr];
        else
            *&dp = NULL;

        return MB_SUCCESS;
    }
    return MB_FAILURE;
}
realType moab::TupleList::get_double ( unsigned int  index,
unsigned int  m 
)

Definition at line 361 of file TupleList.cpp.

References mr, n, and vr.

{
    if( mr > m && n > index ) return vr[index * mr + m];
    return 0;
}
slong moab::TupleList::get_int ( unsigned int  index,
unsigned int  m 
)

Definition at line 349 of file TupleList.cpp.

References ml, n, and vl.

{
    if( ml > m && n > index ) return vl[index * ml + m];
    return 0;
}
sint moab::TupleList::get_sint ( unsigned int  index,
unsigned int  m 
)

get the mth number of return type in the index'th tuple returns 0 if m or index is out of bounds

param index index of the tuple within the TupleList param m index of the value within the tuple return the value at the given position

Definition at line 343 of file TupleList.cpp.

References mi, n, and vi.

{
    if( mi > m && n > index ) return vi[index * mi + m];
    return 0;
}
Ulong moab::TupleList::get_ulong ( unsigned int  index,
unsigned int  m 
)

Definition at line 355 of file TupleList.cpp.

References mul, n, and vul.

{
    if( mul > m && n > index ) return vul[index * mul + m];
    return 0;
}
template<class Value >
void moab::TupleList::index_sort ( const Value *  A,
Index  n,
Index  stride,
Index idx,
SortData< Value > *  work 
) [static, private]

Definition at line 843 of file TupleList.cpp.

References DIGIT_VALUES, merge_index_sort(), and radix_index_sort().

Referenced by sort().

{
    if( n < DIGIT_VALUES )
    {
        if( n == 0 ) return;
        if( n == 1 )
            *idx = 0;
        else
            merge_index_sort( A, n, stride, idx, work );
    }
    else
        radix_index_sort( A, n, stride, idx, work );
}
void moab::TupleList::initialize ( uint  mi,
uint  ml,
uint  mul,
uint  mr,
uint  max 
)

Initializes the starting memory to be used by the TupleList Note: TupleLists must be initialized before they can be used

param mi number of signed ints in each tuple param ml number of long ints in each tuple param mul number of unsigned long ints in each tuple param mr number of reals in each tuple param max starting capacity of max tuples in the TupleList

Examples:
CrystalRouterExample.cpp.

Definition at line 79 of file TupleList.cpp.

References disableWriteAccess(), moab::fail(), max, mi, ml, mr, mul, n, sint, slong, vi, vi_rd, vl, vl_rd, vr, vr_rd, vul, and vul_rd.

Referenced by moab::ParallelComm::augment_default_sets_with_ghosts(), moab::ParCommGraph::compute_partition(), moab::ParallelComm::delete_entities(), moab::ParallelComm::exchange_owned_mesh(), moab::ParCommGraph::form_tuples_to_migrate_mesh(), moab::ParallelComm::get_sent_ents(), iMOAB_SetDoubleTagStorageWithGid(), moab::gs_data::initialize(), moab::DataCoupler::interpolate(), moab::Coupler::interpolate(), moab::Coupler::locate_points(), moab::SpatialLocator::locate_points(), main(), moab::ParallelMergeMesh::PerformMerge(), moab::TempestOnlineMap::ReadParallelMap(), moab::ParallelComm::resolve_shared_ents(), moab::ParallelComm::send_entities(), moab::ScdInterface::tag_shared_vertices(), TupleList(), moab::unpack_tuples(), and moab::TempestOnlineMap::WriteSCRIPMapFile().

{
    this->n   = 0;
    this->max = p_max;
    this->mi  = p_mi;
    this->ml  = p_ml;
    this->mul = p_mul;
    this->mr  = p_mr;
    size_t sz;

    if( max * mi > 0 )
    {
        sz         = max * mi * sizeof( sint );
        void* resi = malloc( sz );
        if( !resi && max * mi > 0 ) fail( "%s: allocation of %d bytes failed\n", __FILE__, (int)sz );
        vi = (sint*)resi;
    }
    else
        vi = NULL;
    if( max * ml > 0 )
    {
        sz         = max * ml * sizeof( slong );
        void* resl = malloc( sz );
        if( !resl && max * ml > 0 ) fail( "%s: allocation of %d bytes failed\n", __FILE__, (int)sz );
        vl = (slong*)resl;
    }
    else
        vl = NULL;
    if( max * mul > 0 )
    {
        sz         = max * mul * sizeof( Ulong );
        void* resu = malloc( sz );
        if( !resu && max * mul > 0 ) fail( "%s: allocation of %d bytes failed\n", __FILE__, (int)sz );
        vul = (Ulong*)resu;
    }
    else
        vul = NULL;
    if( max * mr > 0 )
    {
        sz         = max * mr * sizeof( realType );
        void* resr = malloc( sz );
        if( !resr && max * ml > 0 ) fail( "%s: allocation of %d bytes failed\n", __FILE__, (int)sz );
        vr = (realType*)resr;
    }
    else
        vr = NULL;

    // Begin with write access disabled
    this->disableWriteAccess();

    // Set read variables
    vi_rd  = vi;
    vl_rd  = vl;
    vul_rd = vul;
    vr_rd  = vr;
}
template<class Value >
void moab::TupleList::merge_index_sort ( const Value *  A,
const Index  An,
Index  stride,
Index idx,
SortData< Value > *  work 
) [static, private]

Definition at line 752 of file TupleList.cpp.

References moab::TupleList::SortData< Value >::i, n, and moab::TupleList::SortData< Value >::v.

Referenced by index_sort().

{
    SortData< Value >* const buf[2] = { work + An, work };
    Index n = An, base = -n, odd = 0, c = 0, b = 1;
    Index i = 0;
    for( ;; )
    {
        SortData< Value >* p;
        if( ( c & 1 ) == 0 )
        {
            base += n, n += ( odd & 1 ), c |= 1, b ^= 1;
            while( n > 3 )
                odd <<= 1, odd |= ( n & 1 ), n >>= 1, c <<= 1, b ^= 1;
        }
        else
            base -= n - ( odd & 1 ), n <<= 1, n -= ( odd & 1 ), odd >>= 1, c >>= 1;
        if( c == 0 ) break;
        p = buf[b] + base;
        if( n == 2 )
        {
            Value v[2];
            v[0] = *A, A += stride, v[1] = *A, A += stride;
            if( v[1] < v[0] )
                p[0].v = v[1], p[0].i = i + 1, p[1].v = v[0], p[1].i = i;
            else
                p[0].v = v[0], p[0].i = i, p[1].v = v[1], p[1].i = i + 1;
            i += 2;
        }
        else if( n == 3 )
        {
            Value v[3];
            v[0] = *A, A += stride, v[1] = *A, A += stride, v[2] = *A, A += stride;
            if( v[1] < v[0] )
            {
                if( v[2] < v[1] )
                    p[0].v = v[2], p[1].v = v[1], p[2].v = v[0], p[0].i = i + 2, p[1].i = i + 1, p[2].i = i;
                else
                {
                    if( v[2] < v[0] )
                        p[0].v = v[1], p[1].v = v[2], p[2].v = v[0], p[0].i = i + 1, p[1].i = i + 2, p[2].i = i;
                    else
                        p[0].v = v[1], p[1].v = v[0], p[2].v = v[2], p[0].i = i + 1, p[1].i = i, p[2].i = i + 2;
                }
            }
            else
            {
                if( v[2] < v[0] )
                    p[0].v = v[2], p[1].v = v[0], p[2].v = v[1], p[0].i = i + 2, p[1].i = i, p[2].i = i + 1;
                else
                {
                    if( v[2] < v[1] )
                        p[0].v = v[0], p[1].v = v[2], p[2].v = v[1], p[0].i = i, p[1].i = i + 2, p[2].i = i + 1;
                    else
                        p[0].v = v[0], p[1].v = v[1], p[2].v = v[2], p[0].i = i, p[1].i = i + 1, p[2].i = i + 2;
                }
            }
            i += 3;
        }
        else
        {
            const Index na = n >> 1, nb = ( n + 1 ) >> 1;
            const SortData< Value >*ap = buf[b ^ 1] + base, *ae = ap + na;
            SortData< Value >*bp = p + na, *be = bp + nb;
            for( ;; )
            {
                if( bp->v < ap->v )
                {
                    *p++ = *bp++;
                    if( bp != be ) continue;
                    do
                        *p++ = *ap++;
                    while( ap != ae );
                    break;
                }
                else
                {
                    *p++ = *ap++;
                    if( ap == ae ) break;
                }
            }
        }
    }
    {
        const SortData< Value >*p = buf[0], *pe = p + An;
        do
            *idx++ = ( p++ )->i;
        while( p != pe );
    }
}
void moab::TupleList::permute ( uint perm,
void *  work 
) [private]

Definition at line 510 of file TupleList.cpp.

References mi, ml, mr, mul, n, sint, vi, vl, vr, and vul.

Referenced by sort().

{
    const unsigned int_size = mi * sizeof( sint ), long_size = ml * sizeof( slong ), Ulong_size = mul * sizeof( Ulong ),
                   real_size = mr * sizeof( realType );
    if( mi )
    {
        uint *p = perm, *pe = p + n;
        char* sorted = (char*)work;
        while( p != pe )
            memcpy( (void*)sorted, &vi[mi * ( *p++ )], int_size ), sorted += int_size;
        memcpy( vi, work, int_size * n );
    }
    if( ml )
    {
        uint *p = perm, *pe = p + n;
        char* sorted = (char*)work;
        while( p != pe )
            memcpy( (void*)sorted, &vl[ml * ( *p++ )], long_size ), sorted += long_size;
        memcpy( vl, work, long_size * n );
    }
    if( mul )
    {
        uint *p = perm, *pe = p + n;
        char* sorted = (char*)work;
        while( p != pe )
            memcpy( (void*)sorted, &vul[mul * ( *p++ )], Ulong_size ), sorted += Ulong_size;
        memcpy( vul, work, Ulong_size * n );
    }
    if( mr )
    {
        uint *p = perm, *pe = p + n;
        char* sorted = (char*)work;
        while( p != pe )
            memcpy( (void*)sorted, &vr[mr * ( *p++ )], real_size ), sorted += real_size;
        memcpy( vr, work, real_size * n );
    }
}
void moab::TupleList::print ( const char *  name) const
Examples:
CrystalRouterExample.cpp.

Definition at line 453 of file TupleList.cpp.

References mi, ml, mr, mul, n, vi, vl, vr, and vul.

Referenced by moab::ParallelComm::augment_default_sets_with_ghosts(), and main().

{
    std::cout << "Printing Tuple " << name << "===================" << std::endl;
    unsigned long i = 0, l = 0, ul = 0, r = 0;
    for( uint k = 0; k < n; k++ )
    {
        for( uint j = 0; j < mi; j++ )
        {
            std::cout << vi[i++] << " | ";
        }
        for( uint j = 0; j < ml; j++ )
        {
            std::cout << vl[l++] << " | ";
        }
        for( uint j = 0; j < mul; j++ )
        {
            std::cout << vul[ul++] << " | ";
        }
        for( uint j = 0; j < mr; j++ )
        {
            std::cout << vr[r++] << " | ";
        }
        std::cout << std::endl;
    }
    std::cout << "=======================================" << std::endl << std::endl;
}
void moab::TupleList::print_to_file ( const char *  filename) const

Definition at line 479 of file TupleList.cpp.

References mi, ml, mr, mul, n, vi, vl, vr, and vul.

Referenced by moab::ParCommGraph::compute_partition().

{
    std::ofstream ofs;
    ofs.open( filename, std::ofstream::out | std::ofstream::app );

    ofs << "Printing Tuple " << filename << "===================" << std::endl;
    unsigned long i = 0, l = 0, ul = 0, r = 0;
    for( uint k = 0; k < n; k++ )
    {
        for( uint j = 0; j < mi; j++ )
        {
            ofs << vi[i++] << " | ";
        }
        for( uint j = 0; j < ml; j++ )
        {
            ofs << vl[l++] << " | ";
        }
        for( uint j = 0; j < mul; j++ )
        {
            ofs << vul[ul++] << " | ";
        }
        for( uint j = 0; j < mr; j++ )
        {
            ofs << vr[r++] << " | ";
        }
        ofs << std::endl;
    }
    ofs << "=======================================" << std::endl << std::endl;

    ofs.close();
}
unsigned int moab::TupleList::push_back ( sint sp,
slong ip,
Ulong lp,
realType dp 
)

push back a new tuple on the TupleList;

param *sp pointer to a list of signed ints param *ip pointer to a list of signed longs param *lp pointer to a list of unsigned longs param *dp pointer to a list of reals return index of that tuple

Definition at line 393 of file TupleList.cpp.

References last_sorted, mi, ml, mr, mul, n, reserve(), vi, vl, vr, and vul.

{
    reserve();
    if( mi ) memcpy( &vi[mi * ( n - 1 )], sp, mi * sizeof( sint ) );
    if( ml ) memcpy( &vl[ml * ( n - 1 )], ip, ml * sizeof( long ) );
    if( mul ) memcpy( &vul[mul * ( n - 1 )], lp, mul * sizeof( Ulong ) );
    if( mr ) memcpy( &vr[mr * ( n - 1 )], dp, mr * sizeof( realType ) );

    last_sorted = -1;
    return n - 1;
}
template<class Value >
Value moab::TupleList::radix_count ( const Value *  A,
const Value *  end,
Index  stride,
Index  count[DIGITS][DIGIT_VALUES] 
) [static, private]

Definition at line 614 of file TupleList.cpp.

References COUNT_DIGIT_64, COUNT_SIZE, and DIGITS.

Referenced by radix_index_sort().

{
    Value bitorkey = 0;
    memset( count, 0, COUNT_SIZE * sizeof( Index ) );
    do
    {
        Value val = *A;
        bitorkey |= val;
        COUNT_DIGIT_64( DIGITS, 0 );
        // above macro expands to:
        // if(DIGITS> 0) count[ 0][val&DIGIT_MASK]++, val>>=DIGIT_BITS;
        // if(DIGITS> 1) count[ 1][val&DIGIT_MASK]++, val>>=DIGIT_BITS;
        //  ...
        // if(DIGITS>63) count[63][val&DIGIT_MASK]++, val>>=DIGIT_BITS;

    } while( A += stride, A != end );
    return bitorkey;
}
template<class Value >
void moab::TupleList::radix_index_pass_b ( const Value *  A,
Index  n,
Index  stride,
unsigned  sh,
Index off,
SortData< Value > *  out 
) [static, private]

Definition at line 662 of file TupleList.cpp.

References DIGIT_MASK, moab::TupleList::SortData< Value >::i, and moab::TupleList::SortData< Value >::v.

Referenced by radix_index_sort().

{
    Index i = 0;
    do
    {
        Value v              = *A;
        SortData< Value >* d = &out[off[( v >> sh ) & DIGIT_MASK]++];
        d->v = v, d->i = i++;
    } while( A += stride, i != n );
}
template<class Value >
void moab::TupleList::radix_index_pass_be ( const Value *  A,
Index  n,
Index  stride,
unsigned  sh,
Index off,
Index out 
) [static, private]

Definition at line 705 of file TupleList.cpp.

References DIGIT_MASK.

Referenced by radix_index_sort().

{
    Index i = 0;
    do
        out[off[( *A >> sh ) & DIGIT_MASK]++] = i++;
    while( A += stride, i != n );
}
template<class Value >
void moab::TupleList::radix_index_pass_e ( const SortData< Value > *  src,
const SortData< Value > *  end,
unsigned  sh,
Index off,
Index out 
) [static, private]

Definition at line 693 of file TupleList.cpp.

References DIGIT_MASK, moab::TupleList::SortData< Value >::i, and moab::TupleList::SortData< Value >::v.

Referenced by radix_index_sort().

{
    do
        out[off[( src->v >> sh ) & DIGIT_MASK]++] = src->i;
    while( ++src != end );
}
template<class Value >
void moab::TupleList::radix_index_pass_m ( const SortData< Value > *  src,
const SortData< Value > *  end,
unsigned  sh,
Index off,
SortData< Value > *  out 
) [static, private]

Definition at line 679 of file TupleList.cpp.

References DIGIT_MASK, moab::TupleList::SortData< Value >::i, and moab::TupleList::SortData< Value >::v.

Referenced by radix_index_sort().

{
    do
    {
        SortData< Value >* d = &out[off[( src->v >> sh ) & DIGIT_MASK]++];
        d->v = src->v, d->i = src->i;
    } while( ++src != end );
}
template<class Value >
void moab::TupleList::radix_index_sort ( const Value *  A,
Index  n,
Index  stride,
Index idx,
SortData< Value > *  work 
) [static, private]

Definition at line 714 of file TupleList.cpp.

References DIGIT_VALUES, DIGITS, n, radix_count(), radix_index_pass_b(), radix_index_pass_be(), radix_index_pass_e(), radix_index_pass_m(), radix_zeros(), and t.

Referenced by index_sort().

{
    Index count[DIGITS][DIGIT_VALUES];
    Value bitorkey = radix_count( A, A + n * stride, stride, count );
    unsigned shift[DIGITS];
    Index* offsets[DIGITS];
    unsigned digits = radix_zeros( bitorkey, count, shift, offsets );
    if( digits == 0 )
    {
        Index i = 0;
        do
            *idx++ = i++;
        while( i != n );
    }
    else if( digits == 1 )
    {
        radix_index_pass_be( A, n, stride, shift[0], offsets[0], idx );
    }
    else
    {
        SortData< Value >*src, *dst;
        unsigned d;
        if( ( digits & 1 ) == 0 )
            dst = work, src = dst + n;
        else
            src = work, dst = src + n;
        radix_index_pass_b( A, n, stride, shift[0], offsets[0], src );
        for( d = 1; d != digits - 1; ++d )
        {
            SortData< Value >* t;
            radix_index_pass_m( src, src + n, shift[d], offsets[d], dst );
            t = src, src = dst, dst = t;
        }
        radix_index_pass_e( src, src + n, shift[d], offsets[d], idx );
    }
}
void moab::TupleList::radix_offsets ( Index c) [static, private]

Definition at line 641 of file TupleList.cpp.

References DIGIT_VALUES, moab::sum(), and t.

Referenced by radix_zeros().

{
    Index sum = 0, t, *ce = c + DIGIT_VALUES;
    do
        t = *c, *c++ = sum, sum += t;
    while( c != ce );
}
template<class Value >
unsigned moab::TupleList::radix_zeros ( Value  bitorkey,
Index  count[DIGITS][DIGIT_VALUES],
unsigned *  shift,
Index **  offsets 
) [static, private]

Definition at line 650 of file TupleList.cpp.

References DIGIT_BITS, DIGIT_MASK, radix_offsets(), and VALUE_BITS.

Referenced by radix_index_sort().

{
    unsigned digits = 0, sh = 0;
    Index* c = &count[0][0];
    do
    {
        if( bitorkey & DIGIT_MASK ) *shift++ = sh, *offsets++ = c, ++digits, radix_offsets( c );
    } while( bitorkey >>= DIGIT_BITS, sh += DIGIT_BITS, c += DIGIT_VALUES, sh != VALUE_BITS );
    return digits;
}

Adds one to the total number of in-use tuples and resizes if necessary

Definition at line 226 of file TupleList.cpp.

References last_sorted, max, n, and resize().

Referenced by iMOAB_SetDoubleTagStorageWithGid(), and push_back().

{
    n++;
    while( n > max )
        resize( ( max ? max + max / 2 + 1 : 2 ) );
    last_sorted = -1;
}

Frees all allocated memory in use by the TupleList

Definition at line 205 of file TupleList.cpp.

References disableWriteAccess(), vi, vi_rd, vl, vl_rd, vr, vr_rd, vul, and vul_rd.

Referenced by moab::ParallelMergeMesh::CleanUp(), moab::ParallelComm::exchange_ghost_cells(), moab::ParallelComm::exchange_owned_mesh(), moab::Coupler::get_matching_entities(), moab::gs_data::initialize(), moab::Coupler::locate_points(), moab::ParallelMergeMesh::PerformMerge(), moab::TempestOnlineMap::ReadParallelMap(), moab::ParallelComm::resolve_shared_ents(), moab::ParallelComm::send_entities(), and ~TupleList().

{
    // free up the pointers
    free( vi );
    free( vl );
    free( vul );
    free( vr );
    // Set them all to null
    vr  = NULL;
    vi  = NULL;
    vul = NULL;
    vl  = NULL;
    // Set the read and write pointers to null
    disableWriteAccess();
    vi_rd  = NULL;
    vl_rd  = NULL;
    vul_rd = NULL;
    vr_rd  = NULL;
}

Resizes the TupleList to a new max

param max the new max capacity of the TupleList

Definition at line 137 of file TupleList.cpp.

References moab::fail(), max, MB_MEMORY_ALLOCATION_FAILED, MB_SUCCESS, mi, ml, mr, mul, sint, slong, vi, vi_rd, vi_wr, vl, vl_rd, vl_wr, vr, vr_rd, vr_wr, vul, vul_rd, vul_wr, and writeEnabled.

Referenced by moab::ParallelComm::augment_default_sets_with_ghosts(), moab::Coupler::consolidate_tuples(), moab::gs_data::initialize(), moab::Coupler::locate_points(), moab::ParallelMergeMesh::PopulateMyMatches(), moab::ParallelMergeMesh::PopulateMyTup(), reserve(), and moab::Coupler::test_local_box().

{
    this->max = maxIn;
    size_t sz;

    if( vi || ( max * mi > 0 ) )
    {
        sz         = max * mi * sizeof( sint );
        void* resi = realloc( vi, sz );
        if( !resi && max * mi > 0 )
        {
            fail( "%s: allocation of %d bytes failed\n", __FILE__, (int)sz );
            return moab::MB_MEMORY_ALLOCATION_FAILED;
        }
        vi = (sint*)resi;
    }
    if( vl || ( max * ml > 0 ) )
    {
        sz         = max * ml * sizeof( slong );
        void* resl = realloc( vl, sz );
        if( !resl && max * ml > 0 )
        {
            fail( "%s: allocation of %d bytes failed\n", __FILE__, (int)sz );
            return moab::MB_MEMORY_ALLOCATION_FAILED;
        }
        vl = (slong*)resl;
    }
    if( vul || ( max * mul > 0 ) )
    {
        sz         = max * mul * sizeof( Ulong );
        void* resu = realloc( vul, sz );
        if( !resu && max * mul > 0 )
        {
            fail( "%s: allocation of %d bytes failed\n", __FILE__, (int)sz );
            return moab::MB_MEMORY_ALLOCATION_FAILED;
        }
        vul = (Ulong*)resu;
    }
    if( vr || ( max * mr > 0 ) )
    {
        sz         = max * mr * sizeof( realType );
        void* resr = realloc( vr, sz );
        if( !resr && max * mr > 0 )
        {
            fail( "%s: allocation of %d bytes failed\n", __FILE__, (int)sz );
            return moab::MB_MEMORY_ALLOCATION_FAILED;
        }
        vr = (realType*)resr;
    }

    // Set read variables
    vi_rd  = vi;
    vl_rd  = vl;
    vul_rd = vul;
    vr_rd  = vr;

    // Set the write variables if necessary
    if( writeEnabled )
    {
        vi_wr  = vi;
        vl_wr  = vl;
        vul_wr = vul;
        vr_wr  = vr;
    }
    return moab::MB_SUCCESS;
}

Sorts the TupleList by 'key' if key<mi: key represents the key'th index in vi if mi<key<ml: key represents the (key-mi)'th index in vl if ml<key<mul: key represents the (key-mi-ml)'th index in vul

param key index to be sorted by param *buf buffer space used for sorting

Definition at line 550 of file TupleList.cpp.

References index_sort(), last_sorted, MB_NOT_IMPLEMENTED, MB_SUCCESS, mi, ml, mr, mul, n, permute(), moab::TupleList::buffer::ptr, sint, slong, umax_2, vi, vl, vul, and writeEnabled.

Referenced by moab::Coupler::consolidate_tuples(), moab::ParallelComm::exchange_owned_mesh(), moab::ParallelComm::get_sent_ents(), iMOAB_SetDoubleTagStorageWithGid(), moab::gs_data::initialize(), moab::TempestOnlineMap::ReadParallelMap(), moab::ParallelComm::resolve_shared_ents(), moab::ParallelComm::send_entities(), moab::ParallelMergeMesh::SortMyMatches(), and moab::ScdInterface::tag_shared_vertices().

{
    const unsigned int_size   = mi * sizeof( sint );
    const unsigned long_size  = ml * sizeof( slong );
    const unsigned Ulong_size = mul * sizeof( Ulong );
    const unsigned real_size  = mr * sizeof( realType );
    const unsigned width      = umax_2( umax_2( int_size, long_size ), umax_2( Ulong_size, real_size ) );
    unsigned data_size        = key >= mi ? sizeof( SortData< long > ) : sizeof( SortData< uint > );
#if defined( WIN32 ) || defined( _WIN32 )
    if( key >= mi + ml ) data_size = sizeof( SortData< Ulong > );
#endif

    uint work_min = n * umax_2( 2 * data_size, sizeof( sint ) + width );
    uint* work;
    buf->buffer_reserve( work_min );
    work = (uint*)buf->ptr;
    if( key < mi )
        index_sort( (uint*)&vi[key], n, mi, work, (SortData< uint >*)work );
    else if( key < mi + ml )
        index_sort( (long*)&vl[key - mi], n, ml, work, (SortData< long >*)work );
    else if( key < mi + ml + mul )
        index_sort( (Ulong*)&vul[key - mi - ml], n, mul, work, (SortData< Ulong >*)work );
    else
        return MB_NOT_IMPLEMENTED;

    permute( work, work + n );

    if( !writeEnabled ) last_sorted = key;
    return MB_SUCCESS;
}

Member Data Documentation

Definition at line 309 of file TupleList.hpp.

Referenced by enableWriteAccess(), find(), push_back(), reserve(), and sort().

Definition at line 296 of file TupleList.hpp.

Referenced by get_max(), initialize(), reserve(), and resize().

List of all members.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines