]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Use the preprocessor to avoid writing lots of declarations and instantiations.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 14 Feb 2010 05:37:30 +0000 (05:37 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Sun, 14 Feb 2010 05:37:30 +0000 (05:37 +0000)
git-svn-id: https://svn.dealii.org/trunk@20616 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/source/fe/fe_values.cc
deal.II/deal.II/source/fe/fe_values.decl.1.inst.in [new file with mode: 0644]
deal.II/deal.II/source/fe/fe_values.decl.2.inst.in [new file with mode: 0644]
deal.II/deal.II/source/fe/fe_values.impl.2.inst.in [new file with mode: 0644]
deal.II/deal.II/source/fe/fe_values.impl.inst.in [new file with mode: 0644]

index ab6361add0c27312348e40b73b7a623590a273cd..ae4e039866dfbfe390ce565987ef5211c0273f4c 100644 (file)
@@ -1141,141 +1141,7 @@ class FEValuesBase<dim,spacedim>::CellIteratorBase
     unsigned int
     n_dofs_for_dof_handler () const = 0;
 
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const Vector<double> &in,
-                                Vector<double>       &out) const = 0;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const Vector<float> &in,
-                                Vector<float>       &out) const = 0;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const Vector<long double> &in,
-                                Vector<long double>       &out) const = 0;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const BlockVector<double> &in,
-                                Vector<double>       &out) const = 0;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const BlockVector<float> &in,
-                                Vector<float>       &out) const = 0;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const BlockVector<long double> &in,
-                                Vector<long double>       &out) const = 0;
-
-#ifdef DEAL_II_USE_PETSC
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const PETScWrappers::Vector &in,
-                                Vector<PetscScalar>         &out) const = 0;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const PETScWrappers::BlockVector &in,
-                                Vector<PetscScalar>              &out) const = 0;
-#endif
-
-#ifdef DEAL_II_USE_TRILINOS
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::Vector &in,
-                                Vector<TrilinosScalar>         &out) const = 0;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::BlockVector &in,
-                                Vector<TrilinosScalar>              &out) const = 0;
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::MPI::Vector &in,
-                                Vector<TrilinosScalar>              &out) const = 0;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::MPI::BlockVector &in,
-                                Vector<TrilinosScalar>                   &out) const = 0;
-#endif
+#include "fe_values.decl.1.inst"
 };
 
 
@@ -1335,163 +1201,7 @@ class FEValuesBase<dim,spacedim>::CellIterator : public FEValuesBase<dim,spacedi
     unsigned int
     n_dofs_for_dof_handler () const;
 
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const Vector<double> &in,
-                                Vector<double>       &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const Vector<float> &in,
-                                Vector<float>       &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const Vector<long double> &in,
-                                Vector<long double>       &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const BlockVector<double> &in,
-                                Vector<double>       &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const BlockVector<float> &in,
-                                Vector<float>       &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const BlockVector<long double> &in,
-                                Vector<long double>       &out) const;
-
-#ifdef DEAL_II_USE_PETSC
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const PETScWrappers::Vector &in,
-                                Vector<PetscScalar>         &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const PETScWrappers::BlockVector &in,
-                                Vector<PetscScalar>              &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const PETScWrappers::MPI::Vector &in,
-                                Vector<PetscScalar>         &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const PETScWrappers::MPI::BlockVector &in,
-                                Vector<PetscScalar>              &out) const;
-#endif
-
-#ifdef DEAL_II_USE_TRILINOS
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::Vector &in,
-                                Vector<TrilinosScalar>         &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::BlockVector &in,
-                                Vector<TrilinosScalar>              &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::MPI::Vector &in,
-                                Vector<TrilinosScalar>              &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::MPI::BlockVector &in,
-                                Vector<TrilinosScalar>                   &out) const;
-#endif
+#include "fe_values.decl.2.inst"
 
   private:
                                     /**
@@ -1585,167 +1295,7 @@ class FEValuesBase<dim,spacedim>::TriaCellIterator : public FEValuesBase<dim,spa
     unsigned int
     n_dofs_for_dof_handler () const;
 
-                                    /**
-                                     * Implement the respective
-                                     * function of the base
-                                     * class. Since this is not
-                                     * possible, we just raise an
-                                     * error.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const Vector<double> &in,
-                                Vector<double>       &out) const;
-
-                                    /**
-                                     * Implement the respective
-                                     * function of the base
-                                     * class. Since this is not
-                                     * possible, we just raise an
-                                     * error.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const Vector<float> &in,
-                                Vector<float>       &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const Vector<long double> &in,
-                                Vector<long double>       &out) const;
-
-                                    /**
-                                     * Implement the respective
-                                     * function of the base
-                                     * class. Since this is not
-                                     * possible, we just raise an
-                                     * error.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const BlockVector<double> &in,
-                                Vector<double>       &out) const;
-
-                                    /**
-                                     * Implement the respective
-                                     * function of the base
-                                     * class. Since this is not
-                                     * possible, we just raise an
-                                     * error.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const BlockVector<float> &in,
-                                Vector<float>       &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const BlockVector<long double> &in,
-                                Vector<long double>       &out) const;
-
-#ifdef DEAL_II_USE_PETSC
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const PETScWrappers::Vector &in,
-                                Vector<PetscScalar>         &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const PETScWrappers::BlockVector &in,
-                                Vector<PetscScalar>              &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const PETScWrappers::MPI::Vector &in,
-                                Vector<PetscScalar>         &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const PETScWrappers::MPI::BlockVector &in,
-                                Vector<PetscScalar>              &out) const;
-#endif
-
-#ifdef DEAL_II_USE_TRILINOS
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::Vector &in,
-                                Vector<TrilinosScalar>         &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::BlockVector &in,
-                                Vector<TrilinosScalar>              &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::MPI::Vector &in,
-                                Vector<TrilinosScalar>              &out) const;
-
-                                    /**
-                                     * Call
-                                     * @p get_interpolated_dof_values
-                                     * of the iterator with the
-                                     * given arguments.
-                                     */
-    virtual
-    void
-    get_interpolated_dof_values (const TrilinosWrappers::MPI::BlockVector &in,
-                                Vector<TrilinosScalar>                   &out) const;
-#endif
+#include "fe_values.decl.2.inst"
 
   private:
                                     /**
@@ -1801,171 +1351,7 @@ FEValuesBase<dim,spacedim>::CellIterator<CI>::n_dofs_for_dof_handler () const
 
 
 
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const Vector<double> &in,
-                            Vector<double>       &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const Vector<float> &in,
-                            Vector<float>       &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const Vector<long double> &in,
-                            Vector<long double>       &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const BlockVector<double> &in,
-                            Vector<double>            &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const BlockVector<float> &in,
-                            Vector<float>            &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const BlockVector<long double> &in,
-                            Vector<long double>            &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-
-
-#ifdef DEAL_II_USE_PETSC
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const PETScWrappers::Vector &in,
-                            Vector<PetscScalar>         &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const PETScWrappers::BlockVector &in,
-                            Vector<PetscScalar>              &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const PETScWrappers::MPI::Vector &in,
-                            Vector<PetscScalar>         &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const PETScWrappers::MPI::BlockVector &in,
-                            Vector<PetscScalar>              &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-#endif
-
-#ifdef DEAL_II_USE_TRILINOS
-
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const TrilinosWrappers::Vector &in,
-                            Vector<TrilinosScalar>         &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim,spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const TrilinosWrappers::BlockVector &in,
-                            Vector<TrilinosScalar>              &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim, spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const TrilinosWrappers::MPI::Vector &in,
-                            Vector<TrilinosScalar>              &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-
-
-template <int dim, int spacedim>
-template <typename CI>
-void
-FEValuesBase<dim, spacedim>::CellIterator<CI>::
-get_interpolated_dof_values (const TrilinosWrappers::MPI::BlockVector &in,
-                            Vector<TrilinosScalar>                   &out) const
-{
-  cell->get_interpolated_dof_values (in, out);
-}
-
-#endif
+#include "fe_values.impl.inst"
 
 
 /* ---------------- FEValuesBase<dim,spacedim>::TriaCellIterator --------- */
@@ -2009,158 +1395,7 @@ FEValuesBase<dim,spacedim>::TriaCellIterator::n_dofs_for_dof_handler () const
 }
 
 
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const Vector<double> &,
-                            Vector<double>       &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const Vector<float> &,
-                            Vector<float>       &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const Vector<long double> &,
-                            Vector<long double>       &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const BlockVector<double> &,
-                            Vector<double>            &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const BlockVector<float> &,
-                            Vector<float>            &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const BlockVector<long double> &,
-                            Vector<long double>            &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-
-
-#ifdef DEAL_II_USE_PETSC
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const PETScWrappers::Vector &,
-                            Vector<PetscScalar>         &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const PETScWrappers::BlockVector &,
-                            Vector<PetscScalar>              &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const PETScWrappers::MPI::Vector &,
-                            Vector<PetscScalar>         &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const PETScWrappers::MPI::BlockVector &,
-                            Vector<PetscScalar>              &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-#endif
-
-#ifdef DEAL_II_USE_TRILINOS
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const TrilinosWrappers::Vector &,
-                            Vector<TrilinosScalar>         &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const TrilinosWrappers::BlockVector &,
-                            Vector<TrilinosScalar>              &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim, spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const TrilinosWrappers::MPI::Vector &,
-                            Vector<TrilinosScalar>              &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-
-
-template <int dim, int spacedim>
-void
-FEValuesBase<dim,spacedim>::TriaCellIterator::
-get_interpolated_dof_values (const TrilinosWrappers::MPI::BlockVector &,
-                            Vector<TrilinosScalar>                   &) const
-{
-  Assert (false, ExcMessage (message_string));
-}
-
-#endif
-
+#include "fe_values.impl.2.inst"
 
 
 
diff --git a/deal.II/deal.II/source/fe/fe_values.decl.1.inst.in b/deal.II/deal.II/source/fe/fe_values.decl.1.inst.in
new file mode 100644 (file)
index 0000000..6646733
--- /dev/null
@@ -0,0 +1,30 @@
+//-----------------------------------------------------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2010 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//-----------------------------------------------------------------------
+
+
+// Declarations of member functions of FEValuesBase::CellIteratorBase
+// and derived classes
+
+for (VEC : SERIAL_VECTORS)
+  {
+                                    /**
+                                     * Call
+                                     * @p get_interpolated_dof_values
+                                     * of the iterator with the
+                                     * given arguments.
+                                     */
+    virtual
+    void
+    get_interpolated_dof_values (const VEC &in,
+                                Vector<VEC::value_type> &out) const = 0;
+  }
diff --git a/deal.II/deal.II/source/fe/fe_values.decl.2.inst.in b/deal.II/deal.II/source/fe/fe_values.decl.2.inst.in
new file mode 100644 (file)
index 0000000..6e86304
--- /dev/null
@@ -0,0 +1,30 @@
+//-----------------------------------------------------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2010 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//-----------------------------------------------------------------------
+
+
+// Declarations of member functions of FEValuesBase::CellIteratorBase
+// and derived classes
+
+for (VEC : SERIAL_VECTORS)
+  {
+                                    /**
+                                     * Call
+                                     * @p get_interpolated_dof_values
+                                     * of the iterator with the
+                                     * given arguments.
+                                     */
+    virtual
+    void
+    get_interpolated_dof_values (const VEC &in,
+                                Vector<VEC::value_type> &out) const;
+  }
diff --git a/deal.II/deal.II/source/fe/fe_values.impl.2.inst.in b/deal.II/deal.II/source/fe/fe_values.impl.2.inst.in
new file mode 100644 (file)
index 0000000..0304e04
--- /dev/null
@@ -0,0 +1,24 @@
+//-----------------------------------------------------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2010 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//-----------------------------------------------------------------------
+
+for (VEC : SERIAL_VECTORS)
+  {
+    template <int dim, int spacedim>
+    void
+    FEValuesBase<dim,spacedim>::TriaCellIterator::
+    get_interpolated_dof_values (const VEC &,
+                                Vector<VEC::value_type>       &) const
+    \{
+      Assert (false, ExcMessage (message_string));
+    \}
+  }
diff --git a/deal.II/deal.II/source/fe/fe_values.impl.inst.in b/deal.II/deal.II/source/fe/fe_values.impl.inst.in
new file mode 100644 (file)
index 0000000..ef05984
--- /dev/null
@@ -0,0 +1,25 @@
+//-----------------------------------------------------------------------
+//    $Id$
+//    Version: $Name$
+//
+//    Copyright (C) 2010 by the deal.II authors
+//
+//    This file is subject to QPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//-----------------------------------------------------------------------
+
+for (VEC : SERIAL_VECTORS)
+  {
+    template <int dim, int spacedim>
+    template <typename CI>
+    void
+    FEValuesBase<dim,spacedim>::CellIterator<CI>::
+    get_interpolated_dof_values (const VEC &in,
+                                Vector<VEC::value_type>       &out) const
+    \{
+      cell->get_interpolated_dof_values (in, out);
+    \}
+  }

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.