]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove deprecated functions FunctionParser::initialize().
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Wed, 31 Dec 2014 20:27:15 +0000 (14:27 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 4 Jan 2015 21:13:47 +0000 (15:13 -0600)
doc/news/changes.h
include/deal.II/base/function_parser.h
include/deal.II/base/mg_level_object.h
source/base/function_parser.cc
tests/base/function_parser_02.cc
tests/trilinos/mg_transfer_prebuilt_01.cc

index f721aa59b34737ece4beac576c163252e2ba4831..2e041efdba55287d570b791d747e788098f62f20 100644 (file)
@@ -61,6 +61,12 @@ inconvenience this causes.
   - DataOutBase::create_xdmf_entry with 3 arguments.
   - Algorithms::ThetaTimestepping::operator().
   - Algorithms::Newton::initialize.
+  - MGLevelObject::get_minlevel and MGLevelObject::get_maxlevel.
+  - The versions of FunctionParser::initialize that took a
+    <code>use_degrees</code> or <code>constants</code> argument.
+    The implementation as it is now no longer supports either of
+    these two concepts (since we switched from the FunctionParser
+    library to the muparser library after the deal.II 8.1 release).
 </ol>
 
   <li> Removed: The config.h file no longer exports HAVE_* definitions.
index eede26438e3c28000f2056e163bb8662d834a814..701471d2ac8d7f9aabee2c7180e642721006327c 100644 (file)
@@ -226,7 +226,6 @@ public:
    * the FunctionParser, as declared in the constructor. If this is not the
    * case, an exception is thrown.
    *
-   *
    * <b>constants</b>: a map of constants used to pass any necessary constant
    * that we want to specify in our expressions (in the example above the
    * number pi). An expression is valid if and only if it contains only
@@ -246,36 +245,6 @@ public:
                    const ConstMap                 &constants,
                    const bool time_dependent = false);
 
-  /**
-   * Same as above, but with an additional parameter: <b>use_degrees</b>.
-   * Parameter to decide if the trigonometric functions work in radians or
-   * degrees. The default for this parameter is false, i.e. use radians and
-   * not degrees.
-   *
-   * @note: this function is deprecated. Use the function without this
-   * argument instead (which has the default use_degrees=false).
-   */
-  void initialize (const std::string              &vars,
-                   const std::vector<std::string> &expressions,
-                   const ConstMap                 &constants,
-                   const bool time_dependent,
-                   const bool use_degrees) DEAL_II_DEPRECATED;
-
-
-  /**
-   * Initialize the function. Same as above, but with an additional argument
-   * <b> units </b> - a map of units passed to FunctionParser via AddUnint.
-   *
-   * Can be used as "3cm". Have higher precedence in parsing, i.e. if cm=10
-   * then 3/2cm is 3 /(2*10).
-   */
-  void initialize (const std::string              &vars,
-                   const std::vector<std::string> &expressions,
-                   const ConstMap                 &constants,
-                   const ConstMap                 &units,
-                   const bool time_dependent = false,
-                   const bool use_degrees = false) DEAL_II_DEPRECATED;
-
   /**
    * Initialize the function. Same as above, but accepts a string rather than
    * a vector of strings. If this is a vector valued function, its components
@@ -288,31 +257,6 @@ public:
                    const ConstMap    &constants,
                    const bool time_dependent = false);
 
-  /**
-   * Same as above, but with an additional parameter: <b>use_degrees</b>.
-   * Parameter to decide if the trigonometric functions work in radians or
-   * degrees. The default for this parameter is false, i.e. use radians and
-   * not degrees.
-   *
-   * @note: this function is deprecated. Use the function without this
-   * argument instead (which has the default use_degrees=false).
-   */
-  void initialize (const std::string &vars,
-                   const std::string &expression,
-                   const ConstMap    &constants,
-                   const bool time_dependent,
-                   const bool use_degrees) DEAL_II_DEPRECATED;
-
-  /**
-   * Initialize the function. Same as above, but with <b>units</b>.
-   */
-  void initialize (const std::string &vars,
-                   const std::string &expression,
-                   const ConstMap    &constants,
-                   const ConstMap    &units,
-                   const bool time_dependent = false,
-                   const bool use_degrees = false) DEAL_II_DEPRECATED;
-
   /**
    * A function that returns default names for variables, to be used in the
    * first argument of the initialize() functions: it returns "x" in 1d, "x,y"
index 3b13118cce6a8bf6aff5227762d5d1e9dea13a51..beea1c13f61d02dec5c3278969316556cdc302fa 100644 (file)
@@ -88,16 +88,6 @@ public:
    */
   unsigned int max_level () const;
 
-  /**
-   * @deprecated Replaced by min_level()
-   */
-  unsigned int get_minlevel () const DEAL_II_DEPRECATED;
-
-  /**
-   * @deprecated Replaced by max_level()
-   */
-  unsigned int get_maxlevel () const DEAL_II_DEPRECATED;
-
   /**
    * Memory used by this object.
    */
@@ -188,22 +178,6 @@ MGLevelObject<Object>::clear ()
 }
 
 
-template<class Object>
-unsigned int
-MGLevelObject<Object>::get_minlevel () const
-{
-  return minlevel;
-}
-
-
-template<class Object>
-unsigned int
-MGLevelObject<Object>::get_maxlevel () const
-{
-  return minlevel + objects.size() - 1;
-}
-
-
 template<class Object>
 unsigned int
 MGLevelObject<Object>::min_level () const
index 1182bd31d9da58f76c88eb86b616b81c40d44bb2..bc8e197e89e2add8d8653c7386b0f11f4dfa35c5 100644 (file)
@@ -48,29 +48,57 @@ FunctionParser<dim>::~FunctionParser()
 {}
 
 #ifdef DEAL_II_WITH_MUPARSER
-template <int dim>
-void FunctionParser<dim>::initialize (const std::string                   &variables,
-                                      const std::vector<std::string>      &expressions,
-                                      const std::map<std::string, double> &constants,
-                                      const bool time_dependent,
-                                      const bool use_degrees)
-{
-  initialize (variables,
-              expressions,
-              constants,
-              std::map< std::string, double >(),
-              time_dependent,
-              use_degrees);
-}
-
 
 template <int dim>
-void FunctionParser<dim>::initialize (const std::string              &vars,
+void FunctionParser<dim>::initialize (const std::string              &variables,
                                       const std::vector<std::string> &expressions,
                                       const std::map<std::string, double> &constants,
                                       const bool time_dependent)
 {
-  initialize(vars, expressions, constants, time_dependent, false);
+  this->fp.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) == 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)
+    n_vars = dim+1;
+  else
+    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
+  init_muparser ();
+
+  // finally set the initialization bit
+  initialized = true;
 }
 
 
@@ -265,65 +293,6 @@ void FunctionParser<dim>:: init_muparser() const
 }
 
 
-template <int dim>
-void FunctionParser<dim>::initialize (const std::string   &variables,
-                                      const std::vector<std::string>      &expressions,
-                                      const std::map<std::string, double> &constants,
-                                      const std::map<std::string, double> &units,
-                                      const bool time_dependent,
-                                      const bool use_degrees)
-{
-  this->fp.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) == var_names.size(),
-              ExcMessage("Wrong number of variables"));
-  AssertThrow(!use_degrees, ExcNotImplemented());
-
-  // 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()) );
-
-  // we no longer support units:
-  AssertThrow(units.size()==0, ExcNotImplemented());
-
-  // 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)
-    n_vars = dim+1;
-  else
-    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
-  init_muparser ();
-
-  // finally set the initialization bit
-  initialized = true;
-}
-
-
 
 template <int dim>
 void FunctionParser<dim>::initialize (const std::string &vars,
@@ -331,45 +300,8 @@ void FunctionParser<dim>::initialize (const std::string &vars,
                                       const std::map<std::string, double> &constants,
                                       const bool time_dependent)
 {
-  initialize(vars, expression, constants, time_dependent, false);
-}
-
-
-
-template <int dim>
-void FunctionParser<dim>::initialize (const std::string &variables,
-                                      const std::string &expression,
-                                      const std::map<std::string, double> &constants,
-                                      const bool time_dependent,
-                                      const bool use_degrees)
-{
-  // initialize with the things
-  // we got.
-  initialize (variables,
-              Utilities::split_string_list (expression, ';'),
-              constants,
-              time_dependent,
-              use_degrees);
-}
-
-
-
-template <int dim>
-void FunctionParser<dim>::initialize (const std::string &variables,
-                                      const std::string &expression,
-                                      const std::map<std::string, double> &constants,
-                                      const std::map<std::string, double> &units,
-                                      const bool time_dependent,
-                                      const bool use_degrees)
-{
-  // initialize with the things
-  // we got.
-  initialize (variables,
-              Utilities::split_string_list (expression, ';'),
-              constants,
-              units,
-              time_dependent,
-              use_degrees);
+  initialize(vars, Utilities::split_string_list(expression, ';'),
+             constants, time_dependent);
 }
 
 
@@ -435,55 +367,6 @@ void FunctionParser<dim>::vector_value (const Point<dim> &p,
 #else
 
 
-template <int dim>
-void
-FunctionParser<dim>::initialize(const std::string &,
-                                const std::vector<std::string> &,
-                                const std::map<std::string, double> &,
-                                const bool,
-                                const bool)
-{
-  Assert(false, ExcNeedsFunctionparser());
-}
-
-
-template <int dim>
-void
-FunctionParser<dim>::initialize(const std::string &,
-                                const std::vector<std::string> &,
-                                const std::map<std::string, double> &,
-                                const std::map<std::string, double> &,
-                                const bool,
-                                const bool)
-{
-  Assert(false, ExcNeedsFunctionparser());
-}
-
-
-template <int dim>
-void
-FunctionParser<dim>::initialize(const std::string &,
-                                const std::string &,
-                                const std::map<std::string, double> &,
-                                const bool,
-                                const bool)
-{
-  Assert(false, ExcNeedsFunctionparser());
-}
-
-
-template <int dim>
-void
-FunctionParser<dim>::initialize(const std::string &,
-                                const std::string &,
-                                const std::map<std::string, double> &,
-                                const std::map<std::string, double> &,
-                                const bool,
-                                const bool)
-{
-  Assert(false, ExcNeedsFunctionparser());
-}
-
 template <int dim>
 void
 FunctionParser<dim>::initialize(const std::string &,
index 8dc84db112ed536869f49b63de30ec5459aaa8a4..97cc002948e6fa5b17b2b0b8899520b970b5a8dc 100644 (file)
@@ -41,7 +41,6 @@ int main ()
 
   std::vector<std::string> function(1);
   std::map<std::string, double> constants;
-  std::map<std::string, double> units;
 
   constants["PI"] = 3.141592654;
   constants["cm"] = 10;
@@ -53,7 +52,7 @@ int main ()
   FunctionParser<2> fp;
   function[0] = "x * cm + y * m + PI";
   fp.initialize(FunctionParser<2>::default_variable_names(),
-                function, constants, units);
+                function, constants);
 
   deallog << "Function " << "[" << function[0] << "]" <<
           " @point " << "[" << point << "]" << " is " <<
@@ -64,7 +63,7 @@ int main ()
   //strings
   FunctionParser<2> fp4;
   fp4.initialize(FunctionParser<2>::default_variable_names(),
-                 function[0], constants, units);
+                 function[0], constants);
 
   deallog << "Function " << "[" << function[0] << "]" <<
           " @point " << "[" << point << "]" << " is " <<
index 9db43d904f89599f3d568508b61638fd43900eee..c01191ef4a0cf7764bc434d373b75977d63c6903 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2000 - 2013 by the deal.II authors
+// Copyright (C) 2000 - 2013, 2015 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -48,8 +48,8 @@ reinit_vector (const dealii::MGDoFHandler<dim,spacedim> &mg_dof,
      (&mg_dof.get_tria()));
   AssertThrow(tria!=NULL, ExcMessage("multigrid with Trilinos vectors only works with distributed Triangulation!"));
 
-  for (unsigned int level=v.get_minlevel();
-       level<=v.get_maxlevel(); ++level)
+  for (unsigned int level=v.min_level();
+       level<=v.max_level(); ++level)
     {
       v[level].reinit(mg_dof.locally_owned_mg_dofs(level), tria->get_communicator());
     }

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.