Reference documentation for deal.II version 9.0.0
petsc_parallel_vector.cc
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 #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 <cmath>
22 # include <algorithm>
23 
24 DEAL_II_NAMESPACE_OPEN
25 
26 namespace PETScWrappers
27 {
28  namespace MPI
29  {
30 
32  : communicator (MPI_COMM_SELF)
33  {
34  // virtual functions called in constructors and destructors never use the
35  // override in a derived class
36  // for clarity be explicit on which function is called
38  }
39 
40 
41 
42  Vector::Vector (const MPI_Comm &communicator,
43  const size_type n,
44  const size_type local_size)
45  :
46  communicator (communicator)
47  {
49  }
50 
51 
52 
53  Vector::Vector (const MPI_Comm &communicator,
54  const VectorBase &v,
55  const size_type local_size)
56  :
57  VectorBase (v),
58  communicator (communicator)
59  {
60  // In the past (before it was deprecated) this constructor did a
61  // byte-for-byte copy of v. This choice resulted in two problems:
62  // 1. The created vector will have the same size as v, not local size.
63  // 2. Since both the created vector and v maintain ownership of the same
64  // PETSc Vec, both will try to destroy it: this does not make sense.
65  //
66  // For the sake of backwards compatibility, preserve the behavior of the
67  // copy, but correct the ownership bug. Note that in both this (and the
68  // original) implementation local_size is ultimately unused.
69  (void)local_size;
70  }
71 
72 
73 
74  Vector::Vector (const IndexSet &local,
75  const IndexSet &ghost,
76  const MPI_Comm &communicator)
77  :
78  communicator (communicator)
79  {
81 
82  IndexSet ghost_set = ghost;
83  ghost_set.subtract_set(local);
84 
85  Vector::create_vector(local.size(), local.n_elements(), ghost_set);
86  }
87 
88 
89 
90  Vector::Vector (const IndexSet &local,
91  const MPI_Comm &communicator)
92  :
93  communicator (communicator)
94  {
96  Vector::create_vector(local.size(), local.n_elements());
97  }
98 
99 
100 
101  Vector &
103  {
104  // make sure left- and right-hand side of the assignment are compress()'ed:
106  internal::VectorReference::ExcWrongMode (VectorOperation::unknown,
107  v.last_action));
109  internal::VectorReference::ExcWrongMode (VectorOperation::unknown,
110  last_action));
111 
112  // if the vectors have different sizes,
113  // then first resize the present one
114  if (size() != v.size())
115  {
116  if (v.has_ghost_elements())
118  else
119  reinit (v.communicator, v.size(), v.local_size(), true);
120  }
121 
122  PetscErrorCode ierr = VecCopy (v.vector, vector);
123  AssertThrow (ierr == 0, ExcPETScError(ierr));
124 
125  if (has_ghost_elements())
126  {
127  ierr = VecGhostUpdateBegin(vector, INSERT_VALUES, SCATTER_FORWARD);
128  AssertThrow (ierr == 0, ExcPETScError(ierr));
129  ierr = VecGhostUpdateEnd(vector, INSERT_VALUES, SCATTER_FORWARD);
130  AssertThrow (ierr == 0, ExcPETScError(ierr));
131  }
132  return *this;
133  }
134 
135 
136 
137  void
139  {
140  obtained_ownership = true;
142 
143  create_vector(0, 0);
144  }
145 
146 
147 
148  void
149  Vector::reinit (const MPI_Comm &comm,
150  const size_type n,
151  const size_type local_sz,
152  const bool omit_zeroing_entries)
153  {
154  communicator = comm;
155 
156  // only do something if the sizes
157  // mismatch (may not be true for every proc)
158 
159  int k_global, k = ((size() != n) || (local_size() != local_sz));
160  {
161  const int ierr = MPI_Allreduce (&k, &k_global, 1,
162  MPI_INT, MPI_LOR, communicator);
163  AssertThrowMPI(ierr);
164  }
165 
166  if (k_global || has_ghost_elements())
167  {
168  // FIXME: I'd like to use this here,
169  // but somehow it leads to odd errors
170  // somewhere down the line in some of
171  // the tests:
172 // const PetscErrorCode ierr = VecSetSizes (vector, n, n);
173 // AssertThrow (ierr == 0, ExcPETScError(ierr));
174 
175  // so let's go the slow way:
176 
177  const PetscErrorCode ierr = VecDestroy (&vector);
178  AssertThrow (ierr == 0, ExcPETScError(ierr));
179 
180  create_vector (n, local_sz);
181  }
182 
183  // finally clear the new vector if so
184  // desired
185  if (omit_zeroing_entries == false)
186  *this = 0;
187  }
188 
189 
190 
191  void
193  const bool omit_zeroing_entries)
194  {
195  if (v.has_ghost_elements())
196  {
198  if (!omit_zeroing_entries)
199  {
200  const PetscErrorCode ierr = VecSet(vector, 0.0);
201  AssertThrow (ierr == 0, ExcPETScError(ierr));
202  }
203  }
204  else
205  reinit (v.communicator, v.size(), v.local_size(), omit_zeroing_entries);
206  }
207 
208 
209 
210  void
211  Vector::reinit (const IndexSet &local,
212  const IndexSet &ghost,
213  const MPI_Comm &comm)
214  {
215  const PetscErrorCode ierr = VecDestroy (&vector);
216  AssertThrow (ierr == 0, ExcPETScError(ierr));
217 
218  communicator = comm;
219 
221 
222  IndexSet ghost_set = ghost;
223  ghost_set.subtract_set(local);
224 
225  create_vector(local.size(), local.n_elements(), ghost_set);
226  }
227 
228  void
229  Vector::reinit (const IndexSet &local,
230  const MPI_Comm &comm)
231  {
232  const PetscErrorCode ierr = VecDestroy (&vector);
233  AssertThrow (ierr == 0, ExcPETScError(ierr));
234 
235  communicator = comm;
236 
238  Assert(local.size()>0, ExcMessage("can not create vector of size 0."));
239  create_vector(local.size(), local.n_elements());
240  }
241 
242 
243  void
245  const size_type local_size)
246  {
247  (void)n;
248  Assert (local_size <= n, ExcIndexRange (local_size, 0, n));
249  ghosted = false;
250 
251  const PetscErrorCode ierr = VecCreateMPI (communicator, local_size, PETSC_DETERMINE,
252  &vector);
253  AssertThrow (ierr == 0, ExcPETScError(ierr));
254 
255  Assert (size() == n,
256  ExcDimensionMismatch (size(), n));
257  }
258 
259 
260 
261  void
263  const size_type local_size,
264  const IndexSet &ghostnodes)
265  {
266  (void)n;
267  Assert (local_size <= n, ExcIndexRange (local_size, 0, n));
268  ghosted = true;
269  ghost_indices = ghostnodes;
270 
271  std::vector<size_type> ghostindices;
272  ghostnodes.fill_index_vector(ghostindices);
273 
274  const PetscInt *ptr
275  = (ghostindices.size() > 0
276  ?
277  (const PetscInt *)(&(ghostindices[0]))
278  :
279  nullptr);
280 
281  PetscErrorCode ierr = VecCreateGhost(communicator,
282  local_size,
283  PETSC_DETERMINE,
284  ghostindices.size(),
285  ptr,
286  &vector);
287  AssertThrow (ierr == 0, ExcPETScError(ierr));
288 
289  Assert (size() == n,
290  ExcDimensionMismatch (size(), n));
291 
292 #if DEBUG
293  {
294  // test ghost allocation in debug mode
295  PetscInt begin, end;
296 
297  ierr = VecGetOwnershipRange (vector, &begin, &end);
298  AssertThrow (ierr == 0, ExcPETScError(ierr));
299 
300  Assert(local_size==(size_type)(end-begin), ExcInternalError());
301 
302  Vec l;
303  ierr = VecGhostGetLocalForm(vector, &l);
304  AssertThrow (ierr == 0, ExcPETScError(ierr));
305 
306  PetscInt lsize;
307  ierr = VecGetSize(l, &lsize);
308  AssertThrow (ierr == 0, ExcPETScError(ierr));
309 
310  ierr = VecGhostRestoreLocalForm(vector, &l);
311  AssertThrow (ierr == 0, ExcPETScError(ierr));
312 
313  Assert (lsize==end-begin+(PetscInt)ghost_indices.n_elements(),
314  ExcInternalError());
315  }
316 #endif
317 
318 
319  // in PETSc versions up to 3.5, VecCreateGhost zeroed out the locally
320  // owned vector elements but forgot about the ghost elements. we need to
321  // do this ourselves
322  //
323  // see https://code.google.com/p/dealii/issues/detail?id=233
324 #if DEAL_II_PETSC_VERSION_LT(3,6,0)
326  zero.reinit (communicator, this->size(), local_size);
327  *this = zero;
328 #endif
329 
330  }
331 
332 
333 
334  bool
336  {
337  unsigned int has_nonzero = VectorBase::all_zero()?0:1;
338 #ifdef DEAL_II_WITH_MPI
339  // in parallel, check that the vector
340  // is zero on _all_ processors.
341  unsigned int num_nonzero = Utilities::MPI::sum(has_nonzero, communicator);
342  return num_nonzero == 0;
343 #else
344  return has_nonzero == 0;
345 #endif
346  }
347 
348 
349  void
350  Vector::print (std::ostream &out,
351  const unsigned int precision,
352  const bool scientific,
353  const bool across) const
354  {
355  AssertThrow (out, ExcIO());
356 
357  // get a representation of the vector and
358  // loop over all the elements
359  PetscScalar *val;
360  PetscInt nlocal, istart, iend;
361 
362  PetscErrorCode ierr = VecGetArray (vector, &val);
363  AssertThrow (ierr == 0, ExcPETScError(ierr));
364 
365  ierr = VecGetLocalSize (vector, &nlocal);
366  AssertThrow (ierr == 0, ExcPETScError(ierr));
367 
368  ierr = VecGetOwnershipRange (vector, &istart, &iend);
369  AssertThrow (ierr == 0, ExcPETScError(ierr));
370 
371  // save the state of out stream
372  std::ios::fmtflags old_flags = out.flags();
373  unsigned int old_precision = out.precision (precision);
374 
375  out.precision (precision);
376  if (scientific)
377  out.setf (std::ios::scientific, std::ios::floatfield);
378  else
379  out.setf (std::ios::fixed, std::ios::floatfield);
380 
381  // let each processor produce its output in turn. this requires
382  // synchronizing output between processors using a barrier --
383  // which is clearly slow, but nobody is going to print a whole
384  // matrix this way on a regular basis for production runs, so
385  // the slowness of the barrier doesn't matter
386  for ( unsigned int i = 0;
388  i++)
389  {
390  const int mpi_ierr = MPI_Barrier(communicator);
391  AssertThrowMPI(mpi_ierr);
392 
394  {
395  if (across)
396  {
397  out << "[Proc" << i << " " << istart << "-" << iend-1 << "]" << ' ';
398  for (PetscInt i=0; i<nlocal; ++i)
399  out << val[i] << ' ';
400  }
401  else
402  {
403  out << "[Proc " << i << " " << istart << "-" << iend-1 << "]" << std::endl;
404  for (PetscInt i=0; i<nlocal; ++i)
405  out << val[i] << std::endl;
406  }
407  out << std::endl;
408  }
409  }
410  // reset output format
411  out.flags (old_flags);
412  out.precision(old_precision);
413 
414  // restore the representation of the
415  // vector
416  ierr = VecRestoreArray (vector, &val);
417  AssertThrow (ierr == 0, ExcPETScError(ierr));
418 
419  AssertThrow (out, ExcIO());
420  }
421 
422  }
423 
424 }
425 
426 DEAL_II_NAMESPACE_CLOSE
427 
428 #endif // DEAL_II_WITH_PETSC
types::global_dof_index size_type
static ::ExceptionBase & ExcIO()
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
bool is_ascending_and_one_to_one(const MPI_Comm &communicator) const
Definition: index_set.cc:618
size_type size() const
Definition: index_set.h:1575
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:288
VectorOperation::values last_action
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
unsigned int n_mpi_processes(const MPI_Comm &mpi_communicator)
Definition: mpi.cc:65
void fill_index_vector(std::vector< size_type > &indices) const
Definition: index_set.cc:526
Vector & operator=(const Vector &v)
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1314
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:75
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:1717
static ::ExceptionBase & ExcInternalError()