]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Move the initialize function to the base class.
authorDavid Wells <drwells@email.unc.edu>
Sat, 21 May 2022 19:29:18 +0000 (15:29 -0400)
committerDavid Wells <drwells@email.unc.edu>
Sat, 21 May 2022 21:10:00 +0000 (17:10 -0400)
include/deal.II/base/function_parser.h
include/deal.II/base/mu_parser_internal.h
include/deal.II/base/tensor_function_parser.h
source/base/function_parser.cc
source/base/mu_parser_internal.cc
source/base/tensor_function_parser.cc

index 4564072f3e0814331a1303547951d13cd2d3833c..cf77606c7913a574dbb606d33feb53468aa6f508 100644 (file)
@@ -308,11 +308,11 @@ public:
    * this case is dim+1. The value of this parameter defaults to false, i.e.,
    * do not consider time.
    */
-  void
+  virtual void
   initialize(const std::string &             vars,
              const std::vector<std::string> &expressions,
              const ConstMap &                constants,
-             const bool                      time_dependent = false);
+             const bool                      time_dependent = false) override;
 
   /**
    * Initialize the function. Same as above, but accepts a string rather than
index 562cc6025e879c12822c8d42b43152d63c99ed29..040a49415207fe22a1b84e3c7e83a81ff1a37168 100644 (file)
@@ -164,11 +164,15 @@ namespace internal
       {}
 
       /**
-       * The muParser objects (hidden with the PIMPL idiom) for each thread (and
-       * one for each component).
+       * Initialize the internal state of the object. This is the same as the
+       * inheriting class method - see FunctionParser::initialize() for more
+       * information.
        */
-      mutable Threads::ThreadLocalStorage<internal::FunctionParser::ParserData>
-        parser_data;
+      virtual void
+      initialize(const std::string &                  vars,
+                 const std::vector<std::string> &     expressions,
+                 const std::map<std::string, double> &constants,
+                 const bool                           time_dependent = false);
 
       void
       init_muparser() const;
index 1045c106ff8593445e9f6e82cc5f723bc6903a13..93d9c84e850ac6e572728aa0c19c5f823a916fa5 100644 (file)
@@ -196,11 +196,11 @@ public:
    * method in this case is dim+1. The value of this parameter defaults to
    * false, i.e. do not consider time.
    */
-  void
+  virtual void
   initialize(const std::string &             vars,
              const std::vector<std::string> &expressions,
              const ConstMap &                constants,
-             const bool                      time_dependent = false);
+             const bool                      time_dependent = false) override;
 
   /**
    * Initialize the function. Same as above, but accepts a string rather than
index 13c0812bb00405c8c686f08d73692ffd0dbb5eb8..bec8278c5228983130760fb5d10c15db71b23ca1 100644 (file)
 #include <deal.II/base/function_parser.h>
 #include <deal.II/base/mu_parser_internal.h>
 #include <deal.II/base/patterns.h>
-#include <deal.II/base/thread_management.h>
 #include <deal.II/base/utilities.h>
 
 #include <deal.II/lac/vector.h>
 
-#include <cmath>
 #include <map>
 
 DEAL_II_NAMESPACE_OPEN
@@ -78,38 +76,10 @@ FunctionParser<dim>::initialize(const std::string &             variables,
                                 const std::map<std::string, double> &constants,
                                 const bool time_dependent)
 {
-  this->parser_data.clear(); // this will reset all thread-local objects
-
-  this->constants   = constants;
-  this->var_names   = Utilities::split_string_list(variables, ',');
-  this->expressions = expressions;
-  AssertThrow(((time_dependent) ? dim + 1 : dim) == this->var_names.size(),
-              ExcMessage("Wrong number of variables"));
-
-  // We check that the number of components of this function matches the
-  // number of components passed in as a vector of strings.
   AssertThrow(this->n_components == expressions.size(),
               ExcInvalidExpressionSize(this->n_components, expressions.size()));
-
-  // Now we define how many variables we expect to read in.  We distinguish
-  // between two cases: Time dependent problems, and not time dependent
-  // problems. In the first case the number of variables is given by the
-  // dimension plus one. In the other case, the number of variables is equal
-  // to the dimension. Once we parsed the variables string, if none of this is
-  // the case, then an exception is thrown.
-  if (time_dependent)
-    this->n_vars = dim + 1;
-  else
-    this->n_vars = dim;
-
-  // create a parser object for the current thread we can then query in
-  // value() and vector_value(). this is not strictly necessary because a user
-  // may never call these functions on the current thread, but it gets us
-  // error messages about wrong formulas right away
-  this->init_muparser();
-
-  // finally set the initialization bit
-  this->initialized = true;
+  internal::FunctionParser::ParserImplementation<dim, double>::initialize(
+    variables, expressions, constants, time_dependent);
 }
 
 
index 93c6470bec9d80ca27e4e926ba64c78bd0c383eb..07bf42886c81ff9e10da3fad0965d70af9d73362 100644 (file)
@@ -15,6 +15,7 @@
 
 #include <deal.II/base/mu_parser_internal.h>
 #include <deal.II/base/thread_management.h>
+#include <deal.II/base/utilities.h>
 
 #include <cmath>
 #include <ctime>
@@ -223,6 +224,44 @@ namespace internal
 
 
 
+
+    template <int dim, typename Number>
+    void
+    ParserImplementation<dim, Number>::initialize(
+      const std::string &                  variables,
+      const std::vector<std::string> &     expressions,
+      const std::map<std::string, double> &constants,
+      const bool                           time_dependent)
+    {
+      this->parser_data.clear(); // this will reset all thread-local objects
+
+      this->constants   = constants;
+      this->var_names   = Utilities::split_string_list(variables, ',');
+      this->expressions = expressions;
+      AssertThrow(((time_dependent) ? dim + 1 : dim) == this->var_names.size(),
+                  ExcMessage("Wrong number of variables"));
+
+      // Now we define how many variables we expect to read in. We distinguish
+      // between two cases: Time dependent problems, and not time dependent
+      // problems. In the first case the number of variables is given by the
+      // dimension plus one. In the other case, the number of variables is equal
+      // to the dimension. Once we parsed the variables string, if none of this
+      // is the case, then an exception is thrown.
+      if (time_dependent)
+        this->n_vars = dim + 1;
+      else
+        this->n_vars = dim;
+
+      // create a parser object for the current thread we can then query in
+      // value() and vector_value(). this is not strictly necessary because a
+      // user may never call these functions on the current thread, but it gets
+      // us error messages about wrong formulas right away
+      this->init_muparser();
+      this->initialized = true;
+    }
+
+
+
     template <int dim, typename Number>
     void
     ParserImplementation<dim, Number>::init_muparser() const
index 1584d369ff0225cb0ff18ebe03faba4d58f993ac..08230dca191d042a318281db551a43e73e3be693 100644 (file)
 #include <deal.II/base/tensor.h>
 #include <deal.II/base/tensor_function.h>
 #include <deal.II/base/tensor_function_parser.h>
-#include <deal.II/base/thread_management.h>
 #include <deal.II/base/utilities.h>
 
-#include <cmath>
 #include <map>
 
 DEAL_II_NAMESPACE_OPEN
@@ -79,49 +77,10 @@ TensorFunctionParser<rank, dim, Number>::initialize(
   const std::map<std::string, double> &constants,
   const bool                           time_dependent)
 {
-  this->parser_data.clear(); // this will reset all thread-local objects
-
-  this->constants   = constants;
-  this->var_names   = Utilities::split_string_list(variables, ',');
-  this->expressions = expressions;
-  AssertThrow(((time_dependent) ? dim + 1 : dim) == this->var_names.size(),
-              ExcMessage("Wrong number of variables"));
-
-  // We check that the number of
-  // components of this function
-  // matches the number of components
-  // passed in as a vector of
-  // strings.
   AssertThrow(this->n_components == expressions.size(),
               ExcInvalidExpressionSize(this->n_components, expressions.size()));
-
-  // Now we define how many variables
-  // we expect to read in.  We
-  // distinguish between two cases:
-  // Time dependent problems, and not
-  // time dependent problems. In the
-  // first case the number of
-  // variables is given by the
-  // dimension plus one. In the other
-  // case, the number of variables is
-  // equal to the dimension. Once we
-  // parsed the variables string, if
-  // none of this is the case, then
-  // an exception is thrown.
-  if (time_dependent)
-    this->n_vars = dim + 1;
-  else
-    this->n_vars = dim;
-
-  // create a parser object for the current thread we can then query
-  // in value() and vector_value(). this is not strictly necessary
-  // because a user may never call these functions on the current
-  // thread, but it gets us error messages about wrong formulas right
-  // away
-  this->init_muparser();
-
-  // finally set the initialization bit
-  this->initialized = true;
+  internal::FunctionParser::ParserImplementation<dim, Number>::initialize(
+    variables, expressions, constants, time_dependent);
 }
 
 

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.