]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Restore formatting options in Abaqus_to_UCD
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Mon, 23 Apr 2018 14:06:17 +0000 (16:06 +0200)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Mon, 23 Apr 2018 21:16:47 +0000 (23:16 +0200)
source/grid/grid_in.cc

index e5b389b91e3dc349ffa6a1c3de690050697c9db6..d46cdc11b7b81f4752ec1dffe2df13306b1315ff 100644 (file)
@@ -23,6 +23,8 @@
 #include <deal.II/grid/grid_reordering.h>
 #include <deal.II/grid/grid_tools.h>
 
+#include <boost/io/ios_state.hpp>
+
 #include <map>
 #include <algorithm>
 #include <fstream>
@@ -3533,6 +3535,9 @@ cont:
 
     AssertThrow (output, ExcIO());
 
+    // save old formatting options
+    const boost::io::ios_base_all_saver formatting_saver(output);
+
     // Write out title - Note: No other commented text can be inserted below the
     // title in a UCD file
     output << "# Abaqus to UCD mesh conversion" << std::endl;
@@ -3561,54 +3566,47 @@ cont:
     //    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;
+    output << node_list.size() << "\t"
+           << (cell_list.size() + face_list.size()) << "\t0\t0\t0"
+           << std::endl;
+
+    output.width (16);
+    output.precision (8);
 
     // Write out node numbers
-    for (unsigned int ii = 0; ii < node_list.size(); ++ii) // Loop over all nodes
+    // Loop over all nodes
+    for (unsigned int ii = 0; ii < node_list.size(); ++ii)
       {
-        for (unsigned int jj = 0; jj < dim + 1; ++jj) // Loop over entries to be outputted
+        // Node number
+        output << node_list[ii][0] << "\t";
+
+        // Node coordinates
+        output.setf (std::ios::scientific, std::ios::floatfield);
+        for (unsigned int jj = 1; jj < dim + 1; ++jj)
           {
-            if (jj == 0)        // Node number
-              {
-                output.precision();
-                output << node_list[ii][jj] << "\t";
-              }
-            else                // Node coordinates
-              {
-                output.width (16);
-                output.setf (std::ios::scientific,
-                             std::ios::floatfield);
-                output.precision (8);
-                if (std::abs (node_list[ii][jj]) > tolerance) // invoke tolerance -> set points close to zero equal to zero
-                  output << static_cast<double> (node_list[ii][jj]) << "\t";
-                else
-                  output << 0.0 << "\t";
-              }
+            // invoke tolerance -> set points close to zero equal to zero
+            if (std::abs (node_list[ii][jj]) > tolerance)
+              output << static_cast<double> (node_list[ii][jj]) << "\t";
+            else
+              output << 0.0 << "\t";
           }
         if (dim == 2)
           output << 0.0 << "\t";
 
-        output
-            << std::endl;
+        output << std::endl;
         output.unsetf (std::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"
-            << (dim == 2 ? "quad" : "hex") << "\t";
+        output << ii + 1 << "\t"
+               << cell_list[ii][0] << "\t"
+               << (dim == 2 ? "quad" : "hex") << "\t";
         for (unsigned int jj = 1; jj < GeometryInfo<dim>::vertices_per_cell + 1; ++jj)
-          output
-              << cell_list[ii][jj] << "\t";
+          output << cell_list[ii][jj] << "\t";
 
-        output
-            << std::endl;
+        output << std::endl;
       }
 
     // Write out quad node numbers
@@ -3619,11 +3617,9 @@ cont:
             << face_list[ii][0] << "\t"
             << (dim == 2 ? "line" : "quad") << "\t";
         for (unsigned int jj = 1; jj < GeometryInfo<dim>::vertices_per_face + 1; ++jj)
-          output
-              << face_list[ii][jj] << "\t";
+          output << face_list[ii][jj] << "\t";
 
-        output
-            << std::endl;
+        output << std::endl;
       }
   }
 }

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.