]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Remove a whole bunch of code.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Thu, 3 Sep 2015 15:46:52 +0000 (10:46 -0500)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 6 Sep 2015 18:41:41 +0000 (13:41 -0500)
None of this is needed any more now that the InternalData
class has been largely separated.

include/deal.II/fe/mapping_q.h
include/deal.II/fe/mapping_q1.h
source/fe/mapping_q.cc
source/fe/mapping_q1.cc

index 2d4914eaaf97cb8537c61850df83be1af6e3c2f2..e183478934de699206f0cb465f2d12385f33e8f6 100644 (file)
@@ -441,17 +441,6 @@ protected:
    */
   const unsigned int n_outer;
 
-  /**
-   * Pointer to the @p dim-dimensional tensor product polynomials used as
-   * shape functions for the Qp mapping of cells at the boundary.
-   */
-  const TensorProductPolynomials<dim> *tensor_pols;
-
-  /**
-   * Number of the Qp tensor product shape functions.
-   */
-  const unsigned int n_shape_functions;
-
   /**
    * If this flag is set @p true then @p MappingQ is used on all cells, not
    * only on boundary cells.
index 0bd58ab14721576e70a452bbd06e7121cbe1e3e2..bf07fa14b9ed5ad96dc3b269f8a7d3672378496b 100644 (file)
@@ -538,12 +538,6 @@ protected:
   virtual void compute_mapping_support_points(
     const typename Triangulation<dim,spacedim>::cell_iterator &cell,
     std::vector<Point<spacedim> > &a) const;
-
-  /**
-   * Number of shape functions. Is simply the number of vertices per cell for
-   * the Q1 mapping.
-   */
-  static const unsigned int n_shape_functions = GeometryInfo<dim>::vertices_per_cell;
 };
 
 
index 0a957cc02d09620bab93d4cc33c4ed1c929c046e..dc817420759c8fab6488e6d740146bc8378ad624 100644 (file)
@@ -64,13 +64,13 @@ MappingQ<dim,spacedim>::MappingQ (const unsigned int p,
           ((dim==2) ?
            4+4*(degree-1) :
            8+12*(degree-1)+6*(degree-1)*(degree-1))),
-  n_shape_functions(Utilities::fixed_power<dim>(degree+1)),
   use_mapping_q_on_all_cells (use_mapping_q_on_all_cells
                               || (dim != spacedim)),
   feq(degree),
   line_support_points(degree+1)
 {
-  Assert(n_inner+n_outer==n_shape_functions, ExcInternalError());
+  Assert(n_inner+n_outer==Utilities::fixed_power<dim>(degree+1),
+         ExcInternalError());
 
   // build laplace_on_quad_vector
   if (degree>1)
@@ -90,7 +90,6 @@ MappingQ<dim,spacedim>::MappingQ (const MappingQ<dim,spacedim> &mapping)
   degree(mapping.degree),
   n_inner(mapping.n_inner),
   n_outer(mapping.n_outer),
-  n_shape_functions(mapping.n_shape_functions),
   use_mapping_q_on_all_cells (mapping.use_mapping_q_on_all_cells),
   feq(degree),
   line_support_points(degree+1)
@@ -481,7 +480,8 @@ MappingQ<dim,spacedim>::compute_laplace_vector(Table<2,double> &lvs) const
   const unsigned int n_q_points=quadrature.size();
 
   InternalData quadrature_data(degree);
-  quadrature_data.shape_derivatives.resize(n_shape_functions * n_q_points);
+  quadrature_data.shape_derivatives.resize(quadrature_data.n_shape_functions *
+                                           n_q_points);
   quadrature_data.compute_shape_function_values(quadrature.get_points());
 
   // Compute the stiffness matrix of the inner dofs
@@ -880,7 +880,7 @@ transform (const VectorSlice<const std::vector<Tensor<1,dim> > >   input,
     }
 
   // Now, q1_data should have the right tensors in it and we call the base
-  // classes transform function
+  // class's transform function
   MappingQ1<dim,spacedim>::transform(input, mapping_type, *q1_data, output);
 }
 
@@ -910,7 +910,7 @@ transform (const VectorSlice<const std::vector<DerivativeForm<1, dim ,spacedim>
     }
 
   // Now, q1_data should have the right tensors in it and we call the base
-  // classes transform function
+  // class's transform function
   MappingQ1<dim,spacedim>::transform(input, mapping_type, *q1_data, output);
 }
 
@@ -938,7 +938,7 @@ transform (const VectorSlice<const std::vector<Tensor<2, dim> > >  input,
     }
 
   // Now, q1_data should have the right tensors in it and we call the base
-  // classes transform function
+  // class's transform function
   MappingQ1<dim,spacedim>::transform(input, mapping_type, *q1_data, output);
 }
 
@@ -968,7 +968,7 @@ transform (const VectorSlice<const std::vector<DerivativeForm<2, dim ,spacedim>
     }
 
   // Now, q1_data should have the right tensors in it and we call the base
-  // classes transform function
+  // class's transform function
   MappingQ1<dim,spacedim>::transform(input, mapping_type, *q1_data, output);
 }
 
@@ -996,7 +996,7 @@ transform (const VectorSlice<const std::vector<Tensor<3, dim> > >  input,
     }
 
   // Now, q1_data should have the right tensors in it and we call the base
-  // classes transform function
+  // class's transform function
   MappingQ1<dim,spacedim>::transform(input, mapping_type, *q1_data, output);
 }
 
index 96ab2065de8837a141fd6cb5070edda7765deef2..d4f83839d34ce0781ce72c18d772ae247cbf478d 100644 (file)
@@ -17,6 +17,8 @@
 #include <deal.II/base/derivative_form.h>
 #include <deal.II/base/quadrature.h>
 #include <deal.II/base/qprojector.h>
+#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/base/tensor_product_polynomials.h>
 #include <deal.II/base/memory_consumption.h>
 #include <deal.II/base/std_cxx11/array.h>
 #include <deal.II/base/std_cxx11/unique_ptr.h>
@@ -28,8 +30,6 @@
 #include <deal.II/fe/fe.h>
 #include <deal.II/fe/fe_values.h>
 #include <deal.II/fe/mapping_q1.h>
-#include <deal.II/fe/mapping_q.h>
-#include <deal.II/fe/mapping_q1_eulerian.h>
 
 #include <cmath>
 #include <algorithm>
 DEAL_II_NAMESPACE_OPEN
 
 
-template <int dim, int spacedim>
-const unsigned int MappingQ1<dim,spacedim>::n_shape_functions;
-
-
-
 template<int dim, int spacedim>
 MappingQ1<dim,spacedim>::InternalData::InternalData (const unsigned int polynomial_degree)
   :

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.