18 #ifdef DEAL_II_WITH_PETSC
36 VectorReference::operator PetscScalar()
const
50 VecGetOwnershipRange(vector.vector, &
begin, &
end);
53 Vec locally_stored_elements = PETSC_NULL;
54 ierr = VecGhostGetLocalForm(vector.vector, &locally_stored_elements);
58 ierr = VecGetSize(locally_stored_elements, &lsize);
62 ierr = VecGetArray(locally_stored_elements, &ptr);
77 vector.ghost_indices.index_within_set(index);
83 ierr = VecRestoreArray(locally_stored_elements, &ptr);
87 VecGhostRestoreLocalForm(vector.vector, &locally_stored_elements);
99 PetscErrorCode ierr = VecGetOwnershipRange(vector.vector, &
begin, &
end);
104 ExcAccessToNonlocalElement(index,
begin,
end - 1));
106 PetscInt idx = index;
108 ierr = VecGetValues(vector.vector, 1, &idx, &
value);
120 , obtained_ownership(true)
123 ExcMessage(
"PETSc does not support multi-threaded access, set "
124 "the thread limit to 1 in MPI_InitFinalize()."));
132 , ghost_indices(v.ghost_indices)
134 , obtained_ownership(true)
137 ExcMessage(
"PETSc does not support multi-threaded access, set "
138 "the thread limit to 1 in MPI_InitFinalize()."));
154 , obtained_ownership(false)
157 ExcMessage(
"PETSc does not support multi-threaded access, set "
158 "the thread limit to 1 in MPI_InitFinalize()."));
167 const PetscErrorCode ierr = VecDestroy(&
vector);
180 const PetscErrorCode ierr = VecDestroy(&
vector);
199 PetscErrorCode ierr = VecSet(
vector, s);
204 Vec ghost = PETSC_NULL;
205 ierr = VecGhostGetLocalForm(
vector, &ghost);
208 ierr = VecSet(ghost, s);
211 ierr = VecGhostRestoreLocalForm(
vector, &ghost);
226 const PetscErrorCode ierr = VecEqual(
vector, v.
vector, &flag);
229 return (flag == PETSC_TRUE);
240 const PetscErrorCode ierr = VecEqual(
vector, v.
vector, &flag);
243 return (flag == PETSC_FALSE);
252 const PetscErrorCode ierr = VecGetSize(
vector, &sz);
264 const PetscErrorCode ierr = VecGetLocalSize(
vector, &sz);
272 std::pair<VectorBase::size_type, VectorBase::size_type>
276 const PetscErrorCode ierr =
277 VecGetOwnershipRange(
static_cast<const Vec &
>(
vector), &
begin, &
end);
287 const std::vector<PetscScalar> &values)
289 Assert(indices.size() == values.size(),
290 ExcMessage(
"Function called with arguments of different sizes"));
298 const std::vector<PetscScalar> &values)
300 Assert(indices.size() == values.size(),
301 ExcMessage(
"Function called with arguments of different sizes"));
309 const ::Vector<PetscScalar> &values)
311 Assert(indices.size() == values.size(),
312 ExcMessage(
"Function called with arguments of different sizes"));
321 const PetscScalar *values)
341 const PetscErrorCode ierr = VecDot(vec.
vector,
vector, &result);
365 # ifdef DEAL_II_WITH_MPI
369 int all_int_last_action;
371 const int ierr = MPI_Allreduce(&my_int_last_action,
372 &all_int_last_action,
381 ExcMessage(
"Error: not all processors agree on the last "
382 "VectorOperation before this compress() call."));
391 "Missing compress() or calling with wrong VectorOperation argument."));
407 PetscErrorCode ierr = VecAssemblyBegin(
vector);
409 ierr = VecAssemblyEnd(
vector);
437 const PetscErrorCode ierr = VecSum(
vector, &
sum);
439 return sum /
static_cast<PetscReal
>(
size());
444 PetscScalar * start_ptr;
445 PetscErrorCode ierr = VecGetArray(
vector, &start_ptr);
448 PetscScalar
mean = 0;
450 PetscScalar sum0 = 0, sum1 = 0, sum2 = 0, sum3 = 0;
455 const PetscScalar *ptr = start_ptr;
456 const PetscScalar *eptr = ptr + (
size() / 4) * 4;
465 while (ptr != start_ptr +
size())
468 mean = (sum0 + sum1 + sum2 + sum3) /
static_cast<PetscReal
>(
size());
473 ierr = VecRestoreArray(
vector, &start_ptr);
485 const PetscErrorCode ierr = VecNorm(
vector, NORM_1, &
d);
498 const PetscErrorCode ierr = VecNorm(
vector, NORM_2, &
d);
511 PetscScalar * start_ptr;
512 PetscErrorCode ierr = VecGetArray(
vector, &start_ptr);
517 real_type sum0 = 0, sum1 = 0, sum2 = 0, sum3 = 0;
522 const PetscScalar *ptr = start_ptr;
523 const PetscScalar *eptr = ptr + (
size() / 4) * 4;
532 while (ptr != start_ptr +
size())
535 norm = std::pow(sum0 + sum1 + sum2 + sum3, 1. / p);
540 ierr = VecRestoreArray(
vector, &start_ptr);
553 const PetscErrorCode ierr = VecNorm(
vector, NORM_INFINITY, &
d);
567 const PetscErrorCode ierr = VecMin(
vector, &p, &
d);
580 const PetscErrorCode ierr = VecMax(
vector, &p, &
d);
593 PetscScalar * start_ptr;
594 PetscErrorCode ierr = VecGetArray(
vector, &start_ptr);
597 const PetscScalar *ptr = start_ptr, *eptr = start_ptr +
local_size();
611 ierr = VecRestoreArray(
vector, &start_ptr);
620 template <
typename T>
629 template <
typename T>
635 "whether it is non-negative."))
return true;
646 PetscScalar * start_ptr;
647 PetscErrorCode ierr = VecGetArray(
vector, &start_ptr);
650 const PetscScalar *ptr = start_ptr, *eptr = start_ptr +
local_size();
664 ierr = VecRestoreArray(
vector, &start_ptr);
678 const PetscErrorCode ierr = VecScale(
vector, a);
692 const PetscScalar factor = 1. / a;
695 const PetscErrorCode ierr = VecScale(
vector, factor);
707 const PetscErrorCode ierr = VecAXPY(
vector, 1, v);
719 const PetscErrorCode ierr = VecAXPY(
vector, -1, v);
733 const PetscErrorCode ierr = VecShift(
vector, s);
745 const PetscErrorCode ierr = VecAXPY(
vector, a, v);
761 const PetscScalar weights[2] = {a,
b};
762 Vec addends[2] = {v.
vector,
w.vector};
764 const PetscErrorCode ierr = VecMAXPY(
vector, 2, weights, addends);
776 const PetscErrorCode ierr = VecAYPX(
vector, s, v);
804 const PetscErrorCode ierr = VecPointwiseMult(
vector, factors,
vector);
835 PetscErrorCode ierr =
836 PetscViewerSetFormat(PETSC_VIEWER_STDOUT_WORLD, format);
840 ierr = VecView(
vector, PETSC_VIEWER_STDOUT_WORLD);
848 const unsigned int precision,
849 const bool scientific,
850 const bool across)
const
857 PetscErrorCode ierr = VecGetArray(
vector, &val);
862 const std::ios::fmtflags old_flags = out.flags();
863 const unsigned int old_precision = out.precision(precision);
865 out.precision(precision);
867 out.setf(std::ios::scientific, std::ios::floatfield);
869 out.setf(std::ios::fixed, std::ios::floatfield);
873 out << val[i] <<
' ';
876 out << val[i] << std::endl;
880 out.flags(old_flags);
881 out.precision(old_precision);
885 ierr = VecRestoreArray(
vector, &val);
902 VectorBase::operator
const Vec &()
const
911 std::size_t mem =
sizeof(Vec) +
sizeof(
last_action) +
933 const PetscScalar *values,
934 const bool add_values)
941 internal::VectorReference::ExcWrongMode(action,
last_action));
951 const PetscInt *petsc_indices =
952 reinterpret_cast<const PetscInt *
>(indices);
954 const InsertMode mode = (add_values ? ADD_VALUES : INSERT_VALUES);
955 const PetscErrorCode ierr =
956 VecSetValues(
vector, n_elements, petsc_indices, values, mode);
969 #endif // DEAL_II_WITH_PETSC