24 #include "libmesh/petsc_vector.h"
27 #include "libmesh/petsc_matrix.h"
29 #ifdef LIBMESH_HAVE_PETSC
31 #include "libmesh/dense_subvector.h"
32 #include "libmesh/dense_vector.h"
33 #include "libmesh/int_range.h"
34 #include "libmesh/petsc_macro.h"
50 this->_restore_array();
53 PetscErrorCode
ierr=0;
59 return static_cast<T>(
value);
66 this->_restore_array();
69 PetscErrorCode
ierr=0;
75 return static_cast<Real>(
value);
83 this->_restore_array();
86 PetscErrorCode
ierr=0;
92 return static_cast<Real>(
value);
101 this->_restore_array();
104 PetscErrorCode
ierr=0;
107 ierr = VecNorm (_vec, NORM_INFINITY, &
value);
108 LIBMESH_CHKERR(
ierr);
110 return static_cast<Real>(
value);
116 template <
typename T>
120 this->_restore_array();
130 template <
typename T>
134 this->_restore_array();
144 template <
typename T>
147 this->_restore_array();
148 libmesh_assert_less (i, size());
150 PetscErrorCode
ierr=0;
151 PetscInt i_val = static_cast<PetscInt>(i);
152 PetscScalar petsc_value =
PS(
value);
154 ierr = VecSetValues (_vec, 1, &i_val, &petsc_value, INSERT_VALUES);
155 LIBMESH_CHKERR(
ierr);
157 this->_is_closed =
false;
162 template <
typename T>
165 PetscErrorCode
ierr = 0;
168 ierr = VecReciprocal(_vec);
169 LIBMESH_CHKERR(
ierr);
174 template <
typename T>
177 PetscErrorCode
ierr = 0;
180 ierr = VecConjugate(_vec);
181 LIBMESH_CHKERR(
ierr);
186 template <
typename T>
189 this->_restore_array();
190 libmesh_assert_less (i, size());
192 PetscErrorCode
ierr=0;
193 PetscInt i_val = static_cast<PetscInt>(i);
194 PetscScalar petsc_value =
PS(
value);
196 ierr = VecSetValues (_vec, 1, &i_val, &petsc_value, ADD_VALUES);
197 LIBMESH_CHKERR(
ierr);
199 this->_is_closed =
false;
204 template <
typename T>
206 const std::vector<numeric_index_type> & dof_indices)
209 if (dof_indices.empty())
212 this->_restore_array();
214 PetscErrorCode
ierr=0;
215 const PetscInt * i_val = reinterpret_cast<const PetscInt *>(dof_indices.data());
216 const PetscScalar * petsc_value =
pPS(v);
218 ierr = VecSetValues (_vec, cast_int<PetscInt>(dof_indices.size()),
219 i_val, petsc_value, ADD_VALUES);
220 LIBMESH_CHKERR(
ierr);
222 this->_is_closed =
false;
227 template <
typename T>
231 this->_restore_array();
233 const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
236 PetscErrorCode
ierr=0;
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();
250 LIBMESH_CHKERR(
ierr);
255 template <
typename T>
259 this->_restore_array();
261 const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
264 PetscErrorCode
ierr=0;
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();
278 LIBMESH_CHKERR(
ierr);
282 #if PETSC_VERSION_LESS_THAN(3,1,0)
283 template <
typename T>
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.");
294 template <
typename T>
298 this->_restore_array();
300 const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
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();
314 std::unique_ptr<NumericVector<Number>> this_clone = this->clone();
318 PetscErrorCode
ierr = MatMultHermitianTranspose(
const_cast<PetscMatrix<T> *
>(
A)->mat(), v->
_vec, _vec);
319 LIBMESH_CHKERR(
ierr);
322 this->add(1., *this_clone);
327 template <
typename T>
330 this->_get_array(
false);
333 _values[i] += PetscScalar(v_in);
338 template <
typename T>
346 template <
typename T>
349 this->_restore_array();
351 PetscErrorCode
ierr = 0;
352 PetscScalar a =
PS(a_in);
355 const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
358 libmesh_assert_equal_to (this->size(), v->
size());
363 LIBMESH_CHKERR(
ierr);
369 ierr = VecGhostGetLocalForm (_vec,&loc_vec);
370 LIBMESH_CHKERR(
ierr);
371 ierr = VecGhostGetLocalForm (v->
_vec,&v_loc_vec);
372 LIBMESH_CHKERR(
ierr);
374 ierr = VecAXPY(loc_vec, a, v_loc_vec);
375 LIBMESH_CHKERR(
ierr);
377 ierr = VecGhostRestoreLocalForm (v->
_vec,&v_loc_vec);
378 LIBMESH_CHKERR(
ierr);
379 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
380 LIBMESH_CHKERR(
ierr);
386 template <
typename T>
388 const std::vector<numeric_index_type> & dof_indices)
390 if (dof_indices.empty())
393 this->_restore_array();
395 PetscErrorCode
ierr=0;
397 ierr = VecSetValues (_vec, cast_int<PetscInt>(dof_indices.size()),
398 idx_values,
pPS(v), INSERT_VALUES);
399 LIBMESH_CHKERR(
ierr);
401 this->_is_closed =
false;
406 template <
typename T>
409 this->_restore_array();
411 PetscErrorCode
ierr = 0;
412 PetscScalar factor =
PS(factor_in);
416 ierr = VecScale(_vec, factor);
417 LIBMESH_CHKERR(
ierr);
422 ierr = VecGhostGetLocalForm (_vec,&loc_vec);
423 LIBMESH_CHKERR(
ierr);
425 ierr = VecScale(loc_vec, factor);
426 LIBMESH_CHKERR(
ierr);
428 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
429 LIBMESH_CHKERR(
ierr);
433 template <
typename T>
436 PetscErrorCode
ierr = 0;
438 const PetscVector<T> * v_vec = cast_ptr<const PetscVector<T> *>(&v);
440 ierr = VecPointwiseMult(_vec, _vec, v_vec->
_vec);
441 LIBMESH_CHKERR(
ierr);
446 template <
typename T>
449 PetscErrorCode
ierr = 0;
451 const PetscVector<T> * v_vec = cast_ptr<const PetscVector<T> *>(&v);
453 ierr = VecPointwiseDivide(_vec, _vec, v_vec->
_vec);
454 LIBMESH_CHKERR(
ierr);
459 template <
typename T>
462 this->_restore_array();
464 PetscErrorCode
ierr = 0;
469 LIBMESH_CHKERR(
ierr);
474 ierr = VecGhostGetLocalForm (_vec,&loc_vec);
475 LIBMESH_CHKERR(
ierr);
477 ierr = VecAbs(loc_vec);
478 LIBMESH_CHKERR(
ierr);
480 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
481 LIBMESH_CHKERR(
ierr);
485 template <
typename T>
488 this->_restore_array();
491 PetscErrorCode
ierr = 0;
494 PetscScalar
value=0.;
497 const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
501 LIBMESH_CHKERR(
ierr);
503 return static_cast<T>(
value);
506 template <
typename T>
509 this->_restore_array();
512 PetscErrorCode
ierr = 0;
515 PetscScalar
value=0.;
518 const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
522 LIBMESH_CHKERR(
ierr);
524 return static_cast<T>(
value);
528 template <
typename T>
532 this->_restore_array();
535 PetscErrorCode
ierr = 0;
536 PetscScalar s =
PS(s_in);
538 if (this->size() != 0)
542 ierr = VecSet(_vec, s);
543 LIBMESH_CHKERR(
ierr);
548 ierr = VecGhostGetLocalForm (_vec,&loc_vec);
549 LIBMESH_CHKERR(
ierr);
551 ierr = VecSet(loc_vec, s);
552 LIBMESH_CHKERR(
ierr);
554 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
555 LIBMESH_CHKERR(
ierr);
564 template <
typename T>
569 const PetscVector<T> * v = cast_ptr<const PetscVector<T> *>(&v_in);
578 template <
typename T>
582 this->_restore_array();
585 libmesh_assert_equal_to (this->size(), v.
size());
586 libmesh_assert_equal_to (this->local_size(), v.
local_size());
589 PetscErrorCode
ierr = 0;
599 ierr = VecCopy (v.
_vec, this->_vec);
600 LIBMESH_CHKERR(
ierr);
606 libmesh_assert_equal_to (this->_type, v.
_type);
612 ierr = VecCopy (v.
_vec, this->_vec);
613 LIBMESH_CHKERR(
ierr);
619 ierr = VecGhostGetLocalForm (_vec,&loc_vec);
620 LIBMESH_CHKERR(
ierr);
621 ierr = VecGhostGetLocalForm (v.
_vec,&v_loc_vec);
622 LIBMESH_CHKERR(
ierr);
624 ierr = VecCopy (v_loc_vec, loc_vec);
625 LIBMESH_CHKERR(
ierr);
627 ierr = VecGhostRestoreLocalForm (v.
_vec,&v_loc_vec);
628 LIBMESH_CHKERR(
ierr);
629 ierr = VecGhostRestoreLocalForm (_vec,&loc_vec);
630 LIBMESH_CHKERR(
ierr);
642 template <
typename T>
646 this->_get_array(
false);
652 if (this->size() == v.size())
657 _values[i] =
PS(v[first + i]);
667 _values[i] =
PS(v[i]);
679 template <
typename T>
682 this->_restore_array();
685 PetscVector<T> * v_local = cast_ptr<PetscVector<T> *>(&v_local_in);
688 libmesh_assert_equal_to (v_local->
size(), this->size());
690 PetscErrorCode
ierr = 0;
691 const PetscInt n = this->size();
697 std::vector<PetscInt>
idx(n);
702 LIBMESH_CHKERR(
ierr);
704 ierr = LibMeshVecScatterCreate(_vec,
is,
707 LIBMESH_CHKERR(
ierr);
710 ierr = VecScatterBegin(scatter, _vec, v_local->
_vec,
711 INSERT_VALUES, SCATTER_FORWARD);
712 LIBMESH_CHKERR(
ierr);
714 ierr = VecScatterEnd (scatter, _vec, v_local->
_vec,
715 INSERT_VALUES, SCATTER_FORWARD);
716 LIBMESH_CHKERR(
ierr);
719 ierr = LibMeshISDestroy (&
is);
720 LIBMESH_CHKERR(
ierr);
722 ierr = LibMeshVecScatterDestroy(&scatter);
723 LIBMESH_CHKERR(
ierr);
732 template <
typename T>
734 const std::vector<numeric_index_type> & send_list)
const
748 this->_restore_array();
751 PetscVector<T> * v_local = cast_ptr<PetscVector<T> *>(&v_local_in);
754 libmesh_assert_equal_to (v_local->
size(), this->size());
755 libmesh_assert_less_equal (send_list.size(), v_local->
size());
757 PetscErrorCode
ierr=0;
759 cast_int<numeric_index_type>(send_list.size());
764 std::vector<PetscInt>
idx(n_sl + this->local_size());
767 idx[i] = static_cast<PetscInt>(send_list[i]);
769 idx[n_sl+i] = i + this->first_local_index();
772 PetscInt * idxptr =
idx.empty() ? nullptr :
idx.data();
773 ierr = ISCreateLibMesh(this->comm().
get(), n_sl+this->local_size(),
775 LIBMESH_CHKERR(
ierr);
777 ierr = LibMeshVecScatterCreate(_vec,
is,
780 LIBMESH_CHKERR(
ierr);
784 ierr = VecScatterBegin(scatter, _vec, v_local->
_vec,
785 INSERT_VALUES, SCATTER_FORWARD);
786 LIBMESH_CHKERR(
ierr);
788 ierr = VecScatterEnd (scatter, _vec, v_local->
_vec,
789 INSERT_VALUES, SCATTER_FORWARD);
790 LIBMESH_CHKERR(
ierr);
793 ierr = LibMeshISDestroy (&
is);
794 LIBMESH_CHKERR(
ierr);
796 ierr = LibMeshVecScatterDestroy(&scatter);
797 LIBMESH_CHKERR(
ierr);
806 template <
typename T>
808 const std::vector<numeric_index_type> & indices)
const
811 PetscErrorCode
ierr = 0;
815 ierr = VecCreateSeq(PETSC_COMM_SELF, cast_int<PetscInt>(indices.size()), &dest);
816 LIBMESH_CHKERR(
ierr);
824 ierr = ISCreateLibMesh(this->comm().
get(), cast_int<PetscInt>(indices.size()), idxptr,
826 LIBMESH_CHKERR(
ierr);
830 ierr = LibMeshVecScatterCreate(_vec,
835 LIBMESH_CHKERR(
ierr);
838 ierr = VecScatterBegin(scatter, _vec, dest, INSERT_VALUES, SCATTER_FORWARD);
839 LIBMESH_CHKERR(
ierr);
841 ierr = VecScatterEnd(scatter, _vec, dest, INSERT_VALUES, SCATTER_FORWARD);
842 LIBMESH_CHKERR(
ierr);
845 PetscScalar * values;
846 ierr = VecGetArray (dest, &values);
847 LIBMESH_CHKERR(
ierr);
852 v_local.reserve(indices.size());
854 v_local.insert(v_local.begin(), values, values+indices.size());
857 ierr = VecRestoreArray (dest, &values);
858 LIBMESH_CHKERR(
ierr);
861 ierr = LibMeshVecScatterDestroy(&scatter);
862 LIBMESH_CHKERR(
ierr);
864 ierr = LibMeshISDestroy (&
is);
865 LIBMESH_CHKERR(
ierr);
867 ierr = LibMeshVecDestroy(&dest);
868 LIBMESH_CHKERR(
ierr);
873 template <
typename T>
876 const std::vector<numeric_index_type> & send_list)
878 this->_restore_array();
880 libmesh_assert_less_equal (send_list.size(), this->size());
881 libmesh_assert_less_equal (last_local_idx+1, this->size());
885 PetscErrorCode
ierr=0;
891 if (this->n_processors() == 1)
899 parallel_vec.
init (my_size, my_local_size,
true,
PARALLEL);
908 std::vector<PetscInt>
idx(my_local_size);
912 ierr = ISCreateLibMesh(this->comm().
get(), my_local_size,
914 LIBMESH_CHKERR(
ierr);
916 ierr = LibMeshVecScatterCreate(_vec,
is,
919 LIBMESH_CHKERR(
ierr);
921 ierr = VecScatterBegin(scatter, _vec, parallel_vec.
_vec,
922 INSERT_VALUES, SCATTER_FORWARD);
923 LIBMESH_CHKERR(
ierr);
925 ierr = VecScatterEnd (scatter, _vec, parallel_vec.
_vec,
926 INSERT_VALUES, SCATTER_FORWARD);
927 LIBMESH_CHKERR(
ierr);
930 ierr = LibMeshISDestroy (&
is);
931 LIBMESH_CHKERR(
ierr);
933 ierr = LibMeshVecScatterDestroy(&scatter);
934 LIBMESH_CHKERR(
ierr);
938 parallel_vec.
close();
939 parallel_vec.
localize (*
this, send_list);
945 template <
typename T>
948 this->_restore_array();
951 parallel_object_only();
953 PetscErrorCode
ierr=0;
954 const PetscInt n = this->size();
955 const PetscInt nl = this->local_size();
956 PetscScalar * values;
959 v_local.resize(n, 0.);
961 ierr = VecGetArray (_vec, &values);
962 LIBMESH_CHKERR(
ierr);
966 for (PetscInt i=0; i<nl; i++)
967 v_local[i+ioff] = static_cast<T>(values[i]);
969 ierr = VecRestoreArray (_vec, &values);
970 LIBMESH_CHKERR(
ierr);
972 this->comm().sum(v_local);
978 #ifdef LIBMESH_USE_REAL_NUMBERS
983 timpi_mpi_var(pid))
const
985 this->_restore_array();
987 PetscErrorCode
ierr=0;
988 const PetscInt n = size();
989 PetscScalar * values;
992 if (n_processors() == 1)
996 ierr = VecGetArray (_vec, &values);
997 LIBMESH_CHKERR(
ierr);
999 for (PetscInt i=0; i<n; i++)
1000 v_local[i] = static_cast<Real>(values[i]);
1002 ierr = VecRestoreArray (_vec, &values);
1003 LIBMESH_CHKERR(
ierr);
1007 #ifdef LIBMESH_HAVE_MPI
1015 ierr = VecScatterCreateToZero(_vec, &
ctx, &vout);
1016 LIBMESH_CHKERR(
ierr);
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);
1023 if (processor_id() == 0)
1027 ierr = VecGetArray (vout, &values);
1028 LIBMESH_CHKERR(
ierr);
1030 for (PetscInt i=0; i<n; i++)
1031 v_local[i] = static_cast<Real>(values[i]);
1033 ierr = VecRestoreArray (vout, &values);
1034 LIBMESH_CHKERR(
ierr);
1037 ierr = LibMeshVecScatterDestroy(&
ctx);
1038 LIBMESH_CHKERR(
ierr);
1039 ierr = LibMeshVecDestroy(&vout);
1040 LIBMESH_CHKERR(
ierr);
1048 std::vector<Real> local_values (n, 0.);
1051 ierr = VecGetArray (_vec, &values);
1052 LIBMESH_CHKERR(
ierr);
1054 const PetscInt nl = local_size();
1055 for (PetscInt i=0; i<nl; i++)
1056 local_values[i+ioff] = static_cast<Real>(values[i]);
1058 ierr = VecRestoreArray (_vec, &values);
1059 LIBMESH_CHKERR(
ierr);
1063 MPI_Reduce (local_values.data(), v_local.data(), n,
1064 Parallel::StandardType<Real>(),
1065 Parallel::OpFunction<Real>::sum(), pid,
1066 this->comm().get());
1069 #endif // LIBMESH_HAVE_MPI
1076 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
1082 this->_restore_array();
1084 PetscErrorCode
ierr=0;
1085 const PetscInt n = size();
1086 const PetscInt nl = local_size();
1087 PetscScalar * values;
1093 for (PetscInt i=0; i<n; i++)
1097 if (n_processors() == 1)
1099 ierr = VecGetArray (_vec, &values);
1100 LIBMESH_CHKERR(
ierr);
1102 for (PetscInt i=0; i<n; i++)
1103 v_local[i] = static_cast<Complex>(values[i]);
1105 ierr = VecRestoreArray (_vec, &values);
1106 LIBMESH_CHKERR(
ierr);
1117 std::vector<Real> real_local_values(n, 0.);
1118 std::vector<Real> imag_local_values(n, 0.);
1121 ierr = VecGetArray (_vec, &values);
1122 LIBMESH_CHKERR(
ierr);
1125 for (PetscInt i=0; i<nl; i++)
1127 real_local_values[i+ioff] = static_cast<Complex>(values[i]).real();
1128 imag_local_values[i+ioff] = static_cast<Complex>(values[i]).imag();
1131 ierr = VecRestoreArray (_vec, &values);
1132 LIBMESH_CHKERR(
ierr);
1140 std::vector<Real> real_v_local(n);
1141 std::vector<Real> imag_v_local(n);
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());
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());
1157 for (PetscInt i=0; i<n; i++)
1158 v_local[i] =
Complex(real_v_local[i], imag_v_local[i]);
1166 template <
typename T>
1170 this->_restore_array();
1172 PetscErrorCode
ierr = 0;
1175 const PetscVector<T> * vec1_petsc = cast_ptr<const PetscVector<T> *>(&vec1);
1176 const PetscVector<T> * vec2_petsc = cast_ptr<const PetscVector<T> *>(&vec2);
1182 ierr = VecPointwiseMult(this->vec(),
1185 LIBMESH_CHKERR(
ierr);
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);
1199 ierr = VecPointwiseMult(loc_vec,v1_loc_vec,v2_loc_vec);
1200 LIBMESH_CHKERR(
ierr);
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);
1213 template <
typename T>
1216 this->_restore_array();
1219 PetscErrorCode
ierr=0;
1220 PetscViewer petsc_viewer;
1223 ierr = PetscViewerCreate (this->comm().
get(),
1225 LIBMESH_CHKERR(
ierr);
1233 ierr = PetscViewerASCIIOpen( this->comm().
get(),
1236 LIBMESH_CHKERR(
ierr);
1238 #if PETSC_VERSION_LESS_THAN(3,7,0)
1239 ierr = PetscViewerSetFormat (petsc_viewer,
1240 PETSC_VIEWER_ASCII_MATLAB);
1242 ierr = PetscViewerPushFormat (petsc_viewer,
1243 PETSC_VIEWER_ASCII_MATLAB);
1246 LIBMESH_CHKERR(
ierr);
1248 ierr = VecView (_vec, petsc_viewer);
1249 LIBMESH_CHKERR(
ierr);
1258 #if PETSC_VERSION_LESS_THAN(3,7,0)
1259 ierr = PetscViewerSetFormat (PETSC_VIEWER_STDOUT_WORLD,
1260 PETSC_VIEWER_ASCII_MATLAB);
1262 ierr = PetscViewerPushFormat (PETSC_VIEWER_STDOUT_WORLD,
1263 PETSC_VIEWER_ASCII_MATLAB);
1266 LIBMESH_CHKERR(
ierr);
1268 ierr = VecView (_vec, PETSC_VIEWER_STDOUT_WORLD);
1269 LIBMESH_CHKERR(
ierr);
1276 ierr = LibMeshPetscViewerDestroy (&petsc_viewer);
1277 LIBMESH_CHKERR(
ierr);
1284 template <
typename T>
1286 const std::vector<numeric_index_type> & rows)
const
1288 this->_restore_array();
1291 IS parent_is, subvector_is;
1293 PetscErrorCode
ierr = 0;
1296 PetscVector<T> * petsc_subvector = cast_ptr<PetscVector<T> *>(&subvector);
1309 ierr = VecCreateMPI(this->comm().
get(),
1311 cast_int<PetscInt>(rows.size()),
1312 &(petsc_subvector->
_vec));
1313 LIBMESH_CHKERR(
ierr);
1315 ierr = VecSetFromOptions (petsc_subvector->
_vec);
1316 LIBMESH_CHKERR(
ierr);
1327 std::vector<PetscInt>
idx(rows.size());
1331 ierr = ISCreateLibMesh(this->comm().
get(),
1332 cast_int<PetscInt>(rows.size()),
1336 LIBMESH_CHKERR(
ierr);
1338 ierr = ISCreateLibMesh(this->comm().
get(),
1339 cast_int<PetscInt>(rows.size()),
1343 LIBMESH_CHKERR(
ierr);
1346 ierr = LibMeshVecScatterCreate(this->_vec,
1348 petsc_subvector->
_vec,
1350 &scatter); LIBMESH_CHKERR(
ierr);
1353 ierr = VecScatterBegin(scatter,
1355 petsc_subvector->
_vec,
1357 SCATTER_FORWARD); LIBMESH_CHKERR(
ierr);
1359 ierr = VecScatterEnd(scatter,
1361 petsc_subvector->
_vec,
1363 SCATTER_FORWARD); LIBMESH_CHKERR(
ierr);
1366 ierr = LibMeshISDestroy(&parent_is); LIBMESH_CHKERR(
ierr);
1367 ierr = LibMeshISDestroy(&subvector_is); LIBMESH_CHKERR(
ierr);
1368 ierr = LibMeshVecScatterDestroy(&scatter); LIBMESH_CHKERR(
ierr);
1373 template <
typename T>
1378 #ifdef LIBMESH_HAVE_CXX11_THREAD
1379 std::atomic_thread_fence(std::memory_order_acquire);
1381 Threads::spin_mutex::scoped_lock lock(_petsc_vector_mutex);
1388 if (_array_is_present && !read_only && _values_read_only)
1393 if (_array_is_present && read_only && !_values_read_only)
1394 _read_only_values = _values;
1396 if (!_array_is_present)
1398 #ifdef LIBMESH_HAVE_CXX11_THREAD
1399 std::lock_guard<std::mutex> lock(_petsc_vector_mutex);
1401 if (!_array_is_present)
1403 PetscErrorCode
ierr=0;
1406 #if PETSC_VERSION_LESS_THAN(3,2,0)
1410 ierr = VecGetArray(_vec, const_cast<PetscScalar **>(&_values));
1411 _values_read_only =
false;
1415 ierr = VecGetArrayRead(_vec, &_read_only_values);
1416 _values_read_only =
true;
1420 ierr = VecGetArray(_vec, &_values);
1421 _values_read_only =
false;
1424 LIBMESH_CHKERR(
ierr);
1425 _local_size = this->local_size();
1429 ierr = VecGhostGetLocalForm (_vec,&_local_form);
1430 LIBMESH_CHKERR(
ierr);
1432 #if PETSC_VERSION_LESS_THAN(3,2,0)
1436 ierr = VecGetArray(_local_form, const_cast<PetscScalar **>(&_values));
1437 _values_read_only =
false;
1441 ierr = VecGetArrayRead(_local_form, &_read_only_values);
1442 _values_read_only =
true;
1446 ierr = VecGetArray(_local_form, &_values);
1447 _values_read_only =
false;
1450 LIBMESH_CHKERR(
ierr);
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);
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);
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);
1469 _array_is_present =
true;
1477 template <
typename T>
1480 if (_values_manually_retrieved)
1481 libmesh_error_msg(
"PetscVector values were manually retrieved but are being automatically restored!");
1483 #ifdef LIBMESH_HAVE_CXX11_THREAD
1484 std::atomic_thread_fence(std::memory_order_acquire);
1486 Threads::spin_mutex::scoped_lock lock(_petsc_vector_mutex);
1490 if (_array_is_present)
1492 #ifdef LIBMESH_HAVE_CXX11_THREAD
1493 std::lock_guard<std::mutex> lock(_petsc_vector_mutex);
1496 if (_array_is_present)
1498 PetscErrorCode
ierr=0;
1501 #if PETSC_VERSION_LESS_THAN(3,2,0)
1505 ierr = VecRestoreArray (_vec, const_cast<PetscScalar **>(&_values));
1507 if (_values_read_only)
1508 ierr = VecRestoreArrayRead (_vec, &_read_only_values);
1510 ierr = VecRestoreArray (_vec, &_values);
1513 LIBMESH_CHKERR(
ierr);
1518 #if PETSC_VERSION_LESS_THAN(3,2,0)
1522 ierr = VecRestoreArray (_local_form, const_cast<PetscScalar **>(&_values));
1524 if (_values_read_only)
1525 ierr = VecRestoreArrayRead (_local_form, &_read_only_values);
1527 ierr = VecRestoreArray (_local_form, &_values);
1529 LIBMESH_CHKERR(
ierr);
1531 ierr = VecGhostRestoreLocalForm (_vec,&_local_form);
1532 LIBMESH_CHKERR(
ierr);
1533 _local_form =
nullptr;
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);
1540 _array_is_present =
false;
1557 #endif // #ifdef LIBMESH_HAVE_PETSC