]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
Merge from trunk and some other small fixme patches.
authoryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 13 Feb 2014 11:00:34 +0000 (11:00 +0000)
committeryoung <young@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 13 Feb 2014 11:00:34 +0000 (11:00 +0000)
git-svn-id: https://svn.dealii.org/branches/branch_petscscalar_complex@32477 0785d39b-7218-0410-832d-ea1e28bc413d

17 files changed:
1  2 
deal.II/cmake/config/template-arguments.in
deal.II/include/deal.II/lac/petsc_vector_base.h
deal.II/include/deal.II/numerics/vector_tools.templates.h
deal.II/source/dofs/dof_accessor_get.inst.in
deal.II/source/dofs/dof_accessor_set.inst.in
deal.II/source/fe/fe_tools_interpolate.cc
deal.II/source/fe/fe_values.cc
deal.II/source/fe/fe_values.decl.1.inst.in
deal.II/source/fe/fe_values.decl.2.inst.in
deal.II/source/fe/fe_values.impl.1.inst.in
deal.II/source/fe/fe_values.impl.2.inst.in
deal.II/source/fe/fe_values.inst.in
deal.II/source/meshworker/mesh_worker_vector_selector.inst.in
deal.II/source/numerics/data_out_faces.cc
deal.II/source/numerics/data_out_rotation.cc
deal.II/source/numerics/matrix_tools.cc
tests/tests.h

index 7e97a5eea3338507dd7a07f1a1bbe410d5df0ed1,7e97a5eea3338507dd7a07f1a1bbe410d5df0ed1..c95c0f21d6440d6a93c2463595f9d1d947dc2eb0
@@@ -37,6 -37,6 +37,36 @@@ SERIAL_VECTORS := { Vector<double>
                    @DEAL_II_EXPAND_PETSC_MPI_BLOCKVECTOR@;
                    }
  
++REAL_SERIAL_VECTORS := { Vector<double>;
++                    Vector<float> ;
++                  Vector<long double>;
++
++                  BlockVector<double>;
++                    BlockVector<float>;
++                    BlockVector<long double>;
++
++                  parallel::distributed::Vector<double>;
++                  parallel::distributed::Vector<float> ;
++                  parallel::distributed::Vector<long double>;
++
++                    parallel::distributed::BlockVector<double>;
++                    parallel::distributed::BlockVector<float> ;
++                    parallel::distributed::BlockVector<long double>;
++
++                  @DEAL_II_EXPAND_TRILINOS_VECTOR@;
++                  @DEAL_II_EXPAND_TRILINOS_MPI_VECTOR@;
++
++                  @DEAL_II_EXPAND_TRILINOS_BLOCKVECTOR@;
++                  @DEAL_II_EXPAND_TRILINOS_MPI_BLOCKVECTOR@;
++                  }
++
++COMPLEX_SERIAL_VECTORS := { @DEAL_II_EXPAND_PETSC_VECTOR@;
++                          @DEAL_II_EXPAND_PETSC_MPI_VECTOR@;
++
++                            @DEAL_II_EXPAND_PETSC_BLOCKVECTOR@;
++                          @DEAL_II_EXPAND_PETSC_MPI_BLOCKVECTOR@;
++                          }
++
  EXTERNAL_SEQUENTIAL_VECTORS := { @DEAL_II_EXPAND_TRILINOS_VECTOR@;
                                   @DEAL_II_EXPAND_TRILINOS_BLOCKVECTOR@;
                                   @DEAL_II_EXPAND_PETSC_VECTOR@;
index 00494a4fdd7f63b516dcdf27ea22bf96333b0f7b,163f03d2b7f0350ad49cb53e16aead969bea3dd3..5a99c08bce30ae4f394d6c36447de93a85a940c6
@@@ -41,7 -41,7 +41,7 @@@ DEAL_II_NAMESPACE_OPE
  namespace
  {
    template <class VectorType>
--  double
++  typename VectorType::value_type
    get_vector_element (const VectorType &vector,
                        const types::global_dof_index cell_number)
    {
index e220a5e85bf6df1a00a46eabd0aa8f10f5312142,e220a5e85bf6df1a00a46eabd0aa8f10f5312142..83575b40a9b4984e08853c608698d21e6f396c69
@@@ -19,7 -19,7 +19,7 @@@
  // Declarations of member functions of FEValuesBase::CellIteratorBase
  // and derived classes
  
--for (VEC : SERIAL_VECTORS)
++for (VEC : REAL_SERIAL_VECTORS; S : REAL_SCALARS)
    {
                                      /// Call
                                      /// @p get_interpolated_dof_values
      virtual
      void
      get_interpolated_dof_values (const VEC &in,
--                               Vector<double> &out) const = 0;
++                               Vector<S> &out) const = 0;
++  }
++
++for (VEC : COMPLEX_SERIAL_VECTORS; S : COMPLEX_SCALARS)
++  {
++                                    /// Call
++                                    /// @p get_interpolated_dof_values
++                                    /// of the iterator with the
++                                    /// given arguments.
++    virtual
++    void
++    get_interpolated_dof_values (const VEC &in,
++                               Vector<S> &out) const = 0;
    }
index 567c93f4352e7debd11bec9385f95faedb49bd5e,567c93f4352e7debd11bec9385f95faedb49bd5e..cb183eeda6723d26bf8a076fa131fbde33eba376
@@@ -19,7 -19,7 +19,7 @@@
  // Declarations of member functions of FEValuesBase::CellIteratorBase
  // and derived classes
  
--for (VEC : SERIAL_VECTORS)
++for (VEC : REAL_SERIAL_VECTORS; S : REAL_SCALARS)
    {
                                      /// Call
                                      /// @p get_interpolated_dof_values
      get_interpolated_dof_values (const VEC      &in,
                                 Vector<double> &out) const;
    }
++
++for (VEC : COMPLEX_SERIAL_VECTORS; S : COMPLEX_SCALARS)
++  {
++                                    /// Call
++                                    /// @p get_interpolated_dof_values
++                                    /// of the iterator with the
++                                    /// given arguments.
++    virtual
++    void
++    get_interpolated_dof_values (const VEC      &in,
++                               Vector<S> &out) const;
++  }
index 3dfaaa1a9d03da043795ac3c1d8caf45ad1494d4,3dfaaa1a9d03da043795ac3c1d8caf45ad1494d4..731d6ccb5a9b2a322c7c8ad49ded9c4566a87d0d
@@@ -15,7 -15,7 +15,7 @@@
  // ---------------------------------------------------------------------
  
  
--for (VEC : SERIAL_VECTORS)
++for (VEC : REAL_SERIAL_VECTORS)
    {
      template <int dim, int spacedim>
      template <typename CI>
        cell->get_interpolated_dof_values (in, out);
      \}
    }
++
++for (VEC : COMPLEX_SERIAL_VECTORS)
++  {
++    template <int dim, int spacedim>
++    template <typename CI>
++    void
++    FEValuesBase<dim,spacedim>::CellIterator<CI>::
++    get_interpolated_dof_values (const VEC                     &in,
++                               Vector<std::complex<double> > &out) const
++    \{
++      cell->get_interpolated_dof_values (in, out);
++    \}
++  }
++
index 95883209fd6c7ca205244e2b10c2e7f0d5979ab2,95883209fd6c7ca205244e2b10c2e7f0d5979ab2..4209718bf67d7b1df8231e196284afd5b8bc930d
  // ---------------------------------------------------------------------
  
  
--for (VEC : SERIAL_VECTORS)
++for (VEC : REAL_SERIAL_VECTORS)
    {
      template <int dim, int spacedim>
      void
      FEValuesBase<dim,spacedim>::TriaCellIterator::
--    get_interpolated_dof_values (const VEC &,
++    get_interpolated_dof_values (const VEC      &,
                                 Vector<double> &) const
      \{
        Assert (false, ExcMessage (message_string));
      \}
    }
++
++
++for (VEC : COMPLEX_SERIAL_VECTORS)
++  {
++    template <int dim, int spacedim>
++    void
++    FEValuesBase<dim,spacedim>::TriaCellIterator::
++    get_interpolated_dof_values (const VEC                     &,
++                               Vector<std::complex<double> > &) const
++    \{
++      Assert (false, ExcMessage (message_string));
++    \}
++  }
++
index e63ce3fced7f7cb2c2824daf29e2cb92c2980a70,e63ce3fced7f7cb2c2824daf29e2cb92c2980a70..d0201ca59dccb14b83671793512298d5e1ca398a
@@@ -64,7 -64,7 +64,7 @@@ for (dof_handler : DOFHANDLER_TEMPLATES
  
  
  
--for (VEC : SERIAL_VECTORS; deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS)
++for (VEC : REAL_SERIAL_VECTORS; deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS)
    {
  #if deal_II_dimension <= deal_II_space_dimension
  #if (deal_II_space_dimension == DIM_A) || (deal_II_space_dimension == DIM_B)
  
  
  
--for (VEC : SERIAL_VECTORS; deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS)
++for (VEC : REAL_SERIAL_VECTORS; deal_II_dimension : DIMENSIONS; deal_II_space_dimension :  SPACE_DIMENSIONS)
    {
  #if deal_II_dimension <= deal_II_space_dimension
  #if (deal_II_space_dimension == DIM_A) || (deal_II_space_dimension == DIM_B)
index 4de6186fc85d229041ef1fe3855661451e11a825,4de6186fc85d229041ef1fe3855661451e11a825..7310c234193c0a35a232f4114eb98ac512e02978
@@@ -24,7 -24,7 +24,7 @@@ for (deal_II_dimension : DIMENSIONS
    \}
  }
  
--for (VECTOR : SERIAL_VECTORS ; deal_II_dimension : DIMENSIONS)
++for (VECTOR : REAL_SERIAL_VECTORS ; deal_II_dimension : DIMENSIONS)
  {
    namespace MeshWorker
    \{
index aa5fa1a594820f1b5162eac8e00ab107fe7535fd,a14c34428f21b8434bf1eb21b8a4ea8db7d96b89..2ace9d307cd3af56ee0f84f3f49fda2303b0deb1
@@@ -246,9 -249,8 +249,9 @@@ build_one_patch (const FaceDescriptor *
            Assert (cell_and_face->first->active(), ExcCellNotActiveForCellData());
            const unsigned int cell_number
              = std::distance (this->triangulation->begin_active(),
-                              typename Triangulation<DH::dimension,DH::space_dimension>::active_cell_iterator(cell_and_face->first));
+                              typename Triangulation<dimension,space_dimension>::active_cell_iterator(cell_and_face->first));
  
 +        // @whattodo - this next line does not cause a problem (yet), but probably will one day.
            const double value
              = this->cell_data[dataset]->get_cell_data_value (cell_number);
            for (unsigned int q=0; q<n_q_points; ++q)
index 00238140cab3efc5ea8bd08f1832ed3ab748940b,e53a3f351346b6f5baa31a641324deba5c763e45..d1a9aeda62f58b765a6345f26b806344e79ed6af
@@@ -364,12 -364,10 +364,12 @@@ build_one_patch (const cell_iterator *c
                        ExcMessage("Cell must be active for cell data"));
                const unsigned int cell_number
                  = std::distance (this->triangulation->begin_active(),
-                                  typename Triangulation<DH::dimension,DH::space_dimension>::active_cell_iterator(*cell));
+                                  typename Triangulation<dimension,space_dimension>::active_cell_iterator(*cell));
 +
 +            // @whattodo - this next line does not cause a problem (yet), but probably will one day.
                const double value
                  = this->cell_data[dataset]->get_cell_data_value (cell_number);
-               switch (DH::dimension)
+               switch (dimension)
                  {
                  case 1:
                    for (unsigned int x=0; x<n_points; ++x)
index 34a8f5f962d18607f084a2f79510f63d90a82384,6a726729c127a62bb53e14b7e8c6616667c20e93..0c98dc8b09f2988f84276b9aff5070e207befe34
@@@ -2366,17 -2350,13 +2350,13 @@@ namespace MatrixTool
          // entry from within the part of the
          // matrix that we can see. if we can't
          // find such an entry, take one
-       // @whattodo Note: this is a real mess.
-         PetscScalar average_nonzero_diagonal_entry;
-         // PetscScalar average_nonzero_diagonal_entry = 1;
-         // for (types::global_dof_index i=local_range.first; i<local_range.second; ++i)
-         //   if (matrix.diag_element(i) != 0)
-         //     {
-         //       average_nonzero_diagonal_entry = std::fabs(matrix.diag_element(i));
-         //       break;
-         //     }
-       Assert ((false), ExcMessage ("This function is corrupt: @whattodo"));
 -        PetscScalar average_nonzero_diagonal_entry = 1;
 -        for (types::global_dof_index i=local_range.first; i<local_range.second; ++i)
 -          if (matrix.diag_element(i) != 0)
 -            {
 -              average_nonzero_diagonal_entry = std::fabs(matrix.diag_element(i));
 -              break;
 -            }
++      PetscScalar average_nonzero_diagonal_entry = 1;
++      for (types::global_dof_index i=local_range.first; i<local_range.second; ++i)
++        if (matrix.diag_element(i) != PetscScalar ())
++          {
++            average_nonzero_diagonal_entry = matrix.diag_element(i);
++            break;
++          }
  
          // figure out which rows of the matrix we
          // have to eliminate on this processor
diff --cc tests/tests.h
Simple merge

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.