]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Also instantiate with PETScWrappers::MPI::{Block,}Vector.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 12 Feb 2010 02:41:54 +0000 (02:41 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Fri, 12 Feb 2010 02:41:54 +0000 (02:41 +0000)
git-svn-id: https://svn.dealii.org/trunk@20561 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/include/fe/fe_values.h
deal.II/deal.II/source/fe/fe_values.cc

index c405834f3d72ccbe10326ad03645756cb8e426ab..6594d2092232636936b43998352456de234a7355 100644 (file)
@@ -2,7 +2,7 @@
 //    $Id$
 //    Version: $Name:  $
 //
-//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009 by the deal.II authors
+//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
@@ -58,6 +58,12 @@ namespace PETScWrappers
 {
   class Vector;
   class BlockVector;
+
+  namespace MPI
+  {
+    class Vector;
+    class BlockVector;
+  }
 }
 #endif
 
@@ -2827,6 +2833,28 @@ class FEValuesBase : protected FEValuesData<dim,spacedim>,
         void
         get_interpolated_dof_values (const PETScWrappers::BlockVector &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::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
index 7686ee4d84bd5ac7743fe96597496e8eab9a13bd..0094d67eec094d58eed66e7c8dac0ff8c8e8c4c7 100644 (file)
@@ -2,7 +2,7 @@
 //    $Id$
 //    Version: $Name$
 //
-//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009 by the deal.II authors
+//    Copyright (C) 1998, 1999, 2000, 2001, 2002, 2003, 2004, 2005, 2006, 2007, 2008, 2009, 2010 by the deal.II authors
 //
 //    This file is subject to QPL and may not be  distributed
 //    without copyright and license information. Please refer
@@ -1246,6 +1246,28 @@ class FEValuesBase<dim,spacedim>::CellIterator : public FEValuesBase<dim,spacedi
     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
@@ -1478,6 +1500,28 @@ class FEValuesBase<dim,spacedim>::TriaCellIterator : public FEValuesBase<dim,spa
     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
@@ -1671,6 +1715,28 @@ get_interpolated_dof_values (const PETScWrappers::BlockVector &in,
   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
@@ -1850,6 +1916,26 @@ get_interpolated_dof_values (const PETScWrappers::BlockVector &,
   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

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.