libMesh
petsc_vector.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // C++ includes
21 #include <numeric> // std::iota
22 
23 // Local Includes
24 #include "libmesh/petsc_vector.h"
25 
26 // libMesh includes
27 #include "libmesh/petsc_matrix.h"
28 
29 #ifdef LIBMESH_HAVE_PETSC
30 
31 #include "libmesh/dense_subvector.h"
32 #include "libmesh/dense_vector.h"
33 #include "libmesh/int_range.h"
34 #include "libmesh/petsc_macro.h"
35 
36 // TIMPI includes
37 #include "timpi/op_function.h"
39 #include "timpi/standard_type.h"
40 
41 namespace libMesh
42 {
43 
44 //-----------------------------------------------------------------------
45 // PetscVector members
46 
47 template <typename T>
49 {
50  this->_restore_array();
51  libmesh_assert(this->closed());
52 
53  PetscErrorCode ierr=0;
54  PetscScalar value=0.;
55 
56  ierr = VecSum (_vec, &value);
57  LIBMESH_CHKERR(ierr);
58 
59  return static_cast<T>(value);
60 }
61 
62 
63 template <typename T>
65 {
66  this->_restore_array();
67  libmesh_assert(this->closed());
68 
69  PetscErrorCode ierr=0;
70  PetscReal value=0.;
71 
72  ierr = VecNorm (_vec, NORM_1, &value);
73  LIBMESH_CHKERR(ierr);
74 
75  return static_cast<Real>(value);
76 }
77 
78 
79 
80 template <typename T>
82 {
83  this->_restore_array();
84  libmesh_assert(this->closed());
85 
86  PetscErrorCode ierr=0;
87  PetscReal value=0.;
88 
89  ierr = VecNorm (_vec, NORM_2, &value);
90  LIBMESH_CHKERR(ierr);
91 
92  return static_cast<Real>(value);
93 }
94 
95 
96 
97 
98 template <typename T>
100 {
101  this->_restore_array();
102  libmesh_assert(this->closed());
103 
104  PetscErrorCode ierr=0;
105  PetscReal value=0.;
106 
107  ierr = VecNorm (_vec, NORM_INFINITY, &value);
108  LIBMESH_CHKERR(ierr);
109 
110  return static_cast<Real>(value);
111 }
112 
113 
114 
115 
116 template <typename T>
119 {
120  this->_restore_array();
121  libmesh_assert(this->closed());
122 
123  this->add(1., v);
124 
125  return *this;
126 }
127 
128 
129 
130 template <typename T>
133 {
134  this->_restore_array();
135  libmesh_assert(this->closed());
136 
137  this->add(-1., v);
138 
139  return *this;
140 }
141 
142 
143 
144 template <typename T>
146 {
147  this->_restore_array();
148  libmesh_assert_less (i, size());
149 
150  PetscErrorCode ierr=0;
151  PetscInt i_val = static_cast<PetscInt>(i);
152  PetscScalar petsc_value = PS(value);
153 
154  ierr = VecSetValues (_vec, 1, &i_val, &petsc_value, INSERT_VALUES);
155  LIBMESH_CHKERR(ierr);
156 
157  this->_is_closed = false;
158 }
159 
160 
161 
162 template <typename T>
164 {
165  PetscErrorCode ierr = 0;
166 
167  // VecReciprocal has been in PETSc since at least 2.3.3 days
168  ierr = VecReciprocal(_vec);
169  LIBMESH_CHKERR(ierr);
170 }
171 
172 
173 
174 template <typename T>
176 {
177  PetscErrorCode ierr = 0;
178 
179  // We just call the PETSc VecConjugate
180  ierr = VecConjugate(_vec);
181  LIBMESH_CHKERR(ierr);
182 }
183 
184 
185 
186 template <typename T>
188 {
189  this->_restore_array();
190  libmesh_assert_less (i, size());
191 
192  PetscErrorCode ierr=0;
193  PetscInt i_val = static_cast<PetscInt>(i);
194  PetscScalar petsc_value = PS(value);
195 
196  ierr = VecSetValues (_vec, 1, &i_val, &petsc_value, ADD_VALUES);
197  LIBMESH_CHKERR(ierr);
198 
199  this->_is_closed = false;
200 }
201 
202 
203 
204 template <typename T>
205 void PetscVector<T>::add_vector (const T * v,
206  const std::vector<numeric_index_type> & dof_indices)
207 {
208  // If we aren't adding anything just return
209  if (dof_indices.empty())
210  return;
211 
212  this->_restore_array();
213 
214  PetscErrorCode ierr=0;
215  const PetscInt * i_val = reinterpret_cast<const PetscInt *>(dof_indices.data());
216  const PetscScalar * petsc_value = pPS(v);
217 
218  ierr = VecSetValues (_vec, cast_int<PetscInt>(dof_indices.size()),
219  i_val, petsc_value, ADD_VALUES);
220  LIBMESH_CHKERR(ierr);
221 
222  this->_is_closed = false;
223 }
224 
225 
226 
227 template <typename T>
229  const SparseMatrix<T> & A_in)
230 {
231  this->_restore_array();
232  // Make sure the data passed in are really of Petsc types
233  const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
234  const PetscMatrix<T> * A = cast_ptr<const PetscMatrix<T> *>(&A_in);
235 
236  PetscErrorCode ierr=0;
237 
238  // We shouldn't close() the matrix for you, as that would potentially modify the state of a const object.
239  if (!A->closed())
240  {
241  libmesh_warning("Matrix A must be assembled before calling PetscVector::add_vector(v, A).\n"
242  "Please update your code, as this warning will become an error in a future release.");
243  libmesh_deprecated();
244  const_cast<PetscMatrix<T> *>(A)->close();
245  }
246 
247  // The const_cast<> is not elegant, but it is required since PETSc
248  // expects a non-const Mat.
249  ierr = MatMultAdd(const_cast<PetscMatrix<T> *>(A)->mat(), v->_vec, _vec, _vec);
250  LIBMESH_CHKERR(ierr);
251 }
252 
253 
254 
255 template <typename T>
257  const SparseMatrix<T> & A_in)
258 {
259  this->_restore_array();
260  // Make sure the data passed in are really of Petsc types
261  const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
262  const PetscMatrix<T> * A = cast_ptr<const PetscMatrix<T> *>(&A_in);
263 
264  PetscErrorCode ierr=0;
265 
266  // We shouldn't close() the matrix for you, as that would potentially modify the state of a const object.
267  if (!A->closed())
268  {
269  libmesh_warning("Matrix A must be assembled before calling PetscVector::add_vector_transpose(v, A).\n"
270  "Please update your code, as this warning will become an error in a future release.");
271  libmesh_deprecated();
272  const_cast<PetscMatrix<T> *>(A)->close();
273  }
274 
275  // The const_cast<> is not elegant, but it is required since PETSc
276  // expects a non-const Mat.
277  ierr = MatMultTransposeAdd(const_cast<PetscMatrix<T> *>(A)->mat(), v->_vec, _vec, _vec);
278  LIBMESH_CHKERR(ierr);
279 }
280 
281 
282 #if PETSC_VERSION_LESS_THAN(3,1,0)
283 template <typename T>
285  const SparseMatrix<T> &)
286 {
287  libmesh_error_msg("MatMultHermitianTranspose was introduced in PETSc 3.1.0," \
288  << "No one has made it backwards compatible with older " \
289  << "versions of PETSc so far.");
290 }
291 
292 #else
293 
294 template <typename T>
296  const SparseMatrix<T> & A_in)
297 {
298  this->_restore_array();
299  // Make sure the data passed in are really of Petsc types
300  const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
301  const PetscMatrix<T> * A = cast_ptr<const PetscMatrix<T> *>(&A_in);
302 
303  // We shouldn't close() the matrix for you, as that would potentially modify the state of a const object.
304  if (!A->closed())
305  {
306  libmesh_warning("Matrix A must be assembled before calling PetscVector::add_vector_conjugate_transpose(v, A).\n"
307  "Please update your code, as this warning will become an error in a future release.");
308  libmesh_deprecated();
309  const_cast<PetscMatrix<T> *>(A)->close();
310  }
311 
312  // Store a temporary copy since MatMultHermitianTransposeAdd doesn't seem to work
313  // TODO: Find out why MatMultHermitianTransposeAdd doesn't work, might be a PETSc bug?
314  std::unique_ptr<NumericVector<Number>> this_clone = this->clone();
315 
316  // The const_cast<> is not elegant, but it is required since PETSc
317  // expects a non-const Mat.
318  PetscErrorCode ierr = MatMultHermitianTranspose(const_cast<PetscMatrix<T> *>(A)->mat(), v->_vec, _vec);
319  LIBMESH_CHKERR(ierr);
320 
321  // Add the temporary copy to the matvec result
322  this->add(1., *this_clone);
323 }
324 #endif
325 
326 
327 template <typename T>
328 void PetscVector<T>::add (const T v_in)
329 {
330  this->_get_array(false);
331 
332  for (numeric_index_type i=0; i<_local_size; i++)
333  _values[i] += PetscScalar(v_in);
334 }
335 
336 
337 
338 template <typename T>
340 {
341  this->add (1., v);
342 }
343 
344 
345 
346 template <typename T>
347 void PetscVector<T>::add (const T a_in, const NumericVector<T> & v_in)
348 {
349  this->_restore_array();
350 
351  PetscErrorCode ierr = 0;
352  PetscScalar a = PS(a_in);
353 
354  // Make sure the NumericVector passed in is really a PetscVector
355  const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
356  v->_restore_array();
357 
358  libmesh_assert_equal_to (this->size(), v->size());
359 
360  if (this->type() != GHOSTED)
361  {
362  ierr = VecAXPY(_vec, a, v->_vec);
363  LIBMESH_CHKERR(ierr);
364  }
365  else
366  {
367  Vec loc_vec;
368  Vec v_loc_vec;
369  ierr = VecGhostGetLocalForm (_vec,&loc_vec);
370  LIBMESH_CHKERR(ierr);
371  ierr = VecGhostGetLocalForm (v->_vec,&v_loc_vec);
372  LIBMESH_CHKERR(ierr);
373 
374  ierr = VecAXPY(loc_vec, a, v_loc_vec);
375  LIBMESH_CHKERR(ierr);
376 
377  ierr = VecGhostRestoreLocalForm (v->_vec,&v_loc_vec);
378  LIBMESH_CHKERR(ierr);
379  ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
380  LIBMESH_CHKERR(ierr);
381  }
382 }
383 
384 
385 
386 template <typename T>
387 void PetscVector<T>::insert (const T * v,
388  const std::vector<numeric_index_type> & dof_indices)
389 {
390  if (dof_indices.empty())
391  return;
392 
393  this->_restore_array();
394 
395  PetscErrorCode ierr=0;
396  PetscInt * idx_values = numeric_petsc_cast(dof_indices.data());
397  ierr = VecSetValues (_vec, cast_int<PetscInt>(dof_indices.size()),
398  idx_values, pPS(v), INSERT_VALUES);
399  LIBMESH_CHKERR(ierr);
400 
401  this->_is_closed = false;
402 }
403 
404 
405 
406 template <typename T>
407 void PetscVector<T>::scale (const T factor_in)
408 {
409  this->_restore_array();
410 
411  PetscErrorCode ierr = 0;
412  PetscScalar factor = PS(factor_in);
413 
414  if (this->type() != GHOSTED)
415  {
416  ierr = VecScale(_vec, factor);
417  LIBMESH_CHKERR(ierr);
418  }
419  else
420  {
421  Vec loc_vec;
422  ierr = VecGhostGetLocalForm (_vec,&loc_vec);
423  LIBMESH_CHKERR(ierr);
424 
425  ierr = VecScale(loc_vec, factor);
426  LIBMESH_CHKERR(ierr);
427 
428  ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
429  LIBMESH_CHKERR(ierr);
430  }
431 }
432 
433 template <typename T>
435 {
436  PetscErrorCode ierr = 0;
437 
438  const PetscVector<T> * v_vec = cast_ptr<const PetscVector<T> *>(&v);
439 
440  ierr = VecPointwiseMult(_vec, _vec, v_vec->_vec);
441  LIBMESH_CHKERR(ierr);
442 
443  return *this;
444 }
445 
446 template <typename T>
448 {
449  PetscErrorCode ierr = 0;
450 
451  const PetscVector<T> * v_vec = cast_ptr<const PetscVector<T> *>(&v);
452 
453  ierr = VecPointwiseDivide(_vec, _vec, v_vec->_vec);
454  LIBMESH_CHKERR(ierr);
455 
456  return *this;
457 }
458 
459 template <typename T>
461 {
462  this->_restore_array();
463 
464  PetscErrorCode ierr = 0;
465 
466  if (this->type() != GHOSTED)
467  {
468  ierr = VecAbs(_vec);
469  LIBMESH_CHKERR(ierr);
470  }
471  else
472  {
473  Vec loc_vec;
474  ierr = VecGhostGetLocalForm (_vec,&loc_vec);
475  LIBMESH_CHKERR(ierr);
476 
477  ierr = VecAbs(loc_vec);
478  LIBMESH_CHKERR(ierr);
479 
480  ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
481  LIBMESH_CHKERR(ierr);
482  }
483 }
484 
485 template <typename T>
486 T PetscVector<T>::dot (const NumericVector<T> & v_in) const
487 {
488  this->_restore_array();
489 
490  // Error flag
491  PetscErrorCode ierr = 0;
492 
493  // Return value
494  PetscScalar value=0.;
495 
496  // Make sure the NumericVector passed in is really a PetscVector
497  const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
498 
499  // 2.3.x (at least) style. Untested for previous versions.
500  ierr = VecDot(this->_vec, v->_vec, &value);
501  LIBMESH_CHKERR(ierr);
502 
503  return static_cast<T>(value);
504 }
505 
506 template <typename T>
508 {
509  this->_restore_array();
510 
511  // Error flag
512  PetscErrorCode ierr = 0;
513 
514  // Return value
515  PetscScalar value=0.;
516 
517  // Make sure the NumericVector passed in is really a PetscVector
518  const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
519 
520  // 2.3.x (at least) style. Untested for previous versions.
521  ierr = VecTDot(this->_vec, v->_vec, &value);
522  LIBMESH_CHKERR(ierr);
523 
524  return static_cast<T>(value);
525 }
526 
527 
528 template <typename T>
531 {
532  this->_restore_array();
533  libmesh_assert(this->closed());
534 
535  PetscErrorCode ierr = 0;
536  PetscScalar s = PS(s_in);
537 
538  if (this->size() != 0)
539  {
540  if (this->type() != GHOSTED)
541  {
542  ierr = VecSet(_vec, s);
543  LIBMESH_CHKERR(ierr);
544  }
545  else
546  {
547  Vec loc_vec;
548  ierr = VecGhostGetLocalForm (_vec,&loc_vec);
549  LIBMESH_CHKERR(ierr);
550 
551  ierr = VecSet(loc_vec, s);
552  LIBMESH_CHKERR(ierr);
553 
554  ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
555  LIBMESH_CHKERR(ierr);
556  }
557  }
558 
559  return *this;
560 }
561 
562 
563 
564 template <typename T>
567 {
568  // Make sure the NumericVector passed in is really a PetscVector
569  const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
570 
571  *this = *v;
572 
573  return *this;
574 }
575 
576 
577 
578 template <typename T>
581 {
582  this->_restore_array();
583  v._restore_array();
584 
585  libmesh_assert_equal_to (this->size(), v.size());
586  libmesh_assert_equal_to (this->local_size(), v.local_size());
587  libmesh_assert (v.closed());
588 
589  PetscErrorCode ierr = 0;
590 
591  if (((this->type()==PARALLEL) && (v.type()==GHOSTED)) ||
592  ((this->type()==GHOSTED) && (v.type()==PARALLEL)) ||
593  ((this->type()==GHOSTED) && (v.type()==SERIAL)) ||
594  ((this->type()==SERIAL) && (v.type()==GHOSTED)))
595  {
596  /* Allow assignment of a ghosted to a parallel vector since this
597  causes no difficulty. See discussion in libmesh-devel of
598  June 24, 2010. */
599  ierr = VecCopy (v._vec, this->_vec);
600  LIBMESH_CHKERR(ierr);
601  }
602  else
603  {
604  /* In all other cases, we assert that both vectors are of equal
605  type. */
606  libmesh_assert_equal_to (this->_type, v._type);
607 
608  if (v.size() != 0)
609  {
610  if (this->type() != GHOSTED)
611  {
612  ierr = VecCopy (v._vec, this->_vec);
613  LIBMESH_CHKERR(ierr);
614  }
615  else
616  {
617  Vec loc_vec;
618  Vec v_loc_vec;
619  ierr = VecGhostGetLocalForm (_vec,&loc_vec);
620  LIBMESH_CHKERR(ierr);
621  ierr = VecGhostGetLocalForm (v._vec,&v_loc_vec);
622  LIBMESH_CHKERR(ierr);
623 
624  ierr = VecCopy (v_loc_vec, loc_vec);
625  LIBMESH_CHKERR(ierr);
626 
627  ierr = VecGhostRestoreLocalForm (v._vec,&v_loc_vec);
628  LIBMESH_CHKERR(ierr);
629  ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
630  LIBMESH_CHKERR(ierr);
631  }
632  }
633  }
634 
635  close();
636 
637  return *this;
638 }
639 
640 
641 
642 template <typename T>
644 PetscVector<T>::operator = (const std::vector<T> & v)
645 {
646  this->_get_array(false);
647 
652  if (this->size() == v.size())
653  {
654  numeric_index_type first = first_local_index();
655  numeric_index_type last = last_local_index();
656  for (numeric_index_type i=0; i<last-first; i++)
657  _values[i] = PS(v[first + i]);
658  }
659 
664  else
665  {
666  for (numeric_index_type i=0; i<_local_size; i++)
667  _values[i] = PS(v[i]);
668  }
669 
670  // Make sure ghost dofs are up to date
671  if (this->type() == GHOSTED)
672  this->close();
673 
674  return *this;
675 }
676 
677 
678 
679 template <typename T>
681 {
682  this->_restore_array();
683 
684  // Make sure the NumericVector passed in is really a PetscVector
685  PetscVector<T> * v_local = cast_ptr<PetscVector<T> *>(&v_local_in);
686 
687  libmesh_assert(v_local);
688  libmesh_assert_equal_to (v_local->size(), this->size());
689 
690  PetscErrorCode ierr = 0;
691  const PetscInt n = this->size();
692 
693  IS is;
694  VecScatter scatter;
695 
696  // Create idx, idx[i] = i;
697  std::vector<PetscInt> idx(n);
698  std::iota (idx.begin(), idx.end(), 0);
699 
700  // Create the index set & scatter object
701  ierr = ISCreateLibMesh(this->comm().get(), n, idx.data(), PETSC_USE_POINTER, &is);
702  LIBMESH_CHKERR(ierr);
703 
704  ierr = LibMeshVecScatterCreate(_vec, is,
705  v_local->_vec, is,
706  &scatter);
707  LIBMESH_CHKERR(ierr);
708 
709  // Perform the scatter
710  ierr = VecScatterBegin(scatter, _vec, v_local->_vec,
711  INSERT_VALUES, SCATTER_FORWARD);
712  LIBMESH_CHKERR(ierr);
713 
714  ierr = VecScatterEnd (scatter, _vec, v_local->_vec,
715  INSERT_VALUES, SCATTER_FORWARD);
716  LIBMESH_CHKERR(ierr);
717 
718  // Clean up
719  ierr = LibMeshISDestroy (&is);
720  LIBMESH_CHKERR(ierr);
721 
722  ierr = LibMeshVecScatterDestroy(&scatter);
723  LIBMESH_CHKERR(ierr);
724 
725  // Make sure ghost dofs are up to date
726  if (v_local->type() == GHOSTED)
727  v_local->close();
728 }
729 
730 
731 
732 template <typename T>
734  const std::vector<numeric_index_type> & send_list) const
735 {
736  // FIXME: Workaround for a strange bug at large-scale.
737  // If we have ghosting, PETSc lets us just copy the solution, and
738  // doing so avoids a segfault?
739  if (v_local_in.type() == GHOSTED &&
740  this->type() == PARALLEL)
741  {
742  v_local_in = *this;
743  return;
744  }
745 
746  // Normal code path begins here
747 
748  this->_restore_array();
749 
750  // Make sure the NumericVector passed in is really a PetscVector
751  PetscVector<T> * v_local = cast_ptr<PetscVector<T> *>(&v_local_in);
752 
753  libmesh_assert(v_local);
754  libmesh_assert_equal_to (v_local->size(), this->size());
755  libmesh_assert_less_equal (send_list.size(), v_local->size());
756 
757  PetscErrorCode ierr=0;
758  const numeric_index_type n_sl =
759  cast_int<numeric_index_type>(send_list.size());
760 
761  IS is;
762  VecScatter scatter;
763 
764  std::vector<PetscInt> idx(n_sl + this->local_size());
765 
766  for (numeric_index_type i=0; i<n_sl; i++)
767  idx[i] = static_cast<PetscInt>(send_list[i]);
768  for (auto i : IntRange<numeric_index_type>(0, this->local_size()))
769  idx[n_sl+i] = i + this->first_local_index();
770 
771  // Create the index set & scatter object
772  PetscInt * idxptr = idx.empty() ? nullptr : idx.data();
773  ierr = ISCreateLibMesh(this->comm().get(), n_sl+this->local_size(),
774  idxptr, PETSC_USE_POINTER, &is);
775  LIBMESH_CHKERR(ierr);
776 
777  ierr = LibMeshVecScatterCreate(_vec, is,
778  v_local->_vec, is,
779  &scatter);
780  LIBMESH_CHKERR(ierr);
781 
782 
783  // Perform the scatter
784  ierr = VecScatterBegin(scatter, _vec, v_local->_vec,
785  INSERT_VALUES, SCATTER_FORWARD);
786  LIBMESH_CHKERR(ierr);
787 
788  ierr = VecScatterEnd (scatter, _vec, v_local->_vec,
789  INSERT_VALUES, SCATTER_FORWARD);
790  LIBMESH_CHKERR(ierr);
791 
792  // Clean up
793  ierr = LibMeshISDestroy (&is);
794  LIBMESH_CHKERR(ierr);
795 
796  ierr = LibMeshVecScatterDestroy(&scatter);
797  LIBMESH_CHKERR(ierr);
798 
799  // Make sure ghost dofs are up to date
800  if (v_local->type() == GHOSTED)
801  v_local->close();
802 }
803 
804 
805 
806 template <typename T>
807 void PetscVector<T>::localize (std::vector<T> & v_local,
808  const std::vector<numeric_index_type> & indices) const
809 {
810  // Error code used to check the status of all PETSc function calls.
811  PetscErrorCode ierr = 0;
812 
813  // Create a sequential destination Vec with the right number of entries on each proc.
814  Vec dest;
815  ierr = VecCreateSeq(PETSC_COMM_SELF, cast_int<PetscInt>(indices.size()), &dest);
816  LIBMESH_CHKERR(ierr);
817 
818  // Create an IS using the libmesh routine. PETSc does not own the
819  // IS memory in this case, it is automatically cleaned up by the
820  // std::vector destructor.
821  PetscInt * idxptr =
822  indices.empty() ? nullptr : numeric_petsc_cast(indices.data());
823  IS is;
824  ierr = ISCreateLibMesh(this->comm().get(), cast_int<PetscInt>(indices.size()), idxptr,
826  LIBMESH_CHKERR(ierr);
827 
828  // Create the VecScatter object. "PETSC_NULL" means "use the identity IS".
829  VecScatter scatter;
830  ierr = LibMeshVecScatterCreate(_vec,
831  /*src is=*/is,
832  /*dest vec=*/dest,
833  /*dest is=*/PETSC_NULL,
834  &scatter);
835  LIBMESH_CHKERR(ierr);
836 
837  // Do the scatter
838  ierr = VecScatterBegin(scatter, _vec, dest, INSERT_VALUES, SCATTER_FORWARD);
839  LIBMESH_CHKERR(ierr);
840 
841  ierr = VecScatterEnd(scatter, _vec, dest, INSERT_VALUES, SCATTER_FORWARD);
842  LIBMESH_CHKERR(ierr);
843 
844  // Get access to the values stored in dest.
845  PetscScalar * values;
846  ierr = VecGetArray (dest, &values);
847  LIBMESH_CHKERR(ierr);
848 
849  // Store values into the provided v_local. Make sure there is enough
850  // space reserved and then clear out any existing entries before
851  // inserting.
852  v_local.reserve(indices.size());
853  v_local.clear();
854  v_local.insert(v_local.begin(), values, values+indices.size());
855 
856  // We are done using it, so restore the array.
857  ierr = VecRestoreArray (dest, &values);
858  LIBMESH_CHKERR(ierr);
859 
860  // Clean up PETSc data structures.
861  ierr = LibMeshVecScatterDestroy(&scatter);
862  LIBMESH_CHKERR(ierr);
863 
864  ierr = LibMeshISDestroy (&is);
865  LIBMESH_CHKERR(ierr);
866 
867  ierr = LibMeshVecDestroy(&dest);
868  LIBMESH_CHKERR(ierr);
869 }
870 
871 
872 
873 template <typename T>
874 void PetscVector<T>::localize (const numeric_index_type first_local_idx,
875  const numeric_index_type last_local_idx,
876  const std::vector<numeric_index_type> & send_list)
877 {
878  this->_restore_array();
879 
880  libmesh_assert_less_equal (send_list.size(), this->size());
881  libmesh_assert_less_equal (last_local_idx+1, this->size());
882 
883  const numeric_index_type my_size = this->size();
884  const numeric_index_type my_local_size = (last_local_idx - first_local_idx + 1);
885  PetscErrorCode ierr=0;
886 
887  // Don't bother for serial cases
888  // if ((first_local_idx == 0) &&
889  // (my_local_size == my_size))
890  // But we do need to stay in sync for degenerate cases
891  if (this->n_processors() == 1)
892  return;
893 
894 
895  // Build a parallel vector, initialize it with the local
896  // parts of (*this)
897  PetscVector<T> parallel_vec(this->comm(), PARALLEL);
898 
899  parallel_vec.init (my_size, my_local_size, true, PARALLEL);
900 
901 
902  // Copy part of *this into the parallel_vec
903  {
904  IS is;
905  VecScatter scatter;
906 
907  // Create idx, idx[i] = i+first_local_idx;
908  std::vector<PetscInt> idx(my_local_size);
909  std::iota (idx.begin(), idx.end(), first_local_idx);
910 
911  // Create the index set & scatter object
912  ierr = ISCreateLibMesh(this->comm().get(), my_local_size,
913  my_local_size ? idx.data() : nullptr, PETSC_USE_POINTER, &is);
914  LIBMESH_CHKERR(ierr);
915 
916  ierr = LibMeshVecScatterCreate(_vec, is,
917  parallel_vec._vec, is,
918  &scatter);
919  LIBMESH_CHKERR(ierr);
920 
921  ierr = VecScatterBegin(scatter, _vec, parallel_vec._vec,
922  INSERT_VALUES, SCATTER_FORWARD);
923  LIBMESH_CHKERR(ierr);
924 
925  ierr = VecScatterEnd (scatter, _vec, parallel_vec._vec,
926  INSERT_VALUES, SCATTER_FORWARD);
927  LIBMESH_CHKERR(ierr);
928 
929  // Clean up
930  ierr = LibMeshISDestroy (&is);
931  LIBMESH_CHKERR(ierr);
932 
933  ierr = LibMeshVecScatterDestroy(&scatter);
934  LIBMESH_CHKERR(ierr);
935  }
936 
937  // localize like normal
938  parallel_vec.close();
939  parallel_vec.localize (*this, send_list);
940  this->close();
941 }
942 
943 
944 
945 template <typename T>
946 void PetscVector<T>::localize (std::vector<T> & v_local) const
947 {
948  this->_restore_array();
949 
950  // This function must be run on all processors at once
951  parallel_object_only();
952 
953  PetscErrorCode ierr=0;
954  const PetscInt n = this->size();
955  const PetscInt nl = this->local_size();
956  PetscScalar * values;
957 
958  v_local.clear();
959  v_local.resize(n, 0.);
960 
961  ierr = VecGetArray (_vec, &values);
962  LIBMESH_CHKERR(ierr);
963 
964  numeric_index_type ioff = first_local_index();
965 
966  for (PetscInt i=0; i<nl; i++)
967  v_local[i+ioff] = static_cast<T>(values[i]);
968 
969  ierr = VecRestoreArray (_vec, &values);
970  LIBMESH_CHKERR(ierr);
971 
972  this->comm().sum(v_local);
973 }
974 
975 
976 
977 // Full specialization for Real datatypes
978 #ifdef LIBMESH_USE_REAL_NUMBERS
979 
980 template <>
981 void PetscVector<Real>::localize_to_one (std::vector<Real> & v_local,
982  const processor_id_type
983  timpi_mpi_var(pid)) const
984 {
985  this->_restore_array();
986 
987  PetscErrorCode ierr=0;
988  const PetscInt n = size();
989  PetscScalar * values;
990 
991  // only one processor
992  if (n_processors() == 1)
993  {
994  v_local.resize(n);
995 
996  ierr = VecGetArray (_vec, &values);
997  LIBMESH_CHKERR(ierr);
998 
999  for (PetscInt i=0; i<n; i++)
1000  v_local[i] = static_cast<Real>(values[i]);
1001 
1002  ierr = VecRestoreArray (_vec, &values);
1003  LIBMESH_CHKERR(ierr);
1004  }
1005 
1006  // otherwise multiple processors
1007 #ifdef LIBMESH_HAVE_MPI
1008  else
1009  {
1010  if (pid == 0) // optimized version for localizing to 0
1011  {
1012  Vec vout;
1013  VecScatter ctx;
1014 
1015  ierr = VecScatterCreateToZero(_vec, &ctx, &vout);
1016  LIBMESH_CHKERR(ierr);
1017 
1018  ierr = VecScatterBegin(ctx, _vec, vout, INSERT_VALUES, SCATTER_FORWARD);
1019  LIBMESH_CHKERR(ierr);
1020  ierr = VecScatterEnd(ctx, _vec, vout, INSERT_VALUES, SCATTER_FORWARD);
1021  LIBMESH_CHKERR(ierr);
1022 
1023  if (processor_id() == 0)
1024  {
1025  v_local.resize(n);
1026 
1027  ierr = VecGetArray (vout, &values);
1028  LIBMESH_CHKERR(ierr);
1029 
1030  for (PetscInt i=0; i<n; i++)
1031  v_local[i] = static_cast<Real>(values[i]);
1032 
1033  ierr = VecRestoreArray (vout, &values);
1034  LIBMESH_CHKERR(ierr);
1035  }
1036 
1037  ierr = LibMeshVecScatterDestroy(&ctx);
1038  LIBMESH_CHKERR(ierr);
1039  ierr = LibMeshVecDestroy(&vout);
1040  LIBMESH_CHKERR(ierr);
1041 
1042  }
1043  else
1044  {
1045  v_local.resize(n);
1046 
1047  numeric_index_type ioff = this->first_local_index();
1048  std::vector<Real> local_values (n, 0.);
1049 
1050  {
1051  ierr = VecGetArray (_vec, &values);
1052  LIBMESH_CHKERR(ierr);
1053 
1054  const PetscInt nl = local_size();
1055  for (PetscInt i=0; i<nl; i++)
1056  local_values[i+ioff] = static_cast<Real>(values[i]);
1057 
1058  ierr = VecRestoreArray (_vec, &values);
1059  LIBMESH_CHKERR(ierr);
1060  }
1061 
1062 
1063  MPI_Reduce (local_values.data(), v_local.data(), n,
1064  Parallel::StandardType<Real>(),
1065  Parallel::OpFunction<Real>::sum(), pid,
1066  this->comm().get());
1067  }
1068  }
1069 #endif // LIBMESH_HAVE_MPI
1070 }
1071 
1072 #endif
1073 
1074 
1075 // Full specialization for Complex datatypes
1076 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
1077 
1078 template <>
1079 void PetscVector<Complex>::localize_to_one (std::vector<Complex> & v_local,
1080  const processor_id_type pid) const
1081 {
1082  this->_restore_array();
1083 
1084  PetscErrorCode ierr=0;
1085  const PetscInt n = size();
1086  const PetscInt nl = local_size();
1087  PetscScalar * values;
1088 
1089 
1090  v_local.resize(n);
1091 
1092 
1093  for (PetscInt i=0; i<n; i++)
1094  v_local[i] = 0.;
1095 
1096  // only one processor
1097  if (n_processors() == 1)
1098  {
1099  ierr = VecGetArray (_vec, &values);
1100  LIBMESH_CHKERR(ierr);
1101 
1102  for (PetscInt i=0; i<n; i++)
1103  v_local[i] = static_cast<Complex>(values[i]);
1104 
1105  ierr = VecRestoreArray (_vec, &values);
1106  LIBMESH_CHKERR(ierr);
1107  }
1108 
1109  // otherwise multiple processors
1110  else
1111  {
1112  numeric_index_type ioff = this->first_local_index();
1113 
1114  /* in here the local values are stored, acting as send buffer for MPI
1115  * initialize to zero, since we collect using MPI_SUM
1116  */
1117  std::vector<Real> real_local_values(n, 0.);
1118  std::vector<Real> imag_local_values(n, 0.);
1119 
1120  {
1121  ierr = VecGetArray (_vec, &values);
1122  LIBMESH_CHKERR(ierr);
1123 
1124  // provide my local share to the real and imag buffers
1125  for (PetscInt i=0; i<nl; i++)
1126  {
1127  real_local_values[i+ioff] = static_cast<Complex>(values[i]).real();
1128  imag_local_values[i+ioff] = static_cast<Complex>(values[i]).imag();
1129  }
1130 
1131  ierr = VecRestoreArray (_vec, &values);
1132  LIBMESH_CHKERR(ierr);
1133  }
1134 
1135  /* have buffers of the real and imaginary part of v_local.
1136  * Once MPI_Reduce() collected all the real and imaginary
1137  * parts in these std::vector<Real>, the values can be
1138  * copied to v_local
1139  */
1140  std::vector<Real> real_v_local(n);
1141  std::vector<Real> imag_v_local(n);
1142 
1143  // collect entries from other proc's in real_v_local, imag_v_local
1144  MPI_Reduce (real_local_values.data(),
1145  real_v_local.data(), n,
1146  Parallel::StandardType<Real>(),
1147  Parallel::OpFunction<Real>::sum(),
1148  pid, this->comm().get());
1149 
1150  MPI_Reduce (imag_local_values.data(),
1151  imag_v_local.data(), n,
1152  Parallel::StandardType<Real>(),
1153  Parallel::OpFunction<Real>::sum(),
1154  pid, this->comm().get());
1155 
1156  // copy real_v_local and imag_v_local to v_local
1157  for (PetscInt i=0; i<n; i++)
1158  v_local[i] = Complex(real_v_local[i], imag_v_local[i]);
1159  }
1160 }
1161 
1162 #endif
1163 
1164 
1165 
1166 template <typename T>
1168  const NumericVector<T> & vec2)
1169 {
1170  this->_restore_array();
1171 
1172  PetscErrorCode ierr = 0;
1173 
1174  // Convert arguments to PetscVector*.
1175  const PetscVector<T> * vec1_petsc = cast_ptr<const PetscVector<T> *>(&vec1);
1176  const PetscVector<T> * vec2_petsc = cast_ptr<const PetscVector<T> *>(&vec2);
1177 
1178  // Call PETSc function.
1179 
1180  if (this->type() != GHOSTED)
1181  {
1182  ierr = VecPointwiseMult(this->vec(),
1183  const_cast<PetscVector<T> *>(vec1_petsc)->vec(),
1184  const_cast<PetscVector<T> *>(vec2_petsc)->vec());
1185  LIBMESH_CHKERR(ierr);
1186  }
1187  else
1188  {
1189  Vec loc_vec;
1190  Vec v1_loc_vec;
1191  Vec v2_loc_vec;
1192  ierr = VecGhostGetLocalForm (_vec,&loc_vec);
1193  LIBMESH_CHKERR(ierr);
1194  ierr = VecGhostGetLocalForm (const_cast<PetscVector<T> *>(vec1_petsc)->vec(),&v1_loc_vec);
1195  LIBMESH_CHKERR(ierr);
1196  ierr = VecGhostGetLocalForm (const_cast<PetscVector<T> *>(vec2_petsc)->vec(),&v2_loc_vec);
1197  LIBMESH_CHKERR(ierr);
1198 
1199  ierr = VecPointwiseMult(loc_vec,v1_loc_vec,v2_loc_vec);
1200  LIBMESH_CHKERR(ierr);
1201 
1202  ierr = VecGhostRestoreLocalForm (const_cast<PetscVector<T> *>(vec1_petsc)->vec(),&v1_loc_vec);
1203  LIBMESH_CHKERR(ierr);
1204  ierr = VecGhostRestoreLocalForm (const_cast<PetscVector<T> *>(vec2_petsc)->vec(),&v2_loc_vec);
1205  LIBMESH_CHKERR(ierr);
1206  ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
1207  LIBMESH_CHKERR(ierr);
1208  }
1209 }
1210 
1211 
1212 
1213 template <typename T>
1214 void PetscVector<T>::print_matlab (const std::string & name) const
1215 {
1216  this->_restore_array();
1217  libmesh_assert (this->closed());
1218 
1219  PetscErrorCode ierr=0;
1220  PetscViewer petsc_viewer;
1221 
1222 
1223  ierr = PetscViewerCreate (this->comm().get(),
1224  &petsc_viewer);
1225  LIBMESH_CHKERR(ierr);
1226 
1231  if (name != "")
1232  {
1233  ierr = PetscViewerASCIIOpen( this->comm().get(),
1234  name.c_str(),
1235  &petsc_viewer);
1236  LIBMESH_CHKERR(ierr);
1237 
1238 #if PETSC_VERSION_LESS_THAN(3,7,0)
1239  ierr = PetscViewerSetFormat (petsc_viewer,
1240  PETSC_VIEWER_ASCII_MATLAB);
1241 #else
1242  ierr = PetscViewerPushFormat (petsc_viewer,
1243  PETSC_VIEWER_ASCII_MATLAB);
1244 #endif
1245 
1246  LIBMESH_CHKERR(ierr);
1247 
1248  ierr = VecView (_vec, petsc_viewer);
1249  LIBMESH_CHKERR(ierr);
1250  }
1251 
1255  else
1256  {
1257 
1258 #if PETSC_VERSION_LESS_THAN(3,7,0)
1259  ierr = PetscViewerSetFormat (PETSC_VIEWER_STDOUT_WORLD,
1260  PETSC_VIEWER_ASCII_MATLAB);
1261 #else
1262  ierr = PetscViewerPushFormat (PETSC_VIEWER_STDOUT_WORLD,
1263  PETSC_VIEWER_ASCII_MATLAB);
1264 #endif
1265 
1266  LIBMESH_CHKERR(ierr);
1267 
1268  ierr = VecView (_vec, PETSC_VIEWER_STDOUT_WORLD);
1269  LIBMESH_CHKERR(ierr);
1270  }
1271 
1272 
1276  ierr = LibMeshPetscViewerDestroy (&petsc_viewer);
1277  LIBMESH_CHKERR(ierr);
1278 }
1279 
1280 
1281 
1282 
1283 
1284 template <typename T>
1286  const std::vector<numeric_index_type> & rows) const
1287 {
1288  this->_restore_array();
1289 
1290  // PETSc data structures
1291  IS parent_is, subvector_is;
1292  VecScatter scatter;
1293  PetscErrorCode ierr = 0;
1294 
1295  // Make sure the passed in subvector is really a PetscVector
1296  PetscVector<T> * petsc_subvector = cast_ptr<PetscVector<T> *>(&subvector);
1297 
1298  // If the petsc_subvector is already initialized, we assume that the
1299  // user has already allocated the *correct* amount of space for it.
1300  // If not, we use the appropriate PETSc routines to initialize it.
1301  if (!petsc_subvector->initialized())
1302  {
1303  // Initialize the petsc_subvector to have enough space to hold
1304  // the entries which will be scattered into it. Note: such an
1305  // init() function (where we let PETSc decide the number of local
1306  // entries) is not currently offered by the PetscVector
1307  // class. Should we differentiate here between sequential and
1308  // parallel vector creation based on this->n_processors() ?
1309  ierr = VecCreateMPI(this->comm().get(),
1310  PETSC_DECIDE, // n_local
1311  cast_int<PetscInt>(rows.size()), // n_global
1312  &(petsc_subvector->_vec));
1313  LIBMESH_CHKERR(ierr);
1314 
1315  ierr = VecSetFromOptions (petsc_subvector->_vec);
1316  LIBMESH_CHKERR(ierr);
1317 
1318  // Mark the subvector as initialized
1319  petsc_subvector->_is_initialized = true;
1320  }
1321  else
1322  {
1323  petsc_subvector->_restore_array();
1324  }
1325 
1326  // Use iota to fill an array with entries [0,1,2,3,4,...rows.size()]
1327  std::vector<PetscInt> idx(rows.size());
1328  std::iota (idx.begin(), idx.end(), 0);
1329 
1330  // Construct index sets
1331  ierr = ISCreateLibMesh(this->comm().get(),
1332  cast_int<PetscInt>(rows.size()),
1333  numeric_petsc_cast(rows.data()),
1335  &parent_is);
1336  LIBMESH_CHKERR(ierr);
1337 
1338  ierr = ISCreateLibMesh(this->comm().get(),
1339  cast_int<PetscInt>(rows.size()),
1340  idx.data(),
1342  &subvector_is);
1343  LIBMESH_CHKERR(ierr);
1344 
1345  // Construct the scatter object
1346  ierr = LibMeshVecScatterCreate(this->_vec,
1347  parent_is,
1348  petsc_subvector->_vec,
1349  subvector_is,
1350  &scatter); LIBMESH_CHKERR(ierr);
1351 
1352  // Actually perform the scatter
1353  ierr = VecScatterBegin(scatter,
1354  this->_vec,
1355  petsc_subvector->_vec,
1356  INSERT_VALUES,
1357  SCATTER_FORWARD); LIBMESH_CHKERR(ierr);
1358 
1359  ierr = VecScatterEnd(scatter,
1360  this->_vec,
1361  petsc_subvector->_vec,
1362  INSERT_VALUES,
1363  SCATTER_FORWARD); LIBMESH_CHKERR(ierr);
1364 
1365  // Clean up
1366  ierr = LibMeshISDestroy(&parent_is); LIBMESH_CHKERR(ierr);
1367  ierr = LibMeshISDestroy(&subvector_is); LIBMESH_CHKERR(ierr);
1368  ierr = LibMeshVecScatterDestroy(&scatter); LIBMESH_CHKERR(ierr);
1369 }
1370 
1371 
1372 
1373 template <typename T>
1374 void PetscVector<T>::_get_array(bool read_only) const
1375 {
1376  libmesh_assert (this->initialized());
1377 
1378 #ifdef LIBMESH_HAVE_CXX11_THREAD
1379  std::atomic_thread_fence(std::memory_order_acquire);
1380 #else
1381  Threads::spin_mutex::scoped_lock lock(_petsc_vector_mutex);
1382 #endif
1383 
1384  // If the values have already been retrieved and we're currently
1385  // trying to get a non-read only view (ie read/write) and the
1386  // values are currently read only... then we need to restore
1387  // the array first... and then retrieve it again.
1388  if (_array_is_present && !read_only && _values_read_only)
1389  _restore_array();
1390 
1391  // If we already have a read/write array - and we're trying
1392  // to get a read only array - let's just use the read write
1393  if (_array_is_present && read_only && !_values_read_only)
1394  _read_only_values = _values;
1395 
1396  if (!_array_is_present)
1397  {
1398 #ifdef LIBMESH_HAVE_CXX11_THREAD
1399  std::lock_guard<std::mutex> lock(_petsc_vector_mutex);
1400 #endif
1401  if (!_array_is_present)
1402  {
1403  PetscErrorCode ierr=0;
1404  if (this->type() != GHOSTED)
1405  {
1406 #if PETSC_VERSION_LESS_THAN(3,2,0)
1407  // Vec{Get,Restore}ArrayRead were introduced in PETSc 3.2.0. If you
1408  // have an older PETSc than that, we'll do an ugly
1409  // const_cast and call VecGetArray() instead.
1410  ierr = VecGetArray(_vec, const_cast<PetscScalar **>(&_values));
1411  _values_read_only = false;
1412 #else
1413  if (read_only)
1414  {
1415  ierr = VecGetArrayRead(_vec, &_read_only_values);
1416  _values_read_only = true;
1417  }
1418  else
1419  {
1420  ierr = VecGetArray(_vec, &_values);
1421  _values_read_only = false;
1422  }
1423 #endif
1424  LIBMESH_CHKERR(ierr);
1425  _local_size = this->local_size();
1426  }
1427  else
1428  {
1429  ierr = VecGhostGetLocalForm (_vec,&_local_form);
1430  LIBMESH_CHKERR(ierr);
1431 
1432 #if PETSC_VERSION_LESS_THAN(3,2,0)
1433  // Vec{Get,Restore}ArrayRead were introduced in PETSc 3.2.0. If you
1434  // have an older PETSc than that, we'll do an ugly
1435  // const_cast and call VecGetArray() instead.
1436  ierr = VecGetArray(_local_form, const_cast<PetscScalar **>(&_values));
1437  _values_read_only = false;
1438 #else
1439  if (read_only)
1440  {
1441  ierr = VecGetArrayRead(_local_form, &_read_only_values);
1442  _values_read_only = true;
1443  }
1444  else
1445  {
1446  ierr = VecGetArray(_local_form, &_values);
1447  _values_read_only = false;
1448  }
1449 #endif
1450  LIBMESH_CHKERR(ierr);
1451 
1452  PetscInt my_local_size = 0;
1453  ierr = VecGetLocalSize(_local_form, &my_local_size);
1454  LIBMESH_CHKERR(ierr);
1455  _local_size = static_cast<numeric_index_type>(my_local_size);
1456  }
1457 
1458  { // cache ownership range
1459  PetscInt petsc_first=0, petsc_last=0;
1460  ierr = VecGetOwnershipRange (_vec, &petsc_first, &petsc_last);
1461  LIBMESH_CHKERR(ierr);
1462  _first = static_cast<numeric_index_type>(petsc_first);
1463  _last = static_cast<numeric_index_type>(petsc_last);
1464  }
1465 #ifdef LIBMESH_HAVE_CXX11_THREAD
1466  std::atomic_thread_fence(std::memory_order_release);
1467  _array_is_present.store(true, std::memory_order_relaxed);
1468 #else
1469  _array_is_present = true;
1470 #endif
1471  }
1472  }
1473 }
1474 
1475 
1476 
1477 template <typename T>
1479 {
1480  if (_values_manually_retrieved)
1481  libmesh_error_msg("PetscVector values were manually retrieved but are being automatically restored!");
1482 
1483 #ifdef LIBMESH_HAVE_CXX11_THREAD
1484  std::atomic_thread_fence(std::memory_order_acquire);
1485 #else
1486  Threads::spin_mutex::scoped_lock lock(_petsc_vector_mutex);
1487 #endif
1488 
1489  libmesh_assert (this->initialized());
1490  if (_array_is_present)
1491  {
1492 #ifdef LIBMESH_HAVE_CXX11_THREAD
1493  std::lock_guard<std::mutex> lock(_petsc_vector_mutex);
1494 #endif
1495 
1496  if (_array_is_present)
1497  {
1498  PetscErrorCode ierr=0;
1499  if (this->type() != GHOSTED)
1500  {
1501 #if PETSC_VERSION_LESS_THAN(3,2,0)
1502  // Vec{Get,Restore}ArrayRead were introduced in PETSc 3.2.0. If you
1503  // have an older PETSc than that, we'll do an ugly
1504  // const_cast and call VecRestoreArray() instead.
1505  ierr = VecRestoreArray (_vec, const_cast<PetscScalar **>(&_values));
1506 #else
1507  if (_values_read_only)
1508  ierr = VecRestoreArrayRead (_vec, &_read_only_values);
1509  else
1510  ierr = VecRestoreArray (_vec, &_values);
1511 #endif
1512 
1513  LIBMESH_CHKERR(ierr);
1514  _values = nullptr;
1515  }
1516  else
1517  {
1518 #if PETSC_VERSION_LESS_THAN(3,2,0)
1519  // Vec{Get,Restore}ArrayRead were introduced in PETSc 3.2.0. If you
1520  // have an older PETSc than that, we'll do an ugly
1521  // const_cast and call VecRestoreArray() instead.
1522  ierr = VecRestoreArray (_local_form, const_cast<PetscScalar **>(&_values));
1523 #else
1524  if (_values_read_only)
1525  ierr = VecRestoreArrayRead (_local_form, &_read_only_values);
1526  else
1527  ierr = VecRestoreArray (_local_form, &_values);
1528 #endif
1529  LIBMESH_CHKERR(ierr);
1530  _values = nullptr;
1531  ierr = VecGhostRestoreLocalForm (_vec,&_local_form);
1532  LIBMESH_CHKERR(ierr);
1533  _local_form = nullptr;
1534  _local_size = 0;
1535  }
1536 #ifdef LIBMESH_HAVE_CXX11_THREAD
1537  std::atomic_thread_fence(std::memory_order_release);
1538  _array_is_present.store(false, std::memory_order_relaxed);
1539 #else
1540  _array_is_present = false;
1541 #endif
1542  }
1543  }
1544 }
1545 
1546 
1547 
1548 
1549 //------------------------------------------------------------------
1550 // Explicit instantiations
1551 template class PetscVector<Number>;
1552 
1553 } // namespace libMesh
1554 
1555 
1556 
1557 #endif // #ifdef LIBMESH_HAVE_PETSC
libMesh::PetscVector::pointwise_mult
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2) override
Computes (summation not implied) i.e.
Definition: petsc_vector.C:1167
libMesh::PetscVector::conjugate
virtual void conjugate() override
Negates the imaginary component of each entry in the vector.
Definition: petsc_vector.C:175
libMesh::PetscVector::abs
virtual void abs() override
Sets for each entry in the vector.
Definition: petsc_vector.C:460
libMesh::PetscVector::sum
virtual T sum() const override
Definition: petsc_vector.C:48
libMesh::PetscVector::create_subvector
virtual void create_subvector(NumericVector< T > &subvector, const std::vector< numeric_index_type > &rows) const override
Fills in subvector from this vector using the indices in rows.
Definition: petsc_vector.C:1285
libMesh::PARALLEL
Definition: enum_parallel_type.h:36
libMesh::PetscVector::set
virtual void set(const numeric_index_type i, const T value) override
Sets v(i) = value.
Definition: petsc_vector.C:145
libMesh::SERIAL
Definition: enum_parallel_type.h:35
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::PetscVector::indefinite_dot
T indefinite_dot(const NumericVector< T > &v) const
Definition: petsc_vector.C:507
op_function.h
libMesh::pPS
PetscScalar * pPS(T *ptr)
Definition: petsc_macro.h:174
libMesh::PetscVector::init
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType type=AUTOMATIC) override
Change the dimension of the vector to n.
Definition: petsc_vector.h:640
libMesh::GHOSTED
Definition: enum_parallel_type.h:37
libMesh::PetscVector::add_vector
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices) override
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i].
Definition: petsc_vector.C:205
libMesh::NumericVector::_is_initialized
bool _is_initialized
true once init() has been called.
Definition: numeric_vector.h:732
libMesh::PetscVector::reciprocal
virtual void reciprocal() override
Computes the component-wise reciprocal, .
Definition: petsc_vector.C:163
libMesh::ierr
ierr
Definition: petsc_dm_wrapper.C:72
libMesh::SparseMatrix
Generic sparse matrix.
Definition: vector_fe_ex5.C:45
libMesh::NumericVector::closed
virtual bool closed() const
Definition: numeric_vector.h:171
libMesh::PetscVector::operator*=
virtual NumericVector< T > & operator*=(const NumericVector< T > &v) override
Computes the component-wise multiplication of this vector's entries by another's, .
Definition: petsc_vector.C:434
libMesh::NumericVector
Provides a uniform interface to vector storage schemes for different linear algebra libraries.
Definition: vector_fe_ex5.C:43
libMesh::PetscVector::operator=
PetscVector< T > & operator=(const PetscVector< T > &v)
Copy assignment operator.
Definition: petsc_vector.C:580
PETSC_USE_POINTER
Definition: petsc_macro.h:103
libMesh::NumericVector::_type
ParallelType _type
Type of vector.
Definition: numeric_vector.h:737
libMesh::PetscVector::scale
virtual void scale(const T factor) override
Scale each element of the vector by the given factor.
Definition: petsc_vector.C:407
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::IntRange
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
Definition: int_range.h:53
libMesh::PetscVector::print_matlab
virtual void print_matlab(const std::string &name="") const override
Print the contents of the vector in Matlab's sparse matrix format.
Definition: petsc_vector.C:1214
libMesh::PetscVector::operator-=
virtual NumericVector< T > & operator-=(const NumericVector< T > &v) override
Subtracts v from *this, .
Definition: petsc_vector.C:132
libMesh::PetscVector::dot
virtual T dot(const NumericVector< T > &v) const override
Definition: petsc_vector.C:486
libMesh::PetscVector::close
virtual void close() override
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
Definition: petsc_vector.h:821
libMesh::numeric_petsc_cast
PetscInt * numeric_petsc_cast(const numeric_index_type *p)
Definition: petsc_vector.h:1228
libMesh::PetscVector::localize
virtual void localize(std::vector< T > &v_local) const override
Creates a copy of the global vector in the local vector v_local.
Definition: petsc_vector.C:946
libMesh::PS
PetscScalar PS(T val)
Definition: petsc_macro.h:168
libMesh::NumericVector::initialized
virtual bool initialized() const
Definition: numeric_vector.h:155
libMesh::Utility::iota
void iota(ForwardIter first, ForwardIter last, T value)
Utility::iota is a duplication of the SGI STL extension std::iota.
Definition: utility.h:105
libMesh::PetscVector::add_vector_transpose
virtual void add_vector_transpose(const NumericVector< T > &v, const SparseMatrix< T > &A) override
Computes , i.e.
Definition: petsc_vector.C:256
parallel_implementation.h
libMesh::processor_id_type
uint8_t processor_id_type
Definition: id_types.h:104
libMesh::PetscVector::l2_norm
virtual Real l2_norm() const override
Definition: petsc_vector.C:81
libMesh::PetscVector::operator/=
virtual NumericVector< T > & operator/=(const NumericVector< T > &v) override
Computes the component-wise division of this vector's entries by another's, .
Definition: petsc_vector.C:447
libMesh::MeshTools::Generation::Private::idx
unsigned int idx(const ElemType type, const unsigned int nx, const unsigned int i, const unsigned int j)
A useful inline function which replaces the macros used previously.
Definition: mesh_generation.C:72
libMesh::PetscVector
This class provides a nice interface to PETSc's Vec object.
Definition: petsc_vector.h:72
A
static PetscErrorCode Mat * A
Definition: petscdmlibmeshimpl.C:1026
libMesh::is
PetscErrorCode PetscInt const PetscInt IS * is
Definition: petsc_dm_wrapper.C:60
libMesh::PetscVector::_get_array
void _get_array(bool read_only) const
Queries the array (and the local form if the vector is ghosted) from PETSc.
Definition: petsc_vector.C:1374
libMesh::PetscVector::size
virtual numeric_index_type size() const override
Definition: petsc_vector.h:933
libMesh::numeric_index_type
dof_id_type numeric_index_type
Definition: id_types.h:99
libMesh::Complex
std::complex< Real > Complex
Definition: libmesh_common.h:160
libMesh::initialized
bool initialized()
Checks that library initialization has been done.
Definition: libmesh.C:265
libMesh::PetscVector::add_vector_conjugate_transpose
void add_vector_conjugate_transpose(const NumericVector< T > &v, const SparseMatrix< T > &A)
.
Definition: petsc_vector.C:284
libMesh::PetscVector::local_size
virtual numeric_index_type local_size() const override
Definition: petsc_vector.h:953
libMesh::PetscVector::linfty_norm
virtual Real linfty_norm() const override
Definition: petsc_vector.C:99
libMesh::ReferenceElem::get
const Elem & get(const ElemType type_in)
Definition: reference_elem.C:237
standard_type.h
libMesh::PetscVector::add
virtual void add(const numeric_index_type i, const T value) override
Adds value to each entry of the vector.
Definition: petsc_vector.C:187
value
static const bool value
Definition: xdr_io.C:56
libMesh::PetscMatrix
This class provides a nice interface to the PETSc C-based data structures for parallel,...
Definition: petsc_matrix.h:87
libMesh::PetscVector::_restore_array
void _restore_array() const
Restores the array (and the local form if the vector is ghosted) to PETSc.
Definition: petsc_vector.C:1478
libMesh::ctx
void * ctx
Definition: petsc_dm_wrapper.C:71
libMesh::PetscVector::_vec
Vec _vec
Actual PETSc vector datatype to hold vector entries.
Definition: petsc_vector.h:345
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::PetscVector::l1_norm
virtual Real l1_norm() const override
Definition: petsc_vector.C:64
libMesh::PetscVector::operator+=
virtual NumericVector< T > & operator+=(const NumericVector< T > &v) override
Adds v to *this, .
Definition: petsc_vector.C:118
libMesh::PetscVector::insert
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices) override
Inserts the entries of v in *this at the locations specified by v.
Definition: petsc_vector.C:387
libMesh::NumericVector::type
ParallelType type() const
Definition: numeric_vector.h:160
libMesh::PetscVector::localize_to_one
virtual void localize_to_one(std::vector< T > &v_local, const processor_id_type proc_id=0) const override
Creates a local copy of the global vector in v_local only on processor proc_id.
libMesh::Quality::name
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
libMesh::closed
bool closed()
Checks that the library has been closed.
Definition: libmesh.C:272