protected:
/**
* Create a vector of length @p{n}. For
- * this class, we create a sequential
+ * this class, we create a parallel
* vector. @arg n denotes the total
- * size of the vector to be created.
+ * size of the vector to be
+ * created. @arg local_size denotes how
+ * many of these elements shall be
+ * stored locally. The last argument is
+ * ignored for sequential vectors.
*/
- virtual void create_vector (const unsigned int n);
+ virtual void create_vector (const unsigned int n,
+ const unsigned int local_size = 0);
};
template <typename number>
Vector::Vector (const ::Vector<number> &v)
{
- create_vector (v.size());
+ Vector::create_vector (v.size());
VectorBase::operator = (v);
}
bool operator != (const VectorBase &v) const;
/**
- * Return dimension of the vector.
+ * Return the global dimension of the vector.
*/
unsigned int size () const;
+ /**
+ * Return the local dimension of the
+ * vector, i.e. the number of elements
+ * stored on the present MPI
+ * process. For sequential vectors,
+ * this number is the same as size(),
+ * but for parallel vectors it may be
+ * smaller.
+ */
+ unsigned int local_size () const;
+
/**
* Provide access to a given element,
* both read and write.
mutable LastAction::Values last_action;
/**
- * Create a vector of length
- * @p{n}. Derived classes have to
- * overload this function according to
- * the type of vector they create
- * (e.g. sequential or
- * parallel/distributed vectors).
+ * Create a vector of length @p{n}. For
+ * this class, we create a parallel
+ * vector. @arg n denotes the total
+ * size of the vector to be
+ * created. @arg local_size denotes how
+ * many of these elements shall be
+ * stored locally. The last argument is
+ * ignored for sequential vectors.
+ */
+ virtual void create_vector (const unsigned int n,
+ const unsigned int local_size = 0) = 0;
+
+ /**
+ * Make the reference class a friend.
*/
- virtual void create_vector (const unsigned int n) = 0;
-
friend class internal::VectorReference;
};
Vector::Vector ()
{
- const int n = 0;
- const int ierr
- = VecCreateSeq (PETSC_COMM_SELF, n, &vector);
- AssertThrow (ierr == 0, ExcPETScError(ierr));
+ Vector::create_vector (0);
}
Vector::Vector (const unsigned int n)
{
- const int ierr
- = VecCreateSeq (PETSC_COMM_SELF, n, &vector);
- AssertThrow (ierr == 0, ExcPETScError(ierr));
+ Vector::create_vector (n);
}
Vector::Vector (const VectorBase &v)
{
- int ierr
- = VecCreateSeq (PETSC_COMM_SELF, v.size(), &vector);
- AssertThrow (ierr == 0, ExcPETScError(ierr));
-
+ Vector::create_vector (v.size());
VectorBase::operator = (v);
}
void
- Vector::create_vector (const unsigned int n)
+ Vector::create_vector (const unsigned int n,
+ const unsigned int /*local_size*/)
{
const int ierr
= VecCreateSeq (PETSC_COMM_SELF, n, &vector);
}
+
bool
VectorBase::operator != (const VectorBase &v) const
{
}
+
unsigned int
VectorBase::size () const
{
+ unsigned int
+ VectorBase::local_size () const
+ {
+ int sz;
+ const int ierr = VecGetLocalSize (vector, &sz);
+ AssertThrow (ierr == 0, ExcPETScError(ierr));
+
+ return sz;
+ }
+
+
+
PetscScalar
VectorBase::operator * (const VectorBase &vec) const
{