//---------------------------------------------------------------------------//
// Hex-8 test.
TEUCHOS_UNIT_TEST( MoabEntityIterator, hex_8_test )
{
    // Extract the raw mpi communicator.
    Teuchos::RCP<const Teuchos::Comm<int> > comm = getDefaultComm<int>();
    Teuchos::RCP<const Teuchos::MpiComm<int> > mpi_comm = 
	Teuchos::rcp_dynamic_cast< const Teuchos::MpiComm<int> >( comm );
    Teuchos::RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > opaque_comm = 
	mpi_comm->getRawMpiComm();
    MPI_Comm raw_comm = (*opaque_comm)();

    // Create the mesh.
    int space_dim = 3;
    Teuchos::RCP<moab::Interface> moab_mesh = Teuchos::rcp( new moab::Core() );
    Teuchos::RCP<moab::ParallelComm> parallel_mesh =
	Teuchos::rcp( new moab::ParallelComm(moab_mesh.getRawPtr(),raw_comm) );

    // Create the nodes.
    moab::ErrorCode error = moab::MB_SUCCESS;
    Teuchos::Array<moab::EntityHandle> nodes(8);
    double node_coords[3];
    node_coords[0] = 0.0;
    node_coords[1] = 0.0;
    node_coords[2] = 0.0;
    error = moab_mesh->create_vertex( node_coords, nodes[0] );
    TEST_EQUALITY( error, moab::MB_SUCCESS );

    node_coords[0] = 1.0;
    node_coords[1] = 0.0;
    node_coords[2] = 0.0;
    error = moab_mesh->create_vertex( node_coords, nodes[1] );
    TEST_EQUALITY( error, moab::MB_SUCCESS );

    node_coords[0] = 1.0;
    node_coords[1] = 1.0;
    node_coords[2] = 0.0;
    error = moab_mesh->create_vertex( node_coords, nodes[2] );
    TEST_EQUALITY( error, moab::MB_SUCCESS );

    node_coords[0] = 0.0;
    node_coords[1] = 1.0;
    node_coords[2] = 0.0;
    error = moab_mesh->create_vertex( node_coords, nodes[3] );
    TEST_EQUALITY( error, moab::MB_SUCCESS );

    node_coords[0] = 0.0;
    node_coords[1] = 0.0;
    node_coords[2] = 1.0;
    error = moab_mesh->create_vertex( node_coords, nodes[4] );
    TEST_EQUALITY( error, moab::MB_SUCCESS );

    node_coords[0] = 1.0;
    node_coords[1] = 0.0;
    node_coords[2] = 1.0;
    error = moab_mesh->create_vertex( node_coords, nodes[5] );
    TEST_EQUALITY( error, moab::MB_SUCCESS );

    node_coords[0] = 1.0;
    node_coords[1] = 1.0;
    node_coords[2] = 1.0;
    error = moab_mesh->create_vertex( node_coords, nodes[6] );
    TEST_EQUALITY( error, moab::MB_SUCCESS );

    node_coords[0] = 0.0;
    node_coords[1] = 1.0;
    node_coords[2] = 1.0;
    error = moab_mesh->create_vertex( node_coords, nodes[7] );
    TEST_EQUALITY( error, moab::MB_SUCCESS );

    // Make a hex-8.
    moab::EntityHandle hex_entity;
    error = moab_mesh->create_element( moab::MBHEX,
				       nodes.getRawPtr(),
				       8,
				       hex_entity );
    TEST_EQUALITY( error, moab::MB_SUCCESS );

    // Make 2 entity sets.
    moab::EntityHandle entity_set_1;
    error = moab_mesh->create_meshset( 0, entity_set_1 );
    TEST_EQUALITY( error, moab::MB_SUCCESS );
    moab::EntityHandle entity_set_2;
    error = moab_mesh->create_meshset( 0, entity_set_2 );
    TEST_EQUALITY( error, moab::MB_SUCCESS );

    // Put the hex-8 in the first entity set.
    error = moab_mesh->add_entities( entity_set_1, &hex_entity, 1 );

    // Index the sets in the mesh.
    Teuchos::RCP<DataTransferKit::MoabMeshSetIndexer> set_indexer =
	Teuchos::rcp( new DataTransferKit::MoabMeshSetIndexer(parallel_mesh) );

    // Make a list of hexes.
    unsigned num_hex = 2;
    std::vector<moab::EntityHandle> hex_entities( num_hex, hex_entity );
    
    // Make an iterator for the hex.
    std::function<bool(DataTransferKit::Entity)> all_pred = 
	[=] (DataTransferKit::Entity e){return true;};
    Teuchos::RCP<DataTransferKit::MoabEntityIteratorRange> iterator_range =
	Teuchos::rcp( new DataTransferKit::MoabEntityIteratorRange() );
    iterator_range->d_moab_entities = hex_entities;
    DataTransferKit::EntityIterator entity_iterator = 
	DataTransferKit::MoabEntityIterator(
	    iterator_range, parallel_mesh.ptr(), set_indexer.ptr(), all_pred );

    // Test the entity iterator.
    TEST_EQUALITY( entity_iterator.size(), num_hex );
    TEST_ASSERT( entity_iterator == entity_iterator.begin() );
    TEST_ASSERT( entity_iterator != entity_iterator.end() );

    // Test the first entity under the iterator with a pointer dereference.
    DataTransferKit::EntityId hex_id = 90343;
    DataTransferKit::MoabHelpers::getGlobalIds(
	*parallel_mesh, &hex_entity, 1, &hex_id );
    TEST_EQUALITY( hex_id, entity_iterator->id() );
    TEST_EQUALITY( comm->getRank(), entity_iterator->ownerRank() );
    TEST_EQUALITY( space_dim, entity_iterator->topologicalDimension() );
    TEST_EQUALITY( space_dim, entity_iterator->physicalDimension() );

    int set_1_id = set_indexer->getIndexFromMeshSet( entity_set_1 );
    int set_2_id = set_indexer->getIndexFromMeshSet( entity_set_2 );
    TEST_ASSERT( entity_iterator->inBlock(set_1_id) );
    TEST_ASSERT( !entity_iterator->inBlock(set_2_id) );
    TEST_ASSERT( entity_iterator->onBoundary(set_1_id) );
    TEST_ASSERT( !entity_iterator->onBoundary(set_2_id) );

    Teuchos::RCP<DataTransferKit::EntityExtraData> extra_data_1 =
	entity_iterator->extraData();
    TEST_EQUALITY( hex_entity,
		   Teuchos::rcp_dynamic_cast<DataTransferKit::MoabEntityExtraData>(
		       extra_data_1)->d_moab_entity );

    Teuchos::Tuple<double,6> hex_bounds_1;
    entity_iterator->boundingBox( hex_bounds_1 );
    TEST_EQUALITY( 0.0, hex_bounds_1[0] );
    TEST_EQUALITY( 0.0, hex_bounds_1[1] );
    TEST_EQUALITY( 0.0, hex_bounds_1[2] );
    TEST_EQUALITY( 1.0, hex_bounds_1[3] );
    TEST_EQUALITY( 1.0, hex_bounds_1[4] );
    TEST_EQUALITY( 1.0, hex_bounds_1[5] );

    // Increment the iterator
    ++entity_iterator;

    // Test the second entity under the iterator with a reference dereference.
    TEST_EQUALITY( hex_id, (*entity_iterator).id() );
    TEST_EQUALITY( comm->getRank(), (*entity_iterator).ownerRank() );
    TEST_EQUALITY( space_dim, (*entity_iterator).topologicalDimension() );
    TEST_EQUALITY( space_dim, (*entity_iterator).physicalDimension() );

    TEST_ASSERT( (*entity_iterator).inBlock(set_1_id) );
    TEST_ASSERT( !(*entity_iterator).inBlock(set_2_id) );
    TEST_ASSERT( (*entity_iterator).onBoundary(set_1_id) );
    TEST_ASSERT( !(*entity_iterator).onBoundary(set_2_id) );

    Teuchos::RCP<DataTransferKit::EntityExtraData> extra_data_2 =
	(*entity_iterator).extraData();
    TEST_EQUALITY( hex_entity,
		   Teuchos::rcp_dynamic_cast<DataTransferKit::MoabEntityExtraData>(
		       extra_data_2)->d_moab_entity );

    Teuchos::Tuple<double,6> hex_bounds_2;
    (*entity_iterator).boundingBox( hex_bounds_2 );
    TEST_EQUALITY( 0.0, hex_bounds_2[0] );
    TEST_EQUALITY( 0.0, hex_bounds_2[1] );
    TEST_EQUALITY( 0.0, hex_bounds_2[2] );
    TEST_EQUALITY( 1.0, hex_bounds_2[3] );
    TEST_EQUALITY( 1.0, hex_bounds_2[4] );
    TEST_EQUALITY( 1.0, hex_bounds_2[5] );

    // Test the end of the iterator.
    entity_iterator++;
    TEST_ASSERT( entity_iterator != entity_iterator.begin() );
    TEST_ASSERT( entity_iterator == entity_iterator.end() );

    // Make an iterator with a part 1 predicate.
    DataTransferKit::MoabMeshSetPredicate set_1_pred( entity_set_1, set_indexer );
    DataTransferKit::EntityIterator set_1_iterator =
	DataTransferKit::MoabEntityIterator(
	    iterator_range, parallel_mesh.ptr(), set_indexer.ptr(), set_1_pred.getFunction() );
    TEST_EQUALITY( set_1_iterator.size(), num_hex );

    // Make an iterator with a part 2 predicate.
    DataTransferKit::MoabMeshSetPredicate set_2_pred( entity_set_2, set_indexer );
    DataTransferKit::EntityIterator set_2_iterator =
	DataTransferKit::MoabEntityIterator(
	    iterator_range, parallel_mesh.ptr(), set_indexer.ptr(), set_2_pred.getFunction() );
    TEST_EQUALITY( set_2_iterator.size(), 0 );
}
//---------------------------------------------------------------------------//
// Hex-8 test.
TEUCHOS_UNIT_TEST( STKMeshEntityIterator, hex_8_test )
{
    // Extract the raw mpi communicator.
    Teuchos::RCP<const Teuchos::Comm<int> > comm = getDefaultComm<int>();
    Teuchos::RCP<const Teuchos::MpiComm<int> > mpi_comm = 
	Teuchos::rcp_dynamic_cast< const Teuchos::MpiComm<int> >( comm );
    Teuchos::RCP<const Teuchos::OpaqueWrapper<MPI_Comm> > opaque_comm = 
	mpi_comm->getRawMpiComm();
    MPI_Comm raw_comm = (*opaque_comm)();

    // Create meta data.
    int space_dim = 3;
    stk::mesh::MetaData meta_data( space_dim );

    // Make two parts.
    std::string p1_name = "part_1";
    stk::mesh::Part& part_1 = meta_data.declare_part( p1_name );
    stk::mesh::set_topology( part_1, stk::topology::HEX_8 );
    std::string p2_name = "part_2";
    stk::mesh::Part& part_2 = meta_data.declare_part( p2_name );

    // Make a coordinate field.
    stk::mesh::Field<double, stk::mesh::Cartesian3d>& coord_field =
	meta_data.declare_field<
	stk::mesh::Field<double, stk::mesh::Cartesian3d> >(
	    stk::topology::NODE_RANK, "coordinates");
    meta_data.set_coordinate_field( &coord_field );
    stk::mesh::put_field( coord_field, part_1 );
    meta_data.commit();

    // Create bulk data.
    Teuchos::RCP<stk::mesh::BulkData> bulk_data =
	Teuchos::rcp( new stk::mesh::BulkData(meta_data,raw_comm) );
    bulk_data->modification_begin();

    // Make a hex-8.
    int comm_rank = comm->getRank();
    stk::mesh::EntityId hex_id = 23 + comm_rank;
    stk::mesh::Entity hex_entity = 
	bulk_data->declare_entity( stk::topology::ELEM_RANK, hex_id, part_1 );
    int num_nodes = 8;
    Teuchos::Array<stk::mesh::EntityId> node_ids( num_nodes );
    Teuchos::Array<stk::mesh::Entity> nodes( num_nodes );
    for ( int i = 0; i < num_nodes; ++i )
    {
	node_ids[i] = num_nodes*comm_rank + i + 5;
	nodes[i] = bulk_data->declare_entity( 
	    stk::topology::NODE_RANK, node_ids[i], part_1 );
	bulk_data->declare_relation( hex_entity, nodes[i], i );
    }
    bulk_data->modification_end();

    // Create the node coordinates.
    double* node_coords = 0;
    node_coords = stk::mesh::field_data( coord_field, nodes[0] );
    node_coords[0] = 0.0;
    node_coords[1] = 0.0;
    node_coords[2] = 0.0;

    node_coords = stk::mesh::field_data( coord_field, nodes[1] );
    node_coords[0] = 1.0;
    node_coords[1] = 0.0;
    node_coords[2] = 0.0;

    node_coords = stk::mesh::field_data( coord_field, nodes[2] );
    node_coords[0] = 1.0;
    node_coords[1] = 1.0;
    node_coords[2] = 0.0;

    node_coords = stk::mesh::field_data( coord_field, nodes[3] );
    node_coords[0] = 0.0;
    node_coords[1] = 1.0;
    node_coords[2] = 0.0;

    node_coords = stk::mesh::field_data( coord_field, nodes[4] );
    node_coords[0] = 0.0;
    node_coords[1] = 0.0;
    node_coords[2] = 1.0;

    node_coords = stk::mesh::field_data( coord_field, nodes[5] );
    node_coords[0] = 1.0;
    node_coords[1] = 0.0;
    node_coords[2] = 1.0;

    node_coords = stk::mesh::field_data( coord_field, nodes[6] );
    node_coords[0] = 1.0;
    node_coords[1] = 1.0;
    node_coords[2] = 1.0;

    node_coords = stk::mesh::field_data( coord_field, nodes[7] );
    node_coords[0] = 0.0;
    node_coords[1] = 1.0;
    node_coords[2] = 1.0;

    // Make a list of hexes.
    unsigned num_hex = 2;
    std::vector<stk::mesh::Entity> hex_entities( num_hex, hex_entity );

    // Make a range for the iterators.
    Teuchos::RCP<DataTransferKit::STKMeshEntityIteratorRange> iterator_range =
	Teuchos::rcp( new DataTransferKit::STKMeshEntityIteratorRange() );
    iterator_range->d_stk_entities = hex_entities;
    
    // Test the name predicate for part 1.
    DataTransferKit::STKPartNamePredicate part_1_name_pred( 
	Teuchos::Array<std::string>(1,p1_name), bulk_data );
    DataTransferKit::EntityIterator part_1_name_iterator =
	DataTransferKit::STKMeshEntityIterator(
	    iterator_range, bulk_data, part_1_name_pred.getFunction() );
    TEST_EQUALITY( part_1_name_iterator.size(), num_hex );

    // Test the name predicate for part 2.
    DataTransferKit::STKPartNamePredicate part_2_name_pred( 
	Teuchos::Array<std::string>(2,p2_name), bulk_data );
    DataTransferKit::EntityIterator part_2_name_iterator =
	DataTransferKit::STKMeshEntityIterator(
	    iterator_range, bulk_data, part_2_name_pred.getFunction() );
    TEST_EQUALITY( part_2_name_iterator.size(), 0 );

    // Test the part vector predicate for part 1.
    stk::mesh::PartVector p1_vec( 1, &part_1 );
    DataTransferKit::STKPartVectorPredicate part_1_vec_pred( p1_vec );
    DataTransferKit::EntityIterator part_1_vec_iterator =
	DataTransferKit::STKMeshEntityIterator(
	    iterator_range, bulk_data, part_1_vec_pred.getFunction() );
    TEST_EQUALITY( part_1_vec_iterator.size(), num_hex );

    // Test the part vector predicate for part 2.
    stk::mesh::PartVector p2_vec( 2, &part_2 );
    DataTransferKit::STKPartVectorPredicate part_2_vec_pred( p2_vec );
    DataTransferKit::EntityIterator part_2_vec_iterator =
	DataTransferKit::STKMeshEntityIterator(
	    iterator_range, bulk_data, part_2_vec_pred.getFunction() );
    TEST_EQUALITY( part_2_vec_iterator.size(), 0 );

    // Test a part vector with 2 part 1's.
    stk::mesh::PartVector p11_vec( 2, &part_1 );
    DataTransferKit::STKPartVectorPredicate part_11_vec_pred( p11_vec );
    DataTransferKit::EntityIterator part_11_vec_iterator =
	DataTransferKit::STKMeshEntityIterator(
	    iterator_range, bulk_data, part_11_vec_pred.getFunction() );
    TEST_EQUALITY( part_11_vec_iterator.size(), num_hex );

    // Test a part vector with a part 1 and part 2
    stk::mesh::PartVector p12_vec( 2 );
    p12_vec[0] = &part_1;
    p12_vec[1] = &part_2;
    DataTransferKit::STKPartVectorPredicate part_12_vec_pred( p12_vec );
    DataTransferKit::EntityIterator part_12_vec_iterator =
	DataTransferKit::STKMeshEntityIterator(
	    iterator_range, bulk_data, part_12_vec_pred.getFunction() );
    TEST_EQUALITY( part_12_vec_iterator.size(), 0 );

    // Test the part selector predicate for part 1.
    stk::mesh::Selector p1_sel( part_1 );
    DataTransferKit::STKSelectorPredicate part_1_sel_pred( p1_sel );
    DataTransferKit::EntityIterator part_1_sel_iterator =
	DataTransferKit::STKMeshEntityIterator(
	    iterator_range, bulk_data, part_1_sel_pred.getFunction() );
    TEST_EQUALITY( part_1_sel_iterator.size(), num_hex );

    // Test the part selector predicate for part 2.
    stk::mesh::Selector p2_sel( part_2 );
    DataTransferKit::STKSelectorPredicate part_2_sel_pred( p2_sel );
    DataTransferKit::EntityIterator part_2_sel_iterator =
	DataTransferKit::STKMeshEntityIterator(
	    iterator_range, bulk_data, part_2_sel_pred.getFunction() );
    TEST_EQUALITY( part_2_sel_iterator.size(), 0 );
}