]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use shared_ptr to store MatrixFree data in MatrixFreeOperators::Base 3764/head
authorBruno Turcksin <bruno.turcksin@gmail.com>
Tue, 10 Jan 2017 23:08:04 +0000 (18:08 -0500)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Thu, 12 Jan 2017 14:59:11 +0000 (09:59 -0500)
13 files changed:
examples/step-37/step-37.cc
include/deal.II/matrix_free/operators.h
include/deal.II/numerics/vector_tools.h
include/deal.II/numerics/vector_tools.templates.h
source/numerics/vector_tools_project_qpmf.inst.in
tests/matrix_free/laplace_operator_01.cc
tests/matrix_free/laplace_operator_02.cc
tests/matrix_free/mass_operator_01.cc
tests/matrix_free/mass_operator_02.cc
tests/matrix_free/mass_operator_03.cc
tests/matrix_free/parallel_multigrid_02.cc
tests/matrix_free/parallel_multigrid_adaptive_02.cc
tests/matrix_free/parallel_multigrid_adaptive_04.cc

index 4360a1f0214b2bd6cb5a8f3cd62b67294be46269..9d8d1a0c82079605c944bed41894c9967488a4d6 100644 (file)
@@ -719,7 +719,7 @@ namespace Step37
     DoFHandler<dim>                            dof_handler;
 
     ConstraintMatrix                           constraints;
-    MatrixFree<dim,double>                     system_mf_storage;
+    std_cxx11::shared_ptr<MatrixFree<dim,double> > system_mf_storage;
     typedef LaplaceOperator<dim,degree_finite_element,double> SystemMatrixType;
     SystemMatrixType                           system_matrix;
 
@@ -760,6 +760,7 @@ namespace Step37
 #endif
     fe (degree_finite_element),
     dof_handler (triangulation),
+    system_mf_storage(new MatrixFree<dim,double>()),
     pcout (std::cout, Utilities::MPI::this_mpi_process(MPI_COMM_WORLD) == 0),
     // The LaplaceProblem class holds an additional output stream that
     // collects detailed timings about the setup phase. This stream, called
@@ -840,8 +841,8 @@ namespace Step37
         MatrixFree<dim,double>::AdditionalData::none;
       additional_data.mapping_update_flags = (update_gradients | update_JxW_values |
                                               update_quadrature_points);
-      system_mf_storage.reinit (dof_handler, constraints, QGauss<1>(fe.degree+1),
-                                additional_data);
+      system_mf_storage->reinit (dof_handler, constraints, QGauss<1>(fe.degree+1),
+                                 additional_data);
     }
 
     system_matrix.initialize (system_mf_storage);
@@ -893,7 +894,9 @@ namespace Step37
         mg_mf_storage[level].reinit(dof_handler, level_constraints,
                                     QGauss<1>(fe.degree+1), additional_data);
 
-        mg_matrices[level].initialize(mg_mf_storage[level], mg_constrained_dofs, level);
+        mg_matrices[level].initialize(
+          std_cxx11::make_shared<MatrixFree<dim,float>>(mg_mf_storage[level]),
+          mg_constrained_dofs, level);
         mg_matrices[level].evaluate_coefficient(Coefficient<dim>());
       }
     setup_time += time.wall_time();
@@ -918,8 +921,8 @@ namespace Step37
     Timer time;
 
     system_rhs = 0;
-    FEEvaluation<dim,degree_finite_element> phi(system_mf_storage);
-    for (unsigned int cell=0; cell<system_mf_storage.n_macro_cells(); ++cell)
+    FEEvaluation<dim,degree_finite_element> phi(*system_mf_storage);
+    for (unsigned int cell=0; cell<system_mf_storage->n_macro_cells(); ++cell)
       {
         phi.reinit(cell);
         for (unsigned int q=0; q<phi.n_q_points; ++q)
index c653fc3c7dded069610cdf5006175a3e7110165b..81eb0e137cf71c5aec293ef94a6f5db33d57b0e8 100644 (file)
@@ -79,12 +79,12 @@ namespace MatrixFreeOperators
     /**
      * Initialize operator on fine scale.
      */
-    void initialize (const MatrixFree<dim,Number> &data);
+    void initialize (std_cxx11::shared_ptr<const MatrixFree<dim,Number> > data);
 
     /**
      * Initialize operator on a level @p level.
      */
-    void initialize (const MatrixFree<dim,Number> &data,
+    void initialize (std_cxx11::shared_ptr<const MatrixFree<dim,Number> > data,
                      const MGConstrainedDoFs &mg_constrained_dofs,
                      const unsigned int level);
 
@@ -199,7 +199,7 @@ namespace MatrixFreeOperators
     /**
      * MatrixFree object to be used with this operator.
      */
-    SmartPointer<const MatrixFree<dim,Number>, Base<dim,Number> > data;
+    std_cxx11::shared_ptr<const MatrixFree<dim,Number> > data;
 
     /**
      * A shared pointer to a diagonal matrix that stores the inverse of
@@ -758,9 +758,9 @@ namespace MatrixFreeOperators
   template <int dim, typename Number>
   void
   Base<dim,Number>::
-  initialize (const MatrixFree<dim,Number> &data_)
+  initialize (std_cxx11::shared_ptr<const MatrixFree<dim,Number> > data_)
   {
-    data =  SmartPointer<const MatrixFree<dim,Number>, Base<dim,Number> >(&data_,typeid(*this).name());
+    data = data_;
     edge_constrained_indices.clear();
     have_interface_matrices = false;
   }
@@ -770,19 +770,19 @@ namespace MatrixFreeOperators
   template <int dim, typename Number>
   void
   Base<dim,Number>::
-  initialize (const MatrixFree<dim,Number> &data_,
+  initialize (std_cxx11::shared_ptr<const MatrixFree<dim,Number> > data_,
               const MGConstrainedDoFs      &mg_constrained_dofs,
               const unsigned int            level)
   {
     AssertThrow (level != numbers::invalid_unsigned_int,
                  ExcMessage("level is not set"));
-    if (data_.n_macro_cells() > 0)
+    if (data_->n_macro_cells() > 0)
       {
         AssertDimension(static_cast<int>(level),
-                        data_.get_cell_iterator(0,0)->level());
+                        data_->get_cell_iterator(0,0)->level());
       }
 
-    data = SmartPointer<const MatrixFree<dim,Number>, Base<dim,Number> >(&data_,typeid(*this).name());
+    data = data_;
 
     // setup edge_constrained indices
     std::vector<types::global_dof_index> interface_indices;
@@ -795,7 +795,7 @@ namespace MatrixFreeOperators
       if (locally_owned.is_element(interface_indices[i]))
         edge_constrained_indices.push_back(locally_owned.index_within_set(interface_indices[i]));
     have_interface_matrices = Utilities::MPI::max((unsigned int)edge_constrained_indices.size(),
-                                                  data_.get_vector_partitioner()->get_mpi_communicator()) > 0;
+                                                  data->get_vector_partitioner()->get_mpi_communicator()) > 0;
   }
 
 
index 7e22e5e1c518ee69d5108d923ea4fbbc24d56d06..802a454ee7cfa16d4b7cc6ca7b78c6b451a83232 100644 (file)
@@ -918,7 +918,7 @@ namespace VectorTools
    * which stores quadrature point data.
    */
   template <int dim, typename VectorType>
-  void project (const MatrixFree<dim,typename VectorType::value_type> &data,
+  void project (std_cxx11::shared_ptr<const MatrixFree<dim,typename VectorType::value_type> > data,
                 const ConstraintMatrix &constraints,
                 const unsigned int n_q_points_1d,
                 const std_cxx11::function< VectorizedArray<typename VectorType::value_type> (const unsigned int, const unsigned int)> func,
@@ -928,7 +928,7 @@ namespace VectorTools
    * Same as above but for <code>n_q_points_1d = matrix_free.get_dof_handler().get_fe().degree+1</code>.
    */
   template <int dim, typename VectorType>
-  void project (const MatrixFree<dim,typename VectorType::value_type> &data,
+  void project (std_cxx11::shared_ptr<const MatrixFree<dim,typename VectorType::value_type> > data,
                 const ConstraintMatrix &constraints,
                 const std_cxx11::function< VectorizedArray<typename VectorType::value_type> (const unsigned int, const unsigned int)> func,
                 VectorType &vec_result);
index 0b64ae5dc32bf94aa75662895e9c7f56871f6f62..aed2ef802bb11aecd77973460c8ebfd2a236f56b 100644 (file)
@@ -936,9 +936,10 @@ namespace VectorTools
       additional_data.tasks_parallel_scheme =
         MatrixFree<dim,number>::AdditionalData::partition_color;
       additional_data.mapping_update_flags = (update_values | update_JxW_values);
-      MatrixFree<dim, number> matrix_free;
-      matrix_free.reinit (mapping, dof, constraints,
-                          QGauss<1>(fe_degree+2), additional_data);
+      std_cxx11::shared_ptr<MatrixFree<dim, number> > matrix_free(
+        new MatrixFree<dim, number> ());
+      matrix_free->reinit (mapping, dof, constraints,
+                           QGauss<1>(fe_degree+2), additional_data);
       typedef MatrixFreeOperators::MassOperator<dim, fe_degree, fe_degree+2, components, number> MatrixType;
       MatrixType mass_matrix;
       mass_matrix.initialize(matrix_free);
@@ -946,9 +947,9 @@ namespace VectorTools
 
       typedef LinearAlgebra::distributed::Vector<number> LocalVectorType;
       LocalVectorType vec, rhs, inhomogeneities;
-      matrix_free.initialize_dof_vector(vec);
-      matrix_free.initialize_dof_vector(rhs);
-      matrix_free.initialize_dof_vector(inhomogeneities);
+      matrix_free->initialize_dof_vector(vec);
+      matrix_free->initialize_dof_vector(rhs);
+      matrix_free->initialize_dof_vector(inhomogeneities);
       constraints.distribute(inhomogeneities);
       inhomogeneities*=-1.;
 
@@ -1163,9 +1164,10 @@ namespace VectorTools
       additional_data.tasks_parallel_scheme =
         MatrixFree<dim,Number>::AdditionalData::partition_color;
       additional_data.mapping_update_flags = (update_values | update_JxW_values);
-      MatrixFree<dim, Number> matrix_free;
-      matrix_free.reinit (mapping, dof, constraints,
-                          QGauss<1>(fe_degree+2), additional_data);
+      std_cxx11::shared_ptr<MatrixFree<dim, Number> > matrix_free(
+        new MatrixFree<dim, Number>());
+      matrix_free->reinit (mapping, dof, constraints,
+                           QGauss<1>(fe_degree+2), additional_data);
       typedef MatrixFreeOperators::MassOperator<dim, fe_degree, fe_degree+2, 1, Number> MatrixType;
       MatrixType mass_matrix;
       mass_matrix.initialize(matrix_free);
@@ -1173,9 +1175,9 @@ namespace VectorTools
 
       typedef LinearAlgebra::distributed::Vector<Number> LocalVectorType;
       LocalVectorType vec, rhs, inhomogeneities;
-      matrix_free.initialize_dof_vector(vec);
-      matrix_free.initialize_dof_vector(rhs);
-      matrix_free.initialize_dof_vector(inhomogeneities);
+      matrix_free->initialize_dof_vector(vec);
+      matrix_free->initialize_dof_vector(rhs);
+      matrix_free->initialize_dof_vector(inhomogeneities);
       constraints.distribute(inhomogeneities);
       inhomogeneities*=-1.;
 
@@ -1240,12 +1242,12 @@ namespace VectorTools
 
 
     template <int dim, typename VectorType, int spacedim, int fe_degree, int n_q_points_1d>
-    void project_parallel (const MatrixFree<dim,typename VectorType::value_type> &matrix_free,
+    void project_parallel (std_cxx11::shared_ptr<const MatrixFree<dim,typename VectorType::value_type> > matrix_free,
                            const ConstraintMatrix &constraints,
                            const std_cxx11::function< VectorizedArray<typename VectorType::value_type> (const unsigned int, const unsigned int)> func,
                            VectorType &vec_result)
     {
-      const DoFHandler<dim,spacedim> &dof = matrix_free.get_dof_handler();
+      const DoFHandler<dim,spacedim> &dof = matrix_free->get_dof_handler();
 
       typedef typename VectorType::value_type Number;
       Assert (dof.get_fe().n_components() == 1,
@@ -1263,16 +1265,16 @@ namespace VectorTools
 
       typedef LinearAlgebra::distributed::Vector<Number> LocalVectorType;
       LocalVectorType vec, rhs, inhomogeneities;
-      matrix_free.initialize_dof_vector(vec);
-      matrix_free.initialize_dof_vector(rhs);
-      matrix_free.initialize_dof_vector(inhomogeneities);
+      matrix_free->initialize_dof_vector(vec);
+      matrix_free->initialize_dof_vector(rhs);
+      matrix_free->initialize_dof_vector(inhomogeneities);
       constraints.distribute(inhomogeneities);
       inhomogeneities*=-1.;
 
       // assemble right hand side:
       {
-        FEEvaluation<dim,fe_degree,n_q_points_1d,1,Number> fe_eval(matrix_free);
-        const unsigned int n_cells = matrix_free.n_macro_cells();
+        FEEvaluation<dim,fe_degree,n_q_points_1d,1,Number> fe_eval(*matrix_free);
+        const unsigned int n_cells = matrix_free->n_macro_cells();
         const unsigned int n_q_points = fe_eval.n_q_points;
 
         for (unsigned int cell=0; cell<n_cells; ++cell)
@@ -1356,14 +1358,14 @@ namespace VectorTools
 
 
   template <int dim, typename VectorType>
-  void project (const MatrixFree<dim,typename VectorType::value_type> &matrix_free,
+  void project (std_cxx11::shared_ptr<const MatrixFree<dim,typename VectorType::value_type> > matrix_free,
                 const ConstraintMatrix &constraints,
                 const unsigned int n_q_points_1d,
                 const std_cxx11::function< VectorizedArray<typename VectorType::value_type> (const unsigned int, const unsigned int)> func,
                 VectorType &vec_result)
   {
     (void) n_q_points_1d;
-    const unsigned int fe_degree = matrix_free.get_dof_handler().get_fe().degree;
+    const unsigned int fe_degree = matrix_free->get_dof_handler().get_fe().degree;
 
     Assert (fe_degree+1 == n_q_points_1d,
             ExcNotImplemented());
@@ -1402,14 +1404,14 @@ namespace VectorTools
 
 
   template <int dim, typename VectorType>
-  void project (const MatrixFree<dim,typename VectorType::value_type> &matrix_free,
+  void project (std_cxx11::shared_ptr<const MatrixFree<dim,typename VectorType::value_type> > matrix_free,
                 const ConstraintMatrix &constraints,
                 const std_cxx11::function< VectorizedArray<typename VectorType::value_type> (const unsigned int, const unsigned int)> func,
                 VectorType &vec_result)
   {
     project (matrix_free,
              constraints,
-             matrix_free.get_dof_handler().get_fe().degree+1,
+             matrix_free->get_dof_handler().get_fe().degree+1,
              func,
              vec_result);
   }
index 2434af942c8ed348e4116962c261994bee79747b..e5fd21666a5c470fd19c0c6ba3f6c4a80f878ae3 100644 (file)
@@ -21,14 +21,14 @@ for (VEC: REAL_NONBLOCK_VECTORS; deal_II_dimension : DIMENSIONS)
 
     template
     void project<deal_II_dimension, VEC>(
-        const MatrixFree<deal_II_dimension,VEC::value_type> &matrix_free,
+        std_cxx11::shared_ptr<const MatrixFree<deal_II_dimension,VEC::value_type> > matrix_free,
         const ConstraintMatrix &constraints,
         const std_cxx11::function< VectorizedArray<VEC::value_type> (const unsigned int, const unsigned int)>,
         VEC &);
 
     template
     void project<deal_II_dimension, VEC>(
-        const MatrixFree<deal_II_dimension,VEC::value_type> &matrix_free,
+        std_cxx11::shared_ptr<const MatrixFree<deal_II_dimension,VEC::value_type> > matrix_free,
         const ConstraintMatrix &constraints,
         const unsigned int,
         const std_cxx11::function< VectorizedArray<VEC::value_type> (const unsigned int, const unsigned int)>,
index 6ca7be6bae9d536a6aa04a05ec6df1c09e565ca0..dc46c7b70fd80997d9d9a8b1dadb5d7a34f292e1 100644 (file)
@@ -97,21 +97,21 @@ void test ()
   //std::cout << "Number of degrees of freedom: " << dof.n_dofs() << std::endl;
   //std::cout << "Number of constraints: " << constraints.n_constraints() << std::endl;
 
-  MatrixFree<dim,number> mf_data;
+  std_cxx11::shared_ptr<MatrixFree<dim,number> > mf_data(new MatrixFree<dim,number> ());
   {
     const QGauss<1> quad (fe_degree+1);
     typename MatrixFree<dim,number>::AdditionalData data;
     data.tasks_parallel_scheme =
       MatrixFree<dim,number>::AdditionalData::none;
     data.tasks_block_size = 7;
-    mf_data.reinit (dof, constraints, quad, data);
+    mf_data->reinit (dof, constraints, quad, data);
   }
 
   MatrixFreeOperators::LaplaceOperator<dim,fe_degree,fe_degree+1, 1, number> mf;
   mf.initialize(mf_data);
   mf.compute_diagonal();
   LinearAlgebra::distributed::Vector<number> in, out, ref;
-  mf_data.initialize_dof_vector (in);
+  mf_data->initialize_dof_vector (in);
   out.reinit (in);
   ref.reinit (in);
 
index 0f95334d897125e68fac34c7f669b2f7a67802bb..2bd3657246f728506a13a3f8c55d7169b59a94f3 100644 (file)
@@ -141,7 +141,7 @@ void test ()
   //std::cout << "Number of degrees of freedom: " << dof.n_dofs() << std::endl;
   //std::cout << "Number of constraints: " << constraints.n_constraints() << std::endl;
 
-  MatrixFree<dim,number> mf_data;
+  std_cxx11::shared_ptr<MatrixFree<dim,number> > mf_data(new MatrixFree<dim,number> ());
   {
     const QGauss<1> quad (fe_degree+1);
     typename MatrixFree<dim,number>::AdditionalData data;
@@ -149,15 +149,15 @@ void test ()
       MatrixFree<dim,number>::AdditionalData::none;
     data.tasks_block_size = 7;
     data.mapping_update_flags = update_quadrature_points | update_gradients | update_JxW_values;
-    mf_data.reinit (dof, constraints, quad, data);
+    mf_data->reinit (dof, constraints, quad, data);
   }
 
   std_cxx11::shared_ptr<Table<2, VectorizedArray<number> > > coefficient;
   coefficient = std_cxx11::make_shared<Table<2, VectorizedArray<number> > >();
   {
-    FEEvaluation<dim,fe_degree,fe_degree+1,1,number> fe_eval(mf_data);
+    FEEvaluation<dim,fe_degree,fe_degree+1,1,number> fe_eval(*mf_data);
 
-    const unsigned int n_cells = mf_data.n_macro_cells();
+    const unsigned int n_cells = mf_data->n_macro_cells();
     const unsigned int n_q_points = fe_eval.n_q_points;
 
     coefficient->reinit(n_cells, n_q_points);
@@ -174,7 +174,7 @@ void test ()
   mf.set_coefficient(coefficient);
   mf.compute_diagonal();
   LinearAlgebra::distributed::Vector<number> in, out, ref;
-  mf_data.initialize_dof_vector (in);
+  mf_data->initialize_dof_vector (in);
   out.reinit (in);
   ref.reinit (in);
 
index 43d7d4cd7fffaa635e7c4960a093b1550c2669c5..6d8db2a358b3e1d13dd70c9b5bac3278a9ce55b4 100644 (file)
@@ -105,21 +105,21 @@ void test ()
   //std::cout << "Number of degrees of freedom: " << dof.n_dofs() << std::endl;
   //std::cout << "Number of constraints: " << constraints.n_constraints() << std::endl;
 
-  MatrixFree<dim,number> mf_data;
+  std_cxx11::shared_ptr<MatrixFree<dim,number> > mf_data(new MatrixFree<dim,number> ());
   {
     const QGauss<1> quad (fe_degree+2);
     typename MatrixFree<dim,number>::AdditionalData data;
     data.tasks_parallel_scheme =
       MatrixFree<dim,number>::AdditionalData::none;
     data.tasks_block_size = 7;
-    mf_data.reinit (dof, constraints, quad, data);
+    mf_data->reinit (dof, constraints, quad, data);
   }
 
   MatrixFreeOperators::MassOperator<dim,fe_degree, fe_degree+2, 1, number> mf;
   mf.initialize(mf_data);
   mf.compute_diagonal();
   LinearAlgebra::distributed::Vector<number> in, out, ref;
-  mf_data.initialize_dof_vector (in);
+  mf_data->initialize_dof_vector (in);
   out.reinit (in);
   ref.reinit (in);
 
index 9d356780a12bb493ca5b9e4b1782ebf9537d0002..4504061770fec2a37a83f8c40befb6741beca710 100644 (file)
@@ -69,14 +69,14 @@ void test ()
   //std::cout << "Number of degrees of freedom: " << dof.n_dofs() << std::endl;
   //std::cout << "Number of constraints: " << constraints.n_constraints() << std::endl;
 
-  MatrixFree<dim,number> mf_data;
+  std_cxx11::shared_ptr<MatrixFree<dim,number> > mf_data(new MatrixFree<dim,number> ());
   {
     const QGauss<1> quad (fe_degree+2);
     typename MatrixFree<dim,number>::AdditionalData data;
     data.tasks_parallel_scheme =
       MatrixFree<dim,number>::AdditionalData::none;
     data.tasks_block_size = 7;
-    mf_data.reinit (dof, constraints, quad, data);
+    mf_data->reinit (dof, constraints, quad, data);
   }
 
   MatrixFreeOperators::MassOperator<dim,fe_degree, fe_degree+2, 1, number> mf;
@@ -86,7 +86,7 @@ void test ()
     = mf.get_matrix_diagonal_inverse()->get_vector();
 
   LinearAlgebra::distributed::Vector<number> in, out, ref;
-  mf_data.initialize_dof_vector (in);
+  mf_data->initialize_dof_vector (in);
   out.reinit (in);
   ref.reinit (in);
 
index bd441203139849314cc356ac0ef44c4145c3da65..3de18fe067c657026e06d821cdfea30d42f51b05 100644 (file)
@@ -68,14 +68,14 @@ void test ()
   //std::cout << "Number of degrees of freedom: " << dof.n_dofs() << std::endl;
   //std::cout << "Number of constraints: " << constraints.n_constraints() << std::endl;
 
-  MatrixFree<dim,number> mf_data;
+  std_cxx11::shared_ptr<MatrixFree<dim,number> > mf_data(new MatrixFree<dim,number> ());
   {
     const QGauss<1> quad (fe_degree+2);
     typename MatrixFree<dim,number>::AdditionalData data;
     data.tasks_parallel_scheme =
       MatrixFree<dim,number>::AdditionalData::none;
     data.tasks_block_size = 7;
-    mf_data.reinit (dof, constraints, quad, data);
+    mf_data->reinit (dof, constraints, quad, data);
   }
 
   MatrixFreeOperators::MassOperator<dim,fe_degree, fe_degree+2, 1, number> mf;
@@ -85,7 +85,7 @@ void test ()
     = mf.get_matrix_diagonal_inverse()->get_vector();
 
   LinearAlgebra::distributed::Vector<number> in, out, ref;
-  mf_data.initialize_dof_vector (in);
+  mf_data->initialize_dof_vector (in);
   out.reinit (in);
   ref.reinit (in);
 
index 4742e0fda6bab54ee5b3aa03ce75105caa1a01f4..653fa4918fdeb9f3524ca9210ab248364df55630 100644 (file)
@@ -114,12 +114,12 @@ void do_test (const DoFHandler<dim>  &dof)
   MappingQ<dim> mapping(fe_degree+1);
 
   LaplaceOperator<dim,fe_degree,n_q_points_1d,1,number> fine_matrix;
-  MatrixFree<dim,number> fine_level_data;
+  std_cxx11::shared_ptr<MatrixFree<dim,number> > fine_level_data(new MatrixFree<dim,number> ());
 
   typename MatrixFree<dim,number>::AdditionalData fine_level_additional_data;
   fine_level_additional_data.tasks_parallel_scheme = MatrixFree<dim,number>::AdditionalData::none;
-  fine_level_data.reinit (mapping, dof, constraints, QGauss<1>(n_q_points_1d),
-                          fine_level_additional_data);
+  fine_level_data->reinit (mapping, dof, constraints, QGauss<1>(n_q_points_1d),
+                           fine_level_additional_data);
 
   fine_matrix.initialize(fine_level_data);
   fine_matrix.compute_diagonal();
@@ -155,7 +155,8 @@ void do_test (const DoFHandler<dim>  &dof)
 
       mg_level_data[level].reinit (mapping, dof, level_constraints, QGauss<1>(n_q_points_1d),
                                    mg_additional_data);
-      mg_matrices[level].initialize(mg_level_data[level],
+      mg_matrices[level].initialize(std_cxx11::make_shared<MatrixFree<dim,number> >(
+                                      mg_level_data[level]),
                                     mg_constrained_dofs,
                                     level);
       mg_matrices[level].compute_diagonal();
index bebcc0cc79361981967620df2c43cc895751bf08..4105d03a80c1b47e89cfdda2d71ba72d3daa2fda 100644 (file)
@@ -151,13 +151,13 @@ void do_test (const DoFHandler<dim>  &dof)
   MappingQ<dim> mapping(fe_degree+1);
 
   LaplaceOperator<dim,fe_degree,n_q_points_1d,1,number> fine_matrix;
-  MatrixFree<dim,number> fine_level_data;
+  std_cxx11::shared_ptr<MatrixFree<dim,number> > fine_level_data(new MatrixFree<dim,number> ());
 
   typename MatrixFree<dim,number>::AdditionalData fine_level_additional_data;
   fine_level_additional_data.tasks_parallel_scheme = MatrixFree<dim,number>::AdditionalData::none;
   fine_level_additional_data.tasks_block_size = 3;
-  fine_level_data.reinit (mapping, dof, constraints, QGauss<1>(n_q_points_1d),
-                          fine_level_additional_data);
+  fine_level_data->reinit (mapping, dof, constraints, QGauss<1>(n_q_points_1d),
+                           fine_level_additional_data);
 
   fine_matrix.initialize(fine_level_data);
   fine_matrix.compute_diagonal();
@@ -204,7 +204,8 @@ void do_test (const DoFHandler<dim>  &dof)
 
       mg_level_data[level].reinit (mapping, dof, level_constraints, QGauss<1>(n_q_points_1d),
                                    mg_additional_data);
-      mg_matrices[level].initialize(mg_level_data[level],
+      mg_matrices[level].initialize(std_cxx11::make_shared<MatrixFree<dim,number> >(
+                                      mg_level_data[level]),
                                     mg_constrained_dofs,
                                     level);
       mg_matrices[level].compute_diagonal();
index e5b948975512950e6f4ec9724a1e275a77773789..e11ff631db7e76393aa28168659510af10897600 100644 (file)
@@ -151,13 +151,13 @@ void do_test (const DoFHandler<dim>  &dof)
   MappingQ<dim> mapping(fe_degree+1);
 
   LaplaceOperator<dim,fe_degree,n_q_points_1d,1,number> fine_matrix;
-  MatrixFree<dim,number> fine_level_data;
+  std_cxx11::shared_ptr<MatrixFree<dim,number> > fine_level_data(new MatrixFree<dim,number> ());
 
   typename MatrixFree<dim,number>::AdditionalData fine_level_additional_data;
   fine_level_additional_data.tasks_parallel_scheme = MatrixFree<dim,number>::AdditionalData::none;
   fine_level_additional_data.tasks_block_size = 3;
-  fine_level_data.reinit (mapping, dof, constraints, QGauss<1>(n_q_points_1d),
-                          fine_level_additional_data);
+  fine_level_data->reinit (mapping, dof, constraints, QGauss<1>(n_q_points_1d),
+                           fine_level_additional_data);
 
   fine_matrix.initialize(fine_level_data);
   fine_matrix.compute_diagonal();
@@ -204,7 +204,8 @@ void do_test (const DoFHandler<dim>  &dof)
 
       mg_level_data[level].reinit (mapping, dof, level_constraints, QGauss<1>(n_q_points_1d),
                                    mg_additional_data);
-      mg_matrices[level].initialize(mg_level_data[level],
+      mg_matrices[level].initialize(std_cxx11::make_shared<MatrixFree<dim,number> > (
+                                      mg_level_data[level]),
                                     mg_constrained_dofs,
                                     level);
       mg_matrices[level].compute_diagonal();

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.