<h3>Specific improvements</h3>
<ol>
+ <li>New: There is now a method to copy the content from a
+ PETScWrappers::MPI::Vector to deal.II's parallel distributed vector.
+ <br>
+ (Ben Thompson, 2014/01/31)
+
<li>Fixed: The SolutionTransfer class had all sorts of problems when
used with hp::DoFHandler that made its results at least questionable.
Several parts of this class have been rewritten to make the results
DEAL_II_NAMESPACE_OPEN
+#ifdef DEAL_II_WITH_PETSC
+namespace PETScWrappers
+{
+ namespace MPI
+ {
+ class Vector;
+ }
+}
+#endif
+
+
namespace parallel
{
namespace distributed
Vector<Number> &
operator = (const Vector<Number2> &in_vector);
+#ifdef DEAL_II_WITH_PETSC
+ /**
+ * Copy the content of a PETSc vector into a vector. This function
+ * assumes that the vectors layouts have already been initialized to
+ * match.
+ *
+ * This operator is only available if deal.II was configured with PETSc.
+ */
+ Vector<Number> &
+ operator = (const PETScWrappers::MPI::Vector &petsc_vec);
+#endif
+
/**
* This method copies the local range from another vector with the same
* local range, but possibly different layout of ghost indices.
#include <deal.II/base/config.h>
#include <deal.II/lac/parallel_vector.h>
+#include <deal.II/lac/petsc_parallel_vector.h>
+
DEAL_II_NAMESPACE_OPEN
+#ifdef DEAL_II_WITH_PETSC
+
+ template <typename Number>
+ inline
+ Vector<Number> &
+ Vector<Number>::operator = (const PETScWrappers::MPI::Vector &petsc_vec)
+ {
+ Assert(petsc_vec.locally_owned_elements() == locally_owned_elements(),
+ StandardExceptions::ExcInvalidState());
+
+ // get a representation of the vector and copy it
+ PetscScalar *start_ptr;
+ int ierr = VecGetArray (static_cast<const Vec &>(petsc_vec), &start_ptr);
+ AssertThrow (ierr == 0, ExcPETScError(ierr));
+
+ const size_type vec_size = local_size();
+ std::copy (start_ptr, start_ptr + vec_size, begin());
+
+ // restore the representation of the vector
+ ierr = VecRestoreArray (static_cast<const Vec &>(petsc_vec), &start_ptr);
+ AssertThrow (ierr == 0, ExcPETScError(ierr));
+
+ // spread ghost values between processes?
+ bool must_update_ghost_values = vector_is_ghosted ||
+ petsc_vec.has_ghost_elements();
+ if (must_update_ghost_values)
+ update_ghost_values();
+
+ // return a pointer to this object per normal c++ operator overloading
+ // semantics
+ return *this;
+ }
+
+#endif
+
+
+
template <typename Number>
void
Vector<Number>::copy_from (const Vector<Number> &c,