]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Clean up function interpolate and instantiate correctly. Leave a TODO.
authoryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 27 Mar 2014 09:19:37 +0000 (09:19 +0000)
committeryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 27 Mar 2014 09:19:37 +0000 (09:19 +0000)
git-svn-id: https://svn.dealii.org/branches/branch_petscscalar_complex@32699 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/numerics/vector_tools.templates.h
deal.II/source/numerics/vector_tools_interpolate.inst.in

index e41bd1be6eac16f42e9eafdc0544fe10b9390a6c..7216a3c813a0c8dd873571c4bc41f2ae1c08e3d0 100644 (file)
@@ -280,20 +280,17 @@ namespace VectorTools
 
 
 
-  // @whattodo TODO: To return functionality for complex numbers, this
-  // needs to be templated with Number and instantiated for
-  // REAL_SERIAL_VECTORS+double and
-  // COMPLEX_SERIAL_VECTORS+complex<double>
-  template <int dim, class InVector, class OutVector, int spacedim>
+  // @TODO[TY] also enable instantiation for COMPLEX_SERIAL_VECTORS and std::complex<double>
+  template <int dim, class InVector, class OutVector, int spacedim, typename number>
   void
   interpolate (const DoFHandler<dim,spacedim>  &dof_1,
                const DoFHandler<dim,spacedim>  &dof_2,
-               const FullMatrix<double>        &transfer, // <- FullMatrix<std::complex>
+               const FullMatrix<number>        &transfer,
                const InVector                  &data_1,
                OutVector                       &data_2)
   {
-    Vector<double> cell_data_1(dof_1.get_fe().dofs_per_cell); // <- Vector<std::complex>
-    Vector<double> cell_data_2(dof_2.get_fe().dofs_per_cell); // <- Vector<std::complex>
+    Vector<number> cell_data_1(dof_1.get_fe().dofs_per_cell);
+    Vector<number> cell_data_2(dof_2.get_fe().dofs_per_cell);
 
     std::vector<short unsigned int> touch_count (dof_2.n_dofs(), 0); //TODO: check on datatype... kinda strange (UK)
     std::vector<types::global_dof_index>       local_dof_indices (dof_2.get_fe().dofs_per_cell);
index 9b466454111bf22566edfe6c396b68f6df214c86..00bb3370eec9b768ee4801da73eaf1cf1d41055c 100644 (file)
@@ -33,11 +33,21 @@ for (VEC : REAL_SERIAL_VECTORS ; deal_II_dimension : DIMENSIONS; deal_II_space_d
          const Function<deal_II_space_dimension>&,
          VEC&);
 
+      \}
+#endif
+  }
+
+
+for (VEC : REAL_SERIAL_VECTORS ; deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS; NUMBER : double)
+  {
+#if deal_II_dimension <= deal_II_space_dimension
+    namespace VectorTools \{
+
       template
         void interpolate
         (const DoFHandler<deal_II_dimension,deal_II_space_dimension>&,
          const DoFHandler<deal_II_dimension,deal_II_space_dimension>&,
-         const FullMatrix<double>&,
+         const FullMatrix<NUMBER>&,
          const VEC&,
          VEC&);
 

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.