libMesh::DistributedVector< T > Class Template Referencefinal

#include <distributed_vector.h>

Inheritance diagram for libMesh::DistributedVector< T >:

Public Member Functions

 DistributedVector (const Parallel::Communicator &comm, const ParallelType=AUTOMATIC)
 
 DistributedVector (const Parallel::Communicator &comm, const numeric_index_type n, const ParallelType ptype=AUTOMATIC)
 
 DistributedVector (const Parallel::Communicator &comm, const numeric_index_type n, const numeric_index_type n_local, const ParallelType ptype=AUTOMATIC)
 
 DistributedVector (const Parallel::Communicator &comm, const numeric_index_type N, const numeric_index_type n_local, const std::vector< numeric_index_type > &ghost, const ParallelType ptype=AUTOMATIC)
 
DistributedVectoroperator= (const DistributedVector &)
 
 DistributedVector (DistributedVector &&)=default
 
 DistributedVector (const DistributedVector &)=default
 
DistributedVectoroperator= (DistributedVector &&)=default
 
virtual ~DistributedVector ()=default
 
virtual void close () override
 
virtual void clear () override
 
virtual void zero () override
 
virtual std::unique_ptr< NumericVector< T > > zero_clone () const override
 
virtual std::unique_ptr< NumericVector< T > > clone () const override
 
virtual void init (const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override
 
virtual void init (const numeric_index_type N, const bool fast=false, const ParallelType ptype=AUTOMATIC) override
 
virtual void init (const numeric_index_type N, const numeric_index_type n_local, const std::vector< numeric_index_type > &ghost, const bool fast=false, const ParallelType=AUTOMATIC) override
 
virtual void init (const NumericVector< T > &other, const bool fast=false) override
 
virtual NumericVector< T > & operator= (const T s) override
 
virtual NumericVector< T > & operator= (const NumericVector< T > &v) override
 
virtual NumericVector< T > & operator= (const std::vector< T > &v) override
 
virtual Real min () const override
 
virtual Real max () const override
 
virtual T sum () const override
 
virtual Real l1_norm () const override
 
virtual Real l2_norm () const override
 
virtual Real linfty_norm () const override
 
virtual numeric_index_type size () const override
 
virtual numeric_index_type local_size () const override
 
virtual numeric_index_type first_local_index () const override
 
virtual numeric_index_type last_local_index () const override
 
virtual T operator() (const numeric_index_type i) const override
 
virtual NumericVector< T > & operator+= (const NumericVector< T > &v) override
 
virtual NumericVector< T > & operator-= (const NumericVector< T > &v) override
 
virtual NumericVector< T > & operator/= (const NumericVector< T > &v) override
 
virtual void reciprocal () override
 
virtual void conjugate () override
 
virtual void set (const numeric_index_type i, const T value) override
 
virtual void add (const numeric_index_type i, const T value) override
 
virtual void add (const T s) override
 
virtual void add (const NumericVector< T > &V) override
 
virtual void add (const T a, const NumericVector< T > &v) override
 
virtual void add_vector (const NumericVector< T > &, const SparseMatrix< T > &) override
 
virtual void add_vector_transpose (const NumericVector< T > &, const SparseMatrix< T > &) override
 
virtual void scale (const T factor) override
 
virtual void abs () override
 
virtual T dot (const NumericVector< T > &V) const override
 
virtual void localize (std::vector< T > &v_local) const override
 
virtual void localize (NumericVector< T > &v_local) const override
 
virtual void localize (NumericVector< T > &v_local, const std::vector< numeric_index_type > &send_list) const override
 
virtual void localize (std::vector< T > &v_local, const std::vector< numeric_index_type > &indices) const override
 
virtual void localize (const numeric_index_type first_local_idx, const numeric_index_type last_local_idx, const std::vector< numeric_index_type > &send_list) override
 
virtual void localize_to_one (std::vector< T > &v_local, const processor_id_type proc_id=0) const override
 
virtual void pointwise_mult (const NumericVector< T > &vec1, const NumericVector< T > &vec2) override
 
virtual void swap (NumericVector< T > &v) override
 
virtual bool initialized () const
 
ParallelType type () const
 
ParallelTypetype ()
 
virtual bool closed () const
 
virtual Real subset_l1_norm (const std::set< numeric_index_type > &indices) const
 
virtual Real subset_l2_norm (const std::set< numeric_index_type > &indices) const
 
virtual Real subset_linfty_norm (const std::set< numeric_index_type > &indices) const
 
virtual T el (const numeric_index_type i) const
 
virtual void get (const std::vector< numeric_index_type > &index, T *values) const
 
void get (const std::vector< numeric_index_type > &index, std::vector< T > &values) const
 
NumericVector< T > & operator*= (const T a)
 
NumericVector< T > & operator/= (const T a)
 
virtual void add_vector (const T *v, const std::vector< numeric_index_type > &dof_indices)
 
void add_vector (const std::vector< T > &v, const std::vector< numeric_index_type > &dof_indices)
 
virtual void add_vector (const NumericVector< T > &v, const std::vector< numeric_index_type > &dof_indices)
 
void add_vector (const DenseVector< T > &v, const std::vector< numeric_index_type > &dof_indices)
 
void add_vector (const NumericVector< T > &v, const ShellMatrix< T > &A)
 
virtual void insert (const T *v, const std::vector< numeric_index_type > &dof_indices)
 
void insert (const std::vector< T > &v, const std::vector< numeric_index_type > &dof_indices)
 
virtual void insert (const NumericVector< T > &v, const std::vector< numeric_index_type > &dof_indices)
 
void insert (const DenseVector< T > &v, const std::vector< numeric_index_type > &dof_indices)
 
void insert (const DenseSubVector< T > &v, const std::vector< numeric_index_type > &dof_indices)
 
virtual int compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
 
virtual int local_relative_compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
 
virtual int global_relative_compare (const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
 
virtual void print (std::ostream &os=libMesh::out) const
 
template<>
void print (std::ostream &os) const
 
virtual void print_global (std::ostream &os=libMesh::out) const
 
template<>
void print_global (std::ostream &os) const
 
virtual void print_matlab (const std::string &="") const
 
virtual void create_subvector (NumericVector< T > &, const std::vector< numeric_index_type > &) const
 
const Parallel::Communicatorcomm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static std::unique_ptr< NumericVector< T > > build (const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
 
static std::string get_info ()
 
static void print_info (std::ostream &out=libMesh::out)
 
static unsigned int n_objects ()
 
static void enable_print_counter_info ()
 
static void disable_print_counter_info ()
 

Protected Types

typedef std::map< std::string, std::pair< unsigned int, unsigned int > > Counts
 

Protected Member Functions

void increment_constructor_count (const std::string &name)
 
void increment_destructor_count (const std::string &name)
 

Protected Attributes

bool _is_closed
 
bool _is_initialized
 
ParallelType _type
 
const Parallel::Communicator_communicator
 

Static Protected Attributes

static Counts _counts
 
static Threads::atomic< unsigned int > _n_objects
 
static Threads::spin_mutex _mutex
 
static bool _enable_print_counter = true
 

Private Attributes

std::vector< T > _values
 
numeric_index_type _global_size
 
numeric_index_type _local_size
 
numeric_index_type _first_local_index
 
numeric_index_type _last_local_index
 

Detailed Description

template<typename T>
class libMesh::DistributedVector< T >

This class provides a simple parallel, distributed vector datatype which is specific to libmesh. Offers some collective communication capabilities.

Note
The class will sill function without MPI, but only on one processor. All overridden virtual functions are documented in numeric_vector.h.
Author
Benjamin S. Kirk
Date
2003

Definition at line 52 of file distributed_vector.h.

Member Typedef Documentation

◆ Counts

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts
protectedinherited

Data structure to log the information. The log is identified by the class name.

Definition at line 117 of file reference_counter.h.

Constructor & Destructor Documentation

◆ DistributedVector() [1/6]

template<typename T >
libMesh::DistributedVector< T >::DistributedVector ( const Parallel::Communicator comm,
const ParallelType  ptype = AUTOMATIC 
)
inlineexplicit

Dummy-Constructor. Dimension=0

Definition at line 260 of file distributed_vector.h.

References libMesh::NumericVector< T >::_type.

261  :
262  NumericVector<T>(comm_in, ptype),
263  _global_size (0),
264  _local_size (0),
267 {
268  this->_type = ptype;
269 }
numeric_index_type _last_local_index
numeric_index_type _global_size
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ DistributedVector() [2/6]

template<typename T >
libMesh::DistributedVector< T >::DistributedVector ( const Parallel::Communicator comm,
const numeric_index_type  n,
const ParallelType  ptype = AUTOMATIC 
)
inlineexplicit

Constructor. Set dimension to n and initialize all elements with zero.

Definition at line 275 of file distributed_vector.h.

References libMesh::DistributedVector< T >::init().

278  : NumericVector<T>(comm_in, ptype)
279 {
280  this->init(n, n, false, ptype);
281 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override

◆ DistributedVector() [3/6]

template<typename T >
libMesh::DistributedVector< T >::DistributedVector ( const Parallel::Communicator comm,
const numeric_index_type  n,
const numeric_index_type  n_local,
const ParallelType  ptype = AUTOMATIC 
)
inline

Constructor. Set local dimension to n_local, the global dimension to n, and initialize all elements with zero.

Definition at line 287 of file distributed_vector.h.

References libMesh::DistributedVector< T >::init().

291  : NumericVector<T>(comm_in, ptype)
292 {
293  this->init(n, n_local, false, ptype);
294 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override

◆ DistributedVector() [4/6]

template<typename T >
libMesh::DistributedVector< T >::DistributedVector ( const Parallel::Communicator comm,
const numeric_index_type  N,
const numeric_index_type  n_local,
const std::vector< numeric_index_type > &  ghost,
const ParallelType  ptype = AUTOMATIC 
)
inline

Constructor. Set local dimension to n_local, the global dimension to n, but additionally reserve memory for the indices specified by the ghost argument.

Definition at line 300 of file distributed_vector.h.

References libMesh::DistributedVector< T >::init().

305  : NumericVector<T>(comm_in, ptype)
306 {
307  this->init(n, n_local, ghost, false, ptype);
308 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override

◆ DistributedVector() [5/6]

template<typename T>
libMesh::DistributedVector< T >::DistributedVector ( DistributedVector< T > &&  )
default

The 5 special functions can be defaulted for this class, as it does not manage any memory itself.

◆ DistributedVector() [6/6]

template<typename T>
libMesh::DistributedVector< T >::DistributedVector ( const DistributedVector< T > &  )
default

◆ ~DistributedVector()

template<typename T>
virtual libMesh::DistributedVector< T >::~DistributedVector ( )
virtualdefault

Member Function Documentation

◆ abs()

template<typename T >
void libMesh::DistributedVector< T >::abs ( )
overridevirtual

Sets $ u_i \leftarrow |u_i| $ for each entry in the vector.

Implements libMesh::NumericVector< T >.

Definition at line 260 of file distributed_vector.C.

References std::abs(), and libMesh::initialized().

261 {
262  libmesh_assert (this->initialized());
263  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
264 
265  for (numeric_index_type i=0; i<local_size(); i++)
266  this->set(i,std::abs(_values[i]));
267 }
double abs(double a)
numeric_index_type _last_local_index
virtual bool initialized() const
virtual numeric_index_type local_size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ add() [1/4]

template<typename T >
void libMesh::DistributedVector< T >::add ( const numeric_index_type  i,
const T  value 
)
inlineoverridevirtual

Adds value to each entry of the vector.

Implements libMesh::NumericVector< T >.

Definition at line 579 of file distributed_vector.h.

References libMesh::initialized(), and value.

580 {
581  libmesh_assert (this->initialized());
582  libmesh_assert_equal_to (_values.size(), _local_size);
583  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
584  libmesh_assert_less (i, size());
585  libmesh_assert_less (i-first_local_index(), local_size());
586 
588 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual numeric_index_type local_size() const override
virtual numeric_index_type size() const override
numeric_index_type _local_size
virtual numeric_index_type first_local_index() const override
static const bool value
Definition: xdr_io.C:109
numeric_index_type _first_local_index

◆ add() [2/4]

template<typename T >
void libMesh::DistributedVector< T >::add ( const T  s)
overridevirtual

Adds s to each entry of the vector, $ u_i \leftarrow u_i + s $

Implements libMesh::NumericVector< T >.

Definition at line 212 of file distributed_vector.C.

References libMesh::initialized().

213 {
214  libmesh_assert (this->initialized());
215  libmesh_assert_equal_to (_values.size(), _local_size);
216  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
217 
218  for (numeric_index_type i=0; i<local_size(); i++)
219  _values[i] += v;
220 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual numeric_index_type local_size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ add() [3/4]

template<typename T >
void libMesh::DistributedVector< T >::add ( const NumericVector< T > &  v)
overridevirtual

Adds v to this, $ \vec{u} \leftarrow \vec{u} + \vec{v} $. Equivalent to calling operator+=().

Implements libMesh::NumericVector< T >.

Definition at line 225 of file distributed_vector.C.

References libMesh::initialized().

226 {
227  libmesh_assert (this->initialized());
228  libmesh_assert_equal_to (_values.size(), _local_size);
229  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
230 
231  add (1., v);
232 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual void add(const numeric_index_type i, const T value) override
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ add() [4/4]

template<typename T >
void libMesh::DistributedVector< T >::add ( const T  a,
const NumericVector< T > &  v 
)
overridevirtual

Vector addition with a scalar multiple, $ \vec{u} \leftarrow \vec{u} + a\vec{v} $. Equivalent to calling operator+=().

Implements libMesh::NumericVector< T >.

Definition at line 237 of file distributed_vector.C.

References libMesh::initialized().

238 {
239  libmesh_assert (this->initialized());
240  libmesh_assert_equal_to (_values.size(), _local_size);
241  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
242 
243  add(a, v);
244 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual void add(const numeric_index_type i, const T value) override
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ add_vector() [1/6]

template<typename T>
virtual void libMesh::DistributedVector< T >::add_vector ( const NumericVector< T > &  v,
const SparseMatrix< T > &  A 
)
inlineoverridevirtual

Computes $ \vec{u} \leftarrow \vec{u} + A \vec{v} $, i.e. adds the product of a SparseMatrix A and a NumericVector v to this.

Implements libMesh::NumericVector< T >.

Definition at line 191 of file distributed_vector.h.

193  { libmesh_not_implemented(); }

◆ add_vector() [2/6]

template<typename T>
void libMesh::NumericVector< T >::add_vector ( const T *  v,
const std::vector< numeric_index_type > &  dof_indices 
)
virtualinherited

Computes $ \vec{u} \leftarrow \vec{u} + \vec{v} $, where v is a pointer and each dof_indices[i] specifies where to add value v[i]. This should be overridden in subclasses for efficiency.

Reimplemented in libMesh::PetscVector< T >, and libMesh::EpetraVector< T >.

Definition at line 362 of file numeric_vector.C.

Referenced by libMesh::LinearImplicitSystem::assembly(), libMesh::NewmarkSystem::update_rhs(), and libMesh::SparseMatrix< ValOut >::vector_mult_add().

364 {
365  for (std::size_t i=0, n = dof_indices.size(); i != n; i++)
366  this->add (dof_indices[i], v[i]);
367 }
virtual void add(const numeric_index_type i, const T value)=0

◆ add_vector() [3/6]

template<typename T>
void libMesh::NumericVector< T >::add_vector ( const std::vector< T > &  v,
const std::vector< numeric_index_type > &  dof_indices 
)
inlineinherited

Computes $ \vec{u} \leftarrow \vec{u} + \vec{v} $, where v is a std::vector and each dof_indices[i] specifies where to add value v[i].

Definition at line 835 of file numeric_vector.h.

837 {
838  libmesh_assert(v.size() == dof_indices.size());
839  if (!v.empty())
840  this->add_vector(v.data(), dof_indices);
841 }
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)

◆ add_vector() [4/6]

template<typename T>
void libMesh::NumericVector< T >::add_vector ( const NumericVector< T > &  v,
const std::vector< numeric_index_type > &  dof_indices 
)
virtualinherited

Computes $ \vec{u} \leftarrow \vec{u} + \vec{v} $, where v is a NumericVector and each dof_indices[i] specifies where to add value v(i).

Definition at line 372 of file numeric_vector.C.

374 {
375  const std::size_t n = dof_indices.size();
376  libmesh_assert_equal_to(v.size(), n);
377  for (numeric_index_type i=0; i != n; i++)
378  this->add (dof_indices[i], v(i));
379 }
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual void add(const numeric_index_type i, const T value)=0

◆ add_vector() [5/6]

template<typename T>
void libMesh::NumericVector< T >::add_vector ( const DenseVector< T > &  v,
const std::vector< numeric_index_type > &  dof_indices 
)
inlineinherited

Computes $ \vec{u} \leftarrow \vec{u} + \vec{v} $, where v is a DenseVector and each dof_indices[i] specifies where to add value v(i).

Definition at line 847 of file numeric_vector.h.

849 {
850  libmesh_assert(v.size() == dof_indices.size());
851  if (!v.empty())
852  this->add_vector(&v(0), dof_indices);
853 }
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)

◆ add_vector() [6/6]

template<typename T>
void libMesh::NumericVector< T >::add_vector ( const NumericVector< T > &  v,
const ShellMatrix< T > &  A 
)
inherited

Computes $ \vec{u} \leftarrow \vec{u} + A \vec{v} $, i.e. adds the product of a ShellMatrix A and a NumericVector v to this.

Definition at line 384 of file numeric_vector.C.

386 {
387  a.vector_mult_add(*this,v);
388 }

◆ add_vector_transpose()

template<typename T>
virtual void libMesh::DistributedVector< T >::add_vector_transpose ( const NumericVector< T > &  v,
const SparseMatrix< T > &  A 
)
inlineoverridevirtual

Computes $ \vec{u} \leftarrow \vec{u} + A^T \vec{v} $, i.e. adds the product of the transpose of a SparseMatrix A and a NumericVector v to this.

Implements libMesh::NumericVector< T >.

Definition at line 195 of file distributed_vector.h.

197  { libmesh_not_implemented(); }

◆ build()

template<typename T >
std::unique_ptr< NumericVector< T > > libMesh::NumericVector< T >::build ( const Parallel::Communicator comm,
const SolverPackage  solver_package = libMesh::default_solver_package() 
)
staticinherited

Builds a NumericVector on the processors in communicator comm using the linear solver package specified by solver_package

Definition at line 48 of file numeric_vector.C.

Referenced by libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::System::add_vector(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::System::calculate_norm(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), and libMesh::CondensedEigenSystem::get_eigenpair().

49 {
50  // Build the appropriate vector
51  switch (solver_package)
52  {
53 
54 #ifdef LIBMESH_HAVE_LASPACK
55  case LASPACK_SOLVERS:
56  return libmesh_make_unique<LaspackVector<T>>(comm, AUTOMATIC);
57 #endif
58 
59 #ifdef LIBMESH_HAVE_PETSC
60  case PETSC_SOLVERS:
61  return libmesh_make_unique<PetscVector<T>>(comm, AUTOMATIC);
62 #endif
63 
64 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA
65  case TRILINOS_SOLVERS:
66  return libmesh_make_unique<EpetraVector<T>>(comm, AUTOMATIC);
67 #endif
68 
69 #ifdef LIBMESH_HAVE_EIGEN
70  case EIGEN_SOLVERS:
71  return libmesh_make_unique<EigenSparseVector<T>>(comm, AUTOMATIC);
72 #endif
73 
74  default:
75  return libmesh_make_unique<DistributedVector<T>>(comm, AUTOMATIC);
76  }
77 }
EIGEN_SOLVERS
Definition: libmesh.C:246
TRILINOS_SOLVERS
Definition: libmesh.C:244
const Parallel::Communicator & comm() const
LASPACK_SOLVERS
Definition: libmesh.C:248

◆ clear()

template<typename T >
void libMesh::DistributedVector< T >::clear ( )
inlineoverridevirtual

Restores the NumericVector<T> to a pristine state.

Reimplemented from libMesh::NumericVector< T >.

Definition at line 442 of file distributed_vector.h.

References libMesh::libMeshPrivateData::_is_initialized.

443 {
444  _values.clear();
445 
446  _global_size =
447  _local_size =
449  _last_local_index = 0;
450 
451 
452  this->_is_closed = this->_is_initialized = false;
453 }
numeric_index_type _last_local_index
numeric_index_type _global_size
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ clone()

template<typename T >
std::unique_ptr< NumericVector< T > > libMesh::DistributedVector< T >::clone ( ) const
inlineoverridevirtual
Returns
A copy of this vector wrapped in a smart pointer.
Note
This must be overridden in the derived classes.

Implements libMesh::NumericVector< T >.

Definition at line 485 of file distributed_vector.h.

References libMesh::NumericVector< T >::init().

486 {
487  NumericVector<T> * cloned_vector = new DistributedVector<T>(this->comm());
488  cloned_vector->init(*this, true);
489  *cloned_vector = *this;
490  return std::unique_ptr<NumericVector<T>>(cloned_vector);
491 }
const Parallel::Communicator & comm() const

◆ close()

template<typename T >
void libMesh::DistributedVector< T >::close ( )
inlineoverridevirtual

Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across processors.

Implements libMesh::NumericVector< T >.

Definition at line 431 of file distributed_vector.h.

References libMesh::initialized().

432 {
433  libmesh_assert (this->initialized());
434 
435  this->_is_closed = true;
436 }
virtual bool initialized() const

◆ closed()

template<typename T>
virtual bool libMesh::NumericVector< T >::closed ( ) const
inlinevirtualinherited

◆ comm()

const Parallel::Communicator& libMesh::ParallelObject::comm ( ) const
inlineinherited
Returns
A reference to the Parallel::Communicator object used by this mesh.

Definition at line 89 of file parallel_object.h.

References libMesh::ParallelObject::_communicator.

Referenced by libMesh::__libmesh_petsc_diff_solver_jacobian(), libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_tao_equality_constraints(), libMesh::__libmesh_tao_equality_constraints_jacobian(), libMesh::__libmesh_tao_gradient(), libMesh::__libmesh_tao_hessian(), libMesh::__libmesh_tao_inequality_constraints(), libMesh::__libmesh_tao_inequality_constraints_jacobian(), libMesh::__libmesh_tao_objective(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult_add(), libMesh::EquationSystems::_read_impl(), libMesh::MeshRefinement::_refine_elements(), libMesh::MeshRefinement::_smooth_flags(), libMesh::PetscDMWrapper::add_dofs_helper(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::ImplicitSystem::add_matrix(), libMesh::System::add_vector(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assemble_qoi(), libMesh::MeshCommunication::assign_global_indices(), libMesh::DofMap::attach_matrix(), libMesh::MeshTools::Generation::build_extrusion(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::PetscDMWrapper::build_section(), libMesh::PetscDMWrapper::build_sf(), libMesh::MeshBase::cache_elem_dims(), libMesh::System::calculate_norm(), libMesh::DofMap::check_dirichlet_bcid_consistency(), libMesh::PetscDMWrapper::check_section_n_dofs(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::MeshTools::create_bounding_box(), libMesh::MeshTools::create_nodal_bounding_box(), libMesh::MeshRefinement::create_parent_error_vector(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::MeshTools::create_subdomain_bounding_box(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), DMlibMeshFunction(), DMlibMeshJacobian(), DMlibMeshSetSystem_libMesh(), DMVariableBounds_libMesh(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_elem_fraction(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::MeshRefinement::flag_elements_by_nelem_target(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::DofMap::get_info(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::LocationMap< T >::init(), libMesh::TimeSolver::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::EigenSystem::init_data(), libMesh::EigenSystem::init_matrices(), libMesh::OptimizationSystem::initialize_equality_constraints_storage(), libMesh::OptimizationSystem::initialize_inequality_constraints_storage(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_flags(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_p_levels(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshTools::libmesh_assert_valid_unique_ids(), libMesh::libmesh_petsc_snes_fd_residual(), libMesh::libmesh_petsc_snes_jacobian(), libMesh::libmesh_petsc_snes_mffd_residual(), libMesh::libmesh_petsc_snes_postcheck(), libMesh::libmesh_petsc_snes_residual(), libMesh::libmesh_petsc_snes_residual_helper(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::limit_overrefined_boundary(), libMesh::MeshRefinement::limit_underrefined_boundary(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_new_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_new_nodes_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_unique_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshCommunication::make_p_levels_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::FEMSystem::mesh_position_set(), libMesh::DistributedMesh::n_active_elem(), libMesh::MeshTools::n_active_levels(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::CondensedEigenSystem::n_global_non_condensed_dofs(), libMesh::MeshTools::n_levels(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::MeshTools::n_p_levels(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::DistributedMesh::parallel_max_elem_id(), libMesh::DistributedMesh::parallel_max_node_id(), libMesh::ReplicatedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_n_elem(), libMesh::DistributedMesh::parallel_n_nodes(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::MeshTools::paranoid_n_levels(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::CheckpointIO::select_split_config(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::PetscDiffSolver::setup_petsc_data(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::split_mesh(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

90  { return _communicator; }
const Parallel::Communicator & _communicator

◆ compare()

template<typename T>
int libMesh::NumericVector< T >::compare ( const NumericVector< T > &  other_vector,
const Real  threshold = TOLERANCE 
) const
virtualinherited
Returns
-1 when this is equivalent to other_vector (up to the given threshold), or the first index where abs(a[i]-b[i]) exceeds the threshold.

Definition at line 104 of file numeric_vector.C.

106 {
107  libmesh_assert (this->initialized());
108  libmesh_assert (other_vector.initialized());
109  libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
110  libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
111 
112  int first_different_i = std::numeric_limits<int>::max();
114 
115  do
116  {
117  if (std::abs((*this)(i) - other_vector(i)) > threshold)
118  first_different_i = i;
119  else
120  i++;
121  }
122  while (first_different_i==std::numeric_limits<int>::max()
123  && i<last_local_index());
124 
125  // Find the correct first differing index in parallel
126  this->comm().min(first_different_i);
127 
128  if (first_different_i == std::numeric_limits<int>::max())
129  return -1;
130 
131  return first_different_i;
132 }
double abs(double a)
virtual bool initialized() const
const Parallel::Communicator & comm() const
long double max(long double a, double b)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type first_local_index() const =0
virtual numeric_index_type last_local_index() const =0

◆ conjugate()

template<typename T >
void libMesh::DistributedVector< T >::conjugate ( )
overridevirtual

Negates the imaginary component of each entry in the vector.

Implements libMesh::NumericVector< T >.

Definition at line 198 of file distributed_vector.C.

References libMesh::libmesh_conj().

199 {
200  for (numeric_index_type i=0; i<local_size(); i++)
201  {
202  // Replace values by complex conjugate
203  _values[i] = libmesh_conj(_values[i]);
204  }
205 }
T libmesh_conj(T a)
virtual numeric_index_type local_size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ create_subvector()

template<typename T>
virtual void libMesh::NumericVector< T >::create_subvector ( NumericVector< T > &  ,
const std::vector< numeric_index_type > &   
) const
inlinevirtualinherited

Fills in subvector from this vector using the indices in rows. Similar to the create_submatrix() routine for the SparseMatrix class, it is currently only implemented for PetscVectors.

Reimplemented in libMesh::PetscVector< T >, and libMesh::EpetraVector< T >.

Definition at line 692 of file numeric_vector.h.

694  {
695  libmesh_not_implemented();
696  }

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

Referenced by libMesh::LibMeshInit::LibMeshInit().

107 {
108  _enable_print_counter = false;
109  return;
110 }

◆ dot()

template<typename T >
T libMesh::DistributedVector< T >::dot ( const NumericVector< T > &  v) const
overridevirtual
Returns
$ \vec{u} \cdot \vec{v} $, the dot product of (*this) with the vector v.

Uses the complex-conjugate of v in the complex-valued case.

Implements libMesh::NumericVector< T >.

Definition at line 274 of file distributed_vector.C.

References libMesh::DistributedVector< T >::_values, libMesh::DistributedVector< T >::first_local_index(), and libMesh::DistributedVector< T >::last_local_index().

275 {
276  // This function must be run on all processors at once
277  parallel_object_only();
278 
279  // Make sure the NumericVector passed in is really a DistributedVector
280  const DistributedVector<T> * v = cast_ptr<const DistributedVector<T> *>(&V);
281 
282  // Make sure that the two vectors are distributed in the same way.
283  libmesh_assert_equal_to ( this->first_local_index(), v->first_local_index() );
284  libmesh_assert_equal_to ( this->last_local_index(), v->last_local_index() );
285 
286  // The result of dotting together the local parts of the vector.
287  T local_dot = 0;
288 
289  for (std::size_t i=0; i<this->local_size(); i++)
290  local_dot += this->_values[i] * v->_values[i];
291 
292  // The local dot products are now summed via MPI
293  this->comm().sum(local_dot);
294 
295  return local_dot;
296 }
virtual numeric_index_type last_local_index() const override
virtual numeric_index_type local_size() const override
virtual numeric_index_type first_local_index() const override

◆ el()

template<typename T>
virtual T libMesh::NumericVector< T >::el ( const numeric_index_type  i) const
inlinevirtualinherited
Returns
(*this)(i).

Definition at line 346 of file numeric_vector.h.

346 { return (*this)(i); }

◆ enable_print_counter_info()

void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

Methods to enable/disable the reference counter output from print_info()

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

101 {
102  _enable_print_counter = true;
103  return;
104 }

◆ first_local_index()

template<typename T >
numeric_index_type libMesh::DistributedVector< T >::first_local_index ( ) const
inlineoverridevirtual
Returns
The index of the first vector element actually stored on this processor.
Note
The minimum for this index is 0.

Implements libMesh::NumericVector< T >.

Definition at line 523 of file distributed_vector.h.

References libMesh::initialized().

Referenced by libMesh::DistributedVector< T >::dot().

524 {
525  libmesh_assert (this->initialized());
526  libmesh_assert_equal_to (_values.size(), _local_size);
527  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
528 
529  return _first_local_index;
530 }
numeric_index_type _last_local_index
virtual bool initialized() const
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ get() [1/2]

template<typename T>
void libMesh::NumericVector< T >::get ( const std::vector< numeric_index_type > &  index,
T *  values 
) const
inlinevirtualinherited

Access multiple components at once. values will not be reallocated; it should already have enough space. The default implementation calls operator() for each index, but some implementations may supply faster methods here.

Reimplemented in libMesh::PetscVector< T >.

Definition at line 806 of file numeric_vector.h.

Referenced by libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::DofMap::enforce_adjoint_constraints_exactly(), libMesh::DofMap::enforce_constraints_exactly(), libMesh::AdjointRefinementEstimator::estimate_error(), and libMesh::FEMContext::pre_fe_reinit().

808 {
809  const std::size_t num = index.size();
810  for (std::size_t i=0; i<num; i++)
811  {
812  values[i] = (*this)(index[i]);
813  }
814 }

◆ get() [2/2]

template<typename T>
void libMesh::NumericVector< T >::get ( const std::vector< numeric_index_type > &  index,
std::vector< T > &  values 
) const
inlineinherited

Access multiple components at once. values will be resized, if necessary, and filled. The default implementation calls operator() for each index, but some implementations may supply faster methods here.

Definition at line 820 of file numeric_vector.h.

822 {
823  const std::size_t num = index.size();
824  values.resize(num);
825  if (!num)
826  return;
827 
828  this->get(index, values.data());
829 }

◆ get_info()

std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().

Referenced by libMesh::ReferenceCounter::print_info().

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (const auto & pr : _counts)
59  {
60  const std::string name(pr.first);
61  const unsigned int creations = pr.second.first;
62  const unsigned int destructions = pr.second.second;
63 
64  oss << "| " << name << " reference count information:\n"
65  << "| Creations: " << creations << '\n'
66  << "| Destructions: " << destructions << '\n';
67  }
68 
69  oss << " ---------------------------------------------------------------------------- \n";
70 
71  return oss.str();
72 
73 #else
74 
75  return "";
76 
77 #endif
78 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42

◆ global_relative_compare()

template<typename T>
int libMesh::NumericVector< T >::global_relative_compare ( const NumericVector< T > &  other_vector,
const Real  threshold = TOLERANCE 
) const
virtualinherited
Returns
-1 when this is equivalent to other_vector (up to the given global relative threshold), or the first index where abs(a[i]-b[i])/max_j(a[j],b[j]) exceeds the threshold.

Definition at line 169 of file numeric_vector.C.

171 {
172  libmesh_assert (this->initialized());
173  libmesh_assert (other_vector.initialized());
174  libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
175  libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
176 
177  int first_different_i = std::numeric_limits<int>::max();
179 
180  const Real my_norm = this->linfty_norm();
181  const Real other_norm = other_vector.linfty_norm();
182  const Real abs_threshold = std::max(my_norm, other_norm) * threshold;
183 
184  do
185  {
186  if (std::abs((*this)(i) - other_vector(i) ) > abs_threshold)
187  first_different_i = i;
188  else
189  i++;
190  }
191  while (first_different_i==std::numeric_limits<int>::max()
192  && i<last_local_index());
193 
194  // Find the correct first differing index in parallel
195  this->comm().min(first_different_i);
196 
197  if (first_different_i == std::numeric_limits<int>::max())
198  return -1;
199 
200  return first_different_i;
201 }
double abs(double a)
virtual bool initialized() const
const Parallel::Communicator & comm() const
long double max(long double a, double b)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type first_local_index() const =0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual numeric_index_type last_local_index() const =0
virtual Real linfty_norm() const =0

◆ increment_constructor_count()

void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
inlineprotectedinherited

Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.

Definition at line 181 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

182 {
183  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
184  std::pair<unsigned int, unsigned int> & p = _counts[name];
185 
186  p.first++;
187 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

◆ increment_destructor_count()

void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectedinherited

Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.

Definition at line 194 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

195 {
196  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
197  std::pair<unsigned int, unsigned int> & p = _counts[name];
198 
199  p.second++;
200 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

◆ init() [1/4]

template<typename T >
void libMesh::DistributedVector< T >::init ( const numeric_index_type  n,
const numeric_index_type  n_local,
const bool  fast = false,
const ParallelType  ptype = AUTOMATIC 
)
inlineoverridevirtual

Change the dimension of the vector to n. The reserved memory for this vector remains unchanged if possible. If n==0, all memory is freed. Therefore, if you want to resize the vector and release the memory not needed, you have to first call init(0) and then init(n). This behaviour is analogous to that of the STL containers.

On fast==false, the vector is filled by zeros.

Implements libMesh::NumericVector< T >.

Definition at line 314 of file distributed_vector.h.

References libMesh::libMeshPrivateData::_is_initialized, libMesh::AUTOMATIC, libMesh::initialized(), libMesh::PARALLEL, libMesh::SERIAL, and libMesh::zero.

Referenced by libMesh::DistributedVector< T >::DistributedVector(), and libMesh::DistributedVector< T >::localize().

318 {
319  // This function must be run on all processors at once
320  parallel_object_only();
321 
322  libmesh_assert_less_equal (n_local, n);
323 
324  if (ptype == AUTOMATIC)
325  {
326  if (n == n_local)
327  this->_type = SERIAL;
328  else
329  this->_type = PARALLEL;
330  }
331  else
332  this->_type = ptype;
333 
334  libmesh_assert ((this->_type==SERIAL && n==n_local) ||
335  this->_type==PARALLEL);
336 
337  // Clear the data structures if already initialized
338  if (this->initialized())
339  this->clear();
340 
341  // Initialize data structures
342  _values.resize(n_local);
343  _local_size = n_local;
344  _global_size = n;
345 
346  _first_local_index = 0;
347 
348 #ifdef LIBMESH_HAVE_MPI
349 
350  std::vector<numeric_index_type> local_sizes (this->n_processors(), 0);
351 
352  local_sizes[this->processor_id()] = n_local;
353 
354  this->comm().sum(local_sizes);
355 
356  // _first_local_index is the sum of _local_size
357  // for all processor ids less than ours
358  for (processor_id_type p=0; p!=this->processor_id(); p++)
359  _first_local_index += local_sizes[p];
360 
361 
362 # ifdef DEBUG
363  // Make sure all the local sizes sum up to the global
364  // size, otherwise there is big trouble!
365  numeric_index_type dbg_sum=0;
366 
367  for (processor_id_type p=0; p!=this->n_processors(); p++)
368  dbg_sum += local_sizes[p];
369 
370  libmesh_assert_equal_to (dbg_sum, n);
371 
372 # endif
373 
374 #else
375 
376  // No other options without MPI!
377  if (n != n_local)
378  libmesh_error_msg("ERROR: MPI is required for n != n_local!");
379 
380 #endif
381 
383 
384  // Set the initialized flag
385  this->_is_initialized = true;
386 
387  // Zero the components unless directed otherwise
388  if (!fast)
389  this->zero();
390 }
numeric_index_type _last_local_index
virtual bool initialized() const
uint8_t processor_id_type
Definition: id_types.h:99
const Parallel::Communicator & comm() const
numeric_index_type _global_size
virtual void zero() override
processor_id_type n_processors() const
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
virtual void clear() override
processor_id_type processor_id() const
numeric_index_type _first_local_index

◆ init() [2/4]

template<typename T >
void libMesh::DistributedVector< T >::init ( const numeric_index_type  n,
const bool  fast = false,
const ParallelType  ptype = AUTOMATIC 
)
inlineoverridevirtual

Call init() with n_local = N.

Implements libMesh::NumericVector< T >.

Definition at line 420 of file distributed_vector.h.

References libMesh::TriangleWrapper::init().

423 {
424  this->init(n,n,fast,ptype);
425 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override

◆ init() [3/4]

template<typename T >
void libMesh::DistributedVector< T >::init ( const numeric_index_type  n,
const numeric_index_type  n_local,
const std::vector< numeric_index_type > &  ghost,
const bool  fast = false,
const ParallelType  ptype = AUTOMATIC 
)
inlineoverridevirtual

Create a vector that holds tha local indices plus those specified in the ghost argument.

Implements libMesh::NumericVector< T >.

Definition at line 395 of file distributed_vector.h.

References libMesh::TriangleWrapper::init().

400 {
401  // TODO: we shouldn't ignore the ghost sparsity pattern
402  this->init(n, n_local, fast, ptype);
403 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override

◆ init() [4/4]

template<class T >
void libMesh::DistributedVector< T >::init ( const NumericVector< T > &  other,
const bool  fast = false 
)
overridevirtual

Creates a vector that has the same dimension and storage type as other, including ghost dofs.

Implements libMesh::NumericVector< T >.

Definition at line 410 of file distributed_vector.h.

References libMesh::TriangleWrapper::init(), libMesh::NumericVector< T >::local_size(), libMesh::NumericVector< T >::size(), and libMesh::NumericVector< T >::type().

412 {
413  this->init(other.size(),other.local_size(),fast,other.type());
414 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override

◆ initialized()

template<typename T>
virtual bool libMesh::NumericVector< T >::initialized ( ) const
inlinevirtualinherited

◆ insert() [1/5]

template<typename T>
void libMesh::NumericVector< T >::insert ( const T *  v,
const std::vector< numeric_index_type > &  dof_indices 
)
virtualinherited

Inserts the entries of v in *this at the locations specified by v.

Reimplemented in libMesh::PetscVector< T >, and libMesh::EpetraVector< T >.

Definition at line 82 of file numeric_vector.C.

84 {
85  for (numeric_index_type i=0; i<dof_indices.size(); i++)
86  this->set (dof_indices[i], v[i]);
87 }
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ insert() [2/5]

template<typename T>
void libMesh::NumericVector< T >::insert ( const std::vector< T > &  v,
const std::vector< numeric_index_type > &  dof_indices 
)
inlineinherited

Inserts the entries of v in *this at the locations specified by v.

Definition at line 859 of file numeric_vector.h.

861 {
862  libmesh_assert(v.size() == dof_indices.size());
863  if (!v.empty())
864  this->insert(v.data(), dof_indices);
865 }
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices)

◆ insert() [3/5]

template<typename T>
void libMesh::NumericVector< T >::insert ( const NumericVector< T > &  v,
const std::vector< numeric_index_type > &  dof_indices 
)
virtualinherited

Inserts the entries of v in *this at the locations specified by v.

Definition at line 92 of file numeric_vector.C.

94 {
95  libmesh_assert_equal_to (V.size(), dof_indices.size());
96 
97  for (numeric_index_type i=0; i<dof_indices.size(); i++)
98  this->set (dof_indices[i], V(i));
99 }
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ insert() [4/5]

template<typename T>
void libMesh::NumericVector< T >::insert ( const DenseVector< T > &  v,
const std::vector< numeric_index_type > &  dof_indices 
)
inlineinherited

Inserts the entries of v in *this at the locations specified by v.

Definition at line 871 of file numeric_vector.h.

873 {
874  libmesh_assert(v.size() == dof_indices.size());
875  if (!v.empty())
876  this->insert(&v(0), dof_indices);
877 }
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices)

◆ insert() [5/5]

template<typename T>
void libMesh::NumericVector< T >::insert ( const DenseSubVector< T > &  v,
const std::vector< numeric_index_type > &  dof_indices 
)
inlineinherited

Inserts the entries of v in *this at the locations specified by v.

Definition at line 883 of file numeric_vector.h.

885 {
886  libmesh_assert(v.size() == dof_indices.size());
887  if (!v.empty())
888  this->insert(&v(0), dof_indices);
889 }
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices)

◆ l1_norm()

template<typename T >
Real libMesh::DistributedVector< T >::l1_norm ( ) const
overridevirtual
Returns
The $ \ell_1 $-norm of the vector, i.e. the sum of the absolute values of the entries.

Implements libMesh::NumericVector< T >.

Definition at line 65 of file distributed_vector.C.

References std::abs(), and libMesh::initialized().

66 {
67  // This function must be run on all processors at once
68  parallel_object_only();
69 
70  libmesh_assert (this->initialized());
71  libmesh_assert_equal_to (_values.size(), _local_size);
72  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
73 
74  double local_l1 = 0.;
75 
76  for (numeric_index_type i=0; i<local_size(); i++)
77  local_l1 += std::abs(_values[i]);
78 
79  this->comm().sum(local_l1);
80 
81  return local_l1;
82 }
double abs(double a)
numeric_index_type _last_local_index
virtual bool initialized() const
const Parallel::Communicator & comm() const
virtual numeric_index_type local_size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ l2_norm()

template<typename T >
Real libMesh::DistributedVector< T >::l2_norm ( ) const
overridevirtual
Returns
The $ \ell_2 $-norm of the vector, i.e. the square root of the sum of the squares of the entries.

Implements libMesh::NumericVector< T >.

Definition at line 87 of file distributed_vector.C.

References libMesh::initialized(), and libMesh::TensorTools::norm_sq().

88 {
89  // This function must be run on all processors at once
90  parallel_object_only();
91 
92  libmesh_assert (this->initialized());
93  libmesh_assert_equal_to (_values.size(), _local_size);
94  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
95 
96  double local_l2 = 0.;
97 
98  for (numeric_index_type i=0; i<local_size(); i++)
99  local_l2 += TensorTools::norm_sq(_values[i]);
100 
101  this->comm().sum(local_l2);
102 
103  return std::sqrt(local_l2);
104 }
numeric_index_type _last_local_index
virtual bool initialized() const
const Parallel::Communicator & comm() const
virtual numeric_index_type local_size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ last_local_index()

template<typename T >
numeric_index_type libMesh::DistributedVector< T >::last_local_index ( ) const
inlineoverridevirtual
Returns
The index+1 of the last vector element actually stored on this processor.
Note
The maximum for this index is size().

Implements libMesh::NumericVector< T >.

Definition at line 536 of file distributed_vector.h.

References libMesh::initialized().

Referenced by libMesh::DistributedVector< T >::dot().

537 {
538  libmesh_assert (this->initialized());
539  libmesh_assert_equal_to (_values.size(), _local_size);
540  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
541 
542  return _last_local_index;
543 }
numeric_index_type _last_local_index
virtual bool initialized() const
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ linfty_norm()

template<typename T >
Real libMesh::DistributedVector< T >::linfty_norm ( ) const
overridevirtual
Returns
The $ \ell_{\infty} $-norm of the vector, i.e. the maximum absolute value of the entries of the vector.

Implements libMesh::NumericVector< T >.

Definition at line 109 of file distributed_vector.C.

References std::abs(), libMesh::initialized(), std::max(), and libMesh::Real.

110 {
111  // This function must be run on all processors at once
112  parallel_object_only();
113 
114  libmesh_assert (this->initialized());
115  libmesh_assert_equal_to (_values.size(), _local_size);
116  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
117 
118  Real local_linfty = 0.;
119 
120  for (numeric_index_type i=0; i<local_size(); i++)
121  local_linfty = std::max(local_linfty,
122  static_cast<Real>(std::abs(_values[i]))
123  ); // Note we static_cast so that both
124  // types are the same, as required
125  // by std::max
126 
127  this->comm().max(local_linfty);
128 
129  return local_linfty;
130 }
double abs(double a)
numeric_index_type _last_local_index
virtual bool initialized() const
const Parallel::Communicator & comm() const
long double max(long double a, double b)
virtual numeric_index_type local_size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
numeric_index_type _first_local_index

◆ local_relative_compare()

template<typename T>
int libMesh::NumericVector< T >::local_relative_compare ( const NumericVector< T > &  other_vector,
const Real  threshold = TOLERANCE 
) const
virtualinherited
Returns
-1 when this is equivalent to other_vector, (up to the given local relative threshold), or the first index where abs(a[i]-b[i])/max(a[i],b[i]) exceeds the threshold.

Definition at line 136 of file numeric_vector.C.

138 {
139  libmesh_assert (this->initialized());
140  libmesh_assert (other_vector.initialized());
141  libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
142  libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
143 
144  int first_different_i = std::numeric_limits<int>::max();
146 
147  do
148  {
149  if (std::abs((*this)(i) - other_vector(i)) > threshold *
150  std::max(std::abs((*this)(i)), std::abs(other_vector(i))))
151  first_different_i = i;
152  else
153  i++;
154  }
155  while (first_different_i==std::numeric_limits<int>::max()
156  && i<last_local_index());
157 
158  // Find the correct first differing index in parallel
159  this->comm().min(first_different_i);
160 
161  if (first_different_i == std::numeric_limits<int>::max())
162  return -1;
163 
164  return first_different_i;
165 }
double abs(double a)
virtual bool initialized() const
const Parallel::Communicator & comm() const
long double max(long double a, double b)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type first_local_index() const =0
virtual numeric_index_type last_local_index() const =0

◆ local_size()

template<typename T >
numeric_index_type libMesh::DistributedVector< T >::local_size ( ) const
inlineoverridevirtual
Returns
The local size of the vector, i.e. index_stop - index_start.

Implements libMesh::NumericVector< T >.

Definition at line 510 of file distributed_vector.h.

References libMesh::initialized().

Referenced by libMesh::DistributedVector< T >::operator=().

511 {
512  libmesh_assert (this->initialized());
513  libmesh_assert_equal_to (_values.size(), _local_size);
514  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
515 
516  return _local_size;
517 }
numeric_index_type _last_local_index
virtual bool initialized() const
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ localize() [1/5]

template<typename T >
void libMesh::DistributedVector< T >::localize ( std::vector< T > &  v_local) const
overridevirtual

Creates a copy of the global vector in the local vector v_local.

Implements libMesh::NumericVector< T >.

Definition at line 550 of file distributed_vector.C.

References libMesh::initialized().

Referenced by libMesh::DistributedVector< T >::localize().

551 {
552  // This function must be run on all processors at once
553  parallel_object_only();
554 
555  libmesh_assert (this->initialized());
556  libmesh_assert_equal_to (_values.size(), _local_size);
557  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
558 
559  v_local = this->_values;
560 
561  this->comm().allgather (v_local);
562 
563 #ifndef LIBMESH_HAVE_MPI
564  libmesh_assert_equal_to (local_size(), size());
565 #endif
566 }
numeric_index_type _last_local_index
virtual bool initialized() const
void allgather(const T &send, std::vector< T, A > &recv) const
const Parallel::Communicator & comm() const
virtual numeric_index_type local_size() const override
virtual numeric_index_type size() const override
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ localize() [2/5]

template<typename T >
void libMesh::DistributedVector< T >::localize ( NumericVector< T > &  v_local) const
overridevirtual

Same, but fills a NumericVector<T> instead of a std::vector.

Implements libMesh::NumericVector< T >.

Definition at line 377 of file distributed_vector.C.

References libMesh::DistributedVector< T >::_first_local_index, libMesh::DistributedVector< T >::_global_size, libMesh::NumericVector< T >::_is_closed, libMesh::NumericVector< T >::_is_initialized, libMesh::DistributedVector< T >::_last_local_index, libMesh::DistributedVector< T >::_local_size, libMesh::DistributedVector< T >::_values, and libMesh::initialized().

379 {
380  libmesh_assert (this->initialized());
381  libmesh_assert_equal_to (_values.size(), _local_size);
382  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
383 
384  DistributedVector<T> * v_local = cast_ptr<DistributedVector<T> *>(&v_local_in);
385 
386  v_local->_first_local_index = 0;
387 
388  v_local->_global_size =
389  v_local->_local_size =
390  v_local->_last_local_index = size();
391 
392  v_local->_is_initialized =
393  v_local->_is_closed = true;
394 
395  // Call localize on the vector's values. This will help
396  // prevent code duplication
397  localize (v_local->_values);
398 
399 #ifndef LIBMESH_HAVE_MPI
400 
401  libmesh_assert_equal_to (local_size(), size());
402 
403 #endif
404 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual numeric_index_type local_size() const override
virtual numeric_index_type size() const override
numeric_index_type _local_size
virtual void localize(std::vector< T > &v_local) const override
numeric_index_type _first_local_index

◆ localize() [3/5]

template<typename T >
void libMesh::DistributedVector< T >::localize ( NumericVector< T > &  v_local,
const std::vector< numeric_index_type > &  send_list 
) const
overridevirtual

Creates a local vector v_local containing only information relevant to this processor, as defined by the send_list.

Implements libMesh::NumericVector< T >.

Definition at line 409 of file distributed_vector.C.

References libMesh::initialized().

411 {
412  libmesh_assert (this->initialized());
413  libmesh_assert_equal_to (_values.size(), _local_size);
414  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
415 
416  // TODO: We don't yet support the send list; this is inefficient:
417  localize (v_local_in);
418 }
numeric_index_type _last_local_index
virtual bool initialized() const
numeric_index_type _local_size
virtual void localize(std::vector< T > &v_local) const override
numeric_index_type _first_local_index

◆ localize() [4/5]

template<typename T >
void libMesh::DistributedVector< T >::localize ( std::vector< T > &  v_local,
const std::vector< numeric_index_type > &  indices 
) const
overridevirtual

Fill in the local std::vector "v_local" with the global indices given in "indices".

Note
The indices can be different on every processor, and the same index can be localized to more than one processor. The resulting v_local can be shorter than the original, and the entries will be in the order specified by indices.

Example:

*   On 4 procs *this = {a, b, c, d, e, f, g, h, i} is a parallel vector.
*   On each proc, the indices arrays are set up as:
*   proc0, indices = {1,2,4,5}
*   proc1, indices = {2,5,6,8}
*   proc2, indices = {2,3,6,7}
*   proc3, indices = {0,1,2,3}
*
*   After calling this version of localize, the v_local vectors are:
*   proc0, v_local = {b,c,e,f}
*   proc1, v_local = {c,f,g,i}
*   proc2, v_local = {c,d,g,h}
*   proc3, v_local = {a,b,c,d}
* 

This function is useful in parallel I/O routines, when you have a parallel vector of solution values which you want to write a subset of.

Implements libMesh::NumericVector< T >.

Definition at line 423 of file distributed_vector.C.

References libMesh::Parallel::pull_parallel_vector_data().

425 {
426  // Resize v_local so there is enough room to hold all the local values.
427  v_local.resize(indices.size());
428 
429  // We need to know who has the values we want, so get everyone's _local_size
430  std::vector<numeric_index_type> local_sizes;
431  this->comm().allgather (_local_size, local_sizes);
432 
433  // Make a vector of partial sums of local sizes
434  std::vector<numeric_index_type> local_size_sums(this->n_processors());
435  local_size_sums[0] = local_sizes[0];
436  for (numeric_index_type i=1; i<local_sizes.size(); i++)
437  local_size_sums[i] = local_size_sums[i-1] + local_sizes[i];
438 
439  // We now fill in 'requested_ids' based on the indices. Also keep
440  // track of the local index (in the indices vector) for each of
441  // these, since we need that when unpacking.
442  std::map<processor_id_type, std::vector<numeric_index_type>>
443  requested_ids, local_requested_ids;
444 
445  // We'll use this typedef a couple of times below.
446  typedef typename std::vector<numeric_index_type>::iterator iter_t;
447 
448  // For each index in indices, determine which processor it is on.
449  // This is an O(N*log(p)) algorithm that uses std::upper_bound().
450  // Note: upper_bound() returns an iterator to the first entry which is
451  // greater than the given value.
452  for (numeric_index_type i=0; i<indices.size(); i++)
453  {
454  iter_t ub = std::upper_bound(local_size_sums.begin(),
455  local_size_sums.end(),
456  indices[i]);
457 
458  processor_id_type on_proc = cast_int<processor_id_type>
459  (std::distance(local_size_sums.begin(), ub));
460 
461  requested_ids[on_proc].push_back(indices[i]);
462  local_requested_ids[on_proc].push_back(i);
463  }
464 
465  auto gather_functor =
466  [this]
467  (processor_id_type, const std::vector<dof_id_type> & ids,
468  std::vector<T> & values)
469  {
470  // The first send/receive we did was for indices, the second one will be
471  // for corresponding floating point values, so create storage for that now...
472  const std::size_t ids_size = ids.size();
473  values.resize(ids_size);
474 
475  for (std::size_t i=0; i != ids_size; i++)
476  {
477  // The index of the requested value
478  const numeric_index_type requested_index = ids[i];
479 
480  // Transform into local numbering, and get requested value.
481  values[i] = _values[requested_index - _first_local_index];
482  }
483  };
484 
485  auto action_functor =
486  [& v_local, & local_requested_ids]
487  (processor_id_type pid,
488  const std::vector<dof_id_type> &,
489  const std::vector<T> & values)
490  {
491  // Now write the received values to the appropriate place(s) in v_local
492  for (std::size_t i=0, vsize = values.size(); i != vsize; i++)
493  {
494  libmesh_assert(local_requested_ids.count(pid));
495  libmesh_assert_less(i, local_requested_ids[pid].size());
496 
497  // Get the index in v_local where this value needs to be inserted.
498  const numeric_index_type local_requested_index =
499  local_requested_ids[pid][i];
500 
501  // Actually set the value in v_local
502  v_local[local_requested_index] = values[i];
503  }
504  };
505 
506  const T * ex = nullptr;
508  (this->comm(), requested_ids, gather_functor, action_functor, ex);
509 }
void allgather(const T &send, std::vector< T, A > &recv) const
uint8_t processor_id_type
Definition: id_types.h:99
const Parallel::Communicator & comm() const
processor_id_type n_processors() const
dof_id_type numeric_index_type
Definition: id_types.h:92
void pull_parallel_vector_data(const Communicator &comm, const MapToVectors &queries, RequestContainer &reqs, GatherFunctor &gather_data, ActionFunctor &act_on_data, const datum *example)
virtual numeric_index_type size() const override
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ localize() [5/5]

template<typename T >
void libMesh::DistributedVector< T >::localize ( const numeric_index_type  first_local_idx,
const numeric_index_type  last_local_idx,
const std::vector< numeric_index_type > &  send_list 
)
overridevirtual

Updates a local vector with selected values from neighboring processors, as defined by send_list.

Implements libMesh::NumericVector< T >.

Definition at line 514 of file distributed_vector.C.

References libMesh::DistributedVector< T >::_values, libMesh::DistributedVector< T >::init(), libMesh::DistributedVector< T >::localize(), and libMesh::PARALLEL.

517 {
518  // Only good for serial vectors
519  libmesh_assert_equal_to (this->size(), this->local_size());
520  libmesh_assert_greater (last_local_idx, first_local_idx);
521  libmesh_assert_less_equal (send_list.size(), this->size());
522  libmesh_assert_less (last_local_idx, this->size());
523 
524  const numeric_index_type my_size = this->size();
525  const numeric_index_type my_local_size = (last_local_idx - first_local_idx + 1);
526 
527  // Don't bother for serial cases
528  if ((first_local_idx == 0) &&
529  (my_local_size == my_size))
530  return;
531 
532 
533  // Build a parallel vector, initialize it with the local
534  // parts of (*this)
535  DistributedVector<T> parallel_vec(this->comm());
536 
537  parallel_vec.init (my_size, my_local_size, true, PARALLEL);
538 
539  // Copy part of *this into the parallel_vec
540  for (numeric_index_type i=first_local_idx; i<=last_local_idx; i++)
541  parallel_vec._values[i-first_local_idx] = _values[i];
542 
543  // localize like normal
544  parallel_vec.localize (*this, send_list);
545 }
const Parallel::Communicator & comm() const
virtual numeric_index_type local_size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type size() const override

◆ localize_to_one()

template<typename T >
void libMesh::DistributedVector< T >::localize_to_one ( std::vector< T > &  v_local,
const processor_id_type  proc_id = 0 
) const
overridevirtual

Creates a local copy of the global vector in v_local only on processor proc_id. By default the data is sent to processor 0. This method is useful for outputting data from one processor.

Implements libMesh::NumericVector< T >.

Definition at line 571 of file distributed_vector.C.

References libMesh::initialized().

573 {
574  // This function must be run on all processors at once
575  parallel_object_only();
576 
577  libmesh_assert (this->initialized());
578  libmesh_assert_equal_to (_values.size(), _local_size);
579  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
580 
581  v_local = this->_values;
582 
583  this->comm().gather (pid, v_local);
584 
585 #ifndef LIBMESH_HAVE_MPI
586  libmesh_assert_equal_to (local_size(), size());
587 #endif
588 }
numeric_index_type _last_local_index
virtual bool initialized() const
void gather(const unsigned int root_id, const T &send, std::vector< T, A > &recv) const
const Parallel::Communicator & comm() const
virtual numeric_index_type local_size() const override
virtual numeric_index_type size() const override
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ max()

template<typename T >
Real libMesh::DistributedVector< T >::max ( ) const
inlineoverridevirtual
Returns
The maximum entry in the vector, or the maximum real part in the case of complex numbers.

Implements libMesh::NumericVector< T >.

Definition at line 617 of file distributed_vector.h.

References libMesh::initialized(), libMesh::libmesh_real(), std::max(), and libMesh::Real.

618 {
619  // This function must be run on all processors at once
620  parallel_object_only();
621 
622  libmesh_assert (this->initialized());
623  libmesh_assert_equal_to (_values.size(), _local_size);
624  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
625 
626  Real local_max = _values.size() ?
627  libmesh_real(_values[0]) : -std::numeric_limits<Real>::max();
628  for (numeric_index_type i = 1; i < _values.size(); ++i)
629  local_max = std::max(libmesh_real(_values[i]), local_max);
630 
631  this->comm().max(local_max);
632 
633  return local_max;
634 }
T libmesh_real(T a)
virtual Real max() const override
numeric_index_type _last_local_index
virtual bool initialized() const
const Parallel::Communicator & comm() const
long double max(long double a, double b)
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
numeric_index_type _first_local_index

◆ min()

template<typename T >
Real libMesh::DistributedVector< T >::min ( ) const
inlineoverridevirtual
Returns
The minimum entry in the vector, or the minimum real part in the case of complex numbers.

Implements libMesh::NumericVector< T >.

Definition at line 594 of file distributed_vector.h.

References libMesh::initialized(), libMesh::libmesh_real(), std::max(), std::min(), and libMesh::Real.

595 {
596  // This function must be run on all processors at once
597  parallel_object_only();
598 
599  libmesh_assert (this->initialized());
600  libmesh_assert_equal_to (_values.size(), _local_size);
601  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
602 
603  Real local_min = _values.size() ?
604  libmesh_real(_values[0]) : std::numeric_limits<Real>::max();
605  for (numeric_index_type i = 1; i < _values.size(); ++i)
606  local_min = std::min(libmesh_real(_values[i]), local_min);
607 
608  this->comm().min(local_min);
609 
610  return local_min;
611 }
T libmesh_real(T a)
virtual Real max() const override
numeric_index_type _last_local_index
virtual bool initialized() const
const Parallel::Communicator & comm() const
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
long double min(long double a, double b)
numeric_index_type _first_local_index

◆ n_objects()

static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 83 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

84  { return _n_objects; }
static Threads::atomic< unsigned int > _n_objects

◆ n_processors()

processor_id_type libMesh::ParallelObject::n_processors ( ) const
inlineinherited
Returns
The number of processors in the group.

Definition at line 95 of file parallel_object.h.

References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::size().

Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::DistributedMesh::add_node(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assembly(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshBase::get_info(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::DistributedMesh::insert_elem(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::MeshBase::partition(), libMesh::PetscLinearSolver< T >::PetscLinearSolver(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::CheckpointIO::read_connectivity(), libMesh::XdrIO::read_header(), libMesh::CheckpointIO::read_nodes(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::System::read_serialized_vector(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

96  { return cast_int<processor_id_type>(_communicator.size()); }
processor_id_type size() const
Definition: communicator.h:175
const Parallel::Communicator & _communicator

◆ operator()()

template<typename T >
T libMesh::DistributedVector< T >::operator() ( const numeric_index_type  i) const
inlineoverridevirtual
Returns
A copy of the ith entry of the vector.

Implements libMesh::NumericVector< T >.

Definition at line 549 of file distributed_vector.h.

References libMesh::initialized().

550 {
551  libmesh_assert (this->initialized());
552  libmesh_assert_equal_to (_values.size(), _local_size);
553  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
554  libmesh_assert ( ((i >= first_local_index()) &&
555  (i < last_local_index())) );
556 
557  return _values[i - _first_local_index];
558 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual numeric_index_type last_local_index() const override
numeric_index_type _local_size
virtual numeric_index_type first_local_index() const override
numeric_index_type _first_local_index

◆ operator*=()

template<typename T>
NumericVector<T>& libMesh::NumericVector< T >::operator*= ( const T  a)
inlineinherited

Scales the vector by a, $ \vec{u} \leftarrow a\vec{u} $. Equivalent to u.scale(a)

Returns
A reference to *this.

Definition at line 391 of file numeric_vector.h.

391 { this->scale(a); return *this; }
virtual void scale(const T factor)=0

◆ operator+=()

template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator+= ( const NumericVector< T > &  v)
overridevirtual

Adds v to *this, $ \vec{u} \leftarrow \vec{u} + \vec{v} $. Equivalent to u.add(1, v).

Returns
A reference to *this.

Implements libMesh::NumericVector< T >.

Definition at line 135 of file distributed_vector.C.

References libMesh::closed(), and libMesh::initialized().

136 {
137  libmesh_assert (this->closed());
138  libmesh_assert (this->initialized());
139  libmesh_assert_equal_to (_values.size(), _local_size);
140  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
141 
142  add(1., v);
143 
144  return *this;
145 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual void add(const numeric_index_type i, const T value) override
virtual bool closed() const
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ operator-=()

template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator-= ( const NumericVector< T > &  v)
overridevirtual

Subtracts v from *this, $ \vec{u} \leftarrow \vec{u} - \vec{v} $. Equivalent to u.add(-1, v).

Returns
A reference to *this.

Implements libMesh::NumericVector< T >.

Definition at line 150 of file distributed_vector.C.

References libMesh::closed(), and libMesh::initialized().

151 {
152  libmesh_assert (this->closed());
153  libmesh_assert (this->initialized());
154  libmesh_assert_equal_to (_values.size(), _local_size);
155  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
156 
157  add(-1., v);
158 
159  return *this;
160 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual void add(const numeric_index_type i, const T value) override
virtual bool closed() const
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ operator/=() [1/2]

template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator/= ( const NumericVector< T > &  )
overridevirtual

Computes the pointwise division of this vector's entries by another's, $ u_i \leftarrow \frac{u_i}{v_i} \, \forall i$

Returns
A reference to *this.

Implements libMesh::NumericVector< T >.

Definition at line 165 of file distributed_vector.C.

References libMesh::DistributedVector< T >::size(), and libMesh::NumericVector< T >::size().

166 {
167  libmesh_assert_equal_to(size(), v.size());
168 
169  const DistributedVector<T> & v_vec = cast_ref<const DistributedVector<T> &>(v);
170 
171  std::size_t val_size = _values.size();
172 
173  for (std::size_t i=0; i<val_size; i++)
174  _values[i] = _values[i] / v_vec._values[i];
175 
176  return *this;
177 }
virtual numeric_index_type size() const override

◆ operator/=() [2/2]

template<typename T>
NumericVector<T>& libMesh::NumericVector< T >::operator/= ( const T  a)
inlineinherited

Scales the vector by 1/a, $ \vec{u} \leftarrow \frac{1}{a}\vec{u} $. Equivalent to u.scale(1./a)

Returns
A reference to *this.

Definition at line 400 of file numeric_vector.h.

400 { this->scale(1./a); return *this; }
virtual void scale(const T factor)=0

◆ operator=() [1/5]

template<typename T >
DistributedVector< T > & libMesh::DistributedVector< T >::operator= ( const DistributedVector< T > &  v)

Copy assignment operator. We cannot default this (although it essentially implements the default behavior) because the compiler-generated default attempts to automatically call the base class (NumericVector) copy assignment operator, which we have chosen to make pure virtual for other design reasons.

Definition at line 332 of file distributed_vector.C.

References libMesh::DistributedVector< T >::_first_local_index, libMesh::DistributedVector< T >::_global_size, libMesh::NumericVector< T >::_is_closed, libMesh::libMeshPrivateData::_is_initialized, libMesh::NumericVector< T >::_is_initialized, libMesh::DistributedVector< T >::_last_local_index, libMesh::DistributedVector< T >::_local_size, libMesh::DistributedVector< T >::_values, and libMesh::DistributedVector< T >::local_size().

333 {
334  this->_is_initialized = v._is_initialized;
335  this->_is_closed = v._is_closed;
336 
337  _global_size = v._global_size;
338  _local_size = v._local_size;
339  _first_local_index = v._first_local_index;
340  _last_local_index = v._last_local_index;
341 
342  if (v.local_size() == this->local_size())
343  _values = v._values;
344 
345  else
346  libmesh_error_msg("v.local_size() = " << v.local_size() << " must be equal to this->local_size() = " << this->local_size());
347 
348  return *this;
349 }
numeric_index_type _last_local_index
numeric_index_type _global_size
virtual numeric_index_type local_size() const override
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ operator=() [2/5]

template<typename T>
DistributedVector& libMesh::DistributedVector< T >::operator= ( DistributedVector< T > &&  )
default

◆ operator=() [3/5]

template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator= ( const T  s)
overridevirtual

Sets all entries of the vector to the value s.

Returns
A reference to *this.

Implements libMesh::NumericVector< T >.

Definition at line 302 of file distributed_vector.C.

References libMesh::initialized().

303 {
304  libmesh_assert (this->initialized());
305  libmesh_assert_equal_to (_values.size(), _local_size);
306  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
307 
308  for (std::size_t i=0; i<local_size(); i++)
309  _values[i] = s;
310 
311  return *this;
312 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual numeric_index_type local_size() const override
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ operator=() [4/5]

template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator= ( const NumericVector< T > &  v)
overridevirtual

This looks like a copy assignment operator, but note that, unlike normal copy assignment operators, it is pure virtual. This function should be overridden in derived classes so that they can be copied correctly via references to the base class. This design usually isn't a good idea in general, but in this context it works because we usually don't have a mix of different kinds of NumericVectors active in the library at a single time.

Returns
A reference to *this as the base type.

Implements libMesh::NumericVector< T >.

Definition at line 318 of file distributed_vector.C.

319 {
320  // Make sure the NumericVector passed in is really a DistributedVector
321  const DistributedVector<T> * v = cast_ptr<const DistributedVector<T> *>(&v_in);
322 
323  *this = *v;
324 
325  return *this;
326 }

◆ operator=() [5/5]

template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator= ( const std::vector< T > &  v)
overridevirtual

Sets (*this)(i) = v(i) for each entry of the vector.

Returns
A reference to *this as the base type.

Implements libMesh::NumericVector< T >.

Definition at line 355 of file distributed_vector.C.

References libMesh::initialized().

356 {
357  libmesh_assert (this->initialized());
358  libmesh_assert_equal_to (_values.size(), _local_size);
359  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
360 
361  if (v.size() == local_size())
362  _values = v;
363 
364  else if (v.size() == size())
365  for (std::size_t i=first_local_index(); i<last_local_index(); i++)
366  _values[i-first_local_index()] = v[i];
367 
368  else
369  libmesh_error_msg("Incompatible sizes in DistributedVector::operator=");
370 
371  return *this;
372 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual numeric_index_type last_local_index() const override
virtual numeric_index_type local_size() const override
virtual numeric_index_type size() const override
numeric_index_type _local_size
virtual numeric_index_type first_local_index() const override
numeric_index_type _first_local_index

◆ pointwise_mult()

template<typename T >
void libMesh::DistributedVector< T >::pointwise_mult ( const NumericVector< T > &  vec1,
const NumericVector< T > &  vec2 
)
overridevirtual

Computes $ u_i \leftarrow u_i v_i $ (summation not implied) i.e. the pointwise (component-wise) product of vec1 and vec2, and stores the result in *this.

Implements libMesh::NumericVector< T >.

Definition at line 593 of file distributed_vector.C.

597 {
598  libmesh_not_implemented();
599 }

◆ print() [1/2]

template<typename T >
void libMesh::NumericVector< T >::print ( std::ostream &  os = libMesh::out) const
inlinevirtualinherited

Prints the local contents of the vector, by default to libMesh::out

Definition at line 916 of file numeric_vector.h.

917 {
918  libmesh_assert (this->initialized());
919  os << "Size\tglobal = " << this->size()
920  << "\t\tlocal = " << this->local_size() << std::endl;
921 
922  os << "#\tValue" << std::endl;
923  for (numeric_index_type i=this->first_local_index(); i<this->last_local_index(); i++)
924  os << i << "\t" << (*this)(i) << std::endl;
925 }
virtual bool initialized() const
virtual numeric_index_type size() const =0
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type first_local_index() const =0
virtual numeric_index_type local_size() const =0
virtual numeric_index_type last_local_index() const =0

◆ print() [2/2]

template<>
void libMesh::NumericVector< Complex >::print ( std::ostream &  os) const
inlineinherited

Definition at line 898 of file numeric_vector.h.

899 {
900  libmesh_assert (this->initialized());
901  os << "Size\tglobal = " << this->size()
902  << "\t\tlocal = " << this->local_size() << std::endl;
903 
904  // std::complex<>::operator<<() is defined, but use this form
905  os << "#\tReal part\t\tImaginary part" << std::endl;
906  for (numeric_index_type i=this->first_local_index(); i<this->last_local_index(); i++)
907  os << i << "\t"
908  << (*this)(i).real() << "\t\t"
909  << (*this)(i).imag() << std::endl;
910 }
virtual bool initialized() const
virtual numeric_index_type size() const =0
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type first_local_index() const =0
virtual numeric_index_type local_size() const =0
virtual numeric_index_type last_local_index() const =0

◆ print_global() [1/2]

template<typename T >
void libMesh::NumericVector< T >::print_global ( std::ostream &  os = libMesh::out) const
inlinevirtualinherited

Prints the global contents of the vector, by default to libMesh::out

Definition at line 953 of file numeric_vector.h.

954 {
955  libmesh_assert (this->initialized());
956 
957  std::vector<T> v(this->size());
958  this->localize(v);
959 
960  // Right now we only want one copy of the output
961  if (this->processor_id())
962  return;
963 
964  os << "Size\tglobal = " << this->size() << std::endl;
965  os << "#\tValue" << std::endl;
966  for (numeric_index_type i=0; i!=v.size(); i++)
967  os << i << "\t" << v[i] << std::endl;
968 }
virtual bool initialized() const
virtual numeric_index_type size() const =0
dof_id_type numeric_index_type
Definition: id_types.h:92
processor_id_type processor_id() const
virtual void localize(std::vector< T > &v_local) const =0

◆ print_global() [2/2]

template<>
void libMesh::NumericVector< Complex >::print_global ( std::ostream &  os) const
inlineinherited

Definition at line 931 of file numeric_vector.h.

932 {
933  libmesh_assert (this->initialized());
934 
935  std::vector<Complex> v(this->size());
936  this->localize(v);
937 
938  // Right now we only want one copy of the output
939  if (this->processor_id())
940  return;
941 
942  os << "Size\tglobal = " << this->size() << std::endl;
943  os << "#\tReal part\t\tImaginary part" << std::endl;
944  for (numeric_index_type i=0; i!=v.size(); i++)
945  os << i << "\t"
946  << v[i].real() << "\t\t"
947  << v[i].imag() << std::endl;
948 }
virtual bool initialized() const
virtual numeric_index_type size() const =0
dof_id_type numeric_index_type
Definition: id_types.h:92
processor_id_type processor_id() const
virtual void localize(std::vector< T > &v_local) const =0

◆ print_info()

void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

Prints the reference information, by default to libMesh::out.

Definition at line 87 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().

88 {
90  out_stream << ReferenceCounter::get_info();
91 }
static std::string get_info()

◆ print_matlab()

template<typename T>
virtual void libMesh::NumericVector< T >::print_matlab ( const std::string &  = "") const
inlinevirtualinherited

Print the contents of the vector in Matlab's sparse matrix format. Optionally prints the vector to the file named name. If name is not specified it is dumped to the screen.

Reimplemented in libMesh::PetscVector< T >.

Definition at line 681 of file numeric_vector.h.

682  {
683  libmesh_not_implemented();
684  }

◆ processor_id()

processor_id_type libMesh::ParallelObject::processor_id ( ) const
inlineinherited
Returns
The rank of this processor in the group.

Definition at line 101 of file parallel_object.h.

References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::rank().

Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::EquationSystems::_read_impl(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::BoundaryInfo::add_elements(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::FEMSystem::assembly(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::ExodusII_IO_Helper::close(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::compute_communication_map_parameters(), libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::ExodusII_IO_Helper::create(), libMesh::DistributedMesh::delete_elem(), libMesh::DistributedMesh::delete_node(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DistributedMesh::DistributedMesh(), libMesh::DofMap::end_dof(), libMesh::DofMap::end_old_dof(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshFunction::find_element(), libMesh::MeshFunction::find_elements(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::DofMap::first_dof(), libMesh::DofMap::first_old_dof(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::MeshBase::get_info(), libMesh::DofMap::get_info(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SparsityPattern::Build::handle_vi_vj(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize_element_variables(), libMesh::ExodusII_IO_Helper::initialize_global_variables(), libMesh::ExodusII_IO_Helper::initialize_nodal_variables(), libMesh::DistributedMesh::insert_elem(), libMesh::DofMap::is_evaluable(), libMesh::SparsityPattern::Build::join(), libMesh::DofMap::last_dof(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::DofMap::n_local_dofs(), libMesh::System::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::SparsityPattern::Build::operator()(), libMesh::DistributedMesh::own_node(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::Nemesis_IO_Helper::put_cmap_params(), libMesh::Nemesis_IO_Helper::put_elem_cmap(), libMesh::Nemesis_IO_Helper::put_elem_map(), libMesh::Nemesis_IO_Helper::put_loadbal_param(), libMesh::Nemesis_IO_Helper::put_node_cmap(), libMesh::Nemesis_IO_Helper::put_node_map(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read(), libMesh::ExodusII_IO_Helper::read_elem_num_map(), libMesh::ExodusII_IO_Helper::read_global_values(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::System::read_serialized_data(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::System::read_serialized_vectors(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::CheckpointIO::select_split_config(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::Parallel::Packing< T >::unpack(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::CheckpointIO::write(), libMesh::EquationSystems::write(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), libMesh::ExodusII_IO_Helper::write_element_values(), libMesh::ExodusII_IO_Helper::write_elements(), libMesh::ExodusII_IO::write_global_data(), libMesh::ExodusII_IO_Helper::write_global_values(), libMesh::System::write_header(), libMesh::ExodusII_IO::write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::ExodusII_IO_Helper::write_nodal_values(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::ExodusII_IO_Helper::write_nodesets(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::System::write_serialized_data(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::System::write_serialized_vector(), libMesh::System::write_serialized_vectors(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), libMesh::ExodusII_IO_Helper::write_timestep(), and libMesh::ExodusII_IO::write_timestep_discontinuous().

102  { return cast_int<processor_id_type>(_communicator.rank()); }
const Parallel::Communicator & _communicator
processor_id_type rank() const
Definition: communicator.h:173

◆ reciprocal()

template<typename T >
void libMesh::DistributedVector< T >::reciprocal ( )
overridevirtual

Computes the pointwise reciprocal, $ u_i \leftarrow \frac{1}{u_i} \, \forall i$

Implements libMesh::NumericVector< T >.

Definition at line 183 of file distributed_vector.C.

184 {
185  for (numeric_index_type i=0; i<local_size(); i++)
186  {
187  // Don't divide by zero
188  libmesh_assert_not_equal_to (_values[i], T(0));
189 
190  _values[i] = 1. / _values[i];
191  }
192 }
virtual numeric_index_type local_size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92

◆ scale()

template<typename T >
void libMesh::DistributedVector< T >::scale ( const T  factor)
overridevirtual

Scale each element of the vector by the given factor.

Implements libMesh::NumericVector< T >.

Definition at line 249 of file distributed_vector.C.

References libMesh::initialized().

250 {
251  libmesh_assert (this->initialized());
252  libmesh_assert_equal_to (_values.size(), _local_size);
253  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
254 
255  for (std::size_t i=0; i<local_size(); i++)
256  _values[i] *= factor;
257 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual numeric_index_type local_size() const override
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ set()

template<typename T >
void libMesh::DistributedVector< T >::set ( const numeric_index_type  i,
const T  value 
)
inlineoverridevirtual

Sets v(i) = value.

Implements libMesh::NumericVector< T >.

Definition at line 564 of file distributed_vector.h.

References libMesh::initialized(), and value.

565 {
566  libmesh_assert (this->initialized());
567  libmesh_assert_equal_to (_values.size(), _local_size);
568  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
569  libmesh_assert_less (i, size());
570  libmesh_assert_less (i-first_local_index(), local_size());
571 
573 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual numeric_index_type local_size() const override
virtual numeric_index_type size() const override
numeric_index_type _local_size
virtual numeric_index_type first_local_index() const override
static const bool value
Definition: xdr_io.C:109
numeric_index_type _first_local_index

◆ size()

template<typename T >
numeric_index_type libMesh::DistributedVector< T >::size ( ) const
inlineoverridevirtual
Returns
The size of the vector.

Implements libMesh::NumericVector< T >.

Definition at line 497 of file distributed_vector.h.

References libMesh::initialized().

Referenced by libMesh::DistributedVector< T >::operator/=().

498 {
499  libmesh_assert (this->initialized());
500  libmesh_assert_equal_to (_values.size(), _local_size);
501  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
502 
503  return _global_size;
504 }
numeric_index_type _last_local_index
virtual bool initialized() const
numeric_index_type _global_size
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ subset_l1_norm()

template<class T >
Real libMesh::NumericVector< T >::subset_l1_norm ( const std::set< numeric_index_type > &  indices) const
virtualinherited
Returns
The $ \ell_1 $-norm of the vector, i.e. the sum of the absolute values for the specified entries in the vector.
Note
The indices must necessarily live on this processor.

Definition at line 311 of file numeric_vector.C.

Referenced by libMesh::System::discrete_var_norm().

312 {
313  const NumericVector<T> & v = *this;
314 
315  Real norm = 0;
316 
317  for (const auto & index : indices)
318  norm += std::abs(v(index));
319 
320  this->comm().sum(norm);
321 
322  return norm;
323 }
double abs(double a)
const Parallel::Communicator & comm() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ subset_l2_norm()

template<class T >
Real libMesh::NumericVector< T >::subset_l2_norm ( const std::set< numeric_index_type > &  indices) const
virtualinherited
Returns
The $ \ell_2 $-norm of the vector, i.e. the square root of the sum of the squares of the elements for the specified entries in the vector.
Note
The indices must necessarily live on this processor.

Definition at line 326 of file numeric_vector.C.

Referenced by libMesh::System::discrete_var_norm().

327 {
328  const NumericVector<T> & v = *this;
329 
330  Real norm = 0;
331 
332  for (const auto & index : indices)
333  norm += TensorTools::norm_sq(v(index));
334 
335  this->comm().sum(norm);
336 
337  return std::sqrt(norm);
338 }
const Parallel::Communicator & comm() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ subset_linfty_norm()

template<class T >
Real libMesh::NumericVector< T >::subset_linfty_norm ( const std::set< numeric_index_type > &  indices) const
virtualinherited
Returns
The maximum absolute value of the specified entries of this vector, which is the $ \ell_{\infty} $-norm of a vector.
Note
The indices must necessarily live on this processor.

Definition at line 341 of file numeric_vector.C.

Referenced by libMesh::System::discrete_var_norm().

342 {
343  const NumericVector<T> & v = *this;
344 
345  Real norm = 0;
346 
347  for (const auto & index : indices)
348  {
349  Real value = std::abs(v(index));
350  if (value > norm)
351  norm = value;
352  }
353 
354  this->comm().max(norm);
355 
356  return norm;
357 }
double abs(double a)
const Parallel::Communicator & comm() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const bool value
Definition: xdr_io.C:109

◆ sum()

template<typename T >
T libMesh::DistributedVector< T >::sum ( ) const
overridevirtual
Returns
The sum of all values in the vector.

Implements libMesh::NumericVector< T >.

Definition at line 43 of file distributed_vector.C.

References libMesh::initialized().

44 {
45  // This function must be run on all processors at once
46  parallel_object_only();
47 
48  libmesh_assert (this->initialized());
49  libmesh_assert_equal_to (_values.size(), _local_size);
50  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
51 
52  T local_sum = 0.;
53 
54  for (numeric_index_type i=0; i<local_size(); i++)
55  local_sum += _values[i];
56 
57  this->comm().sum(local_sum);
58 
59  return local_sum;
60 }
numeric_index_type _last_local_index
virtual bool initialized() const
const Parallel::Communicator & comm() const
virtual numeric_index_type local_size() const override
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ swap()

template<typename T >
void libMesh::DistributedVector< T >::swap ( NumericVector< T > &  v)
inlineoverridevirtual

Swaps the contents of this with v. There should be enough indirection in subclasses to make this an O(1) header-swap operation.

Reimplemented from libMesh::NumericVector< T >.

Definition at line 639 of file distributed_vector.h.

References libMesh::DistributedVector< T >::_first_local_index, libMesh::DistributedVector< T >::_global_size, libMesh::DistributedVector< T >::_last_local_index, libMesh::DistributedVector< T >::_local_size, libMesh::DistributedVector< T >::_values, and swap().

640 {
641  DistributedVector<T> & v = cast_ref<DistributedVector<T> &>(other);
642 
643  std::swap(_global_size, v._global_size);
644  std::swap(_local_size, v._local_size);
645  std::swap(_first_local_index, v._first_local_index);
646  std::swap(_last_local_index, v._last_local_index);
647 
648  // This should be O(1) with any reasonable STL implementation
649  std::swap(_values, v._values);
650 }
numeric_index_type _last_local_index
numeric_index_type _global_size
numeric_index_type _local_size
void swap(Iterator &lhs, Iterator &rhs)
numeric_index_type _first_local_index

◆ type() [1/2]

◆ type() [2/2]

template<typename T>
ParallelType& libMesh::NumericVector< T >::type ( )
inlineinherited
Returns
The type (SERIAL, PARALLEL, GHOSTED) of the vector.

Definition at line 159 of file numeric_vector.h.

159 { return _type; }

◆ zero()

template<typename T >
void libMesh::DistributedVector< T >::zero ( )
inlineoverridevirtual

Set all entries to zero. Equivalent to v = 0, but more obvious and faster.

Implements libMesh::NumericVector< T >.

Definition at line 459 of file distributed_vector.h.

References libMesh::initialized().

460 {
461  libmesh_assert (this->initialized());
462  libmesh_assert_equal_to (_values.size(), _local_size);
463  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
464 
465  std::fill (_values.begin(),
466  _values.end(),
467  0.);
468 }
numeric_index_type _last_local_index
virtual bool initialized() const
numeric_index_type _local_size
numeric_index_type _first_local_index

◆ zero_clone()

template<typename T >
std::unique_ptr< NumericVector< T > > libMesh::DistributedVector< T >::zero_clone ( ) const
inlineoverridevirtual
Returns
A smart pointer to a copy of this vector with the same type, size, and partitioning, but with all zero entries.
Note
This must be overridden in the derived classes.

Implements libMesh::NumericVector< T >.

Definition at line 474 of file distributed_vector.h.

References libMesh::NumericVector< T >::init().

475 {
476  NumericVector<T> * cloned_vector = new DistributedVector<T>(this->comm());
477  cloned_vector->init(*this);
478  return std::unique_ptr<NumericVector<T>>(cloned_vector);
479 }
const Parallel::Communicator & comm() const

Member Data Documentation

◆ _communicator

◆ _counts

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited

◆ _enable_print_counter

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

Flag to control whether reference count information is printed when print_info is called.

Definition at line 141 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().

◆ _first_local_index

template<typename T>
numeric_index_type libMesh::DistributedVector< T >::_first_local_index
private

◆ _global_size

template<typename T>
numeric_index_type libMesh::DistributedVector< T >::_global_size
private

◆ _is_closed

template<typename T>
bool libMesh::NumericVector< T >::_is_closed
protectedinherited

Flag which tracks whether the vector's values are consistent on all processors after insertion or addition of values has occurred on some or all processors.

Definition at line 712 of file numeric_vector.h.

Referenced by libMesh::NumericVector< Number >::closed(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::DistributedVector< T >::localize(), libMesh::DistributedVector< T >::operator=(), and libMesh::PetscVector< T >::PetscVector().

◆ _is_initialized

◆ _last_local_index

template<typename T>
numeric_index_type libMesh::DistributedVector< T >::_last_local_index
private

◆ _local_size

◆ _mutex

Threads::spin_mutex libMesh::ReferenceCounter::_mutex
staticprotectedinherited

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 135 of file reference_counter.h.

◆ _n_objects

Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects
staticprotectedinherited

The number of objects. Print the reference count information when the number returns to 0.

Definition at line 130 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().

◆ _type

◆ _values

template<typename T>
std::vector<T> libMesh::DistributedVector< T >::_values
private

The documentation for this class was generated from the following files: