Reference documentation for deal.II version 8.5.1
petsc_parallel_vector.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2017 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/petsc_vector_base.h>
28 # include <deal.II/base/index_set.h>
29 # include <deal.II/lac/vector_type_traits.h>
30 
31 DEAL_II_NAMESPACE_OPEN
32 
33 
37 namespace PETScWrappers
38 {
46  namespace MPI
47  {
48 
157  class Vector : public VectorBase
158  {
159  public:
164 
178  static const bool supports_distributed_data DEAL_II_DEPRECATED = true;
179 
183  Vector ();
184 
201  explicit Vector (const MPI_Comm &communicator,
202  const size_type n,
203  const size_type local_size);
204 
205 
216  template <typename Number>
217  explicit Vector (const MPI_Comm &communicator,
218  const ::Vector<Number> &v,
219  const size_type local_size);
220 
221 
234  explicit Vector (const MPI_Comm &communicator,
235  const VectorBase &v,
236  const size_type local_size) DEAL_II_DEPRECATED;
237 
255  Vector (const IndexSet &local,
256  const IndexSet &ghost,
257  const MPI_Comm &communicator);
258 
263  explicit Vector (const IndexSet &local,
264  const MPI_Comm &communicator);
265 
270  void clear ();
271 
276  Vector &operator= (const Vector &v);
277 
294 
301  Vector &operator= (const PetscScalar s);
302 
312  template <typename number>
313  Vector &operator= (const ::Vector<number> &v);
314 
331  void reinit (const MPI_Comm &communicator,
332  const size_type N,
333  const size_type local_size,
334  const bool omit_zeroing_entries = false);
335 
345  void reinit (const Vector &v,
346  const bool omit_zeroing_entries = false);
347 
355  void reinit (const IndexSet &local,
356  const IndexSet &ghost,
357  const MPI_Comm &communicator);
358 
366  void reinit (const IndexSet &local,
367  const MPI_Comm &communicator);
368 
373  const MPI_Comm &get_mpi_communicator () const;
374 
386  void print (std::ostream &out,
387  const unsigned int precision = 3,
388  const bool scientific = true,
389  const bool across = true) const;
390 
397  bool all_zero () const;
398 
399  protected:
406  virtual void create_vector (const size_type n,
407  const size_type local_size);
408 
409 
410 
416  virtual void create_vector (const size_type n,
417  const size_type local_size,
418  const IndexSet &ghostnodes);
419 
420 
421  private:
425  MPI_Comm communicator;
426  };
427 
428 
429 // ------------------ template and inline functions -------------
430 
431 
440  inline
441  void swap (Vector &u, Vector &v)
442  {
443  u.swap (v);
444  }
445 
446 
447 #ifndef DOXYGEN
448 
449  template <typename number>
450  Vector::Vector (const MPI_Comm &communicator,
451  const ::Vector<number> &v,
452  const size_type local_size)
453  :
454  communicator (communicator)
455  {
456  Vector::create_vector (v.size(), local_size);
457 
458  *this = v;
459  }
460 
461 
462 
463  inline
464  Vector &
465  Vector::operator= (const PetscScalar s)
466  {
468 
469  return *this;
470  }
471 
472 
473 
474  inline
475  Vector &
476  Vector::operator= (const Vector &v)
477  {
478  // make sure left- and right-hand side of the assignment are compress()'ed:
479  Assert(v.last_action == VectorOperation::unknown,
480  internal::VectorReference::ExcWrongMode (VectorOperation::unknown,
481  v.last_action));
483  internal::VectorReference::ExcWrongMode (VectorOperation::unknown,
484  last_action));
485 
486 
487  if (v.size()==0)
488  {
489  // this happens if v has not been initialized to something useful:
490  // Vector x,v;
491  // x=v;
492  // we skip the code below and create a simple serial vector of
493  // length 0
494 
495 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
496  PetscErrorCode ierr = VecDestroy (vector);
497 #else
498  PetscErrorCode ierr = VecDestroy (&vector);
499 #endif
500  AssertThrow (ierr == 0, ExcPETScError(ierr));
501 
502  const int n = 0;
503  ierr = VecCreateSeq (PETSC_COMM_SELF, n, &vector);
504  AssertThrow (ierr == 0, ExcPETScError(ierr));
505  ghosted = false;
507  return *this;
508  }
509 
510  // if the vectors have different sizes,
511  // then first resize the present one
512  if (size() != v.size())
513  {
514  if (v.has_ghost_elements())
515  reinit( v.locally_owned_elements(), v.ghost_indices, v.communicator);
516  else
517  reinit (v.communicator, v.size(), v.local_size(), true);
518  }
519 
520  PetscErrorCode ierr = VecCopy (v.vector, vector);
521  AssertThrow (ierr == 0, ExcPETScError(ierr));
522 
523  if (has_ghost_elements())
524  {
525  ierr = VecGhostUpdateBegin(vector, INSERT_VALUES, SCATTER_FORWARD);
526  AssertThrow (ierr == 0, ExcPETScError(ierr));
527  ierr = VecGhostUpdateEnd(vector, INSERT_VALUES, SCATTER_FORWARD);
528  AssertThrow (ierr == 0, ExcPETScError(ierr));
529  }
530  return *this;
531  }
532 
533 
534 
535  template <typename number>
536  inline
537  Vector &
538  Vector::operator= (const ::Vector<number> &v)
539  {
540  Assert (size() == v.size(),
541  ExcDimensionMismatch (size(), v.size()));
542 
543  // FIXME: the following isn't necessarily fast, but this is due to
544  // the fact that PETSc doesn't offer an inlined access operator.
545  //
546  // if someone wants to contribute some code: to make this code
547  // faster, one could either first convert all values to PetscScalar,
548  // and then set them all at once using VecSetValues. This has the
549  // drawback that it could take quite some memory, if the vector is
550  // large, and it would in addition allocate memory on the heap, which
551  // is expensive. an alternative would be to split the vector into
552  // chunks of, say, 128 elements, convert a chunk at a time and set it
553  // in the output vector using VecSetValues. since 128 elements is
554  // small enough, this could easily be allocated on the stack (as a
555  // local variable) which would make the whole thing much more
556  // efficient.
557  //
558  // a second way to make things faster is for the special case that
559  // number==PetscScalar. we could then declare a specialization of
560  // this template, and omit the conversion. the problem with this is
561  // that the best we can do is to use VecSetValues, but this isn't
562  // very efficient either: it wants to see an array of indices, which
563  // in this case a) again takes up a whole lot of memory on the heap,
564  // and b) is totally dumb since its content would simply be the
565  // sequence 0,1,2,3,...,n. the best of all worlds would probably be a
566  // function in Petsc that would take a pointer to an array of
567  // PetscScalar values and simply copy n elements verbatim into the
568  // vector...
569  for (size_type i=0; i<v.size(); ++i)
570  (*this)(i) = v(i);
571 
573 
574  return *this;
575  }
576 
577 
578 
579  inline
580  const MPI_Comm &
582  {
583  return communicator;
584  }
585 
586 #endif // DOXYGEN
587  }
588 }
589 
590 namespace internal
591 {
592  namespace LinearOperator
593  {
594  template <typename> class ReinitHelper;
595 
600  template<>
601  class ReinitHelper<PETScWrappers::MPI::Vector>
602  {
603  public:
604  template <typename Matrix>
605  static
606  void reinit_range_vector (const Matrix &matrix,
608  bool /*omit_zeroing_entries*/)
609  {
610  v.reinit(matrix.locally_owned_range_indices(), matrix.get_mpi_communicator());
611  }
612 
613  template <typename Matrix>
614  static
615  void reinit_domain_vector(const Matrix &matrix,
617  bool /*omit_zeroing_entries*/)
618  {
619  v.reinit(matrix.locally_owned_domain_indices(), matrix.get_mpi_communicator());
620  }
621  };
622 
623  } /* namespace LinearOperator */
624 } /* namespace internal */
625 
634 template <>
635 struct is_serial_vector< PETScWrappers::MPI::Vector > : std_cxx11::false_type
636 {
637 };
638 
639 
640 DEAL_II_NAMESPACE_CLOSE
641 
642 #endif // DEAL_II_WITH_PETSC
643 
644 /*---------------------------- petsc_parallel_vector.h ---------------------------*/
645 
646 #endif
647 /*---------------------------- petsc_parallel_vector.h ---------------------------*/
types::global_dof_index size_type
Vector & operator=(const Vector &v)
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
VectorBase & operator=(const PetscScalar s)
static void reinit_domain_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
void clear()
Definition: index_set.h:1393
void compress(const VectorOperation::values operation)
unsigned int global_dof_index
Definition: types.h:88
VectorOperation::values last_action
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::size_t size() const
static const bool supports_distributed_data
const MPI_Comm & get_mpi_communicator() const
static void reinit_range_vector(const Matrix &matrix, Vector &v, bool omit_zeroing_entries)
IndexSet locally_owned_elements() const
bool has_ghost_elements() const
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)