From 210a2417c667fd49d58a4e00a576f989283fcedf Mon Sep 17 00:00:00 2001 From: heister Date: Mon, 5 Sep 2011 15:26:36 +0000 Subject: [PATCH] update cubit/abaqus mesh conversion tool and update example files. git-svn-id: https://svn.dealii.org/trunk@24250 0785d39b-7218-0410-832d-ea1e28bc413d --- deal.II/contrib/mesh_conversion/Main.cc | 42 +- .../contrib/mesh_conversion/MeshConversion.cc | 1038 ++++++++--------- .../contrib/mesh_conversion/MeshConversion.h | 16 +- .../mesh/3d/{CC_new.ucd => CC_cubit_new.ucd} | 2 +- .../mesh/3d/{CC_old.ucd => CC_cubit_old.ucd} | 2 +- .../mesh_conversion/mesh/3d/other_simple.inp | 17 + .../mesh_conversion/mesh/3d/other_simple.ucd | 13 + 7 files changed, 517 insertions(+), 613 deletions(-) rename deal.II/contrib/mesh_conversion/mesh/3d/{CC_new.ucd => CC_cubit_new.ucd} (98%) rename deal.II/contrib/mesh_conversion/mesh/3d/{CC_old.ucd => CC_cubit_old.ucd} (98%) create mode 100644 deal.II/contrib/mesh_conversion/mesh/3d/other_simple.inp create mode 100644 deal.II/contrib/mesh_conversion/mesh/3d/other_simple.ucd diff --git a/deal.II/contrib/mesh_conversion/Main.cc b/deal.II/contrib/mesh_conversion/Main.cc index 58fc6e4fb4..b013cfb1fa 100755 --- a/deal.II/contrib/mesh_conversion/Main.cc +++ b/deal.II/contrib/mesh_conversion/Main.cc @@ -9,37 +9,36 @@ using namespace std; void display_help (void) { std::cout << "The first input argument is the spatial dimension of the input file, " << std::endl - << "the second is the type of input file (0 for Abaqus OLD, 1 for Abaqus NEW), " + << "the second is the path to the read-in Abaqus .inp file, " << std::endl - << "the third is the path to the read-in Abaqus .inp file, " - << std::endl - << "and the fourth is the name of the file to which you wish to write the output AVS .ucd file to." + << "and the third is the name of the file to which you wish to write the output AVS .ucd file to." << std::endl << std::endl << "Correct program usage is: " << std::endl - << " './convert_mesh /path/to/input_file.inp /path/to/output_file.ucd'" - << std::endl - << "e.g. './convert_mesh 3 0 mesh/3d/test_in.inp mesh/3d/test_out.ucd'" + << " './convert_mesh /path/to/input_file.inp /path/to/output_file.ucd'" << std::endl - << "NOTE: Old Abaqus files are outputted by Cubit 11.1 and earlier. " + << "e.g. './convert_mesh 3 mesh/3d/test_in.inp mesh/3d/test_out.ucd'" << std::endl - << " New Abaqus files are outputted by Cubit 12.0 and later. They MUST be generated with the 'NOT CUBIT ID's' option selected. " + << "NOTE: New Abaqus files outputted by Cubit 12.0 and later MUST be generated with the 'NOT CUBIT ID's' option selected. " << std::endl; } int main (int argc, char* argv[]) { - if (argc == 5) - { + if (argc != 4) + { + display_help(); + return 0; + } + try { - const unsigned int dimension = atoi (argv[1]); //*argv[1] - '0'; + const unsigned int dimension = atoi (argv[1]); std::cout << "Dimension: " << dimension << std::endl; - const unsigned int input_type = atoi (argv[2]); //*argv[2] - '0'; - std::string input_file = argv[3]; - std::string output_file = argv[4]; + std::string input_file = argv[2]; + std::string output_file = argv[3]; - MeshConversion mesh (dimension, input_type, input_file, output_file); + MeshConversion mesh (dimension, input_file, output_file); const bool success = mesh.convert_mesh (); if (success == false) { @@ -57,17 +56,6 @@ int main (int argc, char* argv[]) { return 1; } - } - else if (argc == 2) { - if (argv[1] == "-h" || argv[1] == "-help") - display_help(); - else - display_help(); - } - else { - std::cerr << "WRONG NUMBER OF INPUT ARGUMENTS:" << std::endl; - display_help(); - } return 0; } diff --git a/deal.II/contrib/mesh_conversion/MeshConversion.cc b/deal.II/contrib/mesh_conversion/MeshConversion.cc index fea24f172f..60060afad6 100755 --- a/deal.II/contrib/mesh_conversion/MeshConversion.cc +++ b/deal.II/contrib/mesh_conversion/MeshConversion.cc @@ -2,60 +2,51 @@ //-------------------------------------------------------------------------------- -MeshConversion::MeshConversion (const unsigned int dimension, const int input_type, const std::string input_file, const std::string output_file) -: -dimension (dimension), -input_file_type (input_type), -input_file_name (input_file), -output_file_name (output_file), -tolerance (5e-16) // Used to offset Cubit tolerance error when outputting value close to zero +MeshConversion::MeshConversion (const unsigned int dimension, const std::string input_file, const std::string output_file) + : + tolerance (5e-16), // Used to offset Cubit tolerance error when outputting value close to zero + dimension (dimension), + input_file_name (input_file), + output_file_name (output_file) { - greeting (); - - if (dimension == 2) - { - node_per_face = 2; // Line edge - node_per_cell = 4; // Quad face - face_per_cell = 2*dimension; // Lines per quad - } - else if (dimension == 3) - { - node_per_face = 4; // Quad face - node_per_cell = 8; // Hex cell - face_per_cell = 2*dimension; // Quads per hex - } - else - { - std::cerr << "ERROR: Chosen spatial dimension is invalid!" << std::endl; - } + greeting (); + + if (dimension == 2) + { + node_per_face = 2; // Line edge + node_per_cell = 4; // Quad face + face_per_cell = 2 * dimension; // Lines per quad + } + else + if (dimension == 3) + { + node_per_face = 4; // Quad face + node_per_cell = 8; // Hex cell + face_per_cell = 2 * dimension; // Quads per hex + } + else + { + std::cerr << "ERROR: Chosen spatial dimension is invalid!" << std::endl; + } } //-------------------------------------------------------------------------------- MeshConversion::~MeshConversion () { - } //-------------------------------------------------------------------------------- -bool MeshConversion::convert_mesh (void) +bool MeshConversion::convert_mesh (void) { - bool readin_successul = false; - if (input_file_type == MeshConversion::abaqus_old) - readin_successul = read_in_abaqus_inp_old (); - else if (input_file_type == MeshConversion::abaqus_new) - readin_successul = read_in_abaqus_inp_new (); - else - std::cerr << "ERROR: File input type invalid." << std::endl; - - bool writeout_successul = false; - if (readin_successul == true) - writeout_successul = write_out_avs_ucd(); - else - std::cerr << "ERROR: Input file not found." << std::endl; - - return readin_successul && writeout_successul; + if (!read_in_abaqus()) + { + std::cerr << "ERROR: could not parse file." << std::endl; + return false; + } + + return write_out_avs_ucd(); } //-------------------------------------------------------------------------------- @@ -64,559 +55,458 @@ bool MeshConversion::convert_mesh (void) // ======================================================== // http://www.codeguru.com/forum/showthread.php?t=231054 // ======================================================== -template bool from_string(T& t, const std::string& s, std::ios_base& (*f)(std::ios_base&)) +template bool from_string (T& t, const std::string& s, std::ios_base& (*f) (std::ios_base&)) { - std::istringstream iss(s); - return !(iss >> f >> t).fail(); + std::istringstream iss (s); + return ! (iss >> f >> t).fail(); } //-------------------------------------------------------------------------------- -void MeshConversion::greeting () { - std::cout << std::endl; std::cout << std::endl; std::cout << std::endl; - std::cout << "************************************************************************************" << std::endl; - std::cout << "*** FEM Mesh conversion tool ***" << std::endl; - std::cout << "************************************************************************************" << std::endl; - std::cout << "*** Author: Jean-Paul Pelteret ***" << std::endl; - std::cout << "*** Date: June 2010 ***" << std::endl; - std::cout << "*** References: ***" << std::endl; - std::cout << "*** http://www.dealii.org/developer/doxygen/deal.II/structGeometryInfo.html ***" << std::endl; - std::cout << "*** http://people.scs.fsu.edu/~burkardt/html/ucd_format.html ***" << std::endl; - std::cout << "*** http://people.scs.fsu.edu/~burkardt/data/ucd/ucd.html ***" << std::endl; - std::cout << "*** http://www.cprogramming.com/tutorial/string.html ***" << std::endl; - std::cout << "*** http://www.codeguru.com/forum/showthread.php?t=231054 ***" << std::endl; - std::cout << "************************************************************************************" << std::endl; - std::cout << "*** FEATURES: ***" << std::endl; - std::cout << "*** Read-in file types: ***" << std::endl; - std::cout << "*** - Abaqus inp (from Cubit) ***" << std::endl; - std::cout << "*** Write-out file types: ***" << std::endl; - std::cout << "*** - AVS UCD ***" << std::endl; - std::cout << "************************************************************************************" << std::endl; - std::cout << std::endl; +void MeshConversion::greeting () +{ + std::cout << std::endl; + std::cout << std::endl; + std::cout << std::endl; + std::cout << "************************************************************************************" << std::endl; + std::cout << "*** FEM Mesh conversion tool ***" << std::endl; + std::cout << "************************************************************************************" << std::endl; + std::cout << "*** Author: Jean-Paul Pelteret, Timo Heister ***" << std::endl; + std::cout << "*** Date: June 2010, September 2011 ***" << std::endl; + std::cout << "*** References: ***" << std::endl; + std::cout << "*** http://www.dealii.org/developer/doxygen/deal.II/structGeometryInfo.html ***" << std::endl; + std::cout << "*** http://people.scs.fsu.edu/~burkardt/html/ucd_format.html ***" << std::endl; + std::cout << "*** http://people.scs.fsu.edu/~burkardt/data/ucd/ucd.html ***" << std::endl; + std::cout << "*** http://www.cprogramming.com/tutorial/string.html ***" << std::endl; + std::cout << "*** http://www.codeguru.com/forum/showthread.php?t=231054 ***" << std::endl; + std::cout << "************************************************************************************" << std::endl; + std::cout << "*** FEATURES: ***" << std::endl; + std::cout << "*** Read-in file types: ***" << std::endl; + std::cout << "*** - Abaqus inp (from Cubit) ***" << std::endl; + std::cout << "*** Write-out file types: ***" << std::endl; + std::cout << "*** - AVS UCD ***" << std::endl; + std::cout << "************************************************************************************" << std::endl; + std::cout << std::endl; } -//-------------------------------------------------------------------------------- -bool MeshConversion::read_in_abaqus_inp_old (void) + +bool MeshConversion::read_in_abaqus () { - bool read_in_successful = true; - std::cout << "Reading in ABAQUS .inp FILE: " << input_file_name << std::endl; - - std::ifstream input_stream; - - input_stream.open(input_file_name.c_str(), std::ifstream::in); - - // Check to see if file exists - if (input_stream.good() == true) - { - std::string buffer; - std::vector < std::string > temp_data; - - // Loop over all the contents of the fibre data file and read it into buffer - while (input_stream >> buffer) - { - for(int i = 0; i < buffer.length(); ++i) - if (buffer[i] == ',') // Get rid of the .inp file's useless comma's - { - buffer.erase(i); - --i; - } - - temp_data.push_back(buffer); - } - - const bool print_input_deck = false; - if (print_input_deck == true) - { - // =========== TEMP =========== - std::string filename_out ("checks/input_stream.txt"); - std::ofstream output; - output.open(filename_out.c_str()); - - for (int ii = 0; ii < temp_data.size(); ii++ ) - output << temp_data[ii] << std::endl; - - output.close(); - // =========== TEMP =========== - } - - for (int k = 0; k < temp_data.size(); ++k) - { - // ================================ NODES =================================== - // ABAQUS formatting - // *NODE - // , , , - // ========================================================================== - - if (temp_data[k] == "*NODE") - { - int j = 0; - float temp_float; - - while (from_string (temp_float, temp_data[k + 1 + j*(dimension+1+(dimension==2 ? 1 : 0))], std::dec) == true) - { - // Initilise storage variables - std::vector node (dimension+1); - - // Convert from string to the variable types - for (int i = 0; i < dimension+1; ++i) - from_string (node[i], temp_data[k + 1 + j*(dimension+1+(dimension==2 ? 1 : 0)) + i], std::dec); - - // Add to the global node number vector - node_list.push_back(node); - - ++j; - } - } - - // ===================== 2D QUADS / 3D HEX CELLS ======================= - // ABAQUS formatting - // *ELEMENT, TYPE=C3D8R, ELSET=EB< block-set number (material-id) > - // , <8 x node number> - // ===================================================================== - - else if (temp_data[k] == "*ELEMENT") - { - int j = 0; - float temp_float; - const int data_per_cell = node_per_cell + 1; - - // Get material id - std::string material_id_line = temp_data[k+ 2]; - std::string material_id_temp; - for (int ll = 8 /*Characters in "ELSET=EB" */; ll < material_id_line.length(); ++ll) - material_id_temp += material_id_line[ll]; - - int material_id = 0; - for (int ll = 0; ll < material_id_temp.length(); ++ll) - material_id += (material_id_temp[material_id_temp.length() - ll - 1] - 48 /* ASCII TRICKS */) * pow(10.0,ll); - - while (from_string (temp_float, temp_data[k + 3 + j*(data_per_cell)], std::dec) == true) - { - // Initilise storage variables - std::vector cell (data_per_cell); - - // Material id - cell[0] = material_id; - - // Convert from string to the variable types - for (int i = 1; i < data_per_cell; ++i) - from_string (cell[i], temp_data[k + 3 + j*(data_per_cell) + i], std::dec); - - // Add to the global node number vector - cell_list.push_back(cell); - - ++j; - } - } - - // =========== 2D LINES / 3D BC_QUADS =========== - // ABAQUS formatting - // SURFACE, NAME=SS - // , - // ============================================== - - else if (temp_data[k] == "*SURFACE" ) - { - float temp_float; - const int data_per_quad = node_per_face + 1; - - // Get sideset id - std::string sideset_id_line = temp_data[k + 1]; - std::string sideset_id_temp; - for (int m = 7 /*Characters in "NAME=SS" */; m < sideset_id_line.length(); ++m) - sideset_id_temp += sideset_id_line[m]; - - int sideset_id = 0; - for (int m = 0; m < sideset_id_temp.length(); ++m) - sideset_id += (sideset_id_temp[sideset_id_temp.length() - m - 1] - 48 /* ASCII TRICKS */) * pow(10.0,m); - - const int data_per_face = 2; - int j = 0; - while (( (k + 2) + j * data_per_face ) < temp_data.size() && from_string (temp_float, temp_data[(k + 2) + j * data_per_face], std::dec) == true) - { - // Get cell to which the face belongs - std::string face_cell_no_line = temp_data[(k + 2) + j * data_per_face]; - int face_cell_no = 0; - for (int m = 0; m < face_cell_no_line.length(); ++m) - face_cell_no += (face_cell_no_line[face_cell_no_line.length() - m - 1] - 48 /* ASCII TRICKS */) * pow(10.0,m); - - // Get ABAQUS cell face number - std::string face_cell_face_no_line = temp_data[(k + 2) + j * data_per_face + 1]; - std::string face_cell_face_no_temp; - for (int m = 1 /*Characters in "S" */; m < face_cell_face_no_line.length(); ++m) - face_cell_face_no_temp += face_cell_face_no_line[m]; - - int face_cell_face_no = 0; - for (int m = 0; m < face_cell_face_no_temp.length(); ++m) - face_cell_face_no += (face_cell_face_no_temp[face_cell_face_no_temp.length() - m - 1] - 48 /* ASCII TRICKS */) * pow(10.0,m); - - // Initilise storage variables - std::vector quad (data_per_quad); - std::vector quad_node_list (node_per_face); - - quad_node_list = get_global_node_numbers(face_cell_no, face_cell_face_no); - - // Sideset id - quad[0] = sideset_id; - - // Global node numbers - for (int m = 0; m < node_per_face; ++m) - quad[m + 1] = quad_node_list[m]; - - // Add to the global quad vector - face_list.push_back(quad); - - ++j; - } - } - } - } - else - { - return false; - } - - return read_in_successful; -} + std::cout << "Reading in ABAQUS .inp FILE: " << input_file_name << std::endl; -//-------------------------------------------------------------------------------- + std::ifstream input_stream; -bool MeshConversion::read_in_abaqus_inp_new (void) -{ - bool read_in_successful = true; - std::cout << "Reading in ABAQUS .inp FILE: " << input_file_name << std::endl; - - std::ifstream input_stream; - - input_stream.open(input_file_name.c_str(), std::ifstream::in); - - // Check to see if file exists - if (input_stream.good() == true) { - std::vector temp_data; - std::string buffer; - - // Read all input data into a temp vector - while (input_stream >> buffer) { - if (buffer.compare(buffer.size()-1,buffer.size()-1, ",") == 0) buffer.erase(buffer.size() -1); - temp_data.push_back(buffer); - } - - // Sort the data into node, element and surface data_per_cell - for (unsigned int k=0; k < temp_data.size(); ++k) { - // Nodes - if (temp_data[k] == "*NODE") { - int j = 0; - float temp_float; - const int header_size = 2; // *NODE, NSET= - - while (from_string (temp_float, temp_data[k + header_size + j*(dimension+1+(dimension==2 ? 1 : 0))], std::dec) == true) - { - // Initilise storage variables - std::vector node (dimension+1); - - // Convert from string to the variable types - for (int i = 0; i < dimension+1; ++i) - from_string (node[i], temp_data[k + header_size + j*(dimension+1+(dimension==2 ? 1 : 0)) + i], std::dec); - - // Add to the global node number vector - node_list.push_back(node); - - ++j; - } - - } - if (temp_data[k] == "*ELEMENT") { - int j = 0; - float temp_float; - const int data_per_cell = node_per_cell + 1; - const int header_size = 3; // *ELEMENT, TYPE=C3D8R, ELSET= - const int mat_id_line = 2; - - // Get material id - std::string material_id_line = temp_data[k+mat_id_line]; - std::string material_id_temp; - for (int ll = 8 /*Characters in "ELSET=EB" */; ll < material_id_line.length(); ++ll) - material_id_temp += material_id_line[ll]; - - int material_id = 0; - for (int ll = 0; ll < material_id_temp.length(); ++ll) - material_id += (material_id_temp[material_id_temp.length() - ll - 1] - 48 /* ASCII TRICKS */) * pow(10.0,ll); - - while (from_string (temp_float, temp_data[k + header_size + j*(data_per_cell)], std::dec) == true) - { - // Initilise storage variables - std::vector cell (data_per_cell); - - // Material id - cell[0] = material_id; - - // Convert from string to the variable types - for (int i = 1; i < data_per_cell; ++i) - from_string (cell[i], temp_data[k + 3 + j*(data_per_cell) + i], std::dec); - - // Add to the global node number vector - cell_list.push_back(cell); - - ++j; - } - - } - if (temp_data[k] == "*ELSET") { - float temp_float; - const int data_per_quad = node_per_face + 1; - const int header_size = 2; // *ELSET, ELSET= - const std::string elset_name = temp_data[k+1].substr(6, temp_data[k+1].size()); // Take out subset name to extract Abaqus face side from later - - // Find the surface to which this element set belongs - this is the next line on which "*SURFACE" appears as multiple elsets can be attached to one surface - int surface_line_offset = 0; - while (temp_data[k + surface_line_offset] != "*SURFACE") - ++surface_line_offset; - - // Get sideset id - const int face_name_from_surface = 1; - std::string sideset_id_line = temp_data[k + surface_line_offset + face_name_from_surface]; // Name on next line - std::string sideset_id_temp; - for (int m = 7 /*Characters in "NAME=SS" */; m < sideset_id_line.length(); ++m) - sideset_id_temp += sideset_id_line[m]; - - // Get ABAQUS cell face number - std::string face_cell_face_no_temp; - for (unsigned int s=k + surface_line_offset; s < temp_data.size(); ++ s) { - if (temp_data[s] == elset_name) { - face_cell_face_no_temp = temp_data[s+1].substr(1, temp_data[s+1].size()); // Next line contains cell side face number - } - } - - int face_cell_face_no = 0; - for (int m = 0; m < face_cell_face_no_temp.length(); ++m) - face_cell_face_no += (face_cell_face_no_temp[face_cell_face_no_temp.length() - m - 1] - 48 /* ASCII TRICKS */) * pow(10.0,m); - - int sideset_id = 0; - for (int m = 0; m < sideset_id_temp.length(); ++m) - sideset_id += (sideset_id_temp[sideset_id_temp.length() - m - 1] - 48 /* ASCII TRICKS */) * pow(10.0,m); - - const int data_per_face = 1; - int j = 0; - while (( (k + header_size) + j * data_per_face ) < temp_data.size() && from_string (temp_float, temp_data[(k + header_size) + j * data_per_face], std::dec) == true) - { - // Get cell to which the face belongs - std::string face_cell_no_line = temp_data[(k + header_size) + j * data_per_face]; - int face_cell_no = 0; - for (int m = 0; m < face_cell_no_line.length(); ++m) - face_cell_no += (face_cell_no_line[face_cell_no_line.length() - m - 1] - 48 /* ASCII TRICKS */) * pow(10.0,m); - - // Initilise storage variables - std::vector quad (data_per_quad); - std::vector quad_node_list (node_per_face); - - quad_node_list = get_global_node_numbers(face_cell_no, face_cell_face_no); - - // Sideset id - quad[0] = sideset_id; - - // Global node numbers - for (int m = 0; m < node_per_face; ++m) - quad[m + 1] = quad_node_list[m]; - - // Add to the global quad vector - face_list.push_back(quad); - - ++j; - } - - } - } - - input_stream.close(); - } - else - { - return false; - } - - return read_in_successful; + input_stream.open (input_file_name.c_str(), std::ifstream::in); + + if (!input_stream.good()) + { + std::cerr << "ERROR: file '" << input_file_name << "' not found." << std::endl; + return false; + } + + std::string line; + + std::getline (input_stream, line); + + while (!input_stream.eof()) + { + if (0 == line.compare ("*HEADING") + || 0 == line.compare (0, 2, "**") + || 0 == line.compare (0, 5, "*PART")) + { + //skip header and comments + while (!input_stream.eof()) + { + std::getline (input_stream, line); + if (line[0] == '*') + goto cont; + } + } + else + if (0 == line.compare (0, 5, "*NODE")) + { + //header might be: + //*NODE, NSET=ALLNODES + //*NODE + + // contains lines in the form: + //index, x, y, z + int c = 0; + while (!input_stream.eof()) + { + std::getline (input_stream, line); + if (line[0] == '*') + goto cont; + + std::vector node (dimension + 1); + + istringstream iss (line); + char comma; + for (unsigned int i = 0; i < dimension + 1; ++i) + iss >> node[i] >> comma; + + node_list.push_back (node); + + ++c; + } + } + else + if (0 == line.compare (0, 8, "*ELEMENT")) + { + //different header formats: + //*ELEMENT, TYPE=S4R, ELSET=EB + //*ELEMENT, TYPE=C3D8R, ELSET=EB + //*ELEMENT, TYPE=C3D8 + + //elements itself (n=4 or n=8): + //index, i[0], ..., i[n] + + int material = 0; + + //scan for material + { + string before_material = "ELSET=EB"; + size_t idx = line.find (before_material); + if (idx != std::string::npos) + { + from_string (material, line.substr (idx + before_material.size()), std::dec); + } + } + + std::getline (input_stream, line); + while (!input_stream.eof()) + { + if (line[0] == '*') + goto cont; + + istringstream iss (line); + char comma; + + const int n_points = 2 << dimension; + std::vector cell (1 + n_points); + for (int i = 0; i < 1 + n_points; ++i) + iss >> cell[i] >> comma; + + //overwrite cell index from file by material + cell[0] = (double) material; + + cell_list.push_back (cell); + + std::getline (input_stream, line); + } + } + else + if (0 == line.compare (0, 8, "*SURFACE")) + { + //old format: + //*SURFACE, NAME=SS + // , S + + int b_indicator = 0; + { + string before = "NAME=SS"; + size_t idx = line.find (before); + if (idx != std::string::npos) + { + from_string (b_indicator, line.substr (idx + before.size()), std::dec); + } + } + + std::getline (input_stream, line); + while (!input_stream.eof()) + { + if (line[0] == '*') + goto cont; + + istringstream iss (line); + char comma; + char temp; + int el_idx; + int face_number; + string face_number_str; + iss >> el_idx >> comma >> temp >> face_number; + + std::vector quad_node_list; + quad_node_list = get_global_node_numbers (el_idx, face_number); + quad_node_list.insert (quad_node_list.begin(), b_indicator); + + face_list.push_back (quad_node_list); + + std::getline (input_stream, line); + } + } + else + if (0 == line.compare (0, 6, "*ELSET")) + { + //new style line/face indicators + //*ELSET, ELSET=SS1_S3 + //,, ... + //*SURFACE, NAME=SS + //SS1_S3, S + + std::vector elements; + std::getline (input_stream, line); + while (!input_stream.eof()) + { + if (line[0] == '*') + break; + + istringstream iss (line); + char comma; + int elid; + while (!iss.eof()) + { + iss >> elid >> comma; + elements.push_back (elid); + } + + std::getline (input_stream, line); + } + + if (0 != line.compare (0, 8, "*SURFACE")) + { + std::cout << "WARNING: expected SURFACE after ELSET." << std::endl; + goto cont; + } + + int b_indicator = 0; + { + string before = "NAME=SS"; + size_t idx = line.find (before); + if (idx != std::string::npos) + { + from_string (b_indicator, line.substr (idx + before.size()), std::dec); + } + } + + std::getline (input_stream, line); + + istringstream iss (line); + char comma; + std::string temp; + int face_number; + iss >> temp >> comma >> face_number; + + for (unsigned int i = 0;i < elements.size();++i) + { + std::vector quad_node_list; + + quad_node_list = get_global_node_numbers (elements[i], face_number); + + quad_node_list.insert (quad_node_list.begin(), b_indicator); + + face_list.push_back (quad_node_list); + } + } + else + if (0 == line.compare (0, 5, "*NSET")) + { + //skip nodesets, we have no use for them + while (!input_stream.eof()) + { + std::getline (input_stream, line); + if (line[0] == '*') + goto cont; + } + } + else + { + std::cout << "Warning, skipping line:" << line << std::endl; + } + + std::getline (input_stream, line); + cont: + (void) 0; + } + return true; } + + //-------------------------------------------------------------------------------- -std::vector MeshConversion::get_global_node_numbers (const int face_cell_no, const int face_cell_face_no) +std::vector MeshConversion::get_global_node_numbers (const int face_cell_no, const int face_cell_face_no) { - std::vector quad_node_list (node_per_face); - - if (dimension == 2) - { - if (face_cell_face_no == 1) - { - quad_node_list[0] = cell_list[face_cell_no - 1][1]; - quad_node_list[1] = cell_list[face_cell_no - 1][2]; - } - else if (face_cell_face_no == 2) - { - quad_node_list[0] = cell_list[face_cell_no - 1][2]; - quad_node_list[1] = cell_list[face_cell_no - 1][3]; - } - else if (face_cell_face_no == 3) - { - quad_node_list[0] = cell_list[face_cell_no - 1][3]; - quad_node_list[1] = cell_list[face_cell_no - 1][4]; - } - else if (face_cell_face_no == 4) - { - quad_node_list[0] = cell_list[face_cell_no - 1][4]; - quad_node_list[1] = cell_list[face_cell_no - 1][1]; - } - else - { - std::cerr << "ERROR! Invalid face number inputted!" << std::endl; - std::cerr << "face_cell_face_no = " << face_cell_face_no << std::endl; - } - } - else // (dimension == 3) - { - if (face_cell_face_no == 1) - { - quad_node_list[0] = cell_list[face_cell_no - 1][1]; - quad_node_list[1] = cell_list[face_cell_no - 1][4]; - quad_node_list[2] = cell_list[face_cell_no - 1][3]; - quad_node_list[3] = cell_list[face_cell_no - 1][2]; - } - else if (face_cell_face_no == 2) - { - quad_node_list[0] = cell_list[face_cell_no - 1][5]; - quad_node_list[1] = cell_list[face_cell_no - 1][8]; - quad_node_list[2] = cell_list[face_cell_no - 1][7]; - quad_node_list[3] = cell_list[face_cell_no - 1][6]; - } - else if (face_cell_face_no == 3) - { - quad_node_list[0] = cell_list[face_cell_no - 1][1]; - quad_node_list[1] = cell_list[face_cell_no - 1][2]; - quad_node_list[2] = cell_list[face_cell_no - 1][6]; - quad_node_list[3] = cell_list[face_cell_no - 1][5]; - } - else if (face_cell_face_no == 4) - { - quad_node_list[0] = cell_list[face_cell_no - 1][2]; - quad_node_list[1] = cell_list[face_cell_no - 1][3]; - quad_node_list[2] = cell_list[face_cell_no - 1][7]; - quad_node_list[3] = cell_list[face_cell_no - 1][6]; - } - else if (face_cell_face_no == 5) - { - quad_node_list[0] = cell_list[face_cell_no - 1][3]; - quad_node_list[1] = cell_list[face_cell_no - 1][4]; - quad_node_list[2] = cell_list[face_cell_no - 1][8]; - quad_node_list[3] = cell_list[face_cell_no - 1][7]; - } - else if (face_cell_face_no == 6) - { - quad_node_list[0] = cell_list[face_cell_no - 1][1]; - quad_node_list[1] = cell_list[face_cell_no - 1][5]; - quad_node_list[2] = cell_list[face_cell_no - 1][8]; - quad_node_list[3] = cell_list[face_cell_no - 1][4]; - } - else - { - std::cerr << "ERROR! Invalid face number inputted!" << std::endl; - std::cerr << "face_cell_face_no = " << face_cell_face_no << std::endl; - } - } - - return quad_node_list; + std::vector quad_node_list (node_per_face); + + if (dimension == 2) + { + if (face_cell_face_no == 1) + { + quad_node_list[0] = cell_list[face_cell_no - 1][1]; + quad_node_list[1] = cell_list[face_cell_no - 1][2]; + } + else + if (face_cell_face_no == 2) + { + quad_node_list[0] = cell_list[face_cell_no - 1][2]; + quad_node_list[1] = cell_list[face_cell_no - 1][3]; + } + else + if (face_cell_face_no == 3) + { + quad_node_list[0] = cell_list[face_cell_no - 1][3]; + quad_node_list[1] = cell_list[face_cell_no - 1][4]; + } + else + if (face_cell_face_no == 4) + { + quad_node_list[0] = cell_list[face_cell_no - 1][4]; + quad_node_list[1] = cell_list[face_cell_no - 1][1]; + } + else + { + std::cerr << "ERROR! Invalid face number inputted!" << std::endl; + std::cerr << "face_cell_face_no = " << face_cell_face_no << std::endl; + } + } + else // (dimension == 3) + { + if (face_cell_face_no == 1) + { + quad_node_list[0] = cell_list[face_cell_no - 1][1]; + quad_node_list[1] = cell_list[face_cell_no - 1][4]; + quad_node_list[2] = cell_list[face_cell_no - 1][3]; + quad_node_list[3] = cell_list[face_cell_no - 1][2]; + } + else + if (face_cell_face_no == 2) + { + quad_node_list[0] = cell_list[face_cell_no - 1][5]; + quad_node_list[1] = cell_list[face_cell_no - 1][8]; + quad_node_list[2] = cell_list[face_cell_no - 1][7]; + quad_node_list[3] = cell_list[face_cell_no - 1][6]; + } + else + if (face_cell_face_no == 3) + { + quad_node_list[0] = cell_list[face_cell_no - 1][1]; + quad_node_list[1] = cell_list[face_cell_no - 1][2]; + quad_node_list[2] = cell_list[face_cell_no - 1][6]; + quad_node_list[3] = cell_list[face_cell_no - 1][5]; + } + else + if (face_cell_face_no == 4) + { + quad_node_list[0] = cell_list[face_cell_no - 1][2]; + quad_node_list[1] = cell_list[face_cell_no - 1][3]; + quad_node_list[2] = cell_list[face_cell_no - 1][7]; + quad_node_list[3] = cell_list[face_cell_no - 1][6]; + } + else + if (face_cell_face_no == 5) + { + quad_node_list[0] = cell_list[face_cell_no - 1][3]; + quad_node_list[1] = cell_list[face_cell_no - 1][4]; + quad_node_list[2] = cell_list[face_cell_no - 1][8]; + quad_node_list[3] = cell_list[face_cell_no - 1][7]; + } + else + if (face_cell_face_no == 6) + { + quad_node_list[0] = cell_list[face_cell_no - 1][1]; + quad_node_list[1] = cell_list[face_cell_no - 1][5]; + quad_node_list[2] = cell_list[face_cell_no - 1][8]; + quad_node_list[3] = cell_list[face_cell_no - 1][4]; + } + else + { + std::cerr << "ERROR! Invalid face number inputted!" << std::endl; + std::cerr << "face_cell_face_no = " << face_cell_face_no << std::endl; + } + } + + return quad_node_list; } //-------------------------------------------------------------------------------- -bool MeshConversion::write_out_avs_ucd (void) +bool MeshConversion::write_out_avs_ucd (void) { - bool write_out_successful = true; - std::cout << "Writing out AVS .ucd FILE: " << output_file_name << std::endl; - - std::ofstream output; - output.open(output_file_name.c_str()); - - // Write out title - Note: No other commented text can be inserted below the title in a UCD file - output << "# FEM Mesh Converter" << std::endl; - output << "# Mesh type: AVS UCD" << std::endl; - output << "# Input file name: " << input_file_name << std::endl; - - // ======================================================== - // http://people.scs.fsu.edu/~burkardt/html/ucd_format.html - // ======================================================== - // ASCII UCD File Format - // The input file cannot contain blank lines or lines with leading blanks. Comments, if present, must precede all data in the file. Comments within the data will cause read errors. The general order of the data is as follows: - // 1. Numbers defining the overall structure, including the number of nodes, the number of cells, and the length of the vector of data associated with the nodes, cells, and the model. - // e.g. 1: - // - // e.g. 2: - // n_elements = n_hex_cells + n_bc_quads + n_quad_cells + n_bc_edges - // outfile.write(str(n_nodes) + " " + str(n_elements) + " 0 0 0\n") - // 2. For each node, its node id and the coordinates of that node in space. Node-ids must be integers, but any number including non sequential numbers can be used. Mid-edge nodes are treated like any other node. - // 3. For each cell: its cell-id, material, cell type (hexahedral, pyramid, etc.), and the list of node-ids that correspond to each of the cell's vertices. The below table specifies the different cell types and the keyword used to represent them in the file. - - // Write out header - output << node_list.size() << "\t" << (cell_list.size() + face_list.size()) << "\t0\t0\t0" << std::endl; - - // Write out node numbers - for (int ii = 0; ii < node_list.size(); ++ii) // Loop over all nodes - { - for (int jj = 0; jj < dimension + 1; ++jj) // Loop over entries to be outputted - { - if (jj == 0) // Node number - { - output.precision(); - output << node_list[ii][jj] << "\t"; - } - else // Node coordinates - { - output.width(16); - output.setf(ios::scientific, ios::floatfield); - output.precision(8); - if (abs(node_list[ii][jj]) > tolerance) // invoke tolerance -> set points close to zero equal to zero - output << static_cast (node_list[ii][jj]) << "\t"; - else - output << 0.0 << "\t"; - } - } - if (dimension == 2) - output << 0.0 << "\t"; - - output << std::endl; - output.unsetf(ios::floatfield); - } - - // Write out cell node numbers - for (int ii = 0; ii < cell_list.size(); ++ii) - { - output << ii + 1 << "\t" << cell_list[ii][0] << "\t" << (dimension==2 ? "quad" : "hex") << "\t"; - for (int jj = 1; jj < node_per_cell + 1; ++jj) - output << cell_list[ii][jj] << "\t"; - - output << std::endl; - } - - // Write out quad node numbers - for (int ii = 0; ii < face_list.size(); ++ii) - { - output << ii + 1 << "\t" << face_list[ii][0] << "\t" << (dimension==2 ? "line" : "quad") << "\t"; - for (int jj = 1; jj < node_per_face + 1; ++jj) - output << face_list[ii][jj] << "\t"; - - output << std::endl; - } - - output.close(); - - if (write_out_successful == true) - { - std::cout << "Output successful!" << std::endl; - std::cout << " Number of nodes: " << node_list.size() << std::endl; - std::cout << " Number of cells: " << cell_list.size() << std::endl; - std::cout << " Number of boundary faces: " << face_list.size() << std::endl << std::endl; - } - - return write_out_successful; + std::cout << "Writing out AVS .ucd FILE: " << output_file_name << std::endl; + + std::ofstream output; + output.open (output_file_name.c_str()); + + // Write out title - Note: No other commented text can be inserted below the title in a UCD file + output << "# FEM Mesh Converter" << std::endl; + output << "# Mesh type: AVS UCD" << std::endl; + output << "# Input file name: " << input_file_name << std::endl; + + // ======================================================== + // http://people.scs.fsu.edu/~burkardt/html/ucd_format.html + // ======================================================== + // ASCII UCD File Format + // The input file cannot contain blank lines or lines with leading blanks. Comments, if present, must precede all data in the file. Comments within the data will cause read errors. The general order of the data is as follows: + // 1. Numbers defining the overall structure, including the number of nodes, the number of cells, and the length of the vector of data associated with the nodes, cells, and the model. + // e.g. 1: + // + // e.g. 2: + // n_elements = n_hex_cells + n_bc_quads + n_quad_cells + n_bc_edges + // outfile.write(str(n_nodes) + " " + str(n_elements) + " 0 0 0\n") + // 2. For each node, its node id and the coordinates of that node in space. Node-ids must be integers, but any number including non sequential numbers can be used. Mid-edge nodes are treated like any other node. + // 3. For each cell: its cell-id, material, cell type (hexahedral, pyramid, etc.), and the list of node-ids that correspond to each of the cell's vertices. The below table specifies the different cell types and the keyword used to represent them in the file. + + // Write out header + output << node_list.size() << "\t" << (cell_list.size() + face_list.size()) << "\t0\t0\t0" << std::endl; + + // Write out node numbers + for (unsigned int ii = 0; ii < node_list.size(); ++ii) // Loop over all nodes + { + for (unsigned int jj = 0; jj < dimension + 1; ++jj) // Loop over entries to be outputted + { + if (jj == 0) // Node number + { + output.precision(); + output << node_list[ii][jj] << "\t"; + } + else // Node coordinates + { + output.width (16); + output.setf (ios::scientific, ios::floatfield); + output.precision (8); + if (abs (node_list[ii][jj]) > tolerance) // invoke tolerance -> set points close to zero equal to zero + output << static_cast (node_list[ii][jj]) << "\t"; + else + output << 0.0 << "\t"; + } + } + if (dimension == 2) + output << 0.0 << "\t"; + + output << std::endl; + output.unsetf (ios::floatfield); + } + + // Write out cell node numbers + for (unsigned int ii = 0; ii < cell_list.size(); ++ii) + { + output << ii + 1 << "\t" << cell_list[ii][0] << "\t" << (dimension == 2 ? "quad" : "hex") << "\t"; + for (unsigned int jj = 1; jj < node_per_cell + 1; ++jj) + output << cell_list[ii][jj] << "\t"; + + output << std::endl; + } + + // Write out quad node numbers + for (unsigned int ii = 0; ii < face_list.size(); ++ii) + { + output << ii + 1 << "\t" << face_list[ii][0] << "\t" << (dimension == 2 ? "line" : "quad") << "\t"; + for (unsigned int jj = 1; jj < node_per_face + 1; ++jj) + output << face_list[ii][jj] << "\t"; + + output << std::endl; + } + + output.close(); + + std::cout << "Output successful!" << std::endl; + std::cout << " Number of nodes: " << node_list.size() << std::endl; + std::cout << " Number of cells: " << cell_list.size() << std::endl; + std::cout << " Number of boundary faces: " << face_list.size() << std::endl << std::endl; + + return true; } //-------------------------------------------------------------------------------- diff --git a/deal.II/contrib/mesh_conversion/MeshConversion.h b/deal.II/contrib/mesh_conversion/MeshConversion.h index a67fbdc9c2..49dd3df4fb 100755 --- a/deal.II/contrib/mesh_conversion/MeshConversion.h +++ b/deal.II/contrib/mesh_conversion/MeshConversion.h @@ -16,24 +16,20 @@ using namespace std; class MeshConversion { public: - MeshConversion (const unsigned int dimension, const int input_type, const std::string input_file, const std::string output_file); - ~MeshConversion (void); - bool convert_mesh (void); + MeshConversion (const unsigned int dimension, const std::string input_file, const std::string output_file); + ~MeshConversion (); + bool convert_mesh (); - enum {abaqus_old = 0, abaqus_new} input_type; - private: void greeting (); - - bool read_in_abaqus_inp_old (void); - bool read_in_abaqus_inp_new (void); - bool write_out_avs_ucd (void); + + bool read_in_abaqus (); + bool write_out_avs_ucd (); std::vector get_global_node_numbers(const int face_cell_no, const int face_cell_face_no); const double tolerance; const unsigned int dimension; - const int input_file_type; unsigned int node_per_face; unsigned int node_per_cell; diff --git a/deal.II/contrib/mesh_conversion/mesh/3d/CC_new.ucd b/deal.II/contrib/mesh_conversion/mesh/3d/CC_cubit_new.ucd similarity index 98% rename from deal.II/contrib/mesh_conversion/mesh/3d/CC_new.ucd rename to deal.II/contrib/mesh_conversion/mesh/3d/CC_cubit_new.ucd index e7f984891e..cc7b82890d 100755 --- a/deal.II/contrib/mesh_conversion/mesh/3d/CC_new.ucd +++ b/deal.II/contrib/mesh_conversion/mesh/3d/CC_cubit_new.ucd @@ -1,6 +1,6 @@ # FEM Mesh Converter # Mesh type: AVS UCD -# Input file name: CC_cubit_new.inp +# Input file name: mesh/3d/CC_cubit_new.inp 54 48 0 0 0 1 -5.00000000e-01 -5.00000000e-01 5.00000000e-01 2 -5.00000000e-01 -5.00000000e-01 0.00000000e+00 diff --git a/deal.II/contrib/mesh_conversion/mesh/3d/CC_old.ucd b/deal.II/contrib/mesh_conversion/mesh/3d/CC_cubit_old.ucd similarity index 98% rename from deal.II/contrib/mesh_conversion/mesh/3d/CC_old.ucd rename to deal.II/contrib/mesh_conversion/mesh/3d/CC_cubit_old.ucd index 1c134305d2..2b5157374f 100755 --- a/deal.II/contrib/mesh_conversion/mesh/3d/CC_old.ucd +++ b/deal.II/contrib/mesh_conversion/mesh/3d/CC_cubit_old.ucd @@ -1,6 +1,6 @@ # FEM Mesh Converter # Mesh type: AVS UCD -# Input file name: CC_cubit_old.inp +# Input file name: mesh/3d/CC_cubit_old.inp 54 48 0 0 0 1 -5.00000000e-01 -5.00000000e-01 5.00000000e-01 2 -5.00000000e-01 -5.00000000e-01 0.00000000e+00 diff --git a/deal.II/contrib/mesh_conversion/mesh/3d/other_simple.inp b/deal.II/contrib/mesh_conversion/mesh/3d/other_simple.inp new file mode 100644 index 0000000000..1250a41207 --- /dev/null +++ b/deal.II/contrib/mesh_conversion/mesh/3d/other_simple.inp @@ -0,0 +1,17 @@ +*NODE +1, 0, 0, 0 +2, 0, 0, 50 +3, 50, 0, 50 +4, 50, 0, 0 +5, 0, 50, 0 +6, 0.705725, 50, 50 +7, 50, 50, 50 +8, 50, 50, 0 +** +** +*ELEMENT, TYPE=C3D8 +1, 1, 2, 3, 4, 5, 6, 7, 8 +*NSET, NSET=n, GENERATE +1, 8 +*ELSET, ELSET=e, GENERATE +1, 1 diff --git a/deal.II/contrib/mesh_conversion/mesh/3d/other_simple.ucd b/deal.II/contrib/mesh_conversion/mesh/3d/other_simple.ucd new file mode 100644 index 0000000000..073dd734a5 --- /dev/null +++ b/deal.II/contrib/mesh_conversion/mesh/3d/other_simple.ucd @@ -0,0 +1,13 @@ +# FEM Mesh Converter +# Mesh type: AVS UCD +# Input file name: mesh/3d/other_simple.inp +8 1 0 0 0 +1 0.00000000e+00 0.00000000e+00 0.00000000e+00 +2 0.00000000e+00 0.00000000e+00 5.00000000e+01 +3 5.00000000e+01 0.00000000e+00 5.00000000e+01 +4 5.00000000e+01 0.00000000e+00 0.00000000e+00 +5 0.00000000e+00 5.00000000e+01 0.00000000e+00 +6 7.05725000e-01 5.00000000e+01 5.00000000e+01 +7 5.00000000e+01 5.00000000e+01 5.00000000e+01 +8 5.00000000e+01 5.00000000e+01 0.00000000e+00 +1 0 hex 1 2 3 4 5 6 7 8 -- 2.39.5