#include <deal.II/fe/block_mask.h>
#include <deal.II/fe/mapping.h>
+#include <memory>
+
+
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim> class FEValuesBase;
virtual ~FiniteElement ();
/**
- * A sort of virtual copy constructor. Some places in the library, for
+ * A sort of virtual copy constructor, this function returns a copy of
+ * the finite element object. Derived classes need to override the function
+ * here in this base class and return an object of the same type as the
+ * derived class.
+ *
+ * Some places in the library, for
* example the constructors of FESystem as well as the hp::FECollection
* class, need to make copies of finite elements without knowing their exact
* type. They do so through this function.
*/
- virtual FiniteElement<dim,spacedim> *clone() const = 0;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim>>
+ clone() const = 0;
/**
* Return a string that uniquely identifies a finite element. The general
std::vector<double> &nodal_values) const;
virtual std::size_t memory_consumption () const;
- virtual FiniteElement<dim> *clone() const;
+
+ virtual
+ std::unique_ptr<FiniteElement<dim,dim> >
+ clone() const;
private:
/**
*/
virtual std::string get_name () const;
- virtual FiniteElement<dim> *clone () const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,dim> >
+ clone() const;
// documentation inherited from the base class
virtual
*/
virtual std::string get_name () const;
-protected:
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
+protected:
/**
* Only for internal use. Its full name is @p get_dofs_per_object_vector
* Constructor for the vector element of degree @p p.
*/
FE_DGVector (const unsigned int p, MappingType m);
-public:
-
- FiniteElement<dim, spacedim> *clone() const;
/**
* Return a string that uniquely identifies a finite element. This class
*/
virtual std::string get_name () const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
/**
* This function returns @p true, if the shape function @p shape_index has
#include <deal.II/fe/fe_tools.h>
#include <deal.II/base/quadrature_lib.h>
+#include <deal.II/base/std_cxx14/memory.h>
+
+
DEAL_II_NAMESPACE_OPEN
template <class PolynomialType, int dim, int spacedim>
-FiniteElement<dim, spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_DGVector<PolynomialType,dim,spacedim>::clone() const
{
- return new FE_DGVector<PolynomialType, dim, spacedim>(*this);
+ return std_cxx14::make_unique<FE_DGVector<PolynomialType, dim, spacedim>>(*this);
}
virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
get_constant_modes () const;
-protected:
-
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
private:
*/
virtual std::size_t memory_consumption () const;
-protected:
-
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim> *clone() const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,dim> >
+ clone() const;
private:
*/
virtual std::string get_name () const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
+
// for documentation, see the FiniteElement base class
virtual
UpdateFlags
protected:
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
-
/**
* Prepare internal data structures and fill in values independent of the
* cell.
*/
virtual std::size_t memory_consumption () const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
protected:
/**
*/
FE_DGQ (const std::vector<Polynomials::Polynomial<double> > &polynomials);
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim, spacedim> *clone() const;
-
private:
/**
* Only for internal use. Its full name is @p get_dofs_per_object_vector
void
convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
std::vector<double> &nodal_values) const;
-
-protected:
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
};
*/
virtual std::string get_name () const;
-protected:
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
};
*/
virtual std::string get_name () const;
-protected:
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
};
const std::vector<std::vector<std::function<const Function<spacedim> *(const typename Triangulation<dim, spacedim>::cell_iterator &) > > > &functions);
public:
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
virtual
UpdateFlags
*/
FE_FaceQ (const unsigned int p);
- virtual FiniteElement<dim,spacedim> *clone() const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
/**
* Return a string that uniquely identifies a finite element. This class
*/
FE_FaceQ (const unsigned int p);
- /**
- * Clone method.
- */
- virtual FiniteElement<1,spacedim> *clone() const;
+ virtual
+ std::unique_ptr<FiniteElement<1,spacedim>>
+ clone() const;
/**
* Return a string that uniquely identifies a finite element. This class
*/
FE_FaceP(unsigned int p);
- /**
- * Clone method.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
/**
* Return a string that uniquely identifies a finite element. This class
get_constant_modes () const;
virtual std::size_t memory_consumption () const;
- virtual FiniteElement<dim> *clone() const;
+
+ virtual
+ std::unique_ptr<FiniteElement<dim,dim> >
+ clone() const;
private:
/**
FE_Nothing (const unsigned int n_components = 1,
const bool dominate = false);
- /**
- * A sort of virtual copy constructor. Some places in the library, for
- * example the constructors of FESystem as well as the hp::FECollection
- * class, need to make copied of finite elements without knowing their exact
- * type. They do so through this function.
- */
virtual
- FiniteElement<dim,spacedim> *
+ std::unique_ptr<FiniteElement<dim,spacedim> >
clone() const;
/**
virtual UpdateFlags requires_update_flags (const UpdateFlags flags) const;
- virtual FiniteElement<2,2> *clone () const;
+ virtual
+ std::unique_ptr<FiniteElement<2,2>>
+ clone() const;
/**
* Destructor.
*/
virtual std::string get_name () const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
+
/**
* Implementation of the corresponding function in the FiniteElement
* class. Since the current element is interpolatory, the nodal
void
convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
std::vector<double> &nodal_values) const;
-
-protected:
-
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
};
virtual bool has_support_on_face (const unsigned int shape_index,
const unsigned int face_index) const;
-protected:
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
private:
virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
get_constant_modes () const;
-protected:
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
private:
*/
virtual std::string get_name () const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,dim> >
+ clone() const;
+
/**
* This function returns @p true, if the shape function @p shape_index has
* non-zero function values somewhere on the face @p face_index.
virtual std::pair<Table<2,bool>, std::vector<unsigned int> >
get_constant_modes () const;
-protected:
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim> *clone() const;
-
private:
/**
*/
virtual std::string get_name () const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
+
/**
* Implementation of the corresponding function in the FiniteElement
* class. Since the current element is interpolatory, the nodal
FiniteElementDomination::Domination
compare_for_face_domination (const FiniteElement<dim,spacedim> &fe_other) const;
//@}
-
-protected:
-
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
};
const unsigned int n_face_support_points = 2);
virtual std::string get_name() const;
- virtual FiniteElement<dim> *clone() const;
+
+ virtual
+ std::unique_ptr<FiniteElement<dim,dim> >
+ clone() const;
// documentation inherited from the base class
virtual
*/
virtual std::string get_name () const;
+ // documentation inherited from the base class
+ virtual
+ std::unique_ptr<FiniteElement<dim,dim> >
+ clone() const;
/**
* This function returns @p true, if the shape function @p shape_index has
get_constant_modes () const;
virtual std::size_t memory_consumption () const;
- virtual FiniteElement<dim> *clone() const;
private:
/**
*/
virtual std::string get_name () const;
- virtual FiniteElement<dim> *clone () const;
-
// documentation inherited from the base class
+ virtual
+ std::unique_ptr<FiniteElement<dim,dim> >
+ clone() const;
+
virtual
void
convert_generalized_support_point_values_to_nodal_values (const std::vector<Vector<double> > &support_point_values,
virtual std::string get_name () const;
// for documentation, see the FiniteElement base class
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
+
virtual
UpdateFlags
requires_update_flags (const UpdateFlags update_flags) const;
protected:
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
-
-
virtual typename FiniteElement<dim,spacedim>::InternalDataBase *
get_data (const UpdateFlags update_flags,
const Mapping<dim,spacedim> &mapping,
*/
FE_TraceQ(unsigned int p);
- /**
- * @p clone function instead of a copy constructor.
- *
- * This function is needed by the constructors of @p FESystem.
- */
- virtual FiniteElement<dim,spacedim> *clone() const;
-
/**
* Return a string that uniquely identifies a finite element. This class
* returns <tt>FE_DGQ<dim>(degree)</tt>, with <tt>dim</tt> and
*/
virtual std::string get_name () const;
+ virtual
+ std::unique_ptr<FiniteElement<dim,spacedim> >
+ clone() const;
+
/**
* This function returns @p true, if the shape function @p shape_index has
* non-zero function values somewhere on the face @p face_index.
#include <sstream>
#include <iostream>
+#include <deal.II/base/std_cxx14/memory.h>
+
//TODO: implement the adjust_quad_dof_index_for_face_orientation_table and
//adjust_line_dof_index_for_line_orientation_table fields, and write tests
template <int dim>
-FiniteElement<dim> *
+std::unique_ptr<FiniteElement<dim,dim> >
FE_ABF<dim>::clone() const
{
- return new FE_ABF<dim>(rt_order);
+ return std_cxx14::make_unique<FE_ABF<dim>>(rt_order);
}
#include <iostream>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim>
-FiniteElement<dim> *
+std::unique_ptr<FiniteElement<dim,dim> >
FE_BDM<dim>::clone() const
{
- return new FE_BDM<dim>(*this);
+ return std_cxx14::make_unique<FE_BDM<dim>>(*this);
}
#include <vector>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
+
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_Bernstein<dim,spacedim>::clone() const
{
- return new FE_Bernstein<dim,spacedim>(*this);
+ return std_cxx14::make_unique<FE_Bernstein<dim,spacedim>>(*this);
}
#include <deal.II/base/polynomials_nedelec.h>
#include <deal.II/base/polynomials_raviart_thomas.h>
+#include <deal.II/base/std_cxx14/memory.h>
+
+
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
#include <deal.II/fe/fe_tools.h>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
+
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_DGP<dim,spacedim>::clone() const
{
- return new FE_DGP<dim,spacedim>(*this);
+ return std_cxx14::make_unique<FE_DGP<dim,spacedim>>(*this);
}
#include <deal.II/fe/fe_tools.h>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
+
DEAL_II_NAMESPACE_OPEN
template <int dim>
-FiniteElement<dim> *
+std::unique_ptr<FiniteElement<dim,dim> >
FE_DGPMonomial<dim>::clone() const
{
- return new FE_DGPMonomial<dim>(*this);
+ return std_cxx14::make_unique<FE_DGPMonomial<dim>>(*this);
}
#include <deal.II/fe/fe_values.h>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
+
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_DGPNonparametric<dim,spacedim>::clone() const
{
- return new FE_DGPNonparametric<dim,spacedim>(*this);
+ return std_cxx14::make_unique<FE_DGPNonparametric<dim,spacedim>>(*this);
}
#include <iostream>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
+
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_DGQ<dim, spacedim>::clone() const
{
- return new FE_DGQ<dim, spacedim>(*this);
+ return std_cxx14::make_unique<FE_DGQ<dim, spacedim>>(*this);
}
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_DGQArbitraryNodes<dim,spacedim>::clone() const
{
// Construct a dummy quadrature formula containing the FE's nodes:
qpoints[i] = Point<1>(this->unit_support_points[lexicographic[i]][0]);
Quadrature<1> pquadrature(qpoints);
- return new FE_DGQArbitraryNodes<dim,spacedim>(pquadrature);
+ return std_cxx14::make_unique<FE_DGQArbitraryNodes<dim,spacedim>>(pquadrature);
}
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_DGQLegendre<dim,spacedim>::clone() const
{
- return new FE_DGQLegendre<dim,spacedim>(this->degree);
+ return std_cxx14::make_unique<FE_DGQLegendre<dim,spacedim>>(this->degree);
}
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_DGQHermite<dim,spacedim>::clone() const
{
- return new FE_DGQHermite<dim,spacedim>(this->degree);
+ return std_cxx14::make_unique<FE_DGQHermite<dim,spacedim>>(this->degree);
}
#include <deal.II/fe/fe_enriched.h>
-
#include <deal.II/fe/fe_tools.h>
+#include <deal.II/base/std_cxx14/memory.h>
+
+
DEAL_II_NAMESPACE_OPEN
namespace
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_Enriched<dim,spacedim>::clone() const
{
std::vector< const FiniteElement< dim, spacedim > * > fes;
multiplicities.push_back(this->element_multiplicity(i) );
}
- return new FE_Enriched<dim,spacedim>(fes, multiplicities, get_enrichments());
+ return std::unique_ptr<FE_Enriched<dim,spacedim>>(new FE_Enriched<dim,spacedim>(fes,
+ multiplicities,
+ get_enrichments()));
}
#include <deal.II/fe/fe_tools.h>
#include <deal.II/base/quadrature_lib.h>
#include <deal.II/lac/householder.h>
+
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_FaceQ<dim,spacedim>::clone() const
{
- return new FE_FaceQ<dim,spacedim>(this->degree);
+ return std_cxx14::make_unique<FE_FaceQ<dim,spacedim>>(this->degree);
}
template <int spacedim>
-FiniteElement<1,spacedim> *
-FE_FaceQ<1,spacedim>::clone() const
+std::unique_ptr<FiniteElement<1,spacedim>>
+ FE_FaceQ<1,spacedim>::clone() const
{
- return new FE_FaceQ<1,spacedim>(this->degree);
+ return std_cxx14::make_unique<FE_FaceQ<1,spacedim>>(this->degree);
}
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_FaceP<dim,spacedim>::clone() const
{
- return new FE_FaceP<dim,spacedim>(this->degree);
+ return std_cxx14::make_unique<FE_FaceP<dim,spacedim>>(this->degree);
}
#include <deal.II/fe/fe_tools.h>
#include <deal.II/lac/full_matrix.h>
#include <deal.II/lac/vector.h>
+
#include <sstream>
#include <iostream>
+#include <deal.II/base/std_cxx14/memory.h>
//TODO: implement the adjust_quad_dof_index_for_face_orientation_table and
//adjust_line_dof_index_for_line_orientation_table fields, and write tests
template <int dim>
-FiniteElement<dim>
-*FE_Nedelec<dim>::clone () const
+std::unique_ptr<FiniteElement<dim,dim> >
+FE_Nedelec<dim>::clone () const
{
- return new FE_Nedelec<dim> (*this);
+ return std_cxx14::make_unique<FE_Nedelec<dim> >(*this);
}
//---------------------------------------------------------------------------
#include <deal.II/fe/fe_nothing.h>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_Nothing<dim,spacedim>::clone() const
{
- return new FE_Nothing<dim,spacedim>(*this);
+ return std_cxx14::make_unique<FE_Nothing<dim,spacedim>>(*this);
}
#include <deal.II/fe/fe_p1nc.h>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
-FiniteElement<2,2> *FE_P1NC::clone () const
+std::unique_ptr<FiniteElement<2,2>>
+ FE_P1NC::clone () const
{
- return new FE_P1NC(*this);
+ return std_cxx14::make_unique<FE_P1NC>(*this);
}
#include <vector>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_Q<dim,spacedim>::clone() const
{
- return new FE_Q<dim,spacedim>(*this);
+ return std_cxx14::make_unique<FE_Q<dim,spacedim>>(*this);
}
#include <vector>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
#include <vector>
#include <sstream>
#include <memory>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_Q_Bubbles<dim,spacedim>::clone() const
{
- return new FE_Q_Bubbles<dim,spacedim>(*this);
+ return std_cxx14::make_unique<FE_Q_Bubbles<dim,spacedim>>(*this);
}
#include <vector>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_Q_DG0<dim,spacedim>::clone() const
{
- return new FE_Q_DG0<dim,spacedim>(*this);
+ return std_cxx14::make_unique<FE_Q_DG0<dim,spacedim>>(*this);
}
#include <cmath>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
//TODO: implement the adjust_quad_dof_index_for_face_orientation_table and
//adjust_line_dof_index_for_line_orientation_table fields, and write tests
template <int dim>
-FiniteElement<dim> *
+std::unique_ptr<FiniteElement<dim,dim> >
FE_Q_Hierarchical<dim>::clone() const
{
- return new FE_Q_Hierarchical<dim>(*this);
+ return std_cxx14::make_unique<FE_Q_Hierarchical<dim>>(*this);
}
#include <vector>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_Q_iso_Q1<dim,spacedim>::clone() const
{
- return new FE_Q_iso_Q1<dim,spacedim>(*this);
+ return std_cxx14::make_unique<FE_Q_iso_Q1<dim,spacedim>>(*this);
}
#include <algorithm>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim>
-FiniteElement<dim> *FE_RannacherTurek<dim>::clone() const
+std::unique_ptr<FiniteElement<dim,dim> >
+FE_RannacherTurek<dim>::clone() const
{
- return new FE_RannacherTurek<dim>(this->order, this->n_face_support_points);
+ return std_cxx14::make_unique<FE_RannacherTurek<dim>>(this->order, this->n_face_support_points);
}
#include <sstream>
#include <iostream>
+#include <deal.II/base/std_cxx14/memory.h>
//TODO: implement the adjust_quad_dof_index_for_face_orientation_table and
//adjust_line_dof_index_for_line_orientation_table fields, and write tests
template <int dim>
-FiniteElement<dim> *
+std::unique_ptr<FiniteElement<dim,dim> >
FE_RaviartThomas<dim>::clone() const
{
- return new FE_RaviartThomas<dim>(*this);
+ return std_cxx14::make_unique<FE_RaviartThomas<dim>>(*this);
}
#include <deal.II/fe/fe_tools.h>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim>
-FiniteElement<dim> *
+std::unique_ptr<FiniteElement<dim,dim> >
FE_RaviartThomasNodal<dim>::clone() const
{
- return new FE_RaviartThomasNodal<dim>(*this);
+ return std_cxx14::make_unique<FE_RaviartThomasNodal<dim>>(*this);
}
#include <deal.II/fe/fe_tools.h>
#include <sstream>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FESystem<dim,spacedim>::clone() const
{
std::vector< const FiniteElement<dim,spacedim>* > fes;
fes.push_back( & base_element(i) );
multiplicities.push_back(this->element_multiplicity(i) );
}
- return new FESystem<dim,spacedim>(fes, multiplicities);
+ return std_cxx14::make_unique<FESystem<dim,spacedim>>(fes, multiplicities);
}
{
clone_base_elements += Threads::new_task ([ &,i,ind] ()
{
- base_elements[ind]
- = std::make_pair (std::unique_ptr<FiniteElement<dim,spacedim>>(fes[i]->clone()),
- multiplicities[i]);
+ base_elements[ind] = { fes[i]->clone(), multiplicities[i] };
});
++ind;
}
#include <sstream>
#include <fstream>
#include <iostream>
+#include <deal.II/base/std_cxx14/memory.h>
DEAL_II_NAMESPACE_OPEN
template <int dim, int spacedim>
-FiniteElement<dim,spacedim> *
+std::unique_ptr<FiniteElement<dim,spacedim> >
FE_TraceQ<dim,spacedim>::clone() const
{
- return new FE_TraceQ<dim,spacedim>(this->degree);
+ return std_cxx14::make_unique<FE_TraceQ<dim,spacedim>>(this->degree);
}