MOAB: Mesh Oriented datABase  (version 5.4.1)
par_spatial_locator_test.cpp File Reference
#include "moab/Core.hpp"
#include "moab/SpatialLocator.hpp"
#include "moab/Tree.hpp"
#include "moab/HomXform.hpp"
#include "moab/ScdInterface.hpp"
#include "moab/CartVect.hpp"
#include "moab/BVHTree.hpp"
#include "moab/ProgOptions.hpp"
#include "moab/CpuTimer.hpp"
#include "moab/ParallelComm.hpp"
#include "moab_mpi.h"
#include "TestUtil.hpp"
#include <cstdlib>
#include <sstream>
#include <string>
+ Include dependency graph for par_spatial_locator_test.cpp:

Go to the source code of this file.

Functions

void test_kd_tree ()
void test_bvh_tree ()
void test_locator (SpatialLocator *sl)
ErrorCode create_hex_mesh (Interface &mb, Range &elems, int n, int dim)
ErrorCode load_file (Interface &mb, std::string &fn, Range &elems)
int main (int argc, char **argv)
bool is_neg (int is_neg)

Variables

int max_depth = 30
int npoints = 1000
int leaf = 6
int tree = -1
bool print_tree = false
int ints = 10
std::string fname

Function Documentation

ErrorCode create_hex_mesh ( Interface mb,
Range elems,
int  n,
int  dim 
)

Definition at line 158 of file par_spatial_locator_test.cpp.

References CHECK_ERR, moab::ScdInterface::construct_box(), ErrorCode, moab::ScdParData::gDims, moab::Interface::get_entities_by_dimension(), n, moab::ScdParData::partMethod, moab::Interface::query_interface(), and moab::ScdParData::SQIJK.

Referenced by main(), test_bvh_tree(), and test_kd_tree().

{
    ScdInterface* scdi;
    ErrorCode rval = mb.query_interface( scdi );CHECK_ERR( rval );
    ScdParData spd;
    spd.gDims[0] = spd.gDims[1] = spd.gDims[2] = 0;
    spd.gDims[3]                               = n;
    spd.partMethod                             = ScdParData::SQIJK;
    if( dim > 1 ) spd.gDims[4] = n;
    if( dim > 2 ) spd.gDims[5] = n;
    ScdBox* new_box;
    rval = scdi->construct_box( HomCoord( 0, 0, 0 ), HomCoord( 0, 0, 0 ), NULL, 0, new_box, NULL, &spd, false, 0 );CHECK_ERR( rval );

    rval = mb.get_entities_by_dimension( 0, dim, elems );CHECK_ERR( rval );

    return rval;
}
bool is_neg ( int  is_neg)

Definition at line 114 of file par_spatial_locator_test.cpp.

Referenced by test_locator().

{
    return ( is_neg < 0 );
}
ErrorCode load_file ( Interface mb,
std::string &  fn,
Range elems 
)

Definition at line 176 of file par_spatial_locator_test.cpp.

References moab::Range::empty(), ErrorCode, moab::Interface::get_entities_by_dimension(), moab::Interface::load_file(), and MB_SUCCESS.

Referenced by moab::Core::load_mesh(), test_bvh_tree(), and test_kd_tree().

{
    std::string options;
    ErrorCode rval;

    options = "PARALLEL=READ_PART;PARTITION=PARALLEL_PARTITION;";
    rval    = mb.load_file( fn.c_str(), NULL, options.c_str() );
    if( MB_SUCCESS != rval ) return rval;

    rval = mb.get_entities_by_dimension( 0, 3, elems );
    if( MB_SUCCESS != rval ) return rval;

    if( elems.empty() )
    {
        rval = mb.get_entities_by_dimension( 0, 2, elems );
        if( MB_SUCCESS != rval ) return rval;
    }

    return MB_SUCCESS;
}
int main ( int  argc,
char **  argv 
)

Definition at line 36 of file par_spatial_locator_test.cpp.

References ProgOptions::addOpt(), moab::fail(), fname, ints, leaf, max_depth, npoints, ProgOptions::parseCommandLine(), print_tree, RUN_TEST, test_bvh_tree(), test_kd_tree(), and tree.

{
    int fail = MPI_Init( &argc, &argv );
    if( fail ) return fail;

    ProgOptions po( "spatial_locator_test options" );
    po.addOpt< std::string >( "file,f", "Input filename", &fname );
    po.addOpt< int >( "ints,i", "Number of intervals on each side of scd mesh", &ints );
    po.addOpt< int >( "leaf,l", "Maximum number of elements per leaf", &leaf );
    po.addOpt< int >( "max_depth,m", "Maximum depth of tree", &max_depth );
    po.addOpt< int >( "npoints,n", "Number of query points", &npoints );
    po.addOpt< void >( "print,p", "Print tree details", &print_tree );
    po.addOpt< int >( "tree,t", "Tree type (-1=all (default), 0=AdaptiveKD, 1=BVH", &tree );
    po.parseCommandLine( argc, argv );

    if( -1 == tree || 0 == tree ) RUN_TEST( test_kd_tree );
    if( -1 == tree || 1 == tree ) RUN_TEST( test_bvh_tree );

    fail = MPI_Finalize();
    if( fail ) return fail;

    return 0;
}
void test_bvh_tree ( )

Definition at line 85 of file par_spatial_locator_test.cpp.

References CHECK_ERR, create_hex_mesh(), ErrorCode, fname, ints, leaf, load_file(), max_depth, mb, moab::BVHTree::parse_options(), and test_locator().

Referenced by main().

{
    ErrorCode rval;
    Core mb;

    // create a simple mesh to test
    Range elems;
    if( fname.empty() )
    {
        rval = create_hex_mesh( mb, elems, ints, 3 );CHECK_ERR( rval );
    }
    else
    {
        rval = load_file( mb, fname, elems );CHECK_ERR( rval );
    }

    // initialize spatial locator with the elements and a BVH tree
    BVHTree bvh( &mb );
    std::ostringstream opts;
    opts << "MAX_DEPTH=" << max_depth << ";MAX_PER_LEAF=" << leaf;
    FileOptions fo( opts.str().c_str() );
    rval               = bvh.parse_options( fo );
    SpatialLocator* sl = new SpatialLocator( &mb, elems, &bvh );
    test_locator( sl );

    // destroy spatial locator, and tree along with it
    delete sl;
}
void test_kd_tree ( )

Definition at line 60 of file par_spatial_locator_test.cpp.

References CHECK_ERR, create_hex_mesh(), ErrorCode, fname, ints, load_file(), mb, and test_locator().

Referenced by main().

{
    ErrorCode rval;
    Core mb;

    // create a simple mesh to test
    Range elems;
    if( fname.empty() )
    {
        rval = create_hex_mesh( mb, elems, ints, 3 );CHECK_ERR( rval );
    }
    else
    {
        rval = load_file( mb, fname, elems );CHECK_ERR( rval );
    }

    // initialize spatial locator with the elements and the default tree type
    SpatialLocator* sl = new SpatialLocator( &mb, elems );

    test_locator( sl );

    // destroy spatial locator, and tree along with it
    delete sl;
}
void test_locator ( SpatialLocator sl)

Definition at line 119 of file par_spatial_locator_test.cpp.

References moab::BoundBox::bMax, moab::BoundBox::bMin, box(), CHECK, CHECK_ERR, ErrorCode, moab::ParallelComm::get_pcomm(), moab::SpatialLocator::get_tree(), is_neg(), moab::SpatialLocator::local_box(), moab::SpatialLocator::moab(), npoints, moab::TreeStats::output_trav_stats(), moab::SpatialLocator::par_loc_table(), moab::Tree::print(), print_tree, moab::ParallelComm::rank(), moab::Tree::tree_stats(), and moab::TupleList::vi_rd.

Referenced by main(), test_bvh_tree(), and test_kd_tree().

{
    CartVect box_del;
    BoundBox box = sl->local_box();
    box_del      = box.bMax - box.bMin;

    double denom = 1.0 / (double)RAND_MAX;
    std::vector< CartVect > test_pts( npoints );

    for( int i = 0; i < npoints; i++ )
    {
        // generate a small number of random point to test
        double rx = (double)rand() * denom, ry = (double)rand() * denom, rz = (double)rand() * denom;
        test_pts[i] = box.bMin + CartVect( rx * box_del[0], ry * box_del[1], rz * box_del[2] );
    }

    // call spatial locator to locate points
    ParallelComm* pc = ParallelComm::get_pcomm( sl->moab(), 0 );
    CHECK( pc != NULL );

    ErrorCode rval = sl->par_locate_points( pc, test_pts[0].array(), npoints );CHECK_ERR( rval );
    if( pc->rank() == 0 )
    {
        int num_out = std::count_if( sl->par_loc_table().vi_rd, sl->par_loc_table().vi_rd + 2 * npoints, is_neg );
        num_out /= 2;

        std::cout << "Number of points inside an element = " << npoints - num_out << "/" << npoints << " ("
                  << 100.0 * ( (double)npoints - num_out ) / npoints << "%)" << std::endl;
        std::cout << "Traversal stats:" << std::endl;
        sl->get_tree()->tree_stats().output_trav_stats();

        if( print_tree )
        {
            std::cout << "Tree information: " << std::endl;
            rval = sl->get_tree()->print();CHECK_ERR( rval );
        }
    }
}

Variable Documentation

int npoints = 1000

Definition at line 29 of file par_spatial_locator_test.cpp.

Referenced by moab::ReadSms::load_file_impl(), main(), and test_locator().

bool print_tree = false

Definition at line 32 of file par_spatial_locator_test.cpp.

Referenced by main(), and test_locator().

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines