Reference documentation for deal.II version 9.0.0
petsc_vector_base.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/lac/petsc_vector_base.h>
17 
18 #ifdef DEAL_II_WITH_PETSC
19 
20 # include <deal.II/base/memory_consumption.h>
21 # include <deal.II/lac/exceptions.h>
22 # include <deal.II/lac/petsc_compatibility.h>
23 # include <deal.II/lac/petsc_parallel_vector.h>
24 # include <cmath>
25 # include <deal.II/base/multithread_info.h>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 namespace PETScWrappers
30 {
31  namespace internal
32  {
33  VectorReference::operator PetscScalar () const
34  {
35  Assert (index < vector.size(),
36  ExcIndexRange (index, 0, vector.size()));
37 
38  // The vector may have ghost entries. In that case, we first need to
39  // figure out which elements we own locally, then get a pointer to the
40  // elements that are stored here (both the ones we own as well as the
41  // ghost elements). In this array, the locally owned elements come first
42  // followed by the ghost elements whose position we can get from an
43  // index set.
44  if (vector.ghosted)
45  {
46  PetscInt begin, end;
47  PetscErrorCode ierr = VecGetOwnershipRange (vector.vector, &begin,
48  &end);
49  AssertThrow (ierr == 0, ExcPETScError(ierr));
50 
51  Vec locally_stored_elements = PETSC_NULL;
52  ierr = VecGhostGetLocalForm(vector.vector, &locally_stored_elements);
53  AssertThrow (ierr == 0, ExcPETScError(ierr));
54 
55  PetscInt lsize;
56  ierr = VecGetSize(locally_stored_elements, &lsize);
57  AssertThrow (ierr == 0, ExcPETScError(ierr));
58 
59  PetscScalar *ptr;
60  ierr = VecGetArray(locally_stored_elements, &ptr);
61  AssertThrow (ierr == 0, ExcPETScError(ierr));
62 
63  PetscScalar value;
64 
65  if ( index>=static_cast<size_type>(begin)
66  && index<static_cast<size_type>(end) )
67  {
68  //local entry
69  value = *(ptr+index-begin);
70  }
71  else
72  {
73  //ghost entry
74  const size_type ghostidx
75  = vector.ghost_indices.index_within_set(index);
76 
77  Assert(ghostidx+end-begin<(size_type)lsize, ExcInternalError());
78  value = *(ptr+ghostidx+end-begin);
79  }
80 
81  ierr = VecRestoreArray(locally_stored_elements, &ptr);
82  AssertThrow (ierr == 0, ExcPETScError(ierr));
83 
84  ierr = VecGhostRestoreLocalForm(vector.vector, &locally_stored_elements);
85  AssertThrow (ierr == 0, ExcPETScError(ierr));
86 
87  return value;
88  }
89 
90 
91  // first verify that the requested
92  // element is actually locally
93  // available
94  PetscInt begin, end;
95 
96  PetscErrorCode ierr = VecGetOwnershipRange (vector.vector, &begin, &end);
97  AssertThrow (ierr == 0, ExcPETScError(ierr));
98 
99  AssertThrow ((index >= static_cast<size_type>(begin)) &&
100  (index < static_cast<size_type>(end)),
101  ExcAccessToNonlocalElement (index, begin, end-1));
102 
103  PetscInt idx = index;
104  PetscScalar value;
105  ierr = VecGetValues(vector.vector, 1, &idx, &value);
106  AssertThrow (ierr == 0, ExcPETScError(ierr));
107 
108  return value;
109  }
110  }
111 
113  :
114  vector (nullptr),
115  ghosted(false),
116  last_action (::VectorOperation::unknown),
117  obtained_ownership(true)
118  {
120  ExcMessage("PETSc does not support multi-threaded access, set "
121  "the thread limit to 1 in MPI_InitFinalize()."));
122  }
123 
124 
125 
127  :
128  Subscriptor (),
129  ghosted(v.ghosted),
130  ghost_indices(v.ghost_indices),
131  last_action (::VectorOperation::unknown),
132  obtained_ownership(true)
133  {
135  ExcMessage("PETSc does not support multi-threaded access, set "
136  "the thread limit to 1 in MPI_InitFinalize()."));
137 
138  PetscErrorCode ierr = VecDuplicate (v.vector, &vector);
139  AssertThrow (ierr == 0, ExcPETScError(ierr));
140 
141  ierr = VecCopy (v.vector, vector);
142  AssertThrow (ierr == 0, ExcPETScError(ierr));
143  }
144 
145 
146 
147  VectorBase::VectorBase (const Vec &v)
148  :
149  Subscriptor (),
150  vector(v),
151  ghosted(false),
152  last_action (::VectorOperation::unknown),
153  obtained_ownership(false)
154  {
156  ExcMessage("PETSc does not support multi-threaded access, set "
157  "the thread limit to 1 in MPI_InitFinalize()."));
158  }
159 
160 
161 
163  {
164  if (obtained_ownership)
165  {
166  const PetscErrorCode ierr = VecDestroy (&vector);
167  AssertNothrow (ierr == 0, ExcPETScError(ierr));
168  (void) ierr;
169  }
170  }
171 
172 
173 
174  void
176  {
177  if (obtained_ownership)
178  {
179  const PetscErrorCode ierr = VecDestroy (&vector);
180  AssertThrow (ierr == 0, ExcPETScError(ierr));
181  }
182 
183  ghosted = false;
184  ghost_indices.clear ();
186  obtained_ownership = true;
187  }
188 
189 
190 
191  VectorBase &
192  VectorBase::operator = (const PetscScalar s)
193  {
194  AssertIsFinite(s);
195 
196  //TODO[TH]: assert(is_compressed())
197 
198  PetscErrorCode ierr = VecSet (vector, s);
199  AssertThrow (ierr == 0, ExcPETScError(ierr));
200 
201  if (has_ghost_elements())
202  {
203  Vec ghost = PETSC_NULL;
204  ierr = VecGhostGetLocalForm(vector, &ghost);
205  AssertThrow (ierr == 0, ExcPETScError(ierr));
206 
207  ierr = VecSet (ghost, s);
208  AssertThrow (ierr == 0, ExcPETScError(ierr));
209 
210  ierr = VecGhostRestoreLocalForm(vector, &ghost);
211  AssertThrow (ierr == 0, ExcPETScError(ierr));
212  }
213 
214  return *this;
215  }
216 
217 
218 
219  bool
221  {
222  Assert (size() == v.size(),
223  ExcDimensionMismatch(size(), v.size()));
224 
225  PetscBool flag;
226  const PetscErrorCode ierr = VecEqual (vector, v.vector, &flag);
227  AssertThrow (ierr == 0, ExcPETScError(ierr));
228 
229  return (flag == PETSC_TRUE);
230  }
231 
232 
233 
234  bool
236  {
237  Assert (size() == v.size(),
238  ExcDimensionMismatch(size(), v.size()));
239 
240  PetscBool flag;
241  const PetscErrorCode ierr = VecEqual (vector, v.vector, &flag);
242  AssertThrow (ierr == 0, ExcPETScError(ierr));
243 
244  return (flag == PETSC_FALSE);
245  }
246 
247 
248 
249  VectorBase::size_type
251  {
252  PetscInt sz;
253  const PetscErrorCode ierr = VecGetSize (vector, &sz);
254  AssertThrow (ierr == 0, ExcPETScError(ierr));
255 
256  return sz;
257  }
258 
259 
260 
261  VectorBase::size_type
263  {
264  PetscInt sz;
265  const PetscErrorCode ierr = VecGetLocalSize (vector, &sz);
266  AssertThrow (ierr == 0, ExcPETScError(ierr));
267 
268  return sz;
269  }
270 
271 
272 
273  std::pair<VectorBase::size_type, VectorBase::size_type>
275  {
276  PetscInt begin, end;
277  const PetscErrorCode ierr = VecGetOwnershipRange (static_cast<const Vec &>(vector),
278  &begin, &end);
279  AssertThrow (ierr == 0, ExcPETScError(ierr));
280 
281  return std::make_pair (begin, end);
282  }
283 
284 
285 
286  void
287  VectorBase::set (const std::vector<size_type> &indices,
288  const std::vector<PetscScalar> &values)
289  {
290  Assert (indices.size() == values.size(),
291  ExcMessage ("Function called with arguments of different sizes"));
292  do_set_add_operation(indices.size(), indices.data(), values.data(), false);
293  }
294 
295 
296 
297  void
298  VectorBase::add (const std::vector<size_type> &indices,
299  const std::vector<PetscScalar> &values)
300  {
301  Assert (indices.size() == values.size(),
302  ExcMessage ("Function called with arguments of different sizes"));
303  do_set_add_operation(indices.size(), indices.data(), values.data(), true);
304  }
305 
306 
307 
308  void
309  VectorBase::add (const std::vector<size_type> &indices,
310  const ::Vector<PetscScalar> &values)
311  {
312  Assert (indices.size() == values.size(),
313  ExcMessage ("Function called with arguments of different sizes"));
314  do_set_add_operation(indices.size(), indices.data(), values.begin(), true);
315  }
316 
317 
318 
319  void
320  VectorBase::add (const size_type n_elements,
321  const size_type *indices,
322  const PetscScalar *values)
323  {
324  do_set_add_operation(n_elements, indices, values, true);
325  }
326 
327 
328 
329  PetscScalar
331  {
332  Assert (size() == vec.size(),
333  ExcDimensionMismatch(size(), vec.size()));
334 
335  PetscScalar result;
336 
337  // For complex vectors, VecDot() computes
338  // val = (x,y) = y^H x,
339  // where y^H denotes the conjugate transpose of y.
340  // Note that this corresponds to the usual "mathematicians'"
341  // complex inner product where the SECOND argument gets the
342  // complex conjugate, which is also how we document this
343  // operation.
344  const PetscErrorCode ierr = VecDot (vec.vector, vector, &result);
345  AssertThrow (ierr == 0, ExcPETScError(ierr));
346 
347  return result;
348  }
349 
350 
351 
352  PetscScalar
353  VectorBase::add_and_dot (const PetscScalar a,
354  const VectorBase &V,
355  const VectorBase &W)
356  {
357  this->add(a, V);
358  return *this * W;
359  }
360 
361 
362 
363  void
365  {
366  {
367 #ifdef DEBUG
368 #ifdef DEAL_II_WITH_MPI
369  // Check that all processors agree that last_action is the same (or none!)
370 
371  int my_int_last_action = last_action;
372  int all_int_last_action;
373 
374  const int ierr = MPI_Allreduce
375  (&my_int_last_action, &all_int_last_action, 1, MPI_INT, MPI_BOR,
377  AssertThrowMPI(ierr);
378 
379  AssertThrow(all_int_last_action !=
381  ExcMessage("Error: not all processors agree on the last "
382  "VectorOperation before this compress() call."));
383 #endif
384 #endif
385  }
386 
388  || last_action == operation,
389  ExcMessage("Missing compress() or calling with wrong VectorOperation argument."));
390 
391  // note that one may think that
392  // we only need to do something
393  // if in fact the state is
394  // anything but
395  // last_action::unknown. but
396  // that's not true: one
397  // frequently gets into
398  // situations where only one
399  // processor (or a subset of
400  // processors) actually writes
401  // something into a vector, but
402  // we still need to call
403  // VecAssemblyBegin/End on all
404  // processors.
405  PetscErrorCode ierr = VecAssemblyBegin(vector);
406  AssertThrow (ierr == 0, ExcPETScError(ierr));
407  ierr = VecAssemblyEnd(vector);
408  AssertThrow (ierr == 0, ExcPETScError(ierr));
409 
410  // reset the last action field to
411  // indicate that we're back to a
412  // pristine state
414  }
415 
416 
417 
418  VectorBase::real_type
420  {
421  const real_type d = l2_norm();
422  return d*d;
423  }
424 
425 
426 
427  PetscScalar
429  {
430  // We can only use our more efficient
431  // routine in the serial case.
432  if (dynamic_cast<const PETScWrappers::MPI::Vector *>(this) != nullptr)
433  {
434  PetscScalar sum;
435  const PetscErrorCode ierr = VecSum(vector, &sum);
436  AssertThrow (ierr == 0, ExcPETScError(ierr));
437  return sum/static_cast<PetscReal>(size());
438  }
439 
440  // get a representation of the vector and
441  // loop over all the elements
442  PetscScalar *start_ptr;
443  PetscErrorCode ierr = VecGetArray (vector, &start_ptr);
444  AssertThrow (ierr == 0, ExcPETScError(ierr));
445 
446  PetscScalar mean = 0;
447  {
448  PetscScalar sum0 = 0,
449  sum1 = 0,
450  sum2 = 0,
451  sum3 = 0;
452 
453  // use modern processors better by
454  // allowing pipelined commands to be
455  // executed in parallel
456  const PetscScalar *ptr = start_ptr;
457  const PetscScalar *eptr = ptr + (size()/4)*4;
458  while (ptr!=eptr)
459  {
460  sum0 += *ptr++;
461  sum1 += *ptr++;
462  sum2 += *ptr++;
463  sum3 += *ptr++;
464  };
465  // add up remaining elements
466  while (ptr != start_ptr+size())
467  sum0 += *ptr++;
468 
469  mean = (sum0+sum1+sum2+sum3)/static_cast<PetscReal>(size());
470  }
471 
472  // restore the representation of the
473  // vector
474  ierr = VecRestoreArray (vector, &start_ptr);
475  AssertThrow (ierr == 0, ExcPETScError(ierr));
476 
477  return mean;
478  }
479 
480 
481  VectorBase::real_type
483  {
484  real_type d;
485 
486  const PetscErrorCode ierr = VecNorm (vector, NORM_1, &d);
487  AssertThrow (ierr == 0, ExcPETScError(ierr));
488 
489  return d;
490  }
491 
492 
493 
494  VectorBase::real_type
496  {
497  real_type d;
498 
499  const PetscErrorCode ierr = VecNorm (vector, NORM_2, &d);
500  AssertThrow (ierr == 0, ExcPETScError(ierr));
501 
502  return d;
503  }
504 
505 
506 
507  VectorBase::real_type
508  VectorBase::lp_norm (const real_type p) const
509  {
510  // get a representation of the vector and
511  // loop over all the elements
512  PetscScalar *start_ptr;
513  PetscErrorCode ierr = VecGetArray (vector, &start_ptr);
514  AssertThrow (ierr == 0, ExcPETScError(ierr));
515 
516  real_type norm = 0;
517  {
518  real_type sum0 = 0,
519  sum1 = 0,
520  sum2 = 0,
521  sum3 = 0;
522 
523  // use modern processors better by
524  // allowing pipelined commands to be
525  // executed in parallel
526  const PetscScalar *ptr = start_ptr;
527  const PetscScalar *eptr = ptr + (size()/4)*4;
528  while (ptr!=eptr)
529  {
530  sum0 += std::pow(numbers::NumberTraits<value_type>::abs(*ptr++), p);
531  sum1 += std::pow(numbers::NumberTraits<value_type>::abs(*ptr++), p);
532  sum2 += std::pow(numbers::NumberTraits<value_type>::abs(*ptr++), p);
533  sum3 += std::pow(numbers::NumberTraits<value_type>::abs(*ptr++), p);
534  }
535  // add up remaining elements
536  while (ptr != start_ptr+size())
537  sum0 += std::pow(numbers::NumberTraits<value_type>::abs(*ptr++), p);
538 
539  norm = std::pow(sum0+sum1+sum2+sum3,
540  1./p);
541  }
542 
543  // restore the representation of the
544  // vector
545  ierr = VecRestoreArray (vector, &start_ptr);
546  AssertThrow (ierr == 0, ExcPETScError(ierr));
547 
548  return norm;
549  }
550 
551 
552 
553  VectorBase::real_type
555  {
556  real_type d;
557 
558  const PetscErrorCode ierr = VecNorm (vector, NORM_INFINITY, &d);
559  AssertThrow (ierr == 0, ExcPETScError(ierr));
560 
561  return d;
562  }
563 
564 
565 
566  VectorBase::real_type
568  {
569  PetscInt p;
570  real_type d;
571 
572  const PetscErrorCode ierr = VecMin (vector, &p, &d);
573  AssertThrow (ierr == 0, ExcPETScError(ierr));
574 
575  return d;
576  }
577 
578 
579  VectorBase::real_type
581  {
582  PetscInt p;
583  real_type d;
584 
585  const PetscErrorCode ierr = VecMax (vector, &p, &d);
586  AssertThrow (ierr == 0, ExcPETScError(ierr));
587 
588  return d;
589  }
590 
591 
592 
593  bool
595  {
596  // get a representation of the vector and
597  // loop over all the elements
598  PetscScalar *start_ptr;
599  PetscErrorCode ierr = VecGetArray (vector, &start_ptr);
600  AssertThrow (ierr == 0, ExcPETScError(ierr));
601 
602  const PetscScalar *ptr = start_ptr,
603  *eptr = start_ptr + local_size();
604  bool flag = true;
605  while (ptr != eptr)
606  {
607  if (*ptr != value_type())
608  {
609  flag = false;
610  break;
611  }
612  ++ptr;
613  }
614 
615  // restore the representation of the
616  // vector
617  ierr = VecRestoreArray (vector, &start_ptr);
618  AssertThrow (ierr == 0, ExcPETScError(ierr));
619 
620  return flag;
621  }
622 
623 
624  namespace internal
625  {
626  template <typename T>
627  bool is_non_negative (const T &t)
628  {
629  return t >= 0;
630  }
631 
632 
633 
634  template <typename T>
635  bool is_non_negative (const std::complex<T> &)
636  {
637  Assert (false,
638  ExcMessage ("You can't ask a complex value "
639  "whether it is non-negative."))
640  return true;
641  }
642  }
643 
644 
645 
646  bool
647  VectorBase::is_non_negative () const
648  {
649  // get a representation of the vector and
650  // loop over all the elements
651  PetscScalar *start_ptr;
652  PetscErrorCode ierr = VecGetArray (vector, &start_ptr);
653  AssertThrow (ierr == 0, ExcPETScError(ierr));
654 
655  const PetscScalar *ptr = start_ptr,
656  *eptr = start_ptr + local_size();
657  bool flag = true;
658  while (ptr != eptr)
659  {
660  if (! internal::is_non_negative(*ptr))
661  {
662  flag = false;
663  break;
664  }
665  ++ptr;
666  }
667 
668  // restore the representation of the
669  // vector
670  ierr = VecRestoreArray (vector, &start_ptr);
671  AssertThrow (ierr == 0, ExcPETScError(ierr));
672 
673  return flag;
674  }
675 
676 
677 
678  VectorBase &
679  VectorBase::operator *= (const PetscScalar a)
680  {
681  Assert (!has_ghost_elements(), ExcGhostsPresent());
682  AssertIsFinite(a);
683 
684  const PetscErrorCode ierr = VecScale (vector, a);
685  AssertThrow (ierr == 0, ExcPETScError(ierr));
686 
687  return *this;
688  }
689 
690 
691 
692  VectorBase &
693  VectorBase::operator /= (const PetscScalar a)
694  {
695  Assert (!has_ghost_elements(), ExcGhostsPresent());
696  AssertIsFinite(a);
697 
698  const PetscScalar factor = 1./a;
699  AssertIsFinite(factor);
700 
701  const PetscErrorCode ierr = VecScale (vector, factor);
702  AssertThrow (ierr == 0, ExcPETScError(ierr));
703 
704  return *this;
705  }
706 
707 
708 
709  VectorBase &
711  {
712  Assert (!has_ghost_elements(), ExcGhostsPresent());
713  const PetscErrorCode ierr = VecAXPY (vector, 1, v);
714  AssertThrow (ierr == 0, ExcPETScError(ierr));
715 
716  return *this;
717  }
718 
719 
720 
721  VectorBase &
723  {
724  Assert (!has_ghost_elements(), ExcGhostsPresent());
725  const PetscErrorCode ierr = VecAXPY (vector, -1, v);
726  AssertThrow (ierr == 0, ExcPETScError(ierr));
727 
728  return *this;
729  }
730 
731 
732 
733  void
734  VectorBase::add (const PetscScalar s)
735  {
736  Assert (!has_ghost_elements(), ExcGhostsPresent());
737  AssertIsFinite(s);
738 
739  const PetscErrorCode ierr = VecShift (vector, s);
740  AssertThrow (ierr == 0, ExcPETScError(ierr));
741  }
742 
743 
744 
745  void
746  VectorBase::add (const PetscScalar a,
747  const VectorBase &v)
748  {
749  Assert (!has_ghost_elements(), ExcGhostsPresent());
750  AssertIsFinite(a);
751 
752  const PetscErrorCode ierr = VecAXPY (vector, a, v);
753  AssertThrow (ierr == 0, ExcPETScError(ierr));
754  }
755 
756 
757 
758  void
759  VectorBase::add (const PetscScalar a,
760  const VectorBase &v,
761  const PetscScalar b,
762  const VectorBase &w)
763  {
764  Assert (!has_ghost_elements(), ExcGhostsPresent());
765  AssertIsFinite(a);
766  AssertIsFinite(b);
767 
768  const PetscScalar weights[2] = {a,b};
769  Vec addends[2] = {v.vector, w.vector};
770 
771  const PetscErrorCode ierr = VecMAXPY (vector, 2, weights, addends);
772  AssertThrow (ierr == 0, ExcPETScError(ierr));
773  }
774 
775 
776 
777  void
778  VectorBase::sadd (const PetscScalar s,
779  const VectorBase &v)
780  {
781  Assert (!has_ghost_elements(), ExcGhostsPresent());
782  AssertIsFinite(s);
783 
784  const PetscErrorCode ierr = VecAYPX (vector, s, v);
785  AssertThrow (ierr == 0, ExcPETScError(ierr));
786  }
787 
788 
789 
790  void
791  VectorBase::sadd (const PetscScalar s,
792  const PetscScalar a,
793  const VectorBase &v)
794  {
795  Assert (!has_ghost_elements(), ExcGhostsPresent());
796  AssertIsFinite(s);
797  AssertIsFinite(a);
798 
799  // there is nothing like a AXPAY
800  // operation in Petsc, so do it in two
801  // steps
802  *this *= s;
803  add (a,v);
804  }
805 
806 
807 
808  void
809  VectorBase::scale (const VectorBase &factors)
810  {
811  Assert (!has_ghost_elements(), ExcGhostsPresent());
812  const PetscErrorCode ierr
813  = VecPointwiseMult (vector, factors, vector);
814  AssertThrow (ierr == 0, ExcPETScError(ierr));
815  }
816 
817 
818 
819  void
820  VectorBase::equ (const PetscScalar a,
821  const VectorBase &v)
822  {
823  Assert (!has_ghost_elements(), ExcGhostsPresent());
824  AssertIsFinite(a);
825 
826  Assert (size() == v.size(),
827  ExcDimensionMismatch (size(), v.size()));
828 
829  // there is no simple operation for this
830  // in PETSc. there are multiple ways to
831  // emulate it, we choose this one:
832  const PetscErrorCode ierr = VecCopy (v.vector, vector);
833  AssertThrow (ierr == 0, ExcPETScError(ierr));
834 
835  *this *= a;
836  }
837 
838 
839 
840  void
842  const VectorBase &b)
843  {
844  Assert (!has_ghost_elements(), ExcGhostsPresent());
845  const PetscErrorCode ierr = VecPointwiseDivide (vector, a, b);
846  AssertThrow (ierr == 0, ExcPETScError(ierr));
847  }
848 
849 
850 
851  void
852  VectorBase::write_ascii (const PetscViewerFormat format)
853  {
854  //TODO[TH]:assert(is_compressed())
855 
856  // Set options
857  PetscErrorCode ierr = PetscViewerSetFormat (PETSC_VIEWER_STDOUT_WORLD,
858  format);
859  AssertThrow (ierr == 0, ExcPETScError(ierr));
860 
861  // Write to screen
862  ierr = VecView (vector, PETSC_VIEWER_STDOUT_WORLD);
863  AssertThrow (ierr == 0, ExcPETScError(ierr));
864  }
865 
866 
867 
868  void
869  VectorBase::print (std::ostream &out,
870  const unsigned int precision,
871  const bool scientific,
872  const bool across) const
873  {
874  AssertThrow (out, ExcIO());
875 
876  // get a representation of the vector and
877  // loop over all the elements
878  PetscScalar *val;
879  PetscErrorCode ierr = VecGetArray (vector, &val);
880 
881  AssertThrow (ierr == 0, ExcPETScError(ierr));
882 
883  // save the state of out stream
884  const std::ios::fmtflags old_flags = out.flags();
885  const unsigned int old_precision = out.precision (precision);
886 
887  out.precision (precision);
888  if (scientific)
889  out.setf (std::ios::scientific, std::ios::floatfield);
890  else
891  out.setf (std::ios::fixed, std::ios::floatfield);
892 
893  if (across)
894  for (size_type i=0; i<local_size(); ++i)
895  out << val[i] << ' ';
896  else
897  for (size_type i=0; i<local_size(); ++i)
898  out << val[i] << std::endl;
899  out << std::endl;
900 
901  // reset output format
902  out.flags (old_flags);
903  out.precision(old_precision);
904 
905  // restore the representation of the
906  // vector
907  ierr = VecRestoreArray (vector, &val);
908  AssertThrow (ierr == 0, ExcPETScError(ierr));
909 
910  AssertThrow (out, ExcIO());
911  }
912 
913 
914 
915  void
917  {
918  const PetscErrorCode ierr = VecSwap (vector, v.vector);
919  AssertThrow (ierr == 0, ExcPETScError(ierr));
920  }
921 
922 
923 
924  VectorBase::operator const Vec &() const
925  {
926  return vector;
927  }
928 
929 
930  std::size_t
932  {
933  std::size_t mem = sizeof(Vec)+sizeof(last_action)
936 
937  // TH: I am relatively sure that PETSc is
938  // storing the local data in a contiguous
939  // block without indices:
940  mem += local_size()*sizeof(PetscScalar);
941  // assume that PETSc is storing one index
942  // and one double per ghost element
943  if (ghosted)
944  mem += ghost_indices.n_elements()*(sizeof(PetscScalar)+sizeof(int));
945 
946  //TODO[TH]: size of constant memory for PETSc?
947  return mem;
948  }
949 
950 
951 
952  void
953  VectorBase::do_set_add_operation (const size_type n_elements,
954  const size_type *indices,
955  const PetscScalar *values,
956  const bool add_values)
957  {
958  ::VectorOperation::values action = (add_values ?
961  Assert ((last_action == action)
962  ||
963  (last_action == ::VectorOperation::unknown),
964  internal::VectorReference::ExcWrongMode (action,
965  last_action));
966  Assert (!has_ghost_elements(), ExcGhostsPresent());
967  // VecSetValues complains if we
968  // come with an empty
969  // vector. however, it is not a
970  // collective operation, so we
971  // can skip the call if necessary
972  // (unlike the above calls)
973  if (n_elements != 0)
974  {
975  const PetscInt *petsc_indices = reinterpret_cast<const PetscInt *>(indices);
976 
977  const InsertMode mode = (add_values ? ADD_VALUES : INSERT_VALUES);
978  const PetscErrorCode ierr = VecSetValues (vector, n_elements, petsc_indices,
979  values, mode);
980  AssertThrow (ierr == 0, ExcPETScError(ierr));
981  }
982 
983  // set the mode here, independent of whether we have actually
984  // written elements or whether the list was empty
985  last_action = action;
986  }
987 
988 }
989 
990 DEAL_II_NAMESPACE_CLOSE
991 
992 #endif // DEAL_II_WITH_PETSC
#define AssertNothrow(cond, exc)
Definition: exceptions.h:1183
VectorBase & operator-=(const VectorBase &V)
static ::ExceptionBase & ExcIO()
void sadd(const PetscScalar s, const VectorBase &V)
void ratio(const VectorBase &a, const VectorBase &b)
VectorBase & operator/=(const PetscScalar factor)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
VectorBase & operator=(const PetscScalar s)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void scale(const VectorBase &scaling_factors)
std::pair< size_type, size_type > local_range() const
real_type lp_norm(const real_type p) const
void clear()
Definition: index_set.h:1549
static ::ExceptionBase & ExcMessage(std::string arg1)
void compress(const VectorOperation::values operation)
void do_set_add_operation(const size_type n_elements, const size_type *indices, const PetscScalar *values, const bool add_values)
VectorOperation::values last_action
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static bool is_running_single_threaded()
PetscScalar add_and_dot(const PetscScalar a, const VectorBase &V, const VectorBase &W)
std::size_t memory_consumption() const
void equ(const PetscScalar a, const VectorBase &V)
VectorBase & operator*=(const PetscScalar factor)
static ::ExceptionBase & ExcGhostsPresent()
VectorBase & operator+=(const VectorBase &V)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
#define AssertThrowMPI(error_code)
Definition: exceptions.h:1314
PetscScalar operator*(const VectorBase &vec) const
bool has_ghost_elements() const
PetscScalar mean_value() const
bool operator==(const VectorBase &v) const
void add(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)
void write_ascii(const PetscViewerFormat format=PETSC_VIEWER_DEFAULT)
virtual const MPI_Comm & get_mpi_communicator() const
#define AssertIsFinite(number)
Definition: exceptions.h:1299
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
static ::ExceptionBase & ExcInternalError()
bool operator!=(const VectorBase &v) const
void set(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)