[MOAB-dev] r1252 - MOAB/trunk

kraftche at mcs.anl.gov kraftche at mcs.anl.gov
Fri Aug 24 15:05:52 CDT 2007


Author: kraftche
Date: 2007-08-24 15:05:52 -0500 (Fri, 24 Aug 2007)
New Revision: 1252

Modified:
   MOAB/trunk/MBCore.cpp
   MOAB/trunk/MBReaderIface.hpp
   MOAB/trunk/ReadGmsh.cpp
   MOAB/trunk/ReadGmsh.hpp
   MOAB/trunk/ReadHDF5.cpp
   MOAB/trunk/ReadHDF5.hpp
   MOAB/trunk/ReadNCDF.cpp
   MOAB/trunk/ReadNCDF.hpp
   MOAB/trunk/ReadSTL.cpp
   MOAB/trunk/ReadSTL.hpp
   MOAB/trunk/ReadVtk.cpp
   MOAB/trunk/ReadVtk.hpp
   MOAB/trunk/Tqdcfr.cpp
   MOAB/trunk/Tqdcfr.hpp
Log:
o Change MBReaderIface::load_file such that readers are required to
  pass back an entity set containing all entities read from file.
o Update all readers to create and pass back an entity set containing
  entities read from file (all but Tqdcfr already maintained such an
  entity set internally.)





Modified: MOAB/trunk/MBCore.cpp
===================================================================
--- MOAB/trunk/MBCore.cpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/MBCore.cpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -359,7 +359,7 @@
   MBReaderIface* reader = set->get_file_extension_reader( file_name );
   if (reader)
   { 
-    rval = reader->load_file( file_name, block_id_list, num_blocks );
+    rval = reader->load_file( file_name, file_set, block_id_list, num_blocks );
     delete reader;
     return rval;
   }
@@ -371,7 +371,7 @@
     MBReaderIface* reader = iter->make_reader( this );
     if (NULL != reader)
     {
-      rval = reader->load_file( file_name, block_id_list, num_blocks );
+      rval = reader->load_file( file_name, file_set, block_id_list, num_blocks );
       delete reader;
       if (MB_SUCCESS == rval)
         return MB_SUCCESS;

Modified: MOAB/trunk/MBReaderIface.hpp
===================================================================
--- MOAB/trunk/MBReaderIface.hpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/MBReaderIface.hpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -38,12 +38,14 @@
      * Method all readers must provide to import a mesh.
      *
      *\param file_name           The file to read.
+     *\param file_set            Output: a new entity set containing all data read from file.
      *\param material_set_list   A list of material sets to read, or NULL
      *                           if the entire file is to be read.
      *\param material_set_list_len The length of <code>material_set_list</code>
      *\author Jason Kraftcheck
      */
     virtual MBErrorCode load_file( const char* file_name,
+                                   MBEntityHandle& file_set,
                                    const int* material_set_list,
                                    const int material_set_list_len ) = 0;
 };

Modified: MOAB/trunk/ReadGmsh.cpp
===================================================================
--- MOAB/trunk/ReadGmsh.cpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/ReadGmsh.cpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -77,6 +77,7 @@
 
 
 MBErrorCode ReadGmsh::load_file( const char* filename, 
+                                 MBEntityHandle& file_set,
                                  const int* blocks,
                                  const int num_blocks )
 {
@@ -94,6 +95,7 @@
     mCurrentMeshHandle = 0;
   }
   
+  file_set = mCurrentMeshHandle;
   return result;
 }
 

Modified: MOAB/trunk/ReadGmsh.hpp
===================================================================
--- MOAB/trunk/ReadGmsh.hpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/ReadGmsh.hpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -39,6 +39,7 @@
   static MBReaderIface* factory( MBInterface* );
 
   MBErrorCode load_file(const char *file_name,
+                        MBEntityHandle& file_set,
                         const int* material_set_list,
                         const int num_material_sets );
   

Modified: MOAB/trunk/ReadHDF5.cpp
===================================================================
--- MOAB/trunk/ReadHDF5.cpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/ReadHDF5.cpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -107,7 +107,7 @@
   H5Tclose( handleType );
 }
 
-MBErrorCode ReadHDF5::load_file( const char* filename, const int*, const int num_blocks )
+MBErrorCode ReadHDF5::load_file( const char* filename, MBEntityHandle& file_set, const int*, const int num_blocks )
 {
   MBErrorCode rval;
   mhdf_Status status;
@@ -116,9 +116,11 @@
   char** tag_names = NULL;
   char** groups = NULL;
   std::list<ElemSet>::iterator el_itor;
+  std::list<ElemSet>::reverse_iterator rel_itor;
   unsigned int i, num_groups;
   MBEntityHandle all;  // meshset of everything in file.
   bool have_nodes = true;
+  file_set = 0;
 
   if (num_blocks)
     return MB_FAILURE;
@@ -211,6 +213,18 @@
 DEBUGOUT("Finishing read.\n");
   if (MB_SUCCESS != read_qa( all ))
     goto read_fail;
+    
+DEBUGOUT("Creating entity set for file contents\n")
+  if (MB_SUCCESS != iFace->create_meshset( MESHSET_SET, file_set ))
+    goto read_fail;
+  if (MB_SUCCESS != iFace->add_entities( file_set, setSet.range ))
+    goto read_fail;
+  for (rel_itor = elemList.rbegin(); rel_itor != elemList.rend(); ++rel_itor)
+    if (MB_SUCCESS != iFace->add_entities( file_set, rel_itor->range))
+      goto read_fail;
+  if (MB_SUCCESS != iFace->add_entities( file_set, nodeSet.range ))
+    goto read_fail;
+  
 
     // Clean up and exit.
   free( dataBuffer );
@@ -223,6 +237,9 @@
   
 read_fail:
   
+  if (file_set)
+    iFace->delete_entities( &file_set, 1 );
+  
   if (dataBuffer)
   {
     free( dataBuffer );

Modified: MOAB/trunk/ReadHDF5.hpp
===================================================================
--- MOAB/trunk/ReadHDF5.hpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/ReadHDF5.hpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -47,6 +47,7 @@
    * \param export_set_count Length of <code>export_sets</code> array.
    */
   MBErrorCode load_file( const char* filename,
+                         MBEntityHandle& file_set,
                          const int* material_set_list,
                          int material_set_count  );
 

Modified: MOAB/trunk/ReadNCDF.cpp
===================================================================
--- MOAB/trunk/ReadNCDF.cpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/ReadNCDF.cpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -406,9 +406,11 @@
 
 
 MBErrorCode ReadNCDF::load_file(const char *exodus_file_name,
+                                  MBEntityHandle& file_set,
                                   const int *blocks_to_load,
                                   const int num_blocks) 
 {
+  file_set = 0;
     // this function directs the reading of an exoii file, but doesn't do any of
     // the actual work
   
@@ -460,6 +462,7 @@
 
     // what about properties???
 
+  file_set = mCurrentMeshHandle;
   return MB_SUCCESS;
 }
 

Modified: MOAB/trunk/ReadNCDF.hpp
===================================================================
--- MOAB/trunk/ReadNCDF.hpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/ReadNCDF.hpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -64,6 +64,7 @@
   
     //! load an ExoII file
   MBErrorCode load_file(const char *exodus_file_name,
+                         MBEntityHandle& file_set,
                          const int* blocks_to_load,
                          const int num_blocks);
   

Modified: MOAB/trunk/ReadSTL.cpp
===================================================================
--- MOAB/trunk/ReadSTL.cpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/ReadSTL.cpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -78,9 +78,10 @@
 }
 
 
-MBErrorCode ReadSTL::load_file( const char* filename, 
-                                 const int* blocks,
-                                 const int num_blocks )
+MBErrorCode ReadSTL::load_file( const char* filename,
+                                MBEntityHandle& file_set, 
+                                const int* blocks, 
+                                const int num_blocks )
 {
   mCurrentMeshHandle = 0;
   const MBErrorCode result = load_file_impl( filename, blocks, num_blocks );
@@ -96,6 +97,7 @@
     mCurrentMeshHandle = 0;
   }
   
+  file_set = mCurrentMeshHandle;
   return result;
 }
 

Modified: MOAB/trunk/ReadSTL.hpp
===================================================================
--- MOAB/trunk/ReadSTL.hpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/ReadSTL.hpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -74,6 +74,7 @@
     //! Generic file loading code for both binary and ASCII readers.
     //! Calls reader-specific read_triangles function to do actual I/O.
   MBErrorCode load_file(const char *file_name,
+                        MBEntityHandle& file_set,
                         const int* material_set_list,
                         const int num_material_sets );
   

Modified: MOAB/trunk/ReadVtk.cpp
===================================================================
--- MOAB/trunk/ReadVtk.cpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/ReadVtk.cpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -53,9 +53,11 @@
 
 
 MBErrorCode ReadVtk::load_file(const char *filename,
+                               MBEntityHandle& file_set,
                                const int*, const int) 
 {
   MBErrorCode result;
+  file_set = 0;
 
   int major, minor;
   char vendor_string[257];
@@ -182,6 +184,7 @@
       return result;
   }
   
+  file_set = mCurrentMeshHandle;
   return MB_SUCCESS;
 }
 

Modified: MOAB/trunk/ReadVtk.hpp
===================================================================
--- MOAB/trunk/ReadVtk.hpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/ReadVtk.hpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -31,6 +31,7 @@
 
     //! load a file
   MBErrorCode load_file(const char *file_name,
+                        MBEntityHandle& file_set,
                         const int* material_set_list,
                         const int num_material_sets );
   

Modified: MOAB/trunk/Tqdcfr.cpp
===================================================================
--- MOAB/trunk/Tqdcfr.cpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/Tqdcfr.cpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -170,9 +170,11 @@
 
   
 MBErrorCode Tqdcfr::load_file(const char *file_name,
+                              MBEntityHandle& file_set,
                               const int*, const int) 
 {
   MBErrorCode result;
+  file_set = mFileSet = 0;
   
     // open file
   cubFile = fopen(file_name, "rb");
@@ -183,6 +185,12 @@
   if (!(char_buf[0] == 'C' && char_buf[1] == 'U' && 
         char_buf[2] == 'B' && char_buf[3] == 'E')) 
     return MB_FAILURE;
+  
+    // create meshset to contain file
+  result = mdbImpl->create_meshset( MESHSET_SET, file_set );
+  if (MB_SUCCESS != result) 
+    return result;
+  mFileSet = file_set;
 
     // ***********************
     // read model header type information...
@@ -567,7 +575,7 @@
   if (!reverse.empty()) {
       // need to make a new set, add a reverse sense tag, then add to the sideset
     MBEntityHandle reverse_set;
-    MBErrorCode tmp_result = mdbImpl->create_meshset(MESHSET_SET, reverse_set);
+    MBErrorCode tmp_result = create_set(reverse_set);
     if (MB_SUCCESS != tmp_result) result = tmp_result;
     tmp_result = mdbImpl->add_entities(reverse_set, &reverse[0], reverse.size());
     if (tmp_result != MB_SUCCESS) result = tmp_result;
@@ -626,7 +634,7 @@
   if (!reverse.empty()) {
       // need to make a new set, add a reverse sense tag, then add to the sideset
     MBEntityHandle reverse_set;
-    MBErrorCode tmp_result = mdbImpl->create_meshset(MESHSET_SET, reverse_set);
+    MBErrorCode tmp_result = create_set(reverse_set);
     if (MB_SUCCESS != tmp_result) result = tmp_result;
     tmp_result = mdbImpl->add_entities(reverse_set, &reverse[0], reverse.size());
     if (tmp_result != MB_SUCCESS) result = tmp_result;
@@ -941,6 +949,8 @@
                     node_handle+entity->nodeCt-1);
   MBErrorCode result = mdbImpl->add_entities(entity->setHandle, dum_range);
   if (MB_SUCCESS != result) return result;
+  result = mdbImpl->add_entities( mFileSet, dum_range );
+  if (MB_SUCCESS != result) return result;
 
     // check for id contiguity
   long unsigned int node_offset = mdbImpl->id_from_handle( node_handle);
@@ -1109,6 +1119,8 @@
     MBRange dum_range(start_handle, start_handle+num_elem-1);
     result = mdbImpl->add_entities(entity->setHandle, dum_range);
     if (MB_SUCCESS != result) return result;
+    result = mdbImpl->add_entities(mFileSet, dum_range);
+    if (MB_SUCCESS != result) return result;
 
       // set cub ids
     result = mdbImpl->tag_set_data(cubIdTag, dum_range, &int_buf[0]);
@@ -1356,7 +1368,7 @@
   for (int i = 0; i < info.numEntities; i++) {
 
       // create an entity set for this entity
-    result = instance->mdbImpl->create_meshset(MESHSET_SET, geom_headers[i].setHandle);
+    result = instance->create_set(geom_headers[i].setHandle);
     if (MB_SUCCESS != result) return result;
     
     instance->FREADI(8);
@@ -1403,7 +1415,7 @@
   for (int i = 0; i < info.numEntities; i++) {
 
       // create an entity set for this entity
-    result = instance->mdbImpl->create_meshset(MESHSET_SET, group_headers[i].setHandle);
+    result = instance->create_set(group_headers[i].setHandle);
     if (MB_SUCCESS != result) return result;
     static const char group_category[CATEGORY_TAG_SIZE] = "Group\0";
     
@@ -1453,7 +1465,7 @@
   for (int i = 0; i < info.numEntities; i++) {
 
       // create an entity set for this entity
-    result = instance->mdbImpl->create_meshset(MESHSET_SET, block_headers[i].setHandle);
+    result = instance->create_set(block_headers[i].setHandle);
     if (MB_SUCCESS != result) return result;
     static const char material_category[CATEGORY_TAG_SIZE] = "Material Set\0";
     
@@ -1533,7 +1545,7 @@
   for (int i = 0; i < info.numEntities; i++) {
 
       // create an entity set for this entity
-    result = instance->mdbImpl->create_meshset(MESHSET_SET, nodeset_headers[i].setHandle);
+    result = instance->create_set(nodeset_headers[i].setHandle);
     if (MB_SUCCESS != result) return result;
     static const char dirichlet_category[CATEGORY_TAG_SIZE] = "Dirichlet Set\0";
     
@@ -1583,7 +1595,7 @@
   for (int i = 0; i < info.numEntities; i++) {
 
       // create an entity set for this entity
-    result = instance->mdbImpl->create_meshset(MESHSET_SET, sideset_headers[i].setHandle);
+    result = instance->create_set(sideset_headers[i].setHandle);
     if (MB_SUCCESS != result) return result;
     static const char neumann_category[CATEGORY_TAG_SIZE] = "Neumann Set\0";
     
@@ -2380,6 +2392,19 @@
             << std::endl;
 }
 
+MBErrorCode Tqdcfr::create_set( MBEntityHandle& h, int flags )
+{
+  MBErrorCode rval;
+  if (!mFileSet)
+    return MB_FAILURE;
+  rval = mdbImpl->create_meshset( flags, h );
+  if (MB_SUCCESS != rval)
+    return rval;
+  rval = mdbImpl->add_entities( mFileSet, &h, 1 );
+  if (MB_SUCCESS != rval)
+    mdbImpl->delete_entities( &h, 1 );
+  return rval;
+}
 
 #ifdef TEST_TQDCFR
 #include "MBCore.hpp"
@@ -2400,8 +2425,9 @@
 
   MBCore my_impl;
   Tqdcfr my_tqd(&my_impl);
+  MBEntityHandle file_set;
 
-  MBErrorCode result = my_tqd.load_file(file, 0, 0);
+  MBErrorCode result = my_tqd.load_file(file, file_set, 0, 0);
 
   if (MB_SUCCESS == result)
     std::cout << "Success." << std::endl;

Modified: MOAB/trunk/Tqdcfr.hpp
===================================================================
--- MOAB/trunk/Tqdcfr.hpp	2007-08-24 17:16:31 UTC (rev 1251)
+++ MOAB/trunk/Tqdcfr.hpp	2007-08-24 20:05:52 UTC (rev 1252)
@@ -285,6 +285,7 @@
   
     // read cub file
   MBErrorCode load_file(const char *file_name,
+                        MBEntityHandle& file_set,
                         const int* block_list,
                         int num_blocks );
   MBErrorCode read_nodeset(ModelEntry *model,
@@ -318,9 +319,13 @@
 
   static std::string BLOCK_NODESET_OFFSET_TAG_NAME;
   static std::string BLOCK_SIDESET_OFFSET_TAG_NAME;
+  
+  MBErrorCode create_set( MBEntityHandle& h, int flags = MESHSET_SET );
 
 private:
 
+  MBEntityHandle mFileSet; // set containing read entities.
+
   MBErrorCode convert_nodesets_sidesets();
 
   MBErrorCode read_acis_records();




More information about the moab-dev mailing list