From: Wolfgang Bangerth Date: Sat, 23 Jan 2016 16:43:59 +0000 (-0600) Subject: Remove an unused argument of FiniteElementData. X-Git-Tag: v8.4.0-rc2~63^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=6a1d944f68e25d701e0375355effe17eaa8d7f78;p=dealii.git Remove an unused argument of FiniteElementData. The value of this argument was unused by this class, but several derived classes passed something down anyway. Fix this by simply removing the argument. --- diff --git a/doc/news/changes.h b/doc/news/changes.h index 5403fc08b8..571f51156b 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -38,6 +38,15 @@ inconvenience this causes.

    +
  1. Changed: The constructor of FiniteElementData had a last argument + n_blocks that was not actually used by the class to + initialize anything. It has been removed. In addition, the default + constructor of FiniteElementData has also been removed given that it + only creates a dysfunctional element. +
    + (Wolfgang Bangerth, 2016/01/23). +
  2. +
  3. Rework: SLEPcWrappers were reworked to allow usage of PETSc solvers and preconditioners inside SLEPc's eigensolvers. To that end extra methods were introduced to PETSc wrappers. Moreover, initialization of the @@ -48,7 +57,7 @@ inconvenience this causes. (Denis Davydov, 2015/12/29).
  4. -
  5. Removed the deprecated Operator class in the Algorithms namespace. +
  6. Removed: The deprecated Operator class in the Algorithms namespace.
    (Timo Heister, 2015/12/21)
  7. diff --git a/include/deal.II/fe/fe_base.h b/include/deal.II/fe/fe_base.h index 5ccbeab06f..1a290d7780 100644 --- a/include/deal.II/fe/fe_base.h +++ b/include/deal.II/fe/fe_base.h @@ -340,14 +340,11 @@ public: * $H_\text{div}$ conforming; finally, completely discontinuous * elements (implemented by the FE_DGQ class) are only $L_2$ * conforming. - * - * @param[in] n_blocks obsolete and ignored. */ FiniteElementData (const std::vector &dofs_per_object, const unsigned int n_components, const unsigned int degree, - const Conformity conformity = unknown, - const unsigned int n_blocks = numbers::invalid_unsigned_int); + const Conformity conformity = unknown); /** * Number of dofs per vertex. diff --git a/include/deal.II/fe/fe_dg_vector.templates.h b/include/deal.II/fe/fe_dg_vector.templates.h index 17e10b14cc..ae1e3e9a43 100644 --- a/include/deal.II/fe/fe_dg_vector.templates.h +++ b/include/deal.II/fe/fe_dg_vector.templates.h @@ -30,7 +30,7 @@ FE_DGVector::FE_DGVector ( FE_PolyTensor( deg, FiniteElementData( - get_dpo_vector(deg), dim, deg+1, FiniteElementData::L2, 1), + get_dpo_vector(deg), dim, deg+1, FiniteElementData::L2), std::vector(PolynomialType::compute_n_pols(deg), true), std::vector(PolynomialType::compute_n_pols(deg), ComponentMask(dim,true))) diff --git a/source/fe/fe_abf.cc b/source/fe/fe_abf.cc index 79d8c06012..064c480a2a 100644 --- a/source/fe/fe_abf.cc +++ b/source/fe/fe_abf.cc @@ -45,7 +45,7 @@ FE_ABF::FE_ABF (const unsigned int deg) FE_PolyTensor, dim> ( deg, FiniteElementData(get_dpo_vector(deg), - dim, deg+1, FiniteElementData::Hdiv, 1), + dim, deg+1, FiniteElementData::Hdiv), std::vector(PolynomialsABF::compute_n_pols(deg), true), std::vector(PolynomialsABF::compute_n_pols(deg), std::vector(dim,true))), diff --git a/source/fe/fe_bdm.cc b/source/fe/fe_bdm.cc index 2fe36b2f36..559709331c 100644 --- a/source/fe/fe_bdm.cc +++ b/source/fe/fe_bdm.cc @@ -39,7 +39,7 @@ FE_BDM::FE_BDM (const unsigned int deg) FE_PolyTensor, dim> ( deg, FiniteElementData(get_dpo_vector(deg), - dim, deg+1, FiniteElementData::Hdiv, 1), + dim, deg+1, FiniteElementData::Hdiv), get_ria_vector (deg), std::vector(PolynomialsBDM::compute_n_pols(deg), std::vector(dim,true))) diff --git a/source/fe/fe_data.cc b/source/fe/fe_data.cc index b04752ba21..ddcc3173f7 100644 --- a/source/fe/fe_data.cc +++ b/source/fe/fe_data.cc @@ -23,8 +23,7 @@ FiniteElementData:: FiniteElementData (const std::vector &dofs_per_object, const unsigned int n_components, const unsigned int degree, - const Conformity conformity, - const unsigned int) + const Conformity conformity) : dofs_per_vertex(dofs_per_object[0]), dofs_per_line(dofs_per_object[1]), diff --git a/source/fe/fe_nedelec.cc b/source/fe/fe_nedelec.cc index 8934aaad24..dcb488a5bf 100644 --- a/source/fe/fe_nedelec.cc +++ b/source/fe/fe_nedelec.cc @@ -46,7 +46,7 @@ FE_Nedelec::FE_Nedelec (const unsigned int p) : FE_PolyTensor, dim> (p, FiniteElementData (get_dpo_vector (p), dim, p + 1, - FiniteElementData::Hcurl, 1), + FiniteElementData::Hcurl), std::vector (PolynomialsNedelec::compute_n_pols (p), true), std::vector (PolynomialsNedelec::compute_n_pols (p), diff --git a/source/fe/fe_raviart_thomas.cc b/source/fe/fe_raviart_thomas.cc index aa3303a57a..37a1a297ab 100644 --- a/source/fe/fe_raviart_thomas.cc +++ b/source/fe/fe_raviart_thomas.cc @@ -45,7 +45,7 @@ FE_RaviartThomas::FE_RaviartThomas (const unsigned int deg) FE_PolyTensor, dim> ( deg, FiniteElementData(get_dpo_vector(deg), - dim, deg+1, FiniteElementData::Hdiv, 1), + dim, deg+1, FiniteElementData::Hdiv), std::vector(PolynomialsRaviartThomas::compute_n_pols(deg), true), std::vector(PolynomialsRaviartThomas::compute_n_pols(deg), std::vector(dim,true))) diff --git a/source/fe/fe_raviart_thomas_nodal.cc b/source/fe/fe_raviart_thomas_nodal.cc index 09f79d33c6..7a81f17cf1 100644 --- a/source/fe/fe_raviart_thomas_nodal.cc +++ b/source/fe/fe_raviart_thomas_nodal.cc @@ -38,7 +38,7 @@ FE_RaviartThomasNodal::FE_RaviartThomasNodal (const unsigned int deg) FE_PolyTensor, dim> ( deg, FiniteElementData(get_dpo_vector(deg), - dim, deg+1, FiniteElementData::Hdiv, 1), + dim, deg+1, FiniteElementData::Hdiv), get_ria_vector (deg), std::vector(PolynomialsRaviartThomas::compute_n_pols(deg), std::vector(dim,true))) diff --git a/source/fe/fe_system.cc b/source/fe/fe_system.cc index f50d3a2f62..2d1265f22f 100644 --- a/source/fe/fe_system.cc +++ b/source/fe/fe_system.cc @@ -2449,8 +2449,6 @@ FESystem::multiply_dof_numbers ( unsigned int degree = 0; // degree is the maximal degree of the components - unsigned int summed_multiplicities = 0; - for (unsigned int i=0; i0) { @@ -2462,8 +2460,6 @@ FESystem::multiply_dof_numbers ( multiplied_n_components+=fes[i]->n_components() * multiplicities[i]; degree = std::max(degree, fes[i]->tensor_degree() ); - - summed_multiplicities += multiplicities[i]; } // assume conformity of the first finite element and then take away @@ -2497,8 +2493,7 @@ FESystem::multiply_dof_numbers ( return FiniteElementData (dpo, multiplied_n_components, degree, - total_conformity, - summed_multiplicities); + total_conformity); }