Reference documentation for deal.II version 8.4.2
petsc_vector_base.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2015 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/petsc_vector.h>
22 # include <deal.II/lac/petsc_parallel_vector.h>
23 # include <cmath>
24 # include <deal.II/base/multithread_info.h>
25 
26 DEAL_II_NAMESPACE_OPEN
27 
28 namespace PETScWrappers
29 {
30  namespace internal
31  {
32  VectorReference::operator PetscScalar () const
33  {
34  Assert (index < vector.size(),
35  ExcIndexRange (index, 0, vector.size()));
36 
37  // if the vector is local, then
38  // simply access the element we
39  // are interested in
40  if (dynamic_cast<const PETScWrappers::Vector *>(&vector) != 0)
41  {
42  PetscInt idx = index;
43  PetscScalar value;
44  int ierr = VecGetValues(vector.vector, 1, &idx, &value);
45  AssertThrow (ierr == 0, ExcPETScError(ierr));
46  return value;
47  }
48  // else see if we are dealing
49  // with a parallel vector
50  else if (dynamic_cast<const PETScWrappers::MPI::Vector *>(&vector) != 0)
51  {
52  int ierr;
53 
54  // there is the possibility
55  // that the vector has
56  // ghost elements. in that
57  // case, we first need to
58  // figure out which
59  // elements we own locally,
60  // then get a pointer to
61  // the elements that are
62  // stored here (both the
63  // ones we own as well as
64  // the ghost elements). in
65  // this array, the locally
66  // owned elements come
67  // first followed by the
68  // ghost elements whose
69  // position we can get from
70  // an index set
71  if (vector.ghosted)
72  {
73  PetscInt begin, end;
74  ierr = VecGetOwnershipRange (vector.vector, &begin, &end);
75  AssertThrow (ierr == 0, ExcPETScError(ierr));
76 
77  Vec locally_stored_elements = PETSC_NULL;
78  ierr = VecGhostGetLocalForm(vector.vector, &locally_stored_elements);
79  AssertThrow (ierr == 0, ExcPETScError(ierr));
80 
81  PetscInt lsize;
82  ierr = VecGetSize(locally_stored_elements, &lsize);
83  AssertThrow (ierr == 0, ExcPETScError(ierr));
84 
85  PetscScalar *ptr;
86  ierr = VecGetArray(locally_stored_elements, &ptr);
87  AssertThrow (ierr == 0, ExcPETScError(ierr));
88 
89  PetscScalar value;
90 
91  if ( index>=static_cast<size_type>(begin)
92  && index<static_cast<size_type>(end) )
93  {
94  //local entry
95  value = *(ptr+index-begin);
96  }
97  else
98  {
99  //ghost entry
100  const size_type ghostidx
101  = vector.ghost_indices.index_within_set(index);
102 
103  Assert(ghostidx+end-begin<(size_type)lsize, ExcInternalError());
104  value = *(ptr+ghostidx+end-begin);
105  }
106 
107  ierr = VecRestoreArray(locally_stored_elements, &ptr);
108  AssertThrow (ierr == 0, ExcPETScError(ierr));
109 
110  ierr = VecGhostRestoreLocalForm(vector.vector, &locally_stored_elements);
111  AssertThrow (ierr == 0, ExcPETScError(ierr));
112 
113  return value;
114  }
115 
116 
117  // first verify that the requested
118  // element is actually locally
119  // available
120  PetscInt begin, end;
121 
122  ierr = VecGetOwnershipRange (vector.vector, &begin, &end);
123  AssertThrow (ierr == 0, ExcPETScError(ierr));
124 
125 
126 
127  AssertThrow ((index >= static_cast<size_type>(begin)) &&
128  (index < static_cast<size_type>(end)),
129  ExcAccessToNonlocalElement (index, begin, end-1));
130 
131  // old version which only work with
132  // VecGetArray()...
133  PetscInt idx = index;
134  PetscScalar value;
135  ierr = VecGetValues(vector.vector, 1, &idx, &value);
136  AssertThrow (ierr == 0, ExcPETScError(ierr));
137 
138  return value;
139  }
140  else
141  // what? what other kind of vector
142  // exists there?
143  Assert (false, ExcInternalError());
144 
145  return -1e20;
146  }
147  }
148 
150  :
151  ghosted(false),
152  last_action (::VectorOperation::unknown),
153  attained_ownership(true)
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  Subscriptor (),
165  ghosted(v.ghosted),
167  last_action (::VectorOperation::unknown),
168  attained_ownership(true)
169  {
171  ExcMessage("PETSc does not support multi-threaded access, set "
172  "the thread limit to 1 in MPI_InitFinalize()."));
173 
174  int ierr = VecDuplicate (v.vector, &vector);
175  AssertThrow (ierr == 0, ExcPETScError(ierr));
176 
177  ierr = VecCopy (v.vector, vector);
178  AssertThrow (ierr == 0, ExcPETScError(ierr));
179  }
180 
181 
182 
183  VectorBase::VectorBase (const Vec &v)
184  :
185  Subscriptor (),
186  vector(v),
187  ghosted(false),
188  last_action (::VectorOperation::unknown),
189  attained_ownership(false)
190  {
192  ExcMessage("PETSc does not support multi-threaded access, set "
193  "the thread limit to 1 in MPI_InitFinalize()."));
194  }
195 
196 
197 
199  {
200  if (attained_ownership)
201  {
202 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
203  const int ierr = VecDestroy (vector);
204 #else
205  const int ierr = VecDestroy (&vector);
206 #endif
207  AssertThrow (ierr == 0, ExcPETScError(ierr));
208  }
209  }
210 
211 
212 
213  void
215  {
216  if (attained_ownership)
217  {
218 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
219  const int ierr = VecDestroy (vector);
220 #else
221  const int ierr = VecDestroy (&vector);
222 #endif
223  AssertThrow (ierr == 0, ExcPETScError(ierr));
224  }
225 
226  ghosted = false;
227  ghost_indices.clear ();
228  last_action = ::VectorOperation::unknown;
229  attained_ownership = true;
230  }
231 
232 
233 
234  VectorBase &
235  VectorBase::operator = (const PetscScalar s)
236  {
237  AssertIsFinite(s);
238 
239  //TODO[TH]: assert(is_compressed())
240 
241  int ierr = VecSet (vector, s);
242  AssertThrow (ierr == 0, ExcPETScError(ierr));
243 
244  if (has_ghost_elements())
245  {
246  Vec ghost = PETSC_NULL;
247  ierr = VecGhostGetLocalForm(vector, &ghost);
248  AssertThrow (ierr == 0, ExcPETScError(ierr));
249 
250  ierr = VecSet (ghost, s);
251  AssertThrow (ierr == 0, ExcPETScError(ierr));
252 
253  ierr = VecGhostRestoreLocalForm(vector, &ghost);
254  AssertThrow (ierr == 0, ExcPETScError(ierr));
255  }
256 
257  return *this;
258  }
259 
260 
261 
262  bool
264  {
265  Assert (size() == v.size(),
266  ExcDimensionMismatch(size(), v.size()));
267 
268 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
269  PetscTruth
270 #else
271  PetscBool
272 #endif
273  flag;
274 
275  const int ierr = VecEqual (vector, v.vector, &flag);
276  AssertThrow (ierr == 0, ExcPETScError(ierr));
277 
278  return (flag == PETSC_TRUE);
279  }
280 
281 
282 
283  bool
285  {
286  Assert (size() == v.size(),
287  ExcDimensionMismatch(size(), v.size()));
288 
289 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
290  PetscTruth
291 #else
292  PetscBool
293 #endif
294  flag;
295 
296  const int ierr = VecEqual (vector, v.vector, &flag);
297  AssertThrow (ierr == 0, ExcPETScError(ierr));
298 
299  return (flag == PETSC_FALSE);
300  }
301 
302 
303 
304  VectorBase::size_type
306  {
307  PetscInt sz;
308  const int ierr = VecGetSize (vector, &sz);
309  AssertThrow (ierr == 0, ExcPETScError(ierr));
310 
311  return sz;
312  }
313 
314 
315 
316  VectorBase::size_type
318  {
319  PetscInt sz;
320  const int ierr = VecGetLocalSize (vector, &sz);
321  AssertThrow (ierr == 0, ExcPETScError(ierr));
322 
323  return sz;
324  }
325 
326 
327 
328  std::pair<VectorBase::size_type, VectorBase::size_type>
330  {
331  PetscInt begin, end;
332  const int ierr = VecGetOwnershipRange (static_cast<const Vec &>(vector),
333  &begin, &end);
334  AssertThrow (ierr == 0, ExcPETScError(ierr));
335 
336  return std::make_pair (begin, end);
337  }
338 
339 
340 
341  void
342  VectorBase::set (const std::vector<size_type> &indices,
343  const std::vector<PetscScalar> &values)
344  {
345  Assert (indices.size() == values.size(),
346  ExcMessage ("Function called with arguments of different sizes"));
347  do_set_add_operation(indices.size(), &indices[0], &values[0], false);
348  }
349 
350 
351 
352  void
353  VectorBase::add (const std::vector<size_type> &indices,
354  const std::vector<PetscScalar> &values)
355  {
356  Assert (indices.size() == values.size(),
357  ExcMessage ("Function called with arguments of different sizes"));
358  do_set_add_operation(indices.size(), &indices[0], &values[0], true);
359  }
360 
361 
362 
363  void
364  VectorBase::add (const std::vector<size_type> &indices,
365  const ::Vector<PetscScalar> &values)
366  {
367  Assert (indices.size() == values.size(),
368  ExcMessage ("Function called with arguments of different sizes"));
369  do_set_add_operation(indices.size(), &indices[0], values.begin(), true);
370  }
371 
372 
373 
374  void
375  VectorBase::add (const size_type n_elements,
376  const size_type *indices,
377  const PetscScalar *values)
378  {
379  do_set_add_operation(n_elements, indices, values, true);
380  }
381 
382 
383 
384  PetscScalar
386  {
387  Assert (size() == vec.size(),
388  ExcDimensionMismatch(size(), vec.size()));
389 
390  PetscScalar result;
391 
392  //For complex vectors, VecDot() computes
393  // val = (x,y) = y^H x,
394  //where y^H denotes the conjugate transpose of y.
395  //Note that this corresponds to the usual "mathematicians" complex inner product where the SECOND argument gets the complex conjugate.
396  const int ierr = VecDot (vec.vector, vector, &result);
397  AssertThrow (ierr == 0, ExcPETScError(ierr));
398 
399  return result;
400  }
401 
402 
403 
404  PetscScalar
405  VectorBase::add_and_dot (const PetscScalar a,
406  const VectorBase &V,
407  const VectorBase &W)
408  {
409  this->add(a, V);
410  return *this * W;
411  }
412 
413 
414 
415  void
416  VectorBase::compress (const VectorOperation::values operation)
417  {
418 #ifdef DEBUG
419 #ifdef DEAL_II_WITH_MPI
420  // Check that all processors agree that last_action is the same (or none!)
421 
422  int my_int_last_action = last_action;
423  int all_int_last_action;
424 
425  MPI_Allreduce(&my_int_last_action, &all_int_last_action, 1, MPI_INT,
426  MPI_BOR, get_mpi_communicator());
427 
428  AssertThrow(all_int_last_action != (::VectorOperation::add | ::VectorOperation::insert),
429  ExcMessage("Error: not all processors agree on the last VectorOperation before this compress() call."));
430 #endif
431 #endif
432 
433  AssertThrow(last_action == ::VectorOperation::unknown
434  || last_action == operation,
435  ExcMessage("Missing compress() or calling with wrong VectorOperation argument."));
436 
437  // note that one may think that
438  // we only need to do something
439  // if in fact the state is
440  // anything but
441  // last_action::unknown. but
442  // that's not true: one
443  // frequently gets into
444  // situations where only one
445  // processor (or a subset of
446  // processors) actually writes
447  // something into a vector, but
448  // we still need to call
449  // VecAssemblyBegin/End on all
450  // processors.
451  int ierr;
452  ierr = VecAssemblyBegin(vector);
453  AssertThrow (ierr == 0, ExcPETScError(ierr));
454  ierr = VecAssemblyEnd(vector);
455  AssertThrow (ierr == 0, ExcPETScError(ierr));
456 
457  // reset the last action field to
458  // indicate that we're back to a
459  // pristine state
460  last_action = ::VectorOperation::unknown;
461  }
462 
463 
464 
465  VectorBase::real_type
467  {
468  const real_type d = l2_norm();
469  return d*d;
470  }
471 
472 
473 
474  PetscScalar
476  {
477  int ierr;
478 
479  // We can only use our more efficient
480  // routine in the serial case.
481  if (dynamic_cast<const PETScWrappers::MPI::Vector *>(this) != 0)
482  {
483  PetscScalar sum;
484  ierr = VecSum(vector, &sum);
485  AssertThrow (ierr == 0, ExcPETScError(ierr));
486  return sum/static_cast<PetscReal>(size());
487  }
488 
489  // get a representation of the vector and
490  // loop over all the elements
491  PetscScalar *start_ptr;
492  ierr = VecGetArray (vector, &start_ptr);
493  AssertThrow (ierr == 0, ExcPETScError(ierr));
494 
495  PetscScalar mean = 0;
496  {
497  PetscScalar sum0 = 0,
498  sum1 = 0,
499  sum2 = 0,
500  sum3 = 0;
501 
502  // use modern processors better by
503  // allowing pipelined commands to be
504  // executed in parallel
505  const PetscScalar *ptr = start_ptr;
506  const PetscScalar *eptr = ptr + (size()/4)*4;
507  while (ptr!=eptr)
508  {
509  sum0 += *ptr++;
510  sum1 += *ptr++;
511  sum2 += *ptr++;
512  sum3 += *ptr++;
513  };
514  // add up remaining elements
515  while (ptr != start_ptr+size())
516  sum0 += *ptr++;
517 
518  mean = (sum0+sum1+sum2+sum3)/static_cast<PetscReal>(size());
519  }
520 
521  // restore the representation of the
522  // vector
523  ierr = VecRestoreArray (vector, &start_ptr);
524  AssertThrow (ierr == 0, ExcPETScError(ierr));
525 
526  return mean;
527  }
528 
529 
530  VectorBase::real_type
532  {
533  real_type d;
534 
535  const int ierr = VecNorm (vector, NORM_1, &d);
536  AssertThrow (ierr == 0, ExcPETScError(ierr));
537 
538  return d;
539  }
540 
541 
542 
543  VectorBase::real_type
545  {
546  real_type d;
547 
548  const int ierr = VecNorm (vector, NORM_2, &d);
549  AssertThrow (ierr == 0, ExcPETScError(ierr));
550 
551  return d;
552  }
553 
554 
555 
556  VectorBase::real_type
557  VectorBase::lp_norm (const real_type p) const
558  {
559  // get a representation of the vector and
560  // loop over all the elements
561  PetscScalar *start_ptr;
562  int ierr = VecGetArray (vector, &start_ptr);
563  AssertThrow (ierr == 0, ExcPETScError(ierr));
564 
565  real_type norm = 0;
566  {
567  real_type sum0 = 0,
568  sum1 = 0,
569  sum2 = 0,
570  sum3 = 0;
571 
572  // use modern processors better by
573  // allowing pipelined commands to be
574  // executed in parallel
575  const PetscScalar *ptr = start_ptr;
576  const PetscScalar *eptr = ptr + (size()/4)*4;
577  while (ptr!=eptr)
578  {
579  sum0 += std::pow(numbers::NumberTraits<value_type>::abs(*ptr++), p);
580  sum1 += std::pow(numbers::NumberTraits<value_type>::abs(*ptr++), p);
581  sum2 += std::pow(numbers::NumberTraits<value_type>::abs(*ptr++), p);
582  sum3 += std::pow(numbers::NumberTraits<value_type>::abs(*ptr++), p);
583  }
584  // add up remaining elements
585  while (ptr != start_ptr+size())
586  sum0 += std::pow(numbers::NumberTraits<value_type>::abs(*ptr++), p);
587 
588  norm = std::pow(sum0+sum1+sum2+sum3,
589  1./p);
590  }
591 
592  // restore the representation of the
593  // vector
594  ierr = VecRestoreArray (vector, &start_ptr);
595  AssertThrow (ierr == 0, ExcPETScError(ierr));
596 
597  return norm;
598  }
599 
600 
601 
602  VectorBase::real_type
604  {
605  real_type d;
606 
607  const int ierr = VecNorm (vector, NORM_INFINITY, &d);
608  AssertThrow (ierr == 0, ExcPETScError(ierr));
609 
610  return d;
611  }
612 
613 
614 
615  VectorBase::real_type
617  {
618  real_type d;
619  Assert (!has_ghost_elements(), ExcGhostsPresent());
620  const int ierr = VecNormalize (vector, &d);
621  AssertThrow (ierr == 0, ExcPETScError(ierr));
622 
623  return d;
624  }
625 
626 
627  VectorBase::real_type
629  {
630  PetscInt p;
631  real_type d;
632 
633  const int ierr = VecMin (vector, &p, &d);
634  AssertThrow (ierr == 0, ExcPETScError(ierr));
635 
636  return d;
637  }
638 
639 
640  VectorBase::real_type
642  {
643  PetscInt p;
644  real_type d;
645 
646  const int ierr = VecMax (vector, &p, &d);
647  AssertThrow (ierr == 0, ExcPETScError(ierr));
648 
649  return d;
650  }
651 
652 
653  VectorBase &
655  {
656  Assert (!has_ghost_elements(), ExcGhostsPresent());
657 
658  const int ierr = VecAbs (vector);
659  AssertThrow (ierr == 0, ExcPETScError(ierr));
660 
661  return *this;
662  }
663 
664 
665 
666  VectorBase &
668  {
669  Assert (!has_ghost_elements(), ExcGhostsPresent());
670 
671  const int ierr = VecConjugate (vector);
672  AssertThrow (ierr == 0, ExcPETScError(ierr));
673 
674  return *this;
675  }
676 
677 
678 
679  VectorBase &
681  {
682  Assert (!has_ghost_elements(), ExcGhostsPresent());
683 
684  const int ierr = VecPointwiseMult (vector,vector,vector);
685  AssertThrow (ierr == 0, ExcPETScError(ierr));
686 
687  return *this;
688  }
689 
690 
691  VectorBase &
693  {
694  Assert (!has_ghost_elements(), ExcGhostsPresent());
695  const int ierr = VecPointwiseMult (vector,vector,v);
696  AssertThrow (ierr == 0, ExcPETScError(ierr));
697 
698  return *this;
699  }
700 
701 
702  VectorBase &
704  const VectorBase &v)
705  {
706  Assert (!has_ghost_elements(), ExcGhostsPresent());
707  const int ierr = VecPointwiseMult (vector,u,v);
708  AssertThrow (ierr == 0, ExcPETScError(ierr));
709 
710  return *this;
711  }
712 
713  bool
715  {
716  // get a representation of the vector and
717  // loop over all the elements
718  PetscScalar *start_ptr;
719  int ierr = VecGetArray (vector, &start_ptr);
720  AssertThrow (ierr == 0, ExcPETScError(ierr));
721 
722  const PetscScalar *ptr = start_ptr,
723  *eptr = start_ptr + local_size();
724  bool flag = true;
725  while (ptr != eptr)
726  {
727  if (*ptr != value_type())
728  {
729  flag = false;
730  break;
731  }
732  ++ptr;
733  }
734 
735  // restore the representation of the
736  // vector
737  ierr = VecRestoreArray (vector, &start_ptr);
738  AssertThrow (ierr == 0, ExcPETScError(ierr));
739 
740  return flag;
741  }
742 
743 
744  namespace internal
745  {
746  template <typename T>
747  bool is_non_negative (const T &t)
748  {
749  return t >= 0;
750  }
751 
752 
753 
754  template <typename T>
755  bool is_non_negative (const std::complex<T> &)
756  {
757  Assert (false,
758  ExcMessage ("You can't ask a complex value "
759  "whether it is non-negative."))
760  return true;
761  }
762  }
763 
764 
765 
766  bool
768  {
769  // get a representation of the vector and
770  // loop over all the elements
771  PetscScalar *start_ptr;
772  int ierr = VecGetArray (vector, &start_ptr);
773  AssertThrow (ierr == 0, ExcPETScError(ierr));
774 
775  const PetscScalar *ptr = start_ptr,
776  *eptr = start_ptr + local_size();
777  bool flag = true;
778  while (ptr != eptr)
779  {
780  if (! internal::is_non_negative(*ptr))
781  {
782  flag = false;
783  break;
784  }
785  ++ptr;
786  }
787 
788  // restore the representation of the
789  // vector
790  ierr = VecRestoreArray (vector, &start_ptr);
791  AssertThrow (ierr == 0, ExcPETScError(ierr));
792 
793  return flag;
794  }
795 
796 
797 
798  VectorBase &
799  VectorBase::operator *= (const PetscScalar a)
800  {
801  Assert (!has_ghost_elements(), ExcGhostsPresent());
802  AssertIsFinite(a);
803 
804  const int ierr = VecScale (vector, a);
805  AssertThrow (ierr == 0, ExcPETScError(ierr));
806 
807  return *this;
808  }
809 
810 
811 
812  VectorBase &
813  VectorBase::operator /= (const PetscScalar a)
814  {
815  Assert (!has_ghost_elements(), ExcGhostsPresent());
816  AssertIsFinite(a);
817 
818  const PetscScalar factor = 1./a;
819  AssertIsFinite(factor);
820 
821  const int ierr = VecScale (vector, factor);
822  AssertThrow (ierr == 0, ExcPETScError(ierr));
823 
824  return *this;
825  }
826 
827 
828 
829  VectorBase &
831  {
832  Assert (!has_ghost_elements(), ExcGhostsPresent());
833  const int ierr = VecAXPY (vector, 1, v);
834  AssertThrow (ierr == 0, ExcPETScError(ierr));
835 
836  return *this;
837  }
838 
839 
840 
841  VectorBase &
843  {
844  Assert (!has_ghost_elements(), ExcGhostsPresent());
845  const int ierr = VecAXPY (vector, -1, v);
846  AssertThrow (ierr == 0, ExcPETScError(ierr));
847 
848  return *this;
849  }
850 
851 
852 
853  void
854  VectorBase::add (const PetscScalar s)
855  {
856  Assert (!has_ghost_elements(), ExcGhostsPresent());
857  AssertIsFinite(s);
858 
859  const int ierr = VecShift (vector, s);
860  AssertThrow (ierr == 0, ExcPETScError(ierr));
861  }
862 
863 
864 
865  void
867  {
868  *this += v;
869  }
870 
871 
872 
873  void
874  VectorBase::add (const PetscScalar a,
875  const VectorBase &v)
876  {
877  Assert (!has_ghost_elements(), ExcGhostsPresent());
878  AssertIsFinite(a);
879 
880  const int ierr = VecAXPY (vector, a, v);
881  AssertThrow (ierr == 0, ExcPETScError(ierr));
882  }
883 
884 
885 
886  void
887  VectorBase::add (const PetscScalar a,
888  const VectorBase &v,
889  const PetscScalar b,
890  const VectorBase &w)
891  {
892  Assert (!has_ghost_elements(), ExcGhostsPresent());
893  AssertIsFinite(a);
894  AssertIsFinite(b);
895 
896  const PetscScalar weights[2] = {a,b};
897  Vec addends[2] = {v.vector, w.vector};
898 
899  const int ierr = VecMAXPY (vector, 2, weights, addends);
900  AssertThrow (ierr == 0, ExcPETScError(ierr));
901  }
902 
903 
904 
905  void
906  VectorBase::sadd (const PetscScalar s,
907  const VectorBase &v)
908  {
909  Assert (!has_ghost_elements(), ExcGhostsPresent());
910  AssertIsFinite(s);
911 
912  const int ierr = VecAYPX (vector, s, v);
913  AssertThrow (ierr == 0, ExcPETScError(ierr));
914  }
915 
916 
917 
918  void
919  VectorBase::sadd (const PetscScalar s,
920  const PetscScalar a,
921  const VectorBase &v)
922  {
923  Assert (!has_ghost_elements(), ExcGhostsPresent());
924  AssertIsFinite(s);
925  AssertIsFinite(a);
926 
927  // there is nothing like a AXPAY
928  // operation in Petsc, so do it in two
929  // steps
930  *this *= s;
931  add (a,v);
932  }
933 
934 
935 
936  void
937  VectorBase::sadd (const PetscScalar s,
938  const PetscScalar a,
939  const VectorBase &v,
940  const PetscScalar b,
941  const VectorBase &w)
942  {
943  Assert (!has_ghost_elements(), ExcGhostsPresent());
944  AssertIsFinite(s);
945  AssertIsFinite(a);
946  AssertIsFinite(b);
947 
948  // there is no operation like MAXPAY, so
949  // do it in two steps
950  *this *= s;
951 
952  const PetscScalar weights[2] = {a,b};
953  Vec addends[2] = {v.vector,w.vector};
954 
955  const int ierr = VecMAXPY (vector, 2, weights, addends);
956  AssertThrow (ierr == 0, ExcPETScError(ierr));
957  }
958 
959 
960 
961  void
962  VectorBase::sadd (const PetscScalar s,
963  const PetscScalar a,
964  const VectorBase &v,
965  const PetscScalar b,
966  const VectorBase &w,
967  const PetscScalar c,
968  const VectorBase &x)
969  {
970  Assert (!has_ghost_elements(), ExcGhostsPresent());
971  AssertIsFinite(s);
972  AssertIsFinite(a);
973  AssertIsFinite(b);
974  AssertIsFinite(c);
975 
976  // there is no operation like MAXPAY, so
977  // do it in two steps
978  *this *= s;
979 
980  const PetscScalar weights[3] = {a,b,c};
981  Vec addends[3] = {v.vector, w.vector, x.vector};
982 
983  const int ierr = VecMAXPY (vector, 3, weights, addends);
984  AssertThrow (ierr == 0, ExcPETScError(ierr));
985  }
986 
987 
988 
989  void
990  VectorBase::scale (const VectorBase &factors)
991  {
992  Assert (!has_ghost_elements(), ExcGhostsPresent());
993  const int ierr
994  = VecPointwiseMult (vector, factors, vector);
995  AssertThrow (ierr == 0, ExcPETScError(ierr));
996  }
997 
998 
999 
1000  void
1001  VectorBase::equ (const PetscScalar a,
1002  const VectorBase &v)
1003  {
1004  Assert (!has_ghost_elements(), ExcGhostsPresent());
1005  AssertIsFinite(a);
1006 
1007  Assert (size() == v.size(),
1008  ExcDimensionMismatch (size(), v.size()));
1009 
1010  // there is no simple operation for this
1011  // in PETSc. there are multiple ways to
1012  // emulate it, we choose this one:
1013  const int ierr = VecCopy (v.vector, vector);
1014  AssertThrow (ierr == 0, ExcPETScError(ierr));
1015 
1016  *this *= a;
1017  }
1018 
1019 
1020 
1021  void
1022  VectorBase::equ (const PetscScalar a,
1023  const VectorBase &v,
1024  const PetscScalar b,
1025  const VectorBase &w)
1026  {
1027  Assert (!has_ghost_elements(), ExcGhostsPresent());
1028  AssertIsFinite(a);
1029  AssertIsFinite(b);
1030 
1031  Assert (size() == v.size(),
1032  ExcDimensionMismatch (size(), v.size()));
1033 
1034  // there is no simple operation for this
1035  // in PETSc. there are multiple ways to
1036  // emulate it, we choose this one:
1037  const int ierr = VecCopy (v.vector, vector);
1038  AssertThrow (ierr == 0, ExcPETScError(ierr));
1039 
1040  sadd (a, b, w);
1041  }
1042 
1043 
1044 
1045  void
1047  const VectorBase &b)
1048  {
1049  Assert (!has_ghost_elements(), ExcGhostsPresent());
1050  const int ierr = VecPointwiseDivide (vector, a, b);
1051  AssertThrow (ierr == 0, ExcPETScError(ierr));
1052  }
1053 
1054 
1055 
1056  void
1057  VectorBase::write_ascii (const PetscViewerFormat format)
1058  {
1059  //TODO[TH]:assert(is_compressed())
1060 
1061  // Set options
1062  PetscViewerSetFormat (PETSC_VIEWER_STDOUT_WORLD,
1063  format);
1064 
1065  // Write to screen
1066  VecView (vector, PETSC_VIEWER_STDOUT_WORLD);
1067  }
1068 
1069 
1070 
1071  void
1072  VectorBase::print (std::ostream &out,
1073  const unsigned int precision,
1074  const bool scientific,
1075  const bool across) const
1076  {
1077  AssertThrow (out, ExcIO());
1078 
1079  // get a representation of the vector and
1080  // loop over all the elements
1081  PetscScalar *val;
1082  int ierr = VecGetArray (vector, &val);
1083 
1084  AssertThrow (ierr == 0, ExcPETScError(ierr));
1085 
1086  // save the state of out stream
1087  const std::ios::fmtflags old_flags = out.flags();
1088  const unsigned int old_precision = out.precision (precision);
1089 
1090  out.precision (precision);
1091  if (scientific)
1092  out.setf (std::ios::scientific, std::ios::floatfield);
1093  else
1094  out.setf (std::ios::fixed, std::ios::floatfield);
1095 
1096  if (across)
1097  for (size_type i=0; i<local_size(); ++i)
1098  out << val[i] << ' ';
1099  else
1100  for (size_type i=0; i<local_size(); ++i)
1101  out << val[i] << std::endl;
1102  out << std::endl;
1103 
1104  // reset output format
1105  out.flags (old_flags);
1106  out.precision(old_precision);
1107 
1108  // restore the representation of the
1109  // vector
1110  ierr = VecRestoreArray (vector, &val);
1111  AssertThrow (ierr == 0, ExcPETScError(ierr));
1112 
1113  AssertThrow (out, ExcIO());
1114  }
1115 
1116 
1117 
1118  void
1120  {
1121  const int ierr = VecSwap (vector, v.vector);
1122  AssertThrow (ierr == 0, ExcPETScError(ierr));
1123  }
1124 
1125 
1126 
1127  VectorBase::operator const Vec &() const
1128  {
1129  return vector;
1130  }
1131 
1132 
1133  std::size_t
1135  {
1136  std::size_t mem = sizeof(Vec)+sizeof(last_action)
1139 
1140  // TH: I am relatively sure that PETSc is
1141  // storing the local data in a contiguous
1142  // block without indices:
1143  mem += local_size()*sizeof(PetscScalar);
1144  // assume that PETSc is storing one index
1145  // and one double per ghost element
1146  if (ghosted)
1147  mem += ghost_indices.n_elements()*(sizeof(PetscScalar)+sizeof(int));
1148 
1149  //TODO[TH]: size of constant memory for PETSc?
1150  return mem;
1151  }
1152 
1153 
1154 
1155  void
1156  VectorBase::do_set_add_operation (const size_type n_elements,
1157  const size_type *indices,
1158  const PetscScalar *values,
1159  const bool add_values)
1160  {
1161  ::VectorOperation::values action = (add_values ?
1162  ::VectorOperation::add :
1163  ::VectorOperation::insert);
1164  Assert ((last_action == action)
1165  ||
1166  (last_action == ::VectorOperation::unknown),
1167  internal::VectorReference::ExcWrongMode (action,
1168  last_action));
1169  Assert (!has_ghost_elements(), ExcGhostsPresent());
1170  // VecSetValues complains if we
1171  // come with an empty
1172  // vector. however, it is not a
1173  // collective operation, so we
1174  // can skip the call if necessary
1175  // (unlike the above calls)
1176  if (n_elements != 0)
1177  {
1178 #ifdef PETSC_USE_64BIT_INDICES
1179  std::vector<PetscInt> petsc_ind (n_elements);
1180  for (size_type i=0; i<n_elements; ++i)
1181  petsc_ind[i] = indices[i];
1182  const PetscInt *petsc_indices = &petsc_ind[0];
1183 #else
1184  const int *petsc_indices = (const int *)indices;
1185 #endif
1186 
1187  const InsertMode mode = (add_values ? ADD_VALUES : INSERT_VALUES);
1188  const int ierr
1189  = VecSetValues (vector, n_elements, petsc_indices, values,
1190  mode);
1191  AssertThrow (ierr == 0, ExcPETScError(ierr));
1192  }
1193 
1194  // set the mode here, independent of whether we have actually
1195  // written elements or whether the list was empty
1196  last_action = action;
1197  }
1198 
1199 }
1200 
1201 DEAL_II_NAMESPACE_CLOSE
1202 
1203 #endif // DEAL_II_WITH_PETSC
VectorBase & operator-=(const VectorBase &V)
void sadd(const PetscScalar s, const VectorBase &V)
::ExceptionBase & ExcMessage(std::string arg1)
VectorBase & mult() DEAL_II_DEPRECATED
VectorBase & operator/=(const PetscScalar factor)
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
VectorBase & operator=(const PetscScalar s)
void scale(const VectorBase &scaling_factors)
std::pair< size_type, size_type > local_range() const
real_type lp_norm(const real_type p) const
VectorBase & abs() DEAL_II_DEPRECATED
void clear()
Definition: index_set.h:1223
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:294
static bool is_running_single_threaded()
PetscScalar add_and_dot(const PetscScalar a, const VectorBase &V, const VectorBase &W)
real_type normalize() const DEAL_II_DEPRECATED
std::size_t memory_consumption() const
void equ(const PetscScalar a, const VectorBase &V)
VectorBase & operator*=(const PetscScalar factor)
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
VectorBase & operator+=(const VectorBase &V)
void print(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const bool across=true) const
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)
void ratio(const VectorBase &a, const VectorBase &b) DEAL_II_DEPRECATED
virtual const MPI_Comm & get_mpi_communicator() const
size_type n_elements() const
Definition: index_set.h:1380
#define AssertIsFinite(number)
Definition: exceptions.h:1096
VectorBase & conjugate() DEAL_II_DEPRECATED
bool operator!=(const VectorBase &v) const
void set(const std::vector< size_type > &indices, const std::vector< PetscScalar > &values)