Reference documentation for deal.II version 9.0.0
petsc_parallel_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_petsc_parallel_vector_h
17 #define dealii_petsc_parallel_vector_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_PETSC
23 
24 # include <deal.II/base/subscriptor.h>
25 # include <deal.II/lac/exceptions.h>
26 # include <deal.II/lac/vector.h>
27 # include <deal.II/lac/vector_operation.h>
28 # include <deal.II/lac/petsc_vector_base.h>
29 # include <deal.II/base/index_set.h>
30 # include <deal.II/lac/vector_type_traits.h>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 
38 namespace PETScWrappers
39 {
47  namespace MPI
48  {
49 
158  class Vector : public VectorBase
159  {
160  public:
165 
169  Vector ();
170 
187  explicit Vector (const MPI_Comm &communicator,
188  const size_type n,
189  const size_type local_size);
190 
191 
202  template <typename Number>
203  explicit Vector (const MPI_Comm &communicator,
204  const ::Vector<Number> &v,
205  const size_type local_size);
206 
207 
220  DEAL_II_DEPRECATED
221  explicit Vector (const MPI_Comm &communicator,
222  const VectorBase &v,
223  const size_type local_size);
224 
248  Vector (const IndexSet &local,
249  const IndexSet &ghost,
250  const MPI_Comm &communicator);
251 
263  explicit Vector (const IndexSet &local,
264  const MPI_Comm &communicator);
265 
270  void clear ();
271 
276  Vector &operator= (const Vector &v);
277 
284  Vector &operator= (const PetscScalar s);
285 
295  template <typename number>
296  Vector &operator= (const ::Vector<number> &v);
297 
314  void reinit (const MPI_Comm &communicator,
315  const size_type N,
316  const size_type local_size,
317  const bool omit_zeroing_entries = false);
318 
328  void reinit (const Vector &v,
329  const bool omit_zeroing_entries = false);
330 
338  void reinit (const IndexSet &local,
339  const IndexSet &ghost,
340  const MPI_Comm &communicator);
341 
349  void reinit (const IndexSet &local,
350  const MPI_Comm &communicator);
351 
356  const MPI_Comm &get_mpi_communicator () const;
357 
369  void print (std::ostream &out,
370  const unsigned int precision = 3,
371  const bool scientific = true,
372  const bool across = true) const;
373 
380  bool all_zero () const;
381 
382  protected:
389  virtual void create_vector (const size_type n,
390  const size_type local_size);
391 
392 
393 
399  virtual void create_vector (const size_type n,
400  const size_type local_size,
401  const IndexSet &ghostnodes);
402 
403 
404  private:
408  MPI_Comm communicator;
409  };
410 
411 
412 // ------------------ template and inline functions -------------
413 
414 
423  inline
424  void swap (Vector &u, Vector &v)
425  {
426  u.swap (v);
427  }
428 
429 
430 #ifndef DOXYGEN
431 
432  template <typename number>
433  Vector::Vector (const MPI_Comm &communicator,
434  const ::Vector<number> &v,
435  const size_type local_size)
436  :
437  communicator (communicator)
438  {
439  Vector::create_vector (v.size(), local_size);
440 
441  *this = v;
442  }
443 
444 
445 
446  inline
447  Vector &
448  Vector::operator= (const PetscScalar s)
449  {
451 
452  return *this;
453  }
454 
455 
456 
457  template <typename number>
458  inline
459  Vector &
460  Vector::operator= (const ::Vector<number> &v)
461  {
462  Assert (size() == v.size(),
463  ExcDimensionMismatch (size(), v.size()));
464 
465  // FIXME: the following isn't necessarily fast, but this is due to
466  // the fact that PETSc doesn't offer an inlined access operator.
467  //
468  // if someone wants to contribute some code: to make this code
469  // faster, one could either first convert all values to PetscScalar,
470  // and then set them all at once using VecSetValues. This has the
471  // drawback that it could take quite some memory, if the vector is
472  // large, and it would in addition allocate memory on the heap, which
473  // is expensive. an alternative would be to split the vector into
474  // chunks of, say, 128 elements, convert a chunk at a time and set it
475  // in the output vector using VecSetValues. since 128 elements is
476  // small enough, this could easily be allocated on the stack (as a
477  // local variable) which would make the whole thing much more
478  // efficient.
479  //
480  // a second way to make things faster is for the special case that
481  // number==PetscScalar. we could then declare a specialization of
482  // this template, and omit the conversion. the problem with this is
483  // that the best we can do is to use VecSetValues, but this isn't
484  // very efficient either: it wants to see an array of indices, which
485  // in this case a) again takes up a whole lot of memory on the heap,
486  // and b) is totally dumb since its content would simply be the
487  // sequence 0,1,2,3,...,n. the best of all worlds would probably be a
488  // function in Petsc that would take a pointer to an array of
489  // PetscScalar values and simply copy n elements verbatim into the
490  // vector...
491  for (size_type i=0; i<v.size(); ++i)
492  (*this)(i) = v(i);
493 
494  compress (::VectorOperation::insert);
495 
496  return *this;
497  }
498 
499 
500 
501  inline
502  const MPI_Comm &
503  Vector::get_mpi_communicator () const
504  {
505  return communicator;
506  }
507 
508 #endif // DOXYGEN
509  }
510 }
511 
512 namespace internal
513 {
514  namespace LinearOperatorImplementation
515  {
516  template <typename> class ReinitHelper;
517 
522  template <>
523  class ReinitHelper<PETScWrappers::MPI::Vector>
524  {
525  public:
526  template <typename Matrix>
527  static
528  void reinit_range_vector (const Matrix &matrix,
530  bool /*omit_zeroing_entries*/)
531  {
532  v.reinit(matrix.locally_owned_range_indices(), matrix.get_mpi_communicator());
533  }
534 
535  template <typename Matrix>
536  static
537  void reinit_domain_vector(const Matrix &matrix,
539  bool /*omit_zeroing_entries*/)
540  {
541  v.reinit(matrix.locally_owned_domain_indices(), matrix.get_mpi_communicator());
542  }
543  };
544 
545  } /* namespace LinearOperator */
546 } /* namespace internal */
547 
556 template <>
557 struct is_serial_vector< PETScWrappers::MPI::Vector > : std::false_type
558 {
559 };
560 
561 
562 DEAL_II_NAMESPACE_CLOSE
563 
564 #endif // DEAL_II_WITH_PETSC
565 
566 /*---------------------------- petsc_parallel_vector.h ---------------------------*/
567 
568 #endif
569 /*---------------------------- petsc_parallel_vector.h ---------------------------*/
types::global_dof_index size_type
VectorBase & operator=(const PetscScalar s)
Vector< Number > & operator=(const Number s)
static void reinit_range_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
static void reinit_domain_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
void swap(BlockVector &u, BlockVector &v)
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const MPI_Comm & get_mpi_communicator() const
Vector & operator=(const Vector &v)
void reinit(const MPI_Comm &communicator, const size_type N, const size_type local_size, const bool omit_zeroing_entries=false)
virtual void create_vector(const size_type n, const size_type local_size)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
void swap(Vector &u, Vector &v)