]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Prefer std::string interfaces over const char * for public interfaces 7436/head
authorDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Fri, 9 Nov 2018 18:00:18 +0000 (19:00 +0100)
committerDaniel Arndt <daniel.arndt@iwr.uni-heidelberg.de>
Sat, 10 Nov 2018 22:51:20 +0000 (23:51 +0100)
13 files changed:
include/deal.II/base/data_out_base.h
include/deal.II/base/event.h
include/deal.II/base/job_identifier.h
include/deal.II/base/parameter_handler.h
include/deal.II/distributed/tria.h
include/deal.II/lac/scalapack.h
source/base/data_out_base.cc
source/base/event.cc
source/base/function_parser.cc
source/base/job_identifier.cc
source/base/parameter_handler.cc
source/distributed/tria.cc
source/lac/scalapack.cc

index 8b9fa93e5278515faaba627239c3e946d9155f62..eea93b567ce69ae6d4b879ef23619ea3e651e86d 100644 (file)
@@ -2963,7 +2963,7 @@ public:
    * DataOutInterface::write_vtu().
    */
   void
-  write_vtu_in_parallel(const char *filename, MPI_Comm comm) const;
+  write_vtu_in_parallel(const std::string &filename, MPI_Comm comm) const;
 
   /**
    * Some visualization programs, such as ParaView, can read several separate
index 43db4140d77c21ff4d12109ed83c75ac3e6dbfe5..d567c163e33e9e48ed9f8f24ca173328f6588f82 100644 (file)
@@ -54,7 +54,7 @@ namespace Algorithms
      * later use.
      */
     static Event
-    assign(const char *name);
+    assign(const std::string &name);
 
     /**
      * If you forgot to store the result of assign, here is how to retrieve it
index 94f3079f368af1e7e8f164caff4a98fd19b4c114..7763c5e110878a396b1ee4d6f89760af68f7a6d3 100644 (file)
@@ -58,7 +58,7 @@ public:
    * identifier for the program being run.
    */
   static std::string
-  base_name(const char *filename);
+  base_name(const std::string &filename);
 
   /**
    * Return the value of <tt>id</tt>.
index 6d088bd13c4c825bf99603f73cdf646f2deb4821..37ecbbf3542352b5604f46f0c303ea8dd8d841af 100644 (file)
@@ -1009,7 +1009,7 @@ public:
    * there for more information.
    */
   virtual void
-  parse_input_from_string(const char *       s,
+  parse_input_from_string(const std::string &s,
                           const std::string &last_line      = "",
                           const bool         skip_undefined = false);
 
index 2c0888764fd7203e428d723109ae8c37503afc61..161fb6443132dfb2938f6333aa2c4793049ee3e3 100644 (file)
@@ -587,7 +587,7 @@ namespace parallel
        * interface between deal.II and p4est.
        */
       void
-      write_mesh_vtk(const char *file_basename) const;
+      write_mesh_vtk(const std::string &file_basename) const;
 
       /**
        * Produce a check sum of the triangulation.  This is a collective
@@ -604,7 +604,7 @@ namespace parallel
        * cell-based data can be saved using register_data_attach().
        */
       void
-      save(const char *filename) const;
+      save(const std::string &filename) const;
 
       /**
        * Load the refinement information saved with save() back in. The mesh
@@ -625,7 +625,7 @@ namespace parallel
        * if a different number of MPI processes is encountered).
        */
       void
-      load(const char *filename, const bool autopartition = true);
+      load(const std::string &filename, const bool autopartition = true);
 
       /**
        * Register a function that can be used to attach data of fixed size
@@ -1043,8 +1043,8 @@ namespace parallel
          */
         void
         save(const typename dealii::internal::p4est::types<dim>::forest
-               *         parallel_forest,
-             const char *filename) const;
+               *                parallel_forest,
+             const std::string &filename) const;
 
         /**
          * Transfer data from file system.
@@ -1067,7 +1067,7 @@ namespace parallel
         void
         load(const typename dealii::internal::p4est::types<dim>::forest
                *                parallel_forest,
-             const char *       filename,
+             const std::string &filename,
              const unsigned int n_attached_deserialize_fixed,
              const unsigned int n_attached_deserialize_variable);
 
@@ -1302,14 +1302,14 @@ namespace parallel
        * compiler.
        */
       void
-      load(const char *filename, const bool autopartition = true);
+      load(const std::string &filename, const bool autopartition = true);
 
       /**
        * This function is not implemented, but needs to be present for the
        * compiler.
        */
       void
-      save(const char *filename) const;
+      save(const std::string &filename) const;
 
       /**
        * This function is not implemented, but needs to be present for the
index 5546f6b0aa6579fdce655f3d38d89b73448fa9eb..93aecaa8e471d50cc07f2bfe6e2b853066a8a931 100644 (file)
@@ -462,17 +462,6 @@ public:
          std::make_pair(numbers::invalid_unsigned_int,
                         numbers::invalid_unsigned_int)) const;
 
-  /**
-   * Same as above but takes a filename as a <code>char</code> pointer.
-   *
-   * @deprecated This function is deprecated. Use the one taking <code>std::string</code> instead.
-   */
-  DEAL_II_DEPRECATED void
-  save(const char *                                 filename,
-       const std::pair<unsigned int, unsigned int> &chunk_size =
-         std::make_pair(numbers::invalid_unsigned_int,
-                        numbers::invalid_unsigned_int)) const;
-
   /**
    * Loads the distributed matrix from file @p filename using HDF5.
    * In case that deal.II was built without HDF5
@@ -487,14 +476,6 @@ public:
   void
   load(const std::string &filename);
 
-  /**
-   * Same as above but takes a filename as a <code>char</code> pointer.
-   *
-   * @deprecated This function is deprecated. Use the one taking <code>std::string</code> instead.
-   */
-  DEAL_II_DEPRECATED void
-  load(const char *filename);
-
   /**
    * Compute the Cholesky factorization of the matrix using ScaLAPACK
    * function <code>pXpotrf</code>. The result of the factorization is stored in
index 8947f951640f248bf813f2c0b13eca29b9c1b595..afee5a9126e438d43e982e3c7f9fc14632f2e184 100644 (file)
@@ -7552,8 +7552,9 @@ DataOutInterface<dim, spacedim>::write_svg(std::ostream &out) const
 
 template <int dim, int spacedim>
 void
-DataOutInterface<dim, spacedim>::write_vtu_in_parallel(const char *filename,
-                                                       MPI_Comm    comm) const
+DataOutInterface<dim, spacedim>::write_vtu_in_parallel(
+  const std::string &filename,
+  MPI_Comm           comm) const
 {
 #ifndef DEAL_II_WITH_MPI
   // without MPI fall back to the normal way to write a vtu file:
@@ -7570,7 +7571,7 @@ DataOutInterface<dim, spacedim>::write_vtu_in_parallel(const char *filename,
   AssertThrowMPI(ierr);
   MPI_File fh;
   ierr = MPI_File_open(comm,
-                       const_cast<char *>(filename),
+                       DEAL_II_MPI_CONST_CAST(filename.c_str()),
                        MPI_MODE_CREATE | MPI_MODE_WRONLY,
                        info,
                        &fh);
@@ -7594,7 +7595,7 @@ DataOutInterface<dim, spacedim>::write_vtu_in_parallel(const char *filename,
       DataOutBase::write_vtu_header(ss, vtk_flags);
       header_size = ss.str().size();
       ierr = MPI_File_write(fh,
-                            const_cast<char *>(ss.str().c_str()),
+                            DEAL_II_MPI_CONST_CAST(ss.str().c_str()),
                             header_size,
                             MPI_CHAR,
                             MPI_STATUS_IGNORE);
@@ -7614,7 +7615,7 @@ DataOutInterface<dim, spacedim>::write_vtu_in_parallel(const char *filename,
                                 vtk_flags,
                                 ss);
     ierr = MPI_File_write_ordered(fh,
-                                  const_cast<char *>(ss.str().c_str()),
+                                  DEAL_II_MPI_CONST_CAST(ss.str().c_str()),
                                   ss.str().size(),
                                   MPI_CHAR,
                                   MPI_STATUS_IGNORE);
@@ -7628,7 +7629,7 @@ DataOutInterface<dim, spacedim>::write_vtu_in_parallel(const char *filename,
       DataOutBase::write_vtu_footer(ss);
       unsigned int footer_size = ss.str().size();
       ierr = MPI_File_write_shared(fh,
-                                   const_cast<char *>(ss.str().c_str()),
+                                   DEAL_II_MPI_CONST_CAST(ss.str().c_str()),
                                    footer_size,
                                    MPI_CHAR,
                                    MPI_STATUS_IGNORE);
index 6187f4b5d04c2ec89edd676054ed4893b5b6f1da..440c703b11ca4cb5269549c881438ba56efd0f0d 100644 (file)
@@ -25,7 +25,7 @@ namespace Algorithms
   std::vector<std::string> Event::names;
 
   Event
-  Event::assign(const char *name)
+  Event::assign(const std::string &name)
   {
     unsigned int index = names.size();
     names.emplace_back(name);
index 769c681851adb60cd2adad0491c7c0a709f6dbdd..f1c30add30df5feddd03a6398f23bedcc2812929 100644 (file)
@@ -296,12 +296,11 @@ FunctionParser<dim>::init_muparser() const
            constant != constants.end();
            ++constant)
         {
-          fp.get()[component]->DefineConst(constant->first.c_str(),
-                                           constant->second);
+          fp.get()[component]->DefineConst(constant->first, constant->second);
         }
 
       for (unsigned int iv = 0; iv < var_names.size(); ++iv)
-        fp.get()[component]->DefineVar(var_names[iv].c_str(), &vars.get()[iv]);
+        fp.get()[component]->DefineVar(var_names[iv], &vars.get()[iv]);
 
       // define some compatibility functions:
       fp.get()[component]->DefineFun("if", internal::mu_if, true);
@@ -411,7 +410,7 @@ FunctionParser<dim>::init_muparser() const
           std::cerr << "Token:    <" << e.GetToken() << ">\n";
           std::cerr << "Position: <" << e.GetPos() << ">\n";
           std::cerr << "Errc:     <" << e.GetCode() << ">" << std::endl;
-          AssertThrow(false, ExcParseError(e.GetCode(), e.GetMsg().c_str()));
+          AssertThrow(false, ExcParseError(e.GetCode(), e.GetMsg()));
         }
     }
 }
@@ -462,7 +461,7 @@ FunctionParser<dim>::value(const Point<dim> & p,
       std::cerr << "Token:    <" << e.GetToken() << ">\n";
       std::cerr << "Position: <" << e.GetPos() << ">\n";
       std::cerr << "Errc:     <" << e.GetCode() << ">" << std::endl;
-      AssertThrow(false, ExcParseError(e.GetCode(), e.GetMsg().c_str()));
+      AssertThrow(false, ExcParseError(e.GetCode(), e.GetMsg()));
       return 0.0;
     }
 }
index e7fc49573f33115b31be3daa03a3d486c7f506f4..ef0f1c2ccff088c450f2aa3abf1029afe33f7c7c 100644 (file)
@@ -52,7 +52,7 @@ JobIdentifier::operator()() const
 
 
 std::string
-JobIdentifier::base_name(const char *filename)
+JobIdentifier::base_name(const std::string &filename)
 {
   std::string            name(filename);
   std::string::size_type pos;
index 925b64073f2e67266e9b2bb3e062dc6d0a87da84..28987ba50cf02b14e0b613e0bbef8bdd6682a9ed 100644 (file)
@@ -461,7 +461,7 @@ ParameterHandler::parse_input(const std::string &filename,
 
 
 void
-ParameterHandler::parse_input_from_string(const char *       s,
+ParameterHandler::parse_input_from_string(const std::string &s,
                                           const std::string &last_line,
                                           const bool         skip_undefined)
 {
index 9d4d9727be6cac0563522a397e0e56ace748be6a..5e0611472a2cf0dff34fcb014f6eca6a48322e04 100644 (file)
@@ -1788,8 +1788,8 @@ namespace parallel
     void
     Triangulation<dim, spacedim>::DataTransfer::save(
       const typename dealii::internal::p4est::types<dim>::forest
-        *         parallel_forest,
-      const char *filename) const
+        *                parallel_forest,
+      const std::string &filename) const
     {
       // Large fractions of this function have been copied from
       // DataOutInterface::write_vtu_in_parallel.
@@ -1964,7 +1964,7 @@ namespace parallel
     Triangulation<dim, spacedim>::DataTransfer::load(
       const typename dealii::internal::p4est::types<dim>::forest
         *                parallel_forest,
-      const char *       filename,
+      const std::string &filename,
       const unsigned int n_attached_deserialize_fixed,
       const unsigned int n_attached_deserialize_variable)
     {
@@ -2689,20 +2689,19 @@ namespace parallel
     template <int dim, int spacedim>
     void
     Triangulation<dim, spacedim>::write_mesh_vtk(
-      const char *file_basename) const
+      const std::string &file_basename) const
     {
       Assert(parallel_forest != nullptr,
              ExcMessage("Can't produce output when no forest is created yet."));
-      dealii::internal::p4est::functions<dim>::vtk_write_file(parallel_forest,
-                                                              nullptr,
-                                                              file_basename);
+      dealii::internal::p4est::functions<dim>::vtk_write_file(
+        parallel_forest, nullptr, file_basename.c_str());
     }
 
 
 
     template <int dim, int spacedim>
     void
-    Triangulation<dim, spacedim>::save(const char *filename) const
+    Triangulation<dim, spacedim>::save(const std::string &filename) const
     {
       Assert(
         cell_attached_data.n_attached_deserialize == 0,
@@ -2757,7 +2756,7 @@ namespace parallel
           tria->data_transfer.clear();
         }
 
-      dealii::internal::p4est::functions<dim>::save(filename,
+      dealii::internal::p4est::functions<dim>::save(filename.c_str(),
                                                     parallel_forest,
                                                     false);
 
@@ -2778,8 +2777,8 @@ namespace parallel
 
     template <int dim, int spacedim>
     void
-    Triangulation<dim, spacedim>::load(const char *filename,
-                                       const bool  autopartition)
+    Triangulation<dim, spacedim>::load(const std::string &filename,
+                                       const bool         autopartition)
     {
       Assert(
         this->n_cells() > 0,
@@ -2827,7 +2826,7 @@ namespace parallel
         attached_count_fixed + attached_count_variable;
 
       parallel_forest = dealii::internal::p4est::functions<dim>::load_ext(
-        filename,
+        filename.c_str(),
         this->mpi_communicator,
         0,
         false,
@@ -5070,7 +5069,7 @@ namespace parallel
 
     template <int spacedim>
     void
-    Triangulation<1, spacedim>::load(const char *, const bool)
+    Triangulation<1, spacedim>::load(const std::string &, const bool)
     {
       Assert(false, ExcNotImplemented());
     }
@@ -5079,7 +5078,7 @@ namespace parallel
 
     template <int spacedim>
     void
-    Triangulation<1, spacedim>::save(const char *) const
+    Triangulation<1, spacedim>::save(const std::string &) const
     {
       Assert(false, ExcNotImplemented());
     }
index be71ea82a014834df794849a0fb23c7bc39bf57d..8a991335137606999917a2999db192c96024a3f3 100644 (file)
@@ -2649,17 +2649,6 @@ ScaLAPACKMatrix<NumberType>::save(
 
 
 
-template <typename NumberType>
-void
-ScaLAPACKMatrix<NumberType>::save(
-  const char *                                 filename,
-  const std::pair<unsigned int, unsigned int> &chunk_size) const
-{
-  save(std::string(filename), chunk_size);
-}
-
-
-
 template <typename NumberType>
 void
 ScaLAPACKMatrix<NumberType>::save_serial(
@@ -3073,15 +3062,6 @@ ScaLAPACKMatrix<NumberType>::load(const std::string &filename)
 
 
 
-template <typename NumberType>
-void
-ScaLAPACKMatrix<NumberType>::load(const char *filename)
-{
-  load(std::string(filename));
-}
-
-
-
 template <typename NumberType>
 void
 ScaLAPACKMatrix<NumberType>::load_serial(const std::string &filename)

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.