Reference documentation for deal.II version 8.5.1
petsc_parallel_vector.cc
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 #include <deal.II/base/mpi.h>
17 #include <deal.II/lac/petsc_parallel_vector.h>
18 
19 #ifdef DEAL_II_WITH_PETSC
20 
21 # include <deal.II/lac/petsc_vector.h>
22 # include <cmath>
23 # include <algorithm>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 namespace PETScWrappers
28 {
29  namespace MPI
30  {
31 
33  : communicator (MPI_COMM_SELF)
34  {
35  // this is an invalid empty vector, so we can just as well create a
36  // sequential one to avoid all the overhead incurred by parallelism
37  const int n = 0;
38  const PetscErrorCode ierr = VecCreateSeq (PETSC_COMM_SELF, n, &vector);
39  AssertThrow (ierr == 0, ExcPETScError(ierr));
40  ghosted = false;
41  }
42 
43 
44 
45  Vector::Vector (const MPI_Comm &communicator,
46  const size_type n,
47  const size_type local_size)
48  :
49  communicator (communicator)
50  {
52  }
53 
54 
55 
56  Vector::Vector (const MPI_Comm &communicator,
57  const VectorBase &v,
58  const size_type local_size)
59  :
60  VectorBase (v),
61  communicator (communicator)
62  {
63  // In the past (before it was deprecated) this constructor did a
64  // byte-for-byte copy of v. This choice resulted in two problems:
65  // 1. The created vector will have the same size as v, not local size.
66  // 2. Since both the created vector and v maintain ownership of the same
67  // PETSc Vec, both will try to destroy it: this does not make sense.
68  //
69  // For the sake of backwards compatibility, preserve the behavior of the
70  // copy, but correct the ownership bug. Note that in both this (and the
71  // original) implementation local_size is ultimately unused.
72  (void)local_size;
73  }
74 
75 
76 
77  Vector::Vector (const IndexSet &local,
78  const IndexSet &ghost,
79  const MPI_Comm &communicator)
80  :
81  communicator (communicator)
82  {
84 
85  IndexSet ghost_set = ghost;
86  ghost_set.subtract_set(local);
87 
88  Vector::create_vector(local.size(), local.n_elements(), ghost_set);
89  }
90 
91 
92 
93  Vector::Vector (const IndexSet &local,
94  const MPI_Comm &communicator)
95  :
96  communicator (communicator)
97  {
99  Vector::create_vector(local.size(), local.n_elements());
100  }
101 
102 
103 
104  void
106  {
107  // destroy the PETSc Vec and create an invalid empty vector,
108  // so we can just as well create a sequential one to avoid
109  // all the overhead incurred by parallelism
110  attained_ownership = true;
112 
113  const int n = 0;
114  const PetscErrorCode ierr = VecCreateSeq (PETSC_COMM_SELF, n, &vector);
115  AssertThrow (ierr == 0, ExcPETScError(ierr));
116  }
117 
118 
119 
120  void
121  Vector::reinit (const MPI_Comm &comm,
122  const size_type n,
123  const size_type local_sz,
124  const bool omit_zeroing_entries)
125  {
126  communicator = comm;
127 
128  // only do something if the sizes
129  // mismatch (may not be true for every proc)
130 
131  int k_global, k = ((size() != n) || (local_size() != local_sz));
132  {
133  const int ierr = MPI_Allreduce (&k, &k_global, 1,
134  MPI_INT, MPI_LOR, communicator);
135  AssertThrowMPI(ierr);
136  }
137 
138  if (k_global || has_ghost_elements())
139  {
140  // FIXME: I'd like to use this here,
141  // but somehow it leads to odd errors
142  // somewhere down the line in some of
143  // the tests:
144 // const PetscErrorCode ierr = VecSetSizes (vector, n, n);
145 // AssertThrow (ierr == 0, ExcPETScError(ierr));
146 
147  // so let's go the slow way:
148 
149 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
150  const PetscErrorCode ierr = VecDestroy (vector);
151 #else
152  const PetscErrorCode ierr = VecDestroy (&vector);
153 #endif
154  AssertThrow (ierr == 0, ExcPETScError(ierr));
155 
156  create_vector (n, local_sz);
157  }
158 
159  // finally clear the new vector if so
160  // desired
161  if (omit_zeroing_entries == false)
162  *this = 0;
163  }
164 
165 
166 
167  void
169  const bool omit_zeroing_entries)
170  {
171  if (v.has_ghost_elements())
172  {
174  if (!omit_zeroing_entries)
175  {
176  const PetscErrorCode ierr = VecSet(vector, 0.0);
177  AssertThrow (ierr == 0, ExcPETScError(ierr));
178  }
179  }
180  else
181  reinit (v.communicator, v.size(), v.local_size(), omit_zeroing_entries);
182  }
183 
184 
185 
186  void
187  Vector::reinit (const IndexSet &local,
188  const IndexSet &ghost,
189  const MPI_Comm &comm)
190  {
191 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
192  const PetscErrorCode ierr = VecDestroy (vector);
193 #else
194  const PetscErrorCode ierr = VecDestroy (&vector);
195 #endif
196  AssertThrow (ierr == 0, ExcPETScError(ierr));
197 
198  communicator = comm;
199 
201 
202  IndexSet ghost_set = ghost;
203  ghost_set.subtract_set(local);
204 
205  create_vector(local.size(), local.n_elements(), ghost_set);
206  }
207 
208  void
209  Vector::reinit (const IndexSet &local,
210  const MPI_Comm &comm)
211  {
212 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
213  const PetscErrorCode ierr = VecDestroy (vector);
214 #else
215  const PetscErrorCode ierr = VecDestroy (&vector);
216 #endif
217  AssertThrow (ierr == 0, ExcPETScError(ierr));
218 
219  communicator = comm;
220 
222  Assert(local.size()>0, ExcMessage("can not create vector of size 0."));
223  create_vector(local.size(), local.n_elements());
224  }
225 
226 
227  Vector &
229  {
231  ExcMessage("Call to compress() required before calling operator=."));
232  //TODO [TH]: can not access v.last_action here. Implement is_compressed()?
233  //Assert(v.last_action==VectorOperation::unknown,
234  // ExcMessage("Call to compress() required before calling operator=."));
235 
236  // get a pointer to the local memory of
237  // this vector
238  PetscScalar *dest_array;
239  PetscErrorCode ierr = VecGetArray (vector, &dest_array);
240  AssertThrow (ierr == 0, ExcPETScError(ierr));
241 
242  // then also a pointer to the source
243  // vector
244  PetscScalar *src_array;
245  ierr = VecGetArray (static_cast<const Vec &>(v), &src_array);
246  AssertThrow (ierr == 0, ExcPETScError(ierr));
247 
248  // then copy:
249  const std::pair<size_type, size_type>
250  local_elements = local_range ();
251  std::copy (src_array + local_elements.first,
252  src_array + local_elements.second,
253  dest_array);
254 
255  // finally restore the arrays
256  ierr = VecRestoreArray (vector, &dest_array);
257  AssertThrow (ierr == 0, ExcPETScError(ierr));
258 
259  ierr = VecRestoreArray (static_cast<const Vec &>(v), &src_array);
260  AssertThrow (ierr == 0, ExcPETScError(ierr));
261 
262  if (has_ghost_elements())
263  {
264  ierr = VecGhostUpdateBegin(vector, INSERT_VALUES, SCATTER_FORWARD);
265  AssertThrow (ierr == 0, ExcPETScError(ierr));
266  ierr = VecGhostUpdateEnd(vector, INSERT_VALUES, SCATTER_FORWARD);
267  AssertThrow (ierr == 0, ExcPETScError(ierr));
268  }
269  return *this;
270  }
271 
272 
273  void
275  const size_type local_size)
276  {
277  (void)n;
278  Assert (local_size <= n, ExcIndexRange (local_size, 0, n));
279  ghosted = false;
280 
281  const PetscErrorCode ierr = VecCreateMPI (communicator, local_size, PETSC_DETERMINE,
282  &vector);
283  AssertThrow (ierr == 0, ExcPETScError(ierr));
284 
285  Assert (size() == n,
286  ExcDimensionMismatch (size(), n));
287  }
288 
289 
290 
291  void
293  const size_type local_size,
294  const IndexSet &ghostnodes)
295  {
296  (void)n;
297  Assert (local_size <= n, ExcIndexRange (local_size, 0, n));
298  ghosted = true;
299  ghost_indices = ghostnodes;
300 
301  std::vector<size_type> ghostindices;
302  ghostnodes.fill_index_vector(ghostindices);
303 
304  const PetscInt *ptr
305  = (ghostindices.size() > 0
306  ?
307  (const PetscInt *)(&(ghostindices[0]))
308  :
309  0);
310 
311  PetscErrorCode ierr = VecCreateGhost(communicator,
312  local_size,
313  PETSC_DETERMINE,
314  ghostindices.size(),
315  ptr,
316  &vector);
317  AssertThrow (ierr == 0, ExcPETScError(ierr));
318 
319  Assert (size() == n,
320  ExcDimensionMismatch (size(), n));
321 
322 #if DEBUG
323  {
324  // test ghost allocation in debug mode
325  PetscInt begin, end;
326 
327  ierr = VecGetOwnershipRange (vector, &begin, &end);
328  AssertThrow (ierr == 0, ExcPETScError(ierr));
329 
330  Assert(local_size==(size_type)(end-begin), ExcInternalError());
331 
332  Vec l;
333  ierr = VecGhostGetLocalForm(vector, &l);
334  AssertThrow (ierr == 0, ExcPETScError(ierr));
335 
336  PetscInt lsize;
337  ierr = VecGetSize(l, &lsize);
338  AssertThrow (ierr == 0, ExcPETScError(ierr));
339 
340  ierr = VecGhostRestoreLocalForm(vector, &l);
341  AssertThrow (ierr == 0, ExcPETScError(ierr));
342 
343  Assert (lsize==end-begin+(PetscInt)ghost_indices.n_elements(),
344  ExcInternalError());
345  }
346 #endif
347 
348 
349  // in PETSc versions up to 3.5, VecCreateGhost zeroed out the locally
350  // owned vector elements but forgot about the ghost elements. we need to
351  // do this ourselves
352  //
353  // see https://code.google.com/p/dealii/issues/detail?id=233
354 #if DEAL_II_PETSC_VERSION_LT(3,6,0)
356  zero.reinit (communicator, this->size(), local_size);
357  *this = zero;
358 #endif
359 
360  }
361 
362 
363 
364  bool
366  {
367  unsigned int has_nonzero = VectorBase::all_zero()?0:1;
368 #ifdef DEAL_II_WITH_MPI
369  // in parallel, check that the vector
370  // is zero on _all_ processors.
371  unsigned int num_nonzero = Utilities::MPI::sum(has_nonzero, communicator);
372  return num_nonzero == 0;
373 #else
374  return has_nonzero == 0;
375 #endif
376  }
377 
378 
379  void
380  Vector::print (std::ostream &out,
381  const unsigned int precision,
382  const bool scientific,
383  const bool across) const
384  {
385  AssertThrow (out, ExcIO());
386 
387  // get a representation of the vector and
388  // loop over all the elements
389  PetscScalar *val;
390  PetscInt nlocal, istart, iend;
391 
392  PetscErrorCode ierr = VecGetArray (vector, &val);
393  AssertThrow (ierr == 0, ExcPETScError(ierr));
394 
395  ierr = VecGetLocalSize (vector, &nlocal);
396  AssertThrow (ierr == 0, ExcPETScError(ierr));
397 
398  ierr = VecGetOwnershipRange (vector, &istart, &iend);
399  AssertThrow (ierr == 0, ExcPETScError(ierr));
400 
401  // save the state of out stream
402  std::ios::fmtflags old_flags = out.flags();
403  unsigned int old_precision = out.precision (precision);
404 
405  out.precision (precision);
406  if (scientific)
407  out.setf (std::ios::scientific, std::ios::floatfield);
408  else
409  out.setf (std::ios::fixed, std::ios::floatfield);
410 
411  for ( unsigned int i = 0;
413  i++)
414  {
415  // This is slow, but most likely only used to debug.
416  const int mpi_ierr = MPI_Barrier(communicator);
417  AssertThrowMPI(mpi_ierr);
419  {
420  if (across)
421  {
422  out << "[Proc" << i << " " << istart << "-" << iend-1 << "]" << ' ';
423  for (PetscInt i=0; i<nlocal; ++i)
424  out << val[i] << ' ';
425  }
426  else
427  {
428  out << "[Proc " << i << " " << istart << "-" << iend-1 << "]" << std::endl;
429  for (PetscInt i=0; i<nlocal; ++i)
430  out << val[i] << std::endl;
431  }
432  out << std::endl;
433  }
434  }
435  // reset output format
436  out.flags (old_flags);
437  out.precision(old_precision);
438 
439  // restore the representation of the
440  // vector
441  ierr = VecRestoreArray (vector, &val);
442  AssertThrow (ierr == 0, ExcPETScError(ierr));
443 
444  AssertThrow (out, ExcIO());
445  }
446 
447  }
448 
449 }
450 
451 DEAL_II_NAMESPACE_CLOSE
452 
453 #endif // DEAL_II_WITH_PETSC
types::global_dof_index size_type
Vector & operator=(const Vector &v)
static ::ExceptionBase & ExcIO()
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
std::pair< size_type, size_type > local_range() const
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
Definition: index_set.cc:604
size_type size() const
Definition: index_set.h:1419
static ::ExceptionBase & ExcMessage(std::string arg1)
T sum(const T &t, const MPI_Comm &mpi_communicator)
void subtract_set(const IndexSet &other)
Definition: index_set.cc:278
VectorOperation::values last_action
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:63
void fill_index_vector(std::vector< size_type > &indices) const
Definition: index_set.cc:512
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1211
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)
IndexSet locally_owned_elements() const
virtual void create_vector(const size_type n, const size_type local_size)
unsigned int this_mpi_process(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:73
static ::ExceptionBase & ExcNotImplemented()
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
size_type n_elements() const
Definition: index_set.h:1560
static ::ExceptionBase & ExcInternalError()