19 #include "libmesh/abaqus_io.h" 20 #include "libmesh/point.h" 21 #include "libmesh/elem.h" 22 #include "libmesh/enum_to_string.h" 23 #include "libmesh/boundary_info.h" 24 #include "libmesh/utility.h" 27 #ifdef LIBMESH_HAVE_GZSTREAM 28 # include "gzstream.h" 32 #include <unordered_map> 47 struct ElementDefinition
50 std::vector<unsigned> abaqus_zero_based_node_id_to_libmesh_node_id;
53 std::vector<unsigned short> abaqus_zero_based_side_id_to_libmesh_side_id;
59 std::map<ElemType, ElementDefinition> eletypes;
64 void add_eletype_entry(
ElemType libmesh_elem_type,
65 const unsigned * node_map,
66 unsigned node_map_size,
67 const unsigned short * side_map,
68 unsigned side_map_size)
71 ElementDefinition & map_entry = eletypes[libmesh_elem_type];
78 (node_map, node_map+node_map_size).swap
79 (map_entry.abaqus_zero_based_node_id_to_libmesh_node_id);
81 std::vector<unsigned short>
82 (side_map, side_map+side_map_size).swap
83 (map_entry.abaqus_zero_based_side_id_to_libmesh_side_id);
100 const unsigned int node_map[] = {0};
101 const unsigned short side_map[] = {0};
102 add_eletype_entry(
NODEELEM, node_map, 1, side_map, 1);
107 const unsigned int node_map[] = {0,1};
108 const unsigned short side_map[] = {0,1};
109 add_eletype_entry(
EDGE2, node_map, 2, side_map, 2);
114 const unsigned int node_map[] = {0,1,2};
115 const unsigned short side_map[] = {0,1,2};
116 add_eletype_entry(
TRI3, node_map, 3, side_map, 3);
121 const unsigned int node_map[] = {0,1,2,3,4,5};
122 const unsigned short side_map[] = {0,1,2};
123 add_eletype_entry(
TRI6, node_map, 6, side_map, 3);
128 const unsigned int node_map[] = {0,1,2,3};
129 const unsigned short side_map[] = {0,1,2,3};
130 add_eletype_entry(
QUAD4, node_map, 4, side_map, 4);
135 const unsigned int node_map[] = {0,1,2,3,4,5,6,7};
136 const unsigned short side_map[] = {0,1,2,3};
137 add_eletype_entry(
QUAD8, node_map, 8, side_map, 4);
142 const unsigned int node_map[] = {0,1,2,3};
143 const unsigned short side_map[] = {0,1,2,3};
144 add_eletype_entry(
TET4, node_map, 4, side_map, 4);
149 const unsigned int node_map[] = {0,1,2,3,4,5,6,7,8,9};
150 const unsigned short side_map[] = {0,1,2,3};
151 add_eletype_entry(
TET10, node_map, 10, side_map, 4);
156 const unsigned int node_map[] = {0,1,2,3,4,5,6,7};
157 const unsigned short side_map[] = {0,5,1,2,3,4};
158 add_eletype_entry(
HEX8, node_map, 8, side_map, 6);
163 const unsigned int node_map[] =
164 {0,1,2,3,4,5,6,7,8,9,10,11,16,17,18,19,12,13,14,15};
165 const unsigned short side_map[] =
167 add_eletype_entry(
HEX20, node_map, 20, side_map, 6);
172 const unsigned int node_map[] =
173 {0,1,2,3,4,5,6,7,8,9,10,11,16,17,18,19,12,13,14,15,26,20,25,21,22,23,24};
174 const unsigned short side_map[] =
176 add_eletype_entry(
HEX27, node_map, 27, side_map, 6);
181 const unsigned int node_map[] = {0,1,2,3,4,5};
182 const unsigned short side_map[] = {0,4,1,2,3};
183 add_eletype_entry(
PRISM6, node_map, 6, side_map, 5);
188 const unsigned int node_map[] =
189 {0,1,2,3,4,5,6,7,8,12,13,14,9,10,11};
190 const unsigned short side_map[] =
192 add_eletype_entry(
PRISM15, node_map, 15, side_map, 5);
197 const unsigned int node_map[] =
198 {0,1,2,3,4,5,6,7,8,12,13,14,9,10,11,15,16,17};
199 const unsigned short side_map[] =
201 add_eletype_entry(
PRISM18, node_map, 18, side_map, 5);
216 build_sidesets_from_nodesets(false),
217 _already_seen_part(false)
236 const bool gzipped_file = (fname.rfind(
".gz") == fname.size() - 3);
240 #ifdef LIBMESH_HAVE_GZSTREAM 241 auto inf = std::make_unique<igzstream>();
243 inf->open(fname.c_str(), std::ios::in);
244 _in = std::move(inf);
246 libmesh_error_msg(
"ERROR: need gzstream to handle .gz files!!!");
251 auto inf = std::make_unique<std::ifstream>();
255 inf->open(new_name.c_str(), std::ios::in);
257 _in = std::move(inf);
272 std::getline(*
_in, s);
285 std::string upper(s);
286 std::transform(upper.begin(), upper.end(), upper.begin(), ::toupper);
289 if (upper.find(
"*PART") ==
static_cast<std::string::size_type
>(0))
293 "We currently don't support reading Abaqus files with multiple PART sections");
299 if (upper.find(
"*NODE") ==
static_cast<std::string::size_type
>(0))
306 bool skip_this_section =
false;
307 std::vector<std::string> keywords_to_ignore = {
"OUTPUT",
"PRINT",
"FILE",
"RESPONSE"};
308 for (
auto & affix : keywords_to_ignore) {
309 std::stringstream keyword;
310 keyword <<
"*NODE " << affix;
311 if (upper.find(keyword.str()) == static_cast<std::string::size_type>(0))
312 skip_this_section =
true;
314 if (skip_this_section)
319 std::string nset_name = this->
parse_label(s,
"nset");
331 else if (upper.find(
"*ELEMENT,") ==
static_cast<std::string::size_type
>(0))
338 bool skip_this_section =
false;
339 std::vector<std::string> keywords_to_ignore = {
"OUTPUT",
"MATRIX",
"PROPERTIES",
"RESPONSE"};
340 for (
auto & affix : keywords_to_ignore) {
341 std::stringstream keyword;
342 keyword <<
"*ELEMENT " << affix;
343 if (upper.find(keyword.str()) == static_cast<std::string::size_type>(0))
344 skip_this_section =
true;
346 if (skip_this_section)
351 std::string elset_name = this->
parse_label(s,
"elset");
363 else if (upper.find(
"*NSET") ==
static_cast<std::string::size_type
>(0))
365 std::string nset_name = this->
parse_label(s,
"nset");
369 libmesh_error_msg_if(nset_name ==
"",
"Unnamed nset encountered!");
388 else if (upper.find(
"*ELSET") ==
static_cast<std::string::size_type
>(0))
390 std::string elset_name = this->
parse_label(s,
"elset");
394 libmesh_error_msg_if(elset_name ==
"",
"Unnamed elset encountered!");
415 else if (upper.find(
"*SURFACE,") ==
static_cast<std::string::size_type
>(0))
418 std::string sideset_name = this->
parse_label(s,
"name");
423 std::string sideset_type = this->
parse_label(s,
"type");
444 libmesh_error_msg(
"Stream is bad! Perhaps the file: " << fname <<
" does not exist?");
472 for (
auto & elem : the_mesh.element_ptr_range())
473 if (elem->dim() < max_dim)
504 std::vector<dof_id_type> * id_storage =
nullptr;
511 while (
_in->peek() !=
'*' &&
_in->peek() != EOF)
515 std::getline(*
_in, line);
524 std::stringstream ss(line);
532 ss >> abaqus_node_id >> c >> x >> c >> y;
536 if (ss.peek() ==
',')
541 id_storage->push_back(abaqus_node_id);
544 libmesh_error_msg_if(abaqus_node_id < 1,
545 "Invalid Abaqus node ID found");
546 const dof_id_type libmesh_node_id = abaqus_node_id-1;
549 "Duplicate Abaqus node ID found");
570 unsigned n_nodes_per_elem = 0;
573 if (upper.find(
"T3D2") != std::string::npos ||
574 upper.find(
"B31") != std::string::npos)
577 n_nodes_per_elem = 2;
580 else if (upper.find(
"B32") != std::string::npos)
583 n_nodes_per_elem = 3;
586 else if (upper.find(
"S3") != std::string::npos ||
587 upper.find(
"CPE3") != std::string::npos ||
588 upper.find(
"2D3") != std::string::npos)
591 n_nodes_per_elem = 3;
594 else if (upper.find(
"CPE4") != std::string::npos ||
595 upper.find(
"S4") != std::string::npos ||
596 upper.find(
"CPEG4") != std::string::npos ||
597 upper.find(
"2D4") != std::string::npos)
600 n_nodes_per_elem = 4;
603 else if (upper.find(
"CPE6") != std::string::npos ||
604 upper.find(
"S6") != std::string::npos ||
605 upper.find(
"CPEG6") != std::string::npos ||
606 upper.find(
"2D6") != std::string::npos)
609 n_nodes_per_elem = 6;
612 else if (upper.find(
"CPE8") != std::string::npos ||
613 upper.find(
"S8") != std::string::npos ||
614 upper.find(
"CPEG8") != std::string::npos ||
615 upper.find(
"2D8") != std::string::npos)
618 n_nodes_per_elem = 8;
621 else if (upper.find(
"3D8") != std::string::npos)
624 n_nodes_per_elem = 8;
627 else if (upper.find(
"3D4") != std::string::npos)
630 n_nodes_per_elem = 4;
633 else if (upper.find(
"3D20") != std::string::npos)
636 n_nodes_per_elem = 20;
639 else if (upper.find(
"3D27") != std::string::npos)
642 n_nodes_per_elem = 27;
645 else if (upper.find(
"3D6") != std::string::npos)
648 n_nodes_per_elem = 6;
651 else if (upper.find(
"3D15") != std::string::npos)
654 n_nodes_per_elem = 15;
657 else if (upper.find(
"3D10") != std::string::npos)
660 n_nodes_per_elem = 10;
663 else if (upper.find(
"MASS") != std::string::npos)
668 n_nodes_per_elem = 1;
672 libmesh_error_msg(
"Unrecognized element type: " << upper);
678 const ElementDefinition & eledef = eletypes[elem_type];
683 (eledef.abaqus_zero_based_node_id_to_libmesh_node_id.size() == 0,
684 "No Abaqus->LibMesh mapping information for ElemType " 689 while (
_in->peek() !=
'*' &&
_in->peek() != EOF)
696 *
_in >> abaqus_elem_id >> c;
700 std::unique_ptr<Elem> new_elem =
Elem::build(elem_type);
701 new_elem->set_id() = abaqus_elem_id-1;
709 while (id_count < n_nodes_per_elem)
712 std::string csv_line;
713 std::getline(*
_in, csv_line);
716 std::stringstream line_stream(csv_line);
720 while (std::getline(line_stream, cell,
','))
728 libmesh_error_msg_if(abaqus_global_node_id < 1,
729 "Invalid Abaqus node ID found");
730 const dof_id_type libmesh_global_node_id = abaqus_global_node_id-1;
733 Node * node = the_mesh.
node_ptr(libmesh_global_node_id);
739 "Error! Mesh::node_ptr() returned nullptr. Either no node exists with ID " 740 << libmesh_global_node_id
741 <<
" or perhaps this input file has *Elements defined before *Nodes?");
745 unsigned libmesh_elem_local_node_id =
746 eledef.abaqus_zero_based_node_id_to_libmesh_node_id[id_count];
749 elem->
set_node(libmesh_elem_local_node_id, node);
759 (id_count != n_nodes_per_elem,
760 "Error: Needed to read " 762 <<
" nodes, but read " 768 if (elset_name !=
"")
793 upper_label_name(label_name);
794 std::transform(upper_line.begin(), upper_line.end(), upper_line.begin(), ::toupper);
795 std::transform(upper_label_name.begin(), upper_label_name.end(), upper_label_name.begin(), ::toupper);
798 size_t label_index = upper_line.find(upper_label_name +
"=");
800 if (label_index != std::string::npos)
803 size_t comma_index = upper_line.find(
",", label_index);
807 std::string::iterator
808 beg = line.begin() + label_name.size() + 1 + label_index,
809 end = (comma_index == std::string::npos) ? line.end() : line.begin() + comma_index;
811 return std::string(beg, end);
815 return std::string(
"");
828 std::stringstream line_stream(upper);
829 while (std::getline(line_stream, cell,
','))
830 if (cell ==
"GENERATE")
841 std::vector<dof_id_type> & id_storage = container[set_name];
844 while (
_in->peek() !=
'*' &&
_in->peek() != EOF)
847 std::string csv_line;
848 std::getline(*
_in, csv_line);
853 std::stringstream line_stream(csv_line);
854 while (std::getline(line_stream, cell,
','))
863 id_storage.push_back(
id);
874 std::vector<dof_id_type> & id_storage = container[set_name];
879 while (
_in->peek() !=
'*' &&
_in->peek() != EOF)
882 std::string csv_line;
883 std::getline(*
_in, csv_line);
892 std::stringstream line_stream(csv_line);
893 line_stream >> start >> c >> end >> c >> stride;
900 for (
dof_id_type current = start; current <= end; current += stride)
901 id_storage.push_back(current);
909 const std::string & sideset_type,
913 std::vector<std::pair<dof_id_type, unsigned>> & id_storage = container[sideset_name];
918 std::string elem_id_or_set, dummy;
921 while (
_in->peek() !=
'*' &&
_in->peek() != EOF)
924 std::getline(*
_in, elem_id_or_set,
',');
931 if (sideset_type ==
"NODE")
936 libMesh::out <<
"Warning: skipped creating a sideset from a " 937 <<
"nodeset, this is not yet supported." 958 id_storage.emplace_back(elem_id,
side_id);
964 const auto & vec = libmesh_map_find(
_elemset_ids, elem_id_or_set);
965 for (
const auto & elem_id_in_elset : vec)
966 id_storage.emplace_back(elem_id_in_elset,
side_id);
972 std::getline(*
_in, dummy);
989 std::map<ElemType, unsigned> elem_types_map;
993 elem_types_map[type] = ctr++;
1003 for (
unsigned elemset_id=0; it != end; ++it, ++elemset_id)
1006 std::vector<dof_id_type> & id_vector = it->second;
1009 for (
const auto &
id : id_vector)
1012 libmesh_error_msg_if(
id < 1,
1013 "Invalid Abaqus element ID found");
1025 if (elem.
dim() < max_dim)
1031 (elemset_id + (elem_types_map[elem.
type()] * n_elemsets));
1063 std::vector<dof_id_type> & nodeset_ids = it->second;
1065 for (
const auto &
id : nodeset_ids)
1068 libmesh_error_msg_if(
id < 1,
1069 "Invalid Abaqus node ID found");
1073 Node * node = the_mesh.
node_ptr(libmesh_global_node_id);
1075 libmesh_error_msg_if(node ==
nullptr,
1076 "Error! Mesh::node_ptr() returned nullptr!");
1105 std::vector<std::pair<dof_id_type,unsigned>> & sideset_ids = it->second;
1107 for (
const auto & [abaqus_elem_id, abaqus_side_number] : sideset_ids)
1110 libmesh_error_msg_if(abaqus_elem_id < 1,
1111 "Invalid Abaqus element ID found");
1112 const dof_id_type libmesh_elem_id = abaqus_elem_id-1;
1118 const ElementDefinition & eledef = eletypes[elem.
type()];
1122 libmesh_error_msg_if(eledef.abaqus_zero_based_side_id_to_libmesh_side_id.size() == 0,
1123 "No Abaqus->LibMesh mapping information for ElemType " 1132 eledef.abaqus_zero_based_side_id_to_libmesh_side_id[abaqus_side_number-1],
1153 std::unordered_multimap<dof_id_type, std::pair<Elem *, boundary_id_type>> provide_bcs;
1162 std::vector<dof_id_type> & id_vector = it->second;
1165 for (
const auto &
id : id_vector)
1168 libmesh_error_msg_if(
id < 1,
1169 "Invalid Abaqus element ID found");
1178 if (elem.
dim() == max_dim)
1184 if (elem.
dim() == 0)
1191 libmesh_error_msg_if(elem.
dim()+1 != max_dim,
1192 "ERROR: Expected boundary element of dimension " << max_dim-1 <<
" but got " << elem.
dim());
1195 provide_bcs.emplace(elem.
key(), std::make_pair(&elem, elemset_id));
1206 for (
auto & elem : the_mesh.active_element_ptr_range())
1207 if (elem->dim() == max_dim)
1208 for (
auto sn : elem->side_index_range())
1216 auto bounds = provide_bcs.equal_range (elem->key(sn));
1219 for (
const auto & pr :
as_range(bounds))
1222 std::unique_ptr<Elem> side (elem->build_side_ptr(sn));
1225 Elem * lower_dim_elem = pr.second.first;
1229 if (*lower_dim_elem == *side)
1247 if (
_in->peek() ==
'*')
1253 if (
_in->peek() ==
'*')
1257 std::getline(*
_in, dummy);
1283 unsigned char max_dim = 0;
1285 unsigned char elem_dimensions_size = cast_int<unsigned char>
1288 for (
unsigned char i=1; i<elem_dimensions_size; ++i)
1298 output = cast_int<dof_id_type>(std::strtol(input.c_str(), &endptr, 10));
1300 return (output != 0 || endptr != input.c_str());
1305 line.erase(std::remove_if(line.begin(), line.end(),
1306 [](
unsigned char const c)
1307 {
return std::isspace(c); }),
virtual void read(const std::string &name) override
This method implements reading a mesh from a specified file.
void process_and_discard_comments()
Any of the various sections can start with some number of lines of comments, which start with "**"...
std::size_t n_boundary_conds() const
ElemType
Defines an enum for geometric element types.
virtual Node *& set_node(const unsigned int i)
A Node is like a Point, but with more information.
std::string & nodeset_name(boundary_id_type id)
bool _already_seen_part
This flag gets set to true after the first "*PART" section we see.
std::map< std::string, std::vector< dof_id_type > > container_t
The type of data structure used to store Node and Elemset IDs.
sideset_container_t _sideset_ids
virtual ~AbaqusIO()
Destructor.
const boundary_id_type side_id
void strip_ws(std::string &line) const
Removes all whitespace characters from "line".
void read_ids(std::string set_name, container_t &container)
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.
This is the base class from which all geometric element types are derived.
bool build_sidesets_from_nodesets
Default false.
void assign_sideset_ids()
Called at the end of the read() function, assigns any sideset IDs found when reading the file to the ...
The libMesh namespace provides an interface to certain functionality in the library.
const BoundaryInfo & get_boundary_info() const
The information about boundary ids on the mesh.
void read_nodes(std::string nset_name)
This function parses a block of nodes in the Abaqus file once such a block has been found...
virtual Node * add_point(const Point &p, const dof_id_type id=DofObject::invalid_id, const processor_id_type proc_id=DofObject::invalid_processor_id)=0
Add a new Node at Point p to the end of the vertex array, with processor_id procid.
virtual void extended_parsing(const std::string &)
Override this callback to implement support for more sections in derived classes. ...
This is the MeshBase class.
void read_elements(std::string upper, std::string elset_name)
This function parses a block of elements in the Abaqus file.
void read_sideset(const std::string &sideset_name, const std::string &sideset_type, sideset_container_t &container)
This function reads a sideset from the input file.
std::map< std::string, std::vector< std::pair< dof_id_type, unsigned > > > sideset_container_t
Type of the data structure for storing the (elem ID, side) pairs defining sidesets.
std::string parse_label(std::string line, std::string label_name) const
This function parses a label of the form foo=bar from a comma-delimited line of the form ...
void add_node(const Node *node, const boundary_id_type id)
Add Node node with boundary id id to the boundary information data structures.
bool detect_generated_set(std::string upper) const
virtual void delete_elem(Elem *e)=0
Removes element e from the mesh.
void build_side_list_from_node_list()
Adds sides to a sideset if every node on that side are in the same sideset.
void generate_ids(std::string set_name, container_t &container)
This function handles "generated" nset and elset sections.
virtual Elem * add_elem(Elem *e)=0
Add elem e to the end of the element array.
static std::unique_ptr< Elem > build(const ElemType type, Elem *p=nullptr)
virtual dof_id_type key(const unsigned int s) const =0
SimpleRange< IndexType > as_range(const std::pair< IndexType, IndexType > &p)
Helper function that allows us to treat a homogenous pair as a range.
void assign_subdomain_ids()
This function is called after all the elements have been read and assigns element subdomain IDs...
virtual const Node * query_node_ptr(const dof_id_type i) const =0
std::string unzip_file(std::string_view name)
Create an unzipped copy of a bz2 or xz file, returning the name of the now-unzipped file that can be ...
void assign_boundary_node_ids()
This function assigns boundary IDs to node sets based on the alphabetical order in which the sets are...
container_t _nodeset_ids
Abaqus writes nodesets and elemsets with labels.
std::string & subdomain_name(subdomain_id_type id)
void set_mesh_dimension(unsigned char d)
Resets the logical dimension of the mesh.
bool string_to_num(const std::string &input, dof_id_type &output) const
Attempts to convert the input string to a numerical value using strtol.
virtual void clear()
Deletes all the element and node data that is currently stored.
std::string & sideset_name(boundary_id_type id)
std::string enum_to_string(const T e)
AbaqusIO(MeshBase &mesh)
Constructor.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::set< ElemType > _elem_types
A set of the different geometric element types detected when reading the mesh.
subdomain_id_type subdomain_id() const
virtual unsigned short dim() const =0
unsigned char max_elem_dimension_seen()
void add_side(const dof_id_type elem, const unsigned short int side, const boundary_id_type id)
Add side side of element number elem with boundary id id to the boundary information data structure...
virtual const Elem & elem_ref(const dof_id_type i) const
std::unique_ptr< std::istream > _in
Stream object used to interact with the file.
virtual const Node * node_ptr(const dof_id_type i) const =0
virtual ElemType type() const =0
A Point defines a location in LIBMESH_DIM dimensional Real space.