#include <deal.II/base/quadrature_lib.h>
#include <deal.II/base/qprojector.h>
+#include <deal.II/base/std_cxx14/memory.h>
#include <deal.II/base/thread_management.h>
#include <deal.II/base/utilities.h>
#include <deal.II/lac/vector.h>
Vector<typename OutVector::value_type> u1_local(DoFTools::max_dofs_per_cell(dof1));
Vector<typename OutVector::value_type> u2_local(DoFTools::max_dofs_per_cell(dof2));
- // have a map for interpolation
- // matrices. shared_ptr make sure
- // that memory is released again
+ // have a map for interpolation matrices.
+ // Using a unique_ptr makes sure that the
+ // memory is released again automatically.
std::map<const FiniteElement<dim,spacedim> *,
std::map<const FiniteElement<dim,spacedim> *,
- std::shared_ptr<FullMatrix<double> > > >
+ std::unique_ptr<FullMatrix<double> > > >
interpolation_matrices;
typename DoFHandlerType1<dim,spacedim>::active_cell_iterator cell1 = dof1.begin_active(),
// there
if (interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()].get() == nullptr)
{
- std::shared_ptr<FullMatrix<double> >
- interpolation_matrix (new FullMatrix<double> (dofs_per_cell2,
- dofs_per_cell1));
- interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()]
- = interpolation_matrix;
+ auto interpolation_matrix = std_cxx14::make_unique<FullMatrix<double> >
+ (dofs_per_cell2, dofs_per_cell1);
get_interpolation_matrix(cell1->get_fe(),
cell2->get_fe(),
*interpolation_matrix);
+
+ interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()]
+ = std::move(interpolation_matrix);
}
cell1->get_dof_values(u1, u1_local);
// dof1 to the back_interpolation
// matrices
std::map<const FiniteElement<dim> *,
- std::shared_ptr<FullMatrix<double> > > interpolation_matrices;
+ std::unique_ptr<FullMatrix<double> > > interpolation_matrices;
for (; cell!=endc; ++cell)
if ((cell->subdomain_id() == subdomain_id)
// matrix is available
if (interpolation_matrices[&cell->get_fe()] == nullptr)
{
- interpolation_matrices[&cell->get_fe()] =
- std::shared_ptr<FullMatrix<double> >
- (new FullMatrix<double>(dofs_per_cell1, dofs_per_cell1));
+ interpolation_matrices[&cell->get_fe()] = std_cxx14::make_unique<FullMatrix<double> >
+ (dofs_per_cell1, dofs_per_cell1);
get_back_interpolation_matrix(cell->get_fe(), fe2,
*interpolation_matrices[&cell->get_fe()]);
}
{
case 1:
if (spacedim==1)
- q_fine.reset(new QGauss<dim> (degree+1));
+ q_fine = std_cxx14::make_unique<QGauss<dim>>(degree+1);
else if (spacedim==2)
- q_fine.reset(new QAnisotropic<dim>(QGauss<1>(degree+1), q_dummy));
+ q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(QGauss<1>(degree+1), q_dummy);
else
- q_fine.reset(new QAnisotropic<dim>(QGauss<1>(degree+1), q_dummy, q_dummy));
+ q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(QGauss<1>(degree+1), q_dummy, q_dummy);
break;
case 2:
if (spacedim==2)
- q_fine.reset(new QGauss<dim> (degree+1));
+ q_fine = std_cxx14::make_unique<QGauss<dim>>(degree+1);
else
- q_fine.reset(new QAnisotropic<dim>(QGauss<1>(degree+1), QGauss<1>(degree+1), q_dummy));
+ q_fine = std_cxx14::make_unique<QAnisotropic<dim>>(QGauss<1>(degree+1), QGauss<1>(degree+1), q_dummy);
break;
case 3:
- q_fine.reset(new QGauss<dim> (degree+1));
+ q_fine = std_cxx14::make_unique<QGauss<dim>>(degree+1);
break;
default:
Assert(false, ExcInternalError());
#include <deal.II/base/numbers.h>
#include <deal.II/base/quadrature.h>
#include <deal.II/base/signaling_nan.h>
+#include <deal.II/base/std_cxx14/memory.h>
#include <deal.II/differentiation/ad.h>
#include <deal.II/lac/vector.h>
#include <deal.II/lac/block_vector.h>
#include <deal.II/fe/fe.h>
#include <iomanip>
-#include <memory>
#include <type_traits>
if (flags & update_mapping)
this->mapping_data.reset (mapping_get_data.return_value());
else
- this->mapping_data.reset (new typename Mapping<dim,spacedim>::InternalDataBase());
+ this->mapping_data = std_cxx14::make_unique<typename Mapping<dim,spacedim>::InternalDataBase> ();
}
}
else
// if the types don't match, there is nothing we can do here
- present_cell.reset (new Type(new_cell));
+ present_cell = std_cxx14::make_unique<Type> (new_cell);
}
}
if (flags & update_mapping)
this->mapping_data.reset (mapping_get_data.return_value());
else
- this->mapping_data.reset (new typename Mapping<dim,spacedim>::InternalDataBase());
+ this->mapping_data = std_cxx14::make_unique<typename Mapping<dim,spacedim>::InternalDataBase> ();
}
if (flags & update_mapping)
this->mapping_data.reset (mapping_get_data.return_value());
else
- this->mapping_data.reset (new typename Mapping<dim,spacedim>::InternalDataBase());
+ this->mapping_data = std_cxx14::make_unique<typename Mapping<dim,spacedim>::InternalDataBase> ();
}
//
// we only need to replace the Qp mapping because that's the one that's
// used on boundary cells where it matters
- this->qp_mapping.reset (new MappingC1<dim,spacedim>::MappingC1Generic());
+ this->qp_mapping = std::make_shared<MappingC1<dim,spacedim>::MappingC1Generic>();
}
(dim != spacedim)),
// create a Q1 mapping for use on interior cells (if necessary)
// or to create a good initial guess in transform_real_to_unit_cell()
- q1_mapping (new MappingQGeneric<dim,spacedim>(1)),
+ q1_mapping (std::make_shared<MappingQGeneric<dim,spacedim> >(1)),
// create a Q_p mapping; if p=1, simply share the Q_1 mapping already
// created via the shared_ptr objects
qp_mapping (this->polynomial_degree>1
?
- std::make_shared<const MappingQGeneric<dim,spacedim>>(degree)
+ std::make_shared<const MappingQGeneric<dim,spacedim> >(degree)
:
q1_mapping)
{}
// reset the q1 mapping we use for interior cells (and previously
// set by the MappingQ constructor) to a MappingQ1Eulerian with the
// current vector
- this->q1_mapping.reset (new MappingQ1Eulerian<dim,VectorType,spacedim>(euler_dof_handler,
- euler_vector));
+ this->q1_mapping = std::make_shared<MappingQ1Eulerian<dim,VectorType,spacedim>>
+ (euler_dof_handler, euler_vector);
// also reset the qp mapping pointer with our own class
- this->qp_mapping.reset (new MappingQEulerianGeneric(degree,*this));
+ this->qp_mapping = std::make_shared<MappingQEulerianGeneric>(degree,*this);
}