]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Initialize the last non-const field of FiniteElementData from constructor arguments. 2112/head
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sat, 23 Jan 2016 18:24:25 +0000 (12:24 -0600)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 24 Jan 2016 19:17:52 +0000 (13:17 -0600)
include/deal.II/fe/fe_base.h
include/deal.II/lac/block_indices.h
source/fe/fe_data.cc
source/fe/fe_system.cc

index 1a290d7780a4a0f0674c2f6a22e888c71e252869..83cacd0061525d89d7f8ce71224c791c65849805 100644 (file)
@@ -308,7 +308,7 @@ public:
    * element. For an element which is not an FESystem, this contains only a
    * single block with length #dofs_per_cell.
    */
-  BlockIndices block_indices_data;
+  const BlockIndices block_indices_data;
 
   /**
    * Constructor, computing all necessary values from the distribution of dofs
@@ -340,11 +340,19 @@ public:
    *   $H_\text{div}$ conforming; finally, completely discontinuous
    *   elements (implemented by the FE_DGQ class) are only $L_2$
    *   conforming.
+   *
+   * @param[in] block_indices An argument that describes how the base elements
+   *   of a finite element are grouped. The default value constructs a single
+   *   block that consists of all @p dofs_per_cell degrees of freedom. This
+   *   is appropriate for all "atomic" elements (including non-primitive ones)
+   *   and these can therefore omit this argument. On the other hand, composed
+   *   elements such as FESystem will want to pass a different value here.
    */
   FiniteElementData (const std::vector<unsigned int> &dofs_per_object,
                      const unsigned int               n_components,
                      const unsigned int               degree,
-                     const Conformity                 conformity = unknown);
+                     const Conformity                 conformity = unknown,
+                     const BlockIndices              &block_indices = BlockIndices());
 
   /**
    * Number of dofs per vertex.
index 7f4f7156044b88523892b7a305721fb3b293cbdd..d6227fd242d4ba6c8e8e55d1a2f0ac5535b63919 100644 (file)
@@ -74,7 +74,8 @@ public:
   /**
    * Specialized constructor for a structure with blocks of equal size.
    */
-  explicit BlockIndices(const unsigned int n_blocks, const size_type block_size = 0);
+  explicit BlockIndices(const unsigned int n_blocks,
+                        const size_type block_size = 0);
 
   /**
    * Reinitialize the number of blocks and assign each block the same number
@@ -259,9 +260,8 @@ BlockIndices::BlockIndices ()
 
 
 inline
-BlockIndices::BlockIndices (
-  const unsigned int n_blocks,
-  const size_type block_size)
+BlockIndices::BlockIndices (const unsigned int n_blocks,
+                            const size_type block_size)
   :
   n_blocks(n_blocks),
   start_indices(n_blocks+1)
index ddcc3173f7e5d6180e14253af111454bb7ffe046..10e2229baa45a7fc36a9ee7a2a32339f989134c0 100644 (file)
@@ -21,9 +21,10 @@ DEAL_II_NAMESPACE_OPEN
 template <int dim>
 FiniteElementData<dim>::
 FiniteElementData (const std::vector<unsigned int> &dofs_per_object,
-                   const unsigned int n_components,
-                   const unsigned int degree,
-                   const Conformity conformity)
+                   const unsigned int               n_components,
+                   const unsigned int               degree,
+                   const Conformity                 conformity,
+                   const BlockIndices              &block_indices)
   :
   dofs_per_vertex(dofs_per_object[0]),
   dofs_per_line(dofs_per_object[1]),
@@ -56,7 +57,11 @@ FiniteElementData (const std::vector<unsigned int> &dofs_per_object,
   components(n_components),
   degree(degree),
   conforming_space(conformity),
-  block_indices_data(1, dofs_per_cell)
+  block_indices_data(block_indices.size() == 0
+                     ?
+                     BlockIndices(1, dofs_per_cell)
+                     :
+                     block_indices)
 {
   Assert(dofs_per_object.size()==dim+1, ExcDimensionMismatch(dofs_per_object.size()-1,dim));
 }
index aefbfa6849bbd548e511c4a4f150d6ebfbd94ead..c9e5c34c04ff515ebc5a085d75e2f28323000184 100644 (file)
@@ -123,7 +123,7 @@ namespace
   multiply_dof_numbers (const std::vector<const FiniteElement<dim,spacedim>*> &fes,
                         const std::vector<unsigned int>                       &multiplicities)
   {
-    Assert (fes.size() == multiplicities.size(), ExcDimensionMismatch (fes.size(), multiplicities.size()));
+    AssertDimension(fes.size(), multiplicities.size());
 
     unsigned int multiplied_dofs_per_vertex = 0;
     unsigned int multiplied_dofs_per_line = 0;
@@ -175,10 +175,17 @@ namespace
     if (dim>1) dpo.push_back(multiplied_dofs_per_quad);
     if (dim>2) dpo.push_back(multiplied_dofs_per_hex);
 
+    BlockIndices block_indices (0,0);
+
+    for (unsigned int base=0; base < fes.size(); ++base)
+      for (unsigned int m = 0; m < multiplicities[base]; ++m)
+        block_indices.push_back(fes[base]->dofs_per_cell);
+
     return FiniteElementData<dim> (dpo,
                                    multiplied_n_components,
                                    degree,
-                                   total_conformity);
+                                   total_conformity,
+                                   block_indices);
   }
 
   /**
@@ -225,7 +232,7 @@ namespace
   compute_restriction_is_additive_flags (const std::vector<const FiniteElement<dim,spacedim>*> &fes,
                                          const std::vector<unsigned int>              &multiplicities)
   {
-    Assert (fes.size() == multiplicities.size(), ExcInternalError());
+    AssertDimension(fes.size(), multiplicities.size());
 
     // first count the number of dofs and components that will emerge from the
     // given FEs
@@ -391,7 +398,7 @@ namespace
   compute_nonzero_components (const std::vector<const FiniteElement<dim,spacedim>*> &fes,
                               const std::vector<unsigned int>              &multiplicities)
   {
-    Assert (fes.size() == multiplicities.size(), ExcInternalError());
+    AssertDimension(fes.size(), multiplicities.size());
 
     // first count the number of dofs and components that will emerge from the
     // given FEs
@@ -666,6 +673,7 @@ FESystem<dim,spacedim>::FESystem (const FiniteElement<dim,spacedim> &fe1,
 }
 
 
+
 template <int dim, int spacedim>
 FESystem<dim,spacedim>::FESystem (const FiniteElement<dim,spacedim> &fe1,
                                   const unsigned int        n1,
@@ -702,6 +710,8 @@ FESystem<dim,spacedim>::FESystem (const FiniteElement<dim,spacedim> &fe1,
   initialize(fes, multiplicities);
 }
 
+
+
 template <int dim, int spacedim>
 FESystem<dim,spacedim>::FESystem (const FiniteElement<dim,spacedim> &fe1,
                                   const unsigned int        n1,
@@ -761,6 +771,7 @@ FESystem<dim,spacedim>::FESystem (
 }
 
 
+
 template <int dim, int spacedim>
 FESystem<dim,spacedim>::~FESystem ()
 {}
@@ -1756,12 +1767,10 @@ FESystem<dim,spacedim>::build_cell_tables()
   this->face_system_to_component_table.resize(this->dofs_per_face);
 
   unsigned int total_index = 0;
-  this->block_indices_data.reinit(0,0);
 
   for (unsigned int base=0; base < this->n_base_elements(); ++base)
     for (unsigned int m = 0; m < this->element_multiplicity(base); ++m)
       {
-        this->block_indices_data.push_back(base_element(base).dofs_per_cell);
         for (unsigned int k=0; k<base_element(base).n_components(); ++k)
           this->component_to_base_table[total_index++]
             = std::make_pair(std::make_pair(base,k), m);

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.