libMesh::AbaqusIO Class Reference

#include <abaqus_io.h>

Inheritance diagram for libMesh::AbaqusIO:

List of all members.

Public Member Functions

 AbaqusIO (MeshBase &mesh)
virtual ~AbaqusIO ()
virtual void read (const std::string &name)

Public Attributes

bool build_sidesets_from_nodesets

Protected Member Functions

MeshBasemesh ()
void skip_comment_lines (std::istream &in, const char comment_start)

Protected Attributes

std::vector< bool > elems_of_dimension

Private Types

typedef std::map< std::string,
std::vector< dof_id_type > > 
container_t
typedef std::map< std::string,
std::vector< std::pair
< dof_id_type, unsigned > > > 
sideset_container_t

Private Member Functions

void read_nodes ()
void read_elements (std::string upper)
std::string parse_label (std::string line, std::string label_name)
void read_ids (std::string set_name, container_t &container)
void assign_subdomain_ids ()
void read_sideset (std::string sideset_name, sideset_container_t &container)
void assign_boundary_node_ids ()
void assign_sideset_ids ()
void process_and_discard_comments ()

Private Attributes

container_t _nodeset_ids
container_t _elemset_ids
sideset_container_t _sideset_ids
std::ifstream _in
std::set< ElemType > _elem_types
std::map< dof_id_type,
dof_id_type
_abaqus_to_libmesh_elem_mapping
std::map< dof_id_type,
dof_id_type
_abaqus_to_libmesh_node_mapping
bool _already_seen_part

Detailed Description

The AbaqusIO class is a preliminary implementation for reading Abaqus mesh files in ASCII format.

Author:
John W. Peterson, 2011.

Definition at line 37 of file abaqus_io.h.


Member Typedef Documentation

typedef std::map<std::string, std::vector<dof_id_type> > libMesh::AbaqusIO::container_t [private]

The type of data structure used to store Node and Elemset IDs.

Definition at line 67 of file abaqus_io.h.

typedef std::map<std::string, std::vector<std::pair<dof_id_type, unsigned> > > libMesh::AbaqusIO::sideset_container_t [private]

Type of the data structure for storing the (elem ID, side) pairs defining sidesets. These come from the *Surface sections of the input file.

Definition at line 74 of file abaqus_io.h.


Constructor & Destructor Documentation

libMesh::AbaqusIO::AbaqusIO ( MeshBase mesh  )  [explicit]

Constructor. Takes a writeable reference to a mesh object.

Definition at line 170 of file abaqus_io.C.

00170                                        :
00171     MeshInput<MeshBase> (mesh_in),
00172     build_sidesets_from_nodesets(false),
00173     _already_seen_part(false)
00174   {
00175   }

libMesh::AbaqusIO::~AbaqusIO (  )  [virtual]

Destructor.

Definition at line 180 of file abaqus_io.C.

00181   {
00182   }


Member Function Documentation

void libMesh::AbaqusIO::assign_boundary_node_ids (  )  [private]

This function assigns boundary IDs to node sets based on the alphabetical order in which the sets are labelled in the Abaqus file. We choose the alphabetical ordering simply because Abaqus does not provide a numerical one within the file.

Definition at line 854 of file abaqus_io.C.

References _abaqus_to_libmesh_node_mapping, _nodeset_ids, libMesh::MeshBase::boundary_info, libMesh::err, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::node_ptr(), and libMesh::out.

Referenced by read().

00855   {
00856     // Get a reference to the mesh we are reading
00857     MeshBase& the_mesh = MeshInput<MeshBase>::mesh();
00858 
00859     // Iterate over the container of nodesets
00860     container_t::iterator it=_nodeset_ids.begin();
00861     for (unsigned current_id=0; it != _nodeset_ids.end(); ++it, ++current_id)
00862       {
00863         libMesh::out << "Assigning node boundary ID " << current_id << " to nodeset '"
00864                      << (*it).first
00865                      << "'." << std::endl;
00866 
00867         // Get a reference to the current vector of nodeset ID values
00868         std::vector<dof_id_type>& nodeset_ids = (*it).second;
00869 
00870         for (std::size_t i=0; i<nodeset_ids.size(); ++i)
00871           {
00872             // Map the Abaqus global node ID to the libmesh node ID
00873             dof_id_type libmesh_global_node_id = _abaqus_to_libmesh_node_mapping[nodeset_ids[i]];
00874 
00875             // Get node pointer from the mesh
00876             Node* node = the_mesh.node_ptr(libmesh_global_node_id);
00877 
00878             if (node == NULL)
00879               {
00880                 libMesh::err << "Error! Mesh returned NULL node pointer!" << std::endl;
00881                 libmesh_error();
00882               }
00883 
00884             // Add this node with the current_id (which is determined by the
00885             // alphabetical ordering of the map) to the BoundaryInfo object
00886             the_mesh.boundary_info->add_node(node, current_id);
00887           }
00888       }
00889 
00890   } // assign_boundary_node_ids()

void libMesh::AbaqusIO::assign_sideset_ids (  )  [private]

Called at the end of the read() function, assigns any sideset IDs found when reading the file to the BoundaryInfo object.

Definition at line 895 of file abaqus_io.C.

References _abaqus_to_libmesh_elem_mapping, _sideset_ids, libMesh::MeshBase::boundary_info, libMesh::MeshBase::elem(), libMesh::Utility::enum_to_string(), libMesh::err, libMesh::MeshInput< MeshBase >::mesh(), libMesh::out, and libMesh::Elem::type().

Referenced by read().

00896   {
00897     // Get a reference to the mesh we are reading
00898     MeshBase& the_mesh = MeshInput<MeshBase>::mesh();
00899 
00900     // initialize the eletypes map (eletypes is a file-global variable)
00901     init_eletypes();
00902 
00903     // Iterate over the container of sidesets
00904     sideset_container_t::iterator it=_sideset_ids.begin();
00905     for (unsigned current_id=0; it != _sideset_ids.end(); ++it, ++current_id)
00906       {
00907         libMesh::out << "Assigning sideset ID " << current_id << " to sideset '"
00908                      << (*it).first
00909                      << "'." << std::endl;
00910 
00911         // Get a reference to the current vector of nodeset ID values
00912         std::vector<std::pair<dof_id_type,unsigned> >& sideset_ids = (*it).second;
00913 
00914         for (std::size_t i=0; i<sideset_ids.size(); ++i)
00915           {
00916             // sideset_ids is a vector of pairs (elem id, side id).  Pull them out
00917             // now to make the code below more readable.
00918             dof_id_type  abaqus_elem_id = sideset_ids[i].first;
00919             unsigned abaqus_side_number = sideset_ids[i].second;
00920 
00921             // Map the Abaqus element ID to LibMesh numbering
00922             dof_id_type libmesh_elem_id = _abaqus_to_libmesh_elem_mapping[ abaqus_elem_id ];
00923 
00924             // Get pointer to that element
00925             Elem* elem = the_mesh.elem(libmesh_elem_id);
00926 
00927             // Check that the pointer returned from the Mesh is non-NULL
00928             if (elem == NULL)
00929               {
00930                 libMesh::err << "Mesh returned NULL pointer for Elem " << libmesh_elem_id << std::endl;
00931                 libmesh_error();
00932               }
00933 
00934             // Grab a reference to the element definition for this element type
00935             const ElementDefinition& eledef = eletypes[elem->type()];
00936 
00937             // If the element definition was not found, the call above would have
00938             // created one with an uninitialized struct.  Check for that here...
00939             if (eledef.abaqus_zero_based_side_id_to_libmesh_side_id.size() == 0)
00940               {
00941                 libMesh::err << "No Abaqus->LibMesh mapping information for ElemType " << Utility::enum_to_string(elem->type()) << "!" << std::endl;
00942                 libmesh_error();
00943               }
00944 
00945             // Add this node with the current_id (which is determined by the
00946             // alphabetical ordering of the map).  Side numbers in Abaqus are 1-based,
00947             // so we subtract 1 here before passing the abaqus side number to the
00948             // mapping array
00949             the_mesh.boundary_info->add_side(elem,
00950                                              eledef.abaqus_zero_based_side_id_to_libmesh_side_id[abaqus_side_number-1],
00951                                              current_id);
00952           }
00953       }
00954   } // assign_sideset_ids()

void libMesh::AbaqusIO::assign_subdomain_ids (  )  [private]

This function is called after all the elements have been read and assigns element subdomain IDs.

The IDs are simply chosen in the order in which the elset labels are stored in the map (roughly alphabetically). To make this easy on people who are planning to use Exodus output, we'll assign different geometric elements to different (but related) subdomains, i.e. assuming there are E elemsets:

Elemset 0, Geometric Type 0: ID 0 Elemset 0, Geometric Type 1: ID 0+E ... Elemset 0, Geometric Type N: ID 0+N*E -------------------------------------- Elemset 1, Geometric Type 0: ID 1 Elemset 1, Geometric Type 1: ID 1+E ... Elemset 1, Geometric Type N: ID 1+N*E etc.

Definition at line 799 of file abaqus_io.C.

References _abaqus_to_libmesh_elem_mapping, _elem_types, _elemset_ids, libMesh::MeshBase::elem(), libMesh::err, libMesh::MeshInput< MeshBase >::mesh(), libMesh::Elem::subdomain_id(), and libMesh::Elem::type().

Referenced by read().

00800   {
00801     // Get a reference to the mesh we are reading
00802     MeshBase& the_mesh = MeshInput<MeshBase>::mesh();
00803 
00804     // The number of elemsets we've found while reading
00805     std::size_t n_elemsets = _elemset_ids.size();
00806 
00807     // Fill in a temporary map with (ElemType, index) pairs based on the _elem_types set.  This
00808     // will allow us to easily look up this index in the loop below.
00809     std::map<ElemType, unsigned> elem_types_map;
00810     {
00811       unsigned ctr=0;
00812       for (std::set<ElemType>::iterator it=_elem_types.begin(); it!=_elem_types.end(); ++it)
00813         elem_types_map[*it] = ctr++;
00814     }
00815 
00816     // Loop over each Elemset and assign subdomain IDs to Mesh elements
00817     {
00818       // The elemset_id counter assigns a logical numbering to the _elemset_ids keys
00819       container_t::iterator it=_elemset_ids.begin();
00820       for (unsigned elemset_id=0; it != _elemset_ids.end(); ++it, ++elemset_id)
00821         {
00822           // Grab a reference to the vector of IDs
00823           std::vector<dof_id_type>& id_vector = (*it).second;
00824 
00825           // Loop over this vector
00826           for (std::size_t i=0; i<id_vector.size(); ++i)
00827             {
00828               // Map the id_vector[i]'th element ID (Abaqus numbering) to LibMesh numbering
00829               dof_id_type libmesh_elem_id = _abaqus_to_libmesh_elem_mapping[ id_vector[i] ];
00830 
00831               // Get pointer to that element
00832               Elem* elem = the_mesh.elem(libmesh_elem_id);
00833 
00834               if (elem == NULL)
00835                 {
00836                   libMesh::err << "Mesh returned NULL pointer for Elem " << libmesh_elem_id << std::endl;
00837                   libmesh_error();
00838                 }
00839 
00840               // Compute the proper subdomain ID, based on the formula in the
00841               // documentation for this function.
00842               subdomain_id_type computed_id = elemset_id + (elem_types_map[elem->type()] * n_elemsets);
00843 
00844               // Assign this ID to the element in question
00845               elem->subdomain_id() = computed_id;
00846             }
00847         }
00848     }
00849   } // assign_subdomain_ids()

MeshBase & libMesh::MeshInput< MeshBase >::mesh (  )  [protected, inherited]

Returns the object as a writeable reference.

Referenced by libMesh::GMVIO::_read_materials(), libMesh::GMVIO::_read_nodes(), libMesh::GMVIO::_read_one_cell(), assign_boundary_node_ids(), assign_sideset_ids(), assign_subdomain_ids(), libMesh::VTKIO::cells_to_vtk(), libMesh::GMVIO::copy_nodal_solution(), libMesh::UNVIO::element_in(), libMesh::TetGenIO::element_in(), libMesh::UNVIO::element_out(), libMesh::UNVIO::node_in(), libMesh::TetGenIO::node_in(), libMesh::UNVIO::node_out(), libMesh::VTKIO::nodes_to_vtk(), libMesh::XdrIO::read(), libMesh::VTKIO::read(), libMesh::TetGenIO::read(), libMesh::Nemesis_IO::read(), libMesh::GMVIO::read(), libMesh::ExodusII_IO::read(), read(), libMesh::LegacyXdrIO::read_ascii(), read_elements(), libMesh::UNVIO::read_implementation(), libMesh::UCDIO::read_implementation(), libMesh::LegacyXdrIO::read_mesh(), libMesh::GmshIO::read_mesh(), read_nodes(), libMesh::XdrIO::read_serialized_bcs(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::OFFIO::read_stream(), libMesh::MatlabIO::read_stream(), libMesh::XdrIO::write(), libMesh::TetGenIO::write(), libMesh::Nemesis_IO::write(), libMesh::GMVIO::write_ascii_new_impl(), libMesh::GMVIO::write_ascii_old_impl(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::UNVIO::write_implementation(), libMesh::UCDIO::write_implementation(), libMesh::LegacyXdrIO::write_mesh(), libMesh::GmshIO::write_mesh(), libMesh::UCDIO::write_nodal_data(), libMesh::Nemesis_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::XdrIO::write_parallel(), libMesh::GmshIO::write_post(), libMesh::XdrIO::write_serialized_bcs(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::LegacyXdrIO::write_soln().

std::string libMesh::AbaqusIO::parse_label ( std::string  line,
std::string  label_name 
) [private]

This function parses a label of the form foo=bar from a comma-delimited line of the form ..., foo=bar, ... The input to the function in this case would be foo, the output would be bar

Definition at line 679 of file abaqus_io.C.

References end.

Referenced by read(), and read_elements().

00680   {
00681     // Do all string comparisons in upper-case
00682     std::string upper_line(line), upper_label_name(label_name);
00683     std::transform(upper_line.begin(), upper_line.end(), upper_line.begin(), ::toupper);
00684     std::transform(upper_label_name.begin(), upper_label_name.end(), upper_label_name.begin(), ::toupper);
00685 
00686     // Get index of start of "label="
00687     size_t label_index = upper_line.find(upper_label_name + "=");
00688 
00689     if (label_index != std::string::npos)
00690       {
00691         // Location of the first comma following "label="
00692         size_t comma_index = upper_line.find(",", label_index);
00693 
00694         // Construct iterators from which to build the sub-string.
00695         // Note the +1 is to skip past the "=" which follows the label name
00696         std::string::iterator
00697           beg = line.begin() + label_name.size() + 1 + label_index,
00698           end = (comma_index == std::string::npos) ? line.end() : line.begin()+comma_index;
00699 
00700         return std::string(beg, end);
00701       }
00702 
00703     // The label index was not found, return the empty string
00704     return std::string("");
00705   } // parse_label()

void libMesh::AbaqusIO::process_and_discard_comments (  )  [private]

Any of the various sections can start with some number of lines of comments, which start with "**". This function discards any lines of comments that it finds from the stream, leaving trailing data intact.

Definition at line 958 of file abaqus_io.C.

References _in.

Referenced by read().

00959   {
00960     std::string dummy;
00961     while (true)
00962       {
00963         // We assume we are at the beginning of a line that may be
00964         // comments or may be data.  We need to only discard the line if
00965         // it begins with **, but we must avoid calling std::getline()
00966         // since there's no way to put that back.
00967         if (_in.peek() == '*')
00968           {
00969             // The first character was a star, so actually read it from the stream.
00970             _in.get();
00971 
00972             // Peek at the next character...
00973             if (_in.peek() == '*')
00974               {
00975                 // OK, second character was star also, by definition this
00976                 // line must be a comment!  Read the rest of the line and discard!
00977                 std::getline(_in, dummy);
00978 
00979                 // Debugging:
00980                 // libMesh::out << "Read comment line: " << dummy << std::endl;
00981               }
00982             else
00983               {
00984                 // The second character was _not_ a star, so put back the first star
00985                 // we pulled out so that the line can be parsed correctly by somebody
00986                 // else!
00987                 _in.unget();
00988 
00989                 // Finally, break out of the while loop, we are done parsing comments
00990                 break;
00991               }
00992           }
00993         else
00994           {
00995             // First character was not *, so this line must be data! Break out of the
00996             // while loop!
00997             break;
00998           }
00999       }
01000 
01001   } // process_and_discard_comments()

void libMesh::AbaqusIO::read ( const std::string &  name  )  [virtual]

This method implements reading a mesh from a specified file.

Implements libMesh::MeshInput< MeshBase >.

Definition at line 187 of file abaqus_io.C.

References _already_seen_part, _elemset_ids, _in, _nodeset_ids, _sideset_ids, assign_boundary_node_ids(), assign_sideset_ids(), assign_subdomain_ids(), libMesh::MeshBase::boundary_info, build_sidesets_from_nodesets, libMesh::MeshBase::clear(), libMesh::err, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::n_elem(), libMesh::MeshBase::n_nodes(), libMesh::out, parse_label(), process_and_discard_comments(), read_elements(), read_ids(), read_nodes(), and read_sideset().

00188   {
00189     // Get a reference to the mesh we are reading
00190     MeshBase& the_mesh = MeshInput<MeshBase>::mesh();
00191 
00192     // Clear any existing mesh data
00193     the_mesh.clear();
00194 
00195     // Open stream for reading
00196     _in.open(fname.c_str());
00197     libmesh_assert(_in.good());
00198 
00199     // Read file line-by-line... this is based on a set of different
00200     // test input files.  I have not looked at the full input file
00201     // specs for Abaqus.
00202     std::string s;
00203     while (true)
00204       {
00205         // Try to read something.  This may set EOF!
00206         std::getline(_in, s);
00207 
00208         if (_in)
00209           {
00210             // Process s...
00211             //
00212             // There are many sections in Abaqus files, we read some
00213             // but others are just ignored...  Some sections may occur
00214             // more than once.  For example for a hybrid grid, you
00215             // will have multiple *Element sections...
00216 
00217             // Some Abaqus files use all upper-case for section names,
00218             // so we will just convert s to uppercase
00219             std::string upper(s);
00220             std::transform(upper.begin(), upper.end(), upper.begin(), ::toupper);
00221 
00222             // 0.) Look for the "*Part" section
00223             if (upper.find("*PART") == 0)
00224               {
00225                 // libMesh::out << "Found parts section!" << std::endl;
00226 
00227                 if (_already_seen_part)
00228                   {
00229                     libMesh::err << "We currently don't support reading Abaqus files with multiple PART sections" << std::endl;
00230                     libmesh_error();
00231                   }
00232 
00233                 _already_seen_part = true;
00234               }
00235 
00236             // 1.) Look for the "*Nodes" section
00237             if (upper.find("*NODE") == 0)
00238               {
00239                 // Process any lines of comments that may be present
00240                 this->process_and_discard_comments();
00241 
00242                 // Read a block of nodes
00243                 this->read_nodes();
00244               }
00245 
00246 
00247 
00248             // 2.) Look for the "*Element" section
00249             else if (upper.find("*ELEMENT,") == 0)
00250               {
00251                 // Process any lines of comments that may be present
00252                 this->process_and_discard_comments();
00253 
00254                 // Read a block of elements
00255                 this->read_elements(upper);
00256               }
00257 
00258 
00259 
00260             // 3.) Look for a Nodeset section
00261             else if (upper.find("*NSET") == 0)
00262               {
00263                 std::string nset_name = this->parse_label(s, "nset");
00264 
00265                 // I haven't seen an unnamed elset yet, but let's detect it
00266                 // just in case...
00267                 if (nset_name == "")
00268                   {
00269                     libMesh::err << "Unnamed nset encountered!" << std::endl;
00270                     libmesh_error();
00271                   }
00272 
00273                 // Process any lines of comments that may be present
00274                 this->process_and_discard_comments();
00275 
00276                 // Read the IDs, storing them in _nodeset_ids
00277                 this->read_ids(nset_name, _nodeset_ids);
00278               } // *Nodeset
00279 
00280 
00281 
00282             // 4.) Look for an Elset section
00283             else if (upper.find("*ELSET") == 0)
00284               {
00285                 std::string elset_name = this->parse_label(s, "elset");
00286 
00287                 // I haven't seen an unnamed elset yet, but let's detect it
00288                 // just in case...
00289                 if (elset_name == "")
00290                   {
00291                     libMesh::err << "Unnamed elset encountered!" << std::endl;
00292                     libmesh_error();
00293                   }
00294 
00295                 // Debugging
00296                 // libMesh::out << "Processing ELSET: " << elset_name << std::endl;
00297 
00298                 // Process any lines of comments that may be present
00299                 this->process_and_discard_comments();
00300 
00301                 // Read the IDs, storing them in _elemset_ids
00302                 this->read_ids(elset_name, _elemset_ids);
00303               } // *Elset
00304 
00305 
00306 
00307             // 5.) Look for a Surface section.  Need to be a little
00308             // careful, since there are also "surface interaction"
00309             // sections we don't want to read here.
00310             else if (upper.find("*SURFACE,") == 0)
00311               {
00312                 // libMesh::out << "Found SURFACE section: " << s << std::endl;
00313 
00314                 // Get the name from the Name=Foo label.  This will be the map key.
00315                 std::string sideset_name = this->parse_label(s, "name");
00316 
00317                 // Print name of section we just found
00318                 // libMesh::out << "Found surface section named: " << sideset_name << std::endl;
00319 
00320                 // Process any lines of comments that may be present
00321                 this->process_and_discard_comments();
00322 
00323                 // Read the sideset IDs
00324                 this->read_sideset(sideset_name, _sideset_ids);
00325 
00326                 // Debugging: print status of most recently read sideset
00327                 // libMesh::out << "Read " << _sideset_ids[sideset_name].size() << " sides in " << sideset_name << std::endl;
00328               }
00329 
00330             continue;
00331           } // if (_in)
00332 
00333         // If !file, check to see if EOF was set.  If so, break out
00334         // of while loop.
00335         if (_in.eof())
00336           break;
00337 
00338         // If !in and !in.eof(), stream is in a bad state!
00339         libMesh::err << "Stream is bad!\n";
00340         libMesh::err << "Perhaps the file: " << fname << " does not exist?" << std::endl;
00341         libmesh_error();
00342       } // while
00343 
00344 
00345     //
00346     // We've read everything we can from the file at this point.  Now
00347     // do some more processing.
00348     //
00349     libMesh::out << "Mesh contains "
00350                  << the_mesh.n_elem()
00351                  << " elements, and "
00352                  << the_mesh.n_nodes()
00353                  << " nodes." << std::endl;
00354 
00355     // TODO: Remove these or write a function to do it?
00356 //    {
00357 //      container_t::iterator it=_nodeset_ids.begin();
00358 //      for (; it != _nodeset_ids.end(); ++it)
00359 //      {
00360 //        libMesh::out << "Node set '" << (*it).first << "' contains " << (*it).second.size() << " ID(s)." << std::endl;
00361 //      }
00362 //    }
00363 //
00364 //    {
00365 //      container_t::iterator it=_elemset_ids.begin();
00366 //      for (; it != _elemset_ids.end(); ++it)
00367 //      {
00368 //        libMesh::out << "Elem set '" << (*it).first << "' contains " << (*it).second.size() << " ID(s)." << std::endl;
00369 //      }
00370 //    }
00371 
00372 
00373     // Set element IDs based on the element sets.
00374     this->assign_subdomain_ids();
00375 
00376     // Assign nodeset values to the BoundaryInfo object
00377     this->assign_boundary_node_ids();
00378 
00379     // Assign sideset values in the BoundaryInfo object
00380     this->assign_sideset_ids();
00381 
00382     // Abaqus files only contain nodesets by default.  To be useful in
00383     // applying most types of BCs in libmesh, we will definitely need
00384     // sidesets.  So we can call the new BoundaryInfo function which
00385     // generates sidesets from nodesets.
00386     if (build_sidesets_from_nodesets)
00387       the_mesh.boundary_info->build_side_list_from_node_list();
00388 
00389   } // read()

void libMesh::AbaqusIO::read_elements ( std::string  upper  )  [private]

This function parses a block of elements in the Abaqus file. You must pass it an upper-cased version of the string declaring this section, which is typically something like: *ELEMENT, TYPE=CPS3 so that it can determine the type of elements to read.

Definition at line 480 of file abaqus_io.C.

References _abaqus_to_libmesh_elem_mapping, _abaqus_to_libmesh_node_mapping, _elem_types, _elemset_ids, _in, libMesh::MeshBase::add_elem(), libMesh::Elem::build(), libMesh::err, libMeshEnums::HEX20, libMeshEnums::HEX8, libMesh::DofObject::id(), libMeshEnums::INVALID_ELEM, libMesh::MeshInput< MeshBase >::mesh(), libMesh::MeshBase::node_ptr(), parse_label(), libMeshEnums::PRISM15, libMeshEnums::PRISM6, libMeshEnums::QUAD4, libMesh::MeshBase::set_mesh_dimension(), libMesh::Elem::set_node(), libMeshEnums::TET10, libMeshEnums::TET4, and libMeshEnums::TRI3.

Referenced by read().

00481   {
00482     // Some *Element sections also specify an Elset name on the same line.
00483     // Look for one here.
00484     std::string elset_name = this->parse_label(upper, "elset");
00485 
00486     // Get a reference to the mesh we are reading
00487     MeshBase& the_mesh = MeshInput<MeshBase>::mesh();
00488 
00489     // initialize the eletypes map (eletypes is a file-global variable)
00490     init_eletypes();
00491 
00492     ElemType elem_type = INVALID_ELEM;
00493     unsigned n_nodes_per_elem = 0;
00494 
00495     // Within s, we should have "type=XXXX"
00496     if (upper.find("CPE4") != std::string::npos ||
00497         upper.find("CPS4") != std::string::npos)
00498       {
00499         elem_type = QUAD4;
00500         n_nodes_per_elem = 4;
00501         the_mesh.set_mesh_dimension(2);
00502       }
00503     else if (upper.find("CPS3") != std::string::npos)
00504       {
00505         elem_type = TRI3;
00506         n_nodes_per_elem = 3;
00507         the_mesh.set_mesh_dimension(2);
00508       }
00509     else if (upper.find("C3D8") != std::string::npos)
00510       {
00511         elem_type = HEX8;
00512         n_nodes_per_elem = 8;
00513         the_mesh.set_mesh_dimension(3);
00514       }
00515     else if (upper.find("C3D4") != std::string::npos)
00516       {
00517         elem_type = TET4;
00518         n_nodes_per_elem = 4;
00519         the_mesh.set_mesh_dimension(3);
00520       }
00521     else if (upper.find("C3D20") != std::string::npos)
00522       {
00523         elem_type = HEX20;
00524         n_nodes_per_elem = 20;
00525         the_mesh.set_mesh_dimension(3);
00526       }
00527     else if (upper.find("C3D6") != std::string::npos)
00528       {
00529         elem_type = PRISM6;
00530         n_nodes_per_elem = 6;
00531         the_mesh.set_mesh_dimension(3);
00532       }
00533     else if (upper.find("C3D15") != std::string::npos)
00534       {
00535         elem_type = PRISM15;
00536         n_nodes_per_elem = 15;
00537         the_mesh.set_mesh_dimension(3);
00538       }
00539     else if (upper.find("C3D10") != std::string::npos)
00540       {
00541         elem_type = TET10;
00542         n_nodes_per_elem = 10;
00543         the_mesh.set_mesh_dimension(3);
00544       }
00545     else
00546       {
00547         libMesh::err << "Unrecognized element type: " << upper << std::endl;
00548         libmesh_error();
00549       }
00550 
00551     // Insert the elem type we detected into the set of all elem types for this mesh
00552     _elem_types.insert(elem_type);
00553 
00554     // For reading in line endings
00555     std::string dummy;
00556 
00557     // Grab a reference to the element definition for this element type
00558     const ElementDefinition& eledef = eletypes[elem_type];
00559 
00560     // If the element definition was not found, the call above would have
00561     // created one with an uninitialized struct.  Check for that here...
00562     if (eledef.abaqus_zero_based_node_id_to_libmesh_node_id.size() == 0)
00563       {
00564         // libMesh::err << "No Abaqus->LibMesh mapping information for ElemType "
00565         //           << Utility::enum_to_string(elem_type)
00566         //           << "!"
00567         //           << std::endl;
00568         libmesh_error();
00569       }
00570 
00571     // We will read elements until the next line begins with *, since that will be the
00572     // next section.
00573     while (_in.peek() != '*' && _in.peek() != EOF)
00574       {
00575         // Read the element ID, it is the first number on each line.  It is
00576         // followed by a comma, so read that also.  We will need this ID later
00577         // when we try to assign subdomain IDs
00578         dof_id_type abaqus_elem_id = 0;
00579         char c;
00580         _in >> abaqus_elem_id >> c;
00581 
00582         // Debugging:
00583         // libMesh::out << "Reading data for element " << abaqus_elem_id << std::endl;
00584 
00585         // Add an element of the appropriate type to the Mesh.
00586         Elem* elem = the_mesh.add_elem(Elem::build(elem_type).release());
00587 
00588         // Associate the ID returned from libmesh with the abaqus element ID
00589         //_libmesh_to_abaqus_elem_mapping[elem->id()] = abaqus_elem_id;
00590         _abaqus_to_libmesh_elem_mapping[abaqus_elem_id] = elem->id();
00591 
00592         // The count of the total number of IDs read for the current element.
00593         unsigned id_count=0;
00594 
00595         // Continue reading line-by-line until we have read enough nodes for this element
00596         while (id_count < n_nodes_per_elem)
00597           {
00598             // Read entire line (up to carriage return) of comma-separated values
00599             std::string csv_line;
00600             std::getline(_in, csv_line);
00601 
00602             // Create a stream object out of the current line
00603             std::stringstream line_stream(csv_line);
00604 
00605             // Process the comma-separated values
00606             std::string cell;
00607             while (std::getline(line_stream, cell, ','))
00608               {
00609                 // FIXME: factor out this strtol stuff into a utility function.
00610                 char* endptr;
00611                 long abaqus_global_node_id = std::strtol(cell.c_str(), &endptr, /*base=*/10);
00612 
00613                 if (abaqus_global_node_id!=0 || cell.c_str() != endptr)
00614                   {
00615                     // Use the global node number mapping to determine the corresponding libmesh global node id
00616                     dof_id_type libmesh_global_node_id = _abaqus_to_libmesh_node_mapping[abaqus_global_node_id];
00617 
00618                     // Grab the node pointer from the mesh for this ID
00619                     Node* node = the_mesh.node_ptr(libmesh_global_node_id);
00620 
00621                     // Debugging:
00622                     // libMesh::out << "Assigning global node id: " << abaqus_global_node_id
00623                     //              << "(Abaqus), " << node->id() << "(LibMesh)" << std::endl;
00624 
00625                     // If node_ptr() returns NULL, it may mean we have not yet read the
00626                     // *Nodes section, though I assumed that always came before the *Elements section...
00627                     if (node == NULL)
00628                       {
00629                         libMesh::err << "Error!  Mesh returned NULL Node pointer.\n";
00630                         libMesh::err << "Either no node exists with ID " << libmesh_global_node_id
00631                                      << " or perhaps this input file has *Elements defined before *Nodes?" << std::endl;
00632                         libmesh_error();
00633                       }
00634 
00635                     // Note: id_count is the zero-based abaqus (elem local) node index.  We therefore map
00636                     // it to a libmesh elem local node index using the element definition map
00637                     unsigned libmesh_elem_local_node_id =
00638                       eledef.abaqus_zero_based_node_id_to_libmesh_node_id[id_count];
00639 
00640                     // Set this node pointer within the element.
00641                     elem->set_node(libmesh_elem_local_node_id) = node;
00642 
00643                     // Debugging:
00644                     // libMesh::out << "Setting elem " << elem->id()
00645                     //              << ", local node " << libmesh_elem_local_node_id
00646                     //              << " to global node " << node->id() << std::endl;
00647 
00648                     // Increment the count of IDs read for this element
00649                     id_count++;
00650                   } // end if strtol success
00651               } // end while getline(',')
00652           } // end while (id_count)
00653 
00654         // Ensure that we read *exactly* as many nodes as we were expecting to, no more.
00655         if (id_count != n_nodes_per_elem)
00656           {
00657             libMesh::err << "Error: Needed to read "
00658                          << n_nodes_per_elem
00659                          << " nodes, but read "
00660                          << id_count
00661                          << " instead!" << std::endl;
00662             libmesh_error();
00663           }
00664 
00665         // If we are recording Elset IDs, add this element to the correct set for later processing.
00666         // Make sure to add it with the Abaqus ID, not the libmesh one!
00667         if (elset_name != "")
00668           {
00669             // Debugging:
00670             // libMesh::out << "Adding Elem " << abaqus_elem_id << " to Elmset " << elset_name << std::endl;
00671             _elemset_ids[elset_name].push_back(abaqus_elem_id);
00672           }
00673       } // end while (peek)
00674   } // read_elements()

void libMesh::AbaqusIO::read_ids ( std::string  set_name,
container_t container 
) [private]

This function reads all the IDs for the current node or element set of the given name, storing them in the passed map using the name as key.

Definition at line 710 of file abaqus_io.C.

References _in.

Referenced by read().

00711   {
00712     // Debugging
00713     // libMesh::out << "Reading ids for set: " << set_name << std::endl;
00714 
00715     // Grab a reference to a vector that will hold all the IDs
00716     std::vector<dof_id_type>& id_storage = container[set_name];
00717 
00718     // Read until the start of another section is detected, or EOF is encountered
00719     while (_in.peek() != '*' && _in.peek() != EOF)
00720       {
00721         // Read entire comma-separated line into a string
00722         std::string csv_line;
00723         std::getline(_in, csv_line);
00724 
00725         // On that line, use std::getline again to parse each
00726         // comma-separated entry.
00727         std::string cell;
00728         std::stringstream line_stream(csv_line);
00729         while (std::getline(line_stream, cell, ','))
00730           {
00731             // If no conversion can be performed by strtol, 0 is returned.
00732             //
00733             // If endptr is not NULL, strtol() stores the address of the
00734             // first invalid character in *endptr.  If there were no
00735             // digits at all, however, strtol() stores the original
00736             // value of str in *endptr.
00737             char* endptr;
00738 
00739             // FIXME - this needs to be updated for 64-bit inputs
00740             long id = std::strtol(cell.c_str(), &endptr, /*base=*/10);
00741 
00742             // Note that lists of comma-separated values in abaqus also
00743             // *end* with a comma, so the last call to getline on a given
00744             // line will get an empty string, which we must detect.
00745             if (id!=0 || cell.c_str() != endptr)
00746               {
00747                 // Debugging
00748                 // libMesh::out << "Read id: " << id << std::endl;
00749 
00750                 // 'cell' is now a string with an integer id in it
00751                 id_storage.push_back( id );
00752               }
00753           }
00754       }
00755 
00756     // Status message
00757     // libMesh::out << "Read " << id_storage.size() << " ID(s) for the set " << set_name << std::endl;
00758   } // read_ids()

void libMesh::AbaqusIO::read_nodes (  )  [private]

This function parses a block of nodes in the Abaqus file once such a block has been found.

Definition at line 397 of file abaqus_io.C.

References _abaqus_to_libmesh_node_mapping, _in, libMesh::MeshBase::add_point(), libMesh::MeshInput< MeshBase >::mesh(), and libMesh::Real.

Referenced by read().

00398   {
00399     // Get a reference to the mesh we are reading
00400     MeshBase& the_mesh = MeshInput<MeshBase>::mesh();
00401 
00402     // Debugging: print node count
00403     // libMesh::out << "Before read_nodes(), mesh contains "
00404     //           << the_mesh.n_elem()
00405     //           << " elements, and "
00406     //           << the_mesh.n_nodes()
00407     //           << " nodes." << std::endl;
00408 
00409     // In the input file I have, Abaqus neither tells what
00410     // the mesh dimension is nor how many nodes it has...
00411 
00412     // The node line format is:
00413     // id, x, y, z
00414     // and you do have to parse out the commas.
00415     // The z-coordinate will only be present for 3D meshes
00416 
00417     // Temporary variables for parsing lines of text
00418     dof_id_type abaqus_node_id=0;
00419     Real x=0, y=0, z=0;
00420     char c;
00421     std::string dummy;
00422 
00423     // Defines the sequential node numbering used by libmesh
00424     dof_id_type libmesh_node_id = 0;
00425 
00426     // We will read nodes until the next line begins with *, since that will be the
00427     // next section.
00428     // TODO: Is Abaqus guaranteed to start the line with '*' or can there be leading white space?
00429     while (_in.peek() != '*' && _in.peek() != EOF)
00430       {
00431         // Re-Initialize variables to be read in from file
00432         abaqus_node_id=0;
00433         x = y = z = 0.;
00434 
00435         // Note: we assume *at least* 2D points here, should we worry about
00436         // trying to read 1D Abaqus meshes?
00437         _in >> abaqus_node_id >> c >> x >> c >> y;
00438 
00439         // Peek at the next character.  If it is a comma, then there is another
00440         // value to read!
00441         if (_in.peek() == ',')
00442           _in >> c >> z;
00443 
00444         // Debugging: Print what we just read in.
00445         // libMesh::out << "node_id=" << node_id
00446         //           << ", x=" << x
00447         //           << ", y=" << y
00448         //           << ", z=" << z
00449         //           << std::endl;
00450 
00451         // Read (and discard) the rest of the line, including the newline.
00452         // This is required so that our 'peek()' at the beginning of this
00453         // loop doesn't read the newline character, for example.
00454         std::getline(_in, dummy);
00455 
00456         // Set up the abaqus -> libmesh node mapping.  This is usually just the
00457         // "off-by-one" map.
00458         _abaqus_to_libmesh_node_mapping[abaqus_node_id] = libmesh_node_id;
00459 
00460         // Add the point to the mesh using libmesh's numbering,
00461         // and post-increment the libmesh node counter.
00462         the_mesh.add_point(Point(x,y,z), libmesh_node_id++);
00463       } // while
00464 
00465     // Debugging: print node count.  Note: in serial mesh, this count may
00466     // be off by one, since Abaqus uses one-based numbering, and libmesh
00467     // just reports the length of its _nodes vector for the number of nodes.
00468     // libMesh::out << "After read_nodes(), mesh contains "
00469     //              << the_mesh.n_elem()
00470     //              << " elements, and "
00471     //              << the_mesh.n_nodes()
00472     //              << " nodes." << std::endl;
00473 
00474   } // read_nodes()

void libMesh::AbaqusIO::read_sideset ( std::string  sideset_name,
sideset_container_t container 
) [private]

This function reads a sideset from the input file. This is defined by a "*Surface" section in the file, and then a list of element ID and side IDs for the set.

Definition at line 763 of file abaqus_io.C.

References _in, and elem_id.

Referenced by read().

00764   {
00765     // Grab a reference to a vector that will hold all the IDs
00766     std::vector<std::pair<dof_id_type, unsigned> >& id_storage = container[sideset_name];
00767 
00768     // Variables for storing values read in from file
00769     dof_id_type elem_id=0; 
00770     unsigned side_id=0;
00771     char c;
00772     std::string dummy;
00773 
00774     // Read until the start of another section is detected, or EOF is encountered
00775     while (_in.peek() != '*' && _in.peek() != EOF)
00776       {
00777         // The strings are of the form: "391, S2"
00778 
00779         // Read the element ID and the leading comma
00780         _in >> elem_id >> c;
00781 
00782         // Read another character (the 'S') and finally the side ID
00783         _in >> c >> side_id;
00784 
00785         // Debugging: print status
00786         // std::cout << "Read elem_id=" << elem_id << ", side_id=" << side_id << std::endl;
00787 
00788         // Store this pair of data in the vector
00789         id_storage.push_back( std::make_pair(elem_id, side_id) );
00790 
00791         // Extract remaining characters on line including newline
00792         std::getline(_in, dummy);
00793       } // while
00794   }

void libMesh::MeshInput< MeshBase >::skip_comment_lines ( std::istream &  in,
const char  comment_start 
) [protected, inherited]

Reads input from in, skipping all the lines that start with the character comment_start.

Referenced by libMesh::TetGenIO::read(), and libMesh::UCDIO::read_implementation().


Member Data Documentation

Map from libmesh element number -> abaqus element number, and the converse.

Definition at line 184 of file abaqus_io.h.

Referenced by assign_sideset_ids(), assign_subdomain_ids(), and read_elements().

Map from abaqus node number -> sequential, 0-based libmesh node numbering. Note that in every Abaqus file I've ever seen the node numbers were 1-based, sequential, and all in order, so that this map is probably overkill. Nevertheless, it is the most general solution in case we come across a weird Abaqus file some day.

Definition at line 193 of file abaqus_io.h.

Referenced by assign_boundary_node_ids(), read_elements(), and read_nodes().

This flag gets set to true after the first "*PART" section we see. If it is still true when we see a second PART section, we will print an error message... we don't currently handle input files with multiple parts.

Definition at line 201 of file abaqus_io.h.

Referenced by read().

std::set<ElemType> libMesh::AbaqusIO::_elem_types [private]

A set of the different geometric element types detected when reading the mesh.

Definition at line 177 of file abaqus_io.h.

Referenced by assign_subdomain_ids(), and read_elements().

Definition at line 165 of file abaqus_io.h.

Referenced by assign_subdomain_ids(), read(), and read_elements().

std::ifstream libMesh::AbaqusIO::_in [private]

Stream object used to interact with the file

Definition at line 171 of file abaqus_io.h.

Referenced by process_and_discard_comments(), read(), read_elements(), read_ids(), read_nodes(), and read_sideset().

Abaqus writes nodesets and elemsets with labels. As we read them in, we'll use these maps to provide a natural ordering for them.

Definition at line 164 of file abaqus_io.h.

Referenced by assign_boundary_node_ids(), and read().

Definition at line 166 of file abaqus_io.h.

Referenced by assign_sideset_ids(), and read().

Default false. Abaqus files have only nodesets in them by default. Set this flag to true if you want libmesh to automatically generate sidesets from Abaqus' nodesets.

Definition at line 61 of file abaqus_io.h.

Referenced by read().


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

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:06 UTC

Hosted By:
SourceForge.net Logo