]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove FunctionMap 10597/head
authorDaniel Arndt <arndtd@ornl.gov>
Wed, 24 Jun 2020 19:51:31 +0000 (15:51 -0400)
committerDaniel Arndt <arndtd@ornl.gov>
Wed, 24 Jun 2020 19:57:51 +0000 (15:57 -0400)
include/deal.II/dofs/deprecated_function_map.h [deleted file]
include/deal.II/dofs/dof_handler.h
include/deal.II/dofs/function_map.h [deleted file]
include/deal.II/numerics/error_estimator.h
include/deal.II/numerics/matrix_tools.h
tests/fe/fe_enriched_color_07.cc

diff --git a/include/deal.II/dofs/deprecated_function_map.h b/include/deal.II/dofs/deprecated_function_map.h
deleted file mode 100644 (file)
index e64a2ad..0000000
+++ /dev/null
@@ -1,91 +0,0 @@
-// ---------------------------------------------------------------------
-//
-// Copyright (C) 2018 - 2020 by the deal.II authors
-//
-// This file is part of the deal.II library.
-//
-// The deal.II library is free software; you can use it, redistribute
-// it, and/or modify it under the terms of the GNU Lesser General
-// Public License as published by the Free Software Foundation; either
-// version 2.1 of the License, or (at your option) any later version.
-// The full text of the license can be found in the file LICENSE.md at
-// the top level directory of deal.II.
-//
-// ---------------------------------------------------------------------
-
-#ifndef dealii_deprecated_function_map_h
-#define dealii_deprecated_function_map_h
-
-#include <deal.II/base/config.h>
-
-#include <map>
-
-DEAL_II_NAMESPACE_OPEN
-
-// Forward declaration
-#ifndef DOXYGEN
-template <int spacedim, typename RangeNumberType>
-class Function;
-#endif
-
-
-/**
- * This class declares a local alias that denotes a mapping between a boundary
- * indicator (see
- * @ref GlossBoundaryIndicator)
- * that is used to describe what
- * kind of boundary condition holds on a particular piece of the boundary, and
- * the function describing the actual function that provides the boundary values
- * on this part of the boundary. This type is required in many functions in the
- * library where, for example, we need to know about the functions
- * $h_i(\mathbf x)$ used in boundary conditions
- * @f{align*}{
- * \mathbf n \cdot \nabla u = h_i \qquad \qquad
- * \text{on}\ \Gamma_i\subset\partial\Omega.
- * @f}
- * An example is the function KellyErrorEstimator::estimate() that allows us to
- * provide a set of functions $h_i$ for all those boundary indicators $i$ for
- * which the boundary condition is supposed to be of Neumann type. Of course,
- * the same kind of principle can be applied to cases where we care about
- * Dirichlet values, where one needs to provide a map from boundary indicator
- * $i$ to Dirichlet function $h_i$ if the boundary conditions are given as
- * @f{align*}{
- * u = h_i \qquad \qquad \text{on}\ \Gamma_i\subset\partial\Omega.
- * @f}
- * This is, for example, the case for the VectorTools::interpolate() functions.
- *
- * Tutorial programs step-6, step-7 and step-8 show examples of how to use
- * function arguments of this type in situations where we actually have an empty
- * map (i.e., we want to describe that <i>no</i> part of the boundary is a
- * Neumann boundary). step-16 actually uses it in a case where one of the parts
- * of the boundary uses a boundary indicator for which we want to use a function
- * object.
- *
- * It seems odd at first to declare this alias inside a class, rather than
- * declaring an alias at global scope. The reason is that C++ does not allow to
- * define templated alias, where here in fact we want an alias that depends on
- * the space dimension. (Defining templated alias is something that is possible
- * starting with the C++11 standard, but that wasn't possible within the C++98
- * standard in place when this programming pattern was conceived.)
- *
- * @ingroup functions
- */
-template <int dim, typename Number = double>
-struct DEAL_II_DEPRECATED FunctionMap
-{
-  /**
-   * Declare the type as discussed above. Since we can't name it FunctionMap
-   * (as that would ambiguate a possible constructor of this class), name it
-   * in the fashion of the standard container local alias.
-   *
-   * @deprecated Directly use the type
-   * <tt>std::map<types::boundary_id, const Function<dim, Number> *></tt>
-   * in your code instead of this alias.
-   */
-  using type DEAL_II_DEPRECATED =
-    std::map<types::boundary_id, const Function<dim, Number> *>;
-};
-
-DEAL_II_NAMESPACE_CLOSE
-
-#endif
index fe4934b5a2dd5a97741399c3ede9930fb3c6cebe..ace14db908f5b1327924bda4c1b09fa4577f11f6 100644 (file)
@@ -29,7 +29,6 @@
 #include <deal.II/distributed/tria_base.h>
 
 #include <deal.II/dofs/block_info.h>
-#include <deal.II/dofs/deprecated_function_map.h>
 #include <deal.II/dofs/dof_accessor.h>
 #include <deal.II/dofs/dof_faces.h>
 #include <deal.II/dofs/dof_iterator_selector.h>
diff --git a/include/deal.II/dofs/function_map.h b/include/deal.II/dofs/function_map.h
deleted file mode 100644 (file)
index b263aaf..0000000
+++ /dev/null
@@ -1,25 +0,0 @@
-// ---------------------------------------------------------------------
-//
-// Copyright (C) 2001 - 2018 by the deal.II authors
-//
-// This file is part of the deal.II library.
-//
-// The deal.II library is free software; you can use it, redistribute
-// it, and/or modify it under the terms of the GNU Lesser General
-// Public License as published by the Free Software Foundation; either
-// version 2.1 of the License, or (at your option) any later version.
-// The full text of the license can be found in the file LICENSE.md at
-// the top level directory of deal.II.
-//
-// ---------------------------------------------------------------------
-
-#ifndef dealii_function_map_h
-#define dealii_function_map_h
-
-#include <deal.II/base/config.h>
-
-#include <deal.II/dofs/deprecated_function_map.h>
-
-DEAL_II_WARNING("This file is deprecated.")
-
-#endif
index a5a88a65858af64b96812b69095efdcc32bfd202..220b3fc69238f28e76f90cbd23bc88bbf43ed15f 100644 (file)
@@ -22,8 +22,6 @@
 #include <deal.II/base/exceptions.h>
 #include <deal.II/base/function.h>
 
-#include <deal.II/dofs/deprecated_function_map.h>
-
 #include <deal.II/fe/component_mask.h>
 
 #include <map>
index cdfcc5e61021356770109ab09a3cc29193c3b708..9cc50320515c72b86db7fcad490925fb9c2bb081 100644 (file)
@@ -22,8 +22,6 @@
 #include <deal.II/base/exceptions.h>
 #include <deal.II/base/function.h>
 
-#include <deal.II/dofs/deprecated_function_map.h>
-
 #include <deal.II/lac/affine_constraints.h>
 
 #include <map>
index 16846b9ffed6b77f868a2c58f313e725e3b7991c..438a682bc36f4df0e2781214f16df26899d797be 100644 (file)
@@ -813,11 +813,12 @@ void
 EstimateEnrichmentFunction::refine_grid()
 {
   Vector<float> estimated_error_per_cell(triangulation.n_active_cells());
-  KellyErrorEstimator<1>::estimate(dof_handler,
-                                   QGauss<1 - 1>(3),
-                                   typename FunctionMap<1>::type(),
-                                   solution,
-                                   estimated_error_per_cell);
+  KellyErrorEstimator<1>::estimate(
+    dof_handler,
+    QGauss<1 - 1>(3),
+    std::map<types::boundary_id, const Function<1> *>{},
+    solution,
+    estimated_error_per_cell);
   GridRefinement::refine_and_coarsen_fixed_number(triangulation,
                                                   estimated_error_per_cell,
                                                   0.2,
@@ -1599,15 +1600,16 @@ LaplaceProblem<dim>::refine_grid()
   for (unsigned int i = 0; i < q_collection.size(); ++i)
     q_collection_face.push_back(QGauss<dim - 1>(1));
 
-  KellyErrorEstimator<dim>::estimate(dof_handler,
-                                     q_collection_face,
-                                     typename FunctionMap<dim>::type(),
-                                     localized_solution,
-                                     local_error_per_cell,
-                                     ComponentMask(),
-                                     nullptr,
-                                     n_mpi_processes,
-                                     this_mpi_process);
+  KellyErrorEstimator<dim>::estimate(
+    dof_handler,
+    q_collection_face,
+    std::map<types::boundary_id, const Function<dim> *>{},
+    localized_solution,
+    local_error_per_cell,
+    ComponentMask(),
+    nullptr,
+    n_mpi_processes,
+    this_mpi_process);
   const unsigned int n_local_cells =
     GridTools::count_cells_with_subdomain_association(triangulation,
                                                       this_mpi_process);

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.