#include <petscvec.h>
+#include <vector>
#include <utility>
*/
PetscScalar
operator () (const unsigned int index) const;
+
+ /**
+ * A collective set operation: instead
+ * of setting individual elements of a
+ * vector, this function allows to set
+ * a whole set of elements at once. The
+ * indices of the elements to be set
+ * are stated in the first argument,
+ * the corresponding values in the
+ * second.
+ */
+ void set (const std::vector<unsigned int> &indices,
+ const std::vector<PetscScalar> &values);
/**
* Return the scalar product of two
// $Id$
// Version: $Name$
//
-// Copyright (C) 2004 by the deal.II authors
+// Copyright (C) 2004, 2005 by the deal.II authors
//
// This file is subject to QPL and may not be distributed
// without copyright and license information. Please refer
+ void
+ VectorBase::set (const std::vector<unsigned int> &indices,
+ const std::vector<PetscScalar> &values)
+ {
+ Assert (indices.size() == values.size(),
+ ExcMessage ("Function called with arguments of different sizes"));
+
+ if (last_action != VectorBase::LastAction::insert)
+ {
+ int ierr;
+ ierr = VecAssemblyBegin (vector);
+ AssertThrow (ierr == 0, ExcPETScError(ierr));
+
+ ierr = VecAssemblyEnd (vector);
+ AssertThrow (ierr == 0, ExcPETScError(ierr));
+ }
+
+#if (PETSC_VERSION_MAJOR <= 2) && \
+ ((PETSC_VERSION_MINOR < 2) || \
+ ((PETSC_VERSION_MINOR == 2) && (PETSC_VERSION_SUBMINOR == 0)))
+ const std::vector<int> petsc_indices (indices.begin(),
+ indices.end());
+#else
+ const std::vector<PetscInt> petsc_indices (indices.begin(),
+ indices.end());
+#endif
+
+ const int ierr
+ = VecSetValues (vector, indices.size(),
+ &petsc_indices[0], &values[0],
+ INSERT_VALUES);
+ AssertThrow (ierr == 0, ExcPETScError(ierr));
+
+
+ last_action = LastAction::insert;
+ }
+
+
+
PetscScalar
VectorBase::operator * (const VectorBase &vec) const
{