libMesh::DistributedVector< T > Class Template Reference

#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)
 
 ~DistributedVector ()
 
virtual void close () libmesh_override
 
virtual void clear () libmesh_override
 
virtual void zero () libmesh_override
 
virtual UniquePtr< NumericVector< T > > zero_clone () const libmesh_override
 
virtual UniquePtr< NumericVector< T > > clone () const libmesh_override
 
virtual void init (const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
 
virtual void init (const numeric_index_type N, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
 
virtual void init (const numeric_index_type, const numeric_index_type, const std::vector< numeric_index_type > &, const bool=false, const ParallelType=AUTOMATIC) libmesh_override
 
virtual void init (const NumericVector< T > &other, const bool fast=false) libmesh_override
 
virtual NumericVector< T > & operator= (const T s) libmesh_override
 
virtual NumericVector< T > & operator= (const NumericVector< T > &v) libmesh_override
 
DistributedVector< T > & operator= (const DistributedVector< T > &v)
 
virtual NumericVector< T > & operator= (const std::vector< T > &v) libmesh_override
 
virtual Real min () const libmesh_override
 
virtual Real max () const libmesh_override
 
virtual T sum () const libmesh_override
 
virtual Real l1_norm () const libmesh_override
 
virtual Real l2_norm () const libmesh_override
 
virtual Real linfty_norm () const libmesh_override
 
virtual numeric_index_type size () const libmesh_override
 
virtual numeric_index_type local_size () const libmesh_override
 
virtual numeric_index_type first_local_index () const libmesh_override
 
virtual numeric_index_type last_local_index () const libmesh_override
 
virtual T operator() (const numeric_index_type i) const libmesh_override
 
virtual NumericVector< T > & operator+= (const NumericVector< T > &v) libmesh_override
 
virtual NumericVector< T > & operator-= (const NumericVector< T > &v) libmesh_override
 
virtual NumericVector< T > & operator/= (NumericVector< T > &v) libmesh_override
 
virtual void reciprocal () libmesh_override
 
virtual void conjugate () libmesh_override
 
virtual void set (const numeric_index_type i, const T value) libmesh_override
 
virtual void add (const numeric_index_type i, const T value) libmesh_override
 
virtual void add (const T s) libmesh_override
 
virtual void add (const NumericVector< T > &V) libmesh_override
 
virtual void add (const T a, const NumericVector< T > &v) libmesh_override
 
virtual void add_vector (const NumericVector< T > &, const SparseMatrix< T > &) libmesh_override
 
virtual void add_vector_transpose (const NumericVector< T > &, const SparseMatrix< T > &) libmesh_override
 
virtual void scale (const T factor) libmesh_override
 
virtual void abs () libmesh_override
 
virtual T dot (const NumericVector< T > &V) const libmesh_override
 
virtual void localize (std::vector< T > &v_local) const libmesh_override
 
virtual void localize (NumericVector< T > &v_local) const libmesh_override
 
virtual void localize (NumericVector< T > &v_local, const std::vector< numeric_index_type > &send_list) const libmesh_override
 
virtual void localize (std::vector< T > &v_local, const std::vector< numeric_index_type > &indices) const libmesh_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) libmesh_override
 
virtual void localize_to_one (std::vector< T > &v_local, const processor_id_type proc_id=0) const libmesh_override
 
virtual void pointwise_mult (const NumericVector< T > &vec1, const NumericVector< T > &vec2) libmesh_override
 
virtual void swap (NumericVector< T > &v) libmesh_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 UniquePtr< NumericVector< T > > build (const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
 
static UniquePtr< NumericVector< T > > build (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 >

Distributed vector. Provides an interface for simple parallel, distributed vectors. Offers some collective communication capabilities. Note that the class will sill function without MPI, but only on one processor. This lets us keep the parallel details behind the scenes.

Author
Benjamin S. Kirk
Date
2003

Definition at line 50 of file distributed_vector.h.

Member Typedef Documentation

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 110 of file reference_counter.h.

Constructor & Destructor Documentation

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

Dummy-Constructor. Dimension=0

Definition at line 447 of file distributed_vector.h.

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

448  :
449  NumericVector<T>(comm_in, ptype),
450  _global_size (0),
451  _local_size (0),
454 {
455  this->_type = ptype;
456 }
numeric_index_type _last_local_index
numeric_index_type _global_size
numeric_index_type _local_size
numeric_index_type _first_local_index
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 462 of file distributed_vector.h.

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

465  : NumericVector<T>(comm_in, ptype)
466 {
467  this->init(n, n, false, ptype);
468 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
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 474 of file distributed_vector.h.

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

478  : NumericVector<T>(comm_in, ptype)
479 {
480  this->init(n, n_local, false, ptype);
481 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
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 487 of file distributed_vector.h.

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

492  : NumericVector<T>(comm_in, ptype)
493 {
494  this->init(n, n_local, ghost, false, ptype);
495 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
template<typename T >
libMesh::DistributedVector< T >::~DistributedVector ( )
inline

Destructor, deallocates memory. Made virtual to allow for derived classes to behave properly.

Definition at line 501 of file distributed_vector.h.

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

502 {
503  this->clear ();
504 }
virtual void clear() libmesh_override

Member Function Documentation

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

v = abs(v)... that is, each entry in v is replaced by its absolute value.

Implements libMesh::NumericVector< T >.

Definition at line 259 of file distributed_vector.C.

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

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

260 {
261  libmesh_assert (this->initialized());
262  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
263 
264  for (numeric_index_type i=0; i<local_size(); i++)
265  this->set(i,std::abs(_values[i]));
266 }
double abs(double a)
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type local_size() const libmesh_override
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T >
void libMesh::DistributedVector< T >::add ( const numeric_index_type  i,
const T  value 
)
inlinevirtual

v(i) += value

Implements libMesh::NumericVector< T >.

Definition at line 775 of file distributed_vector.h.

References libMesh::DistributedVector< T >::_first_local_index, libMesh::DistributedVector< T >::_last_local_index, libMesh::DistributedVector< T >::_local_size, libMesh::DistributedVector< T >::_values, libMesh::DistributedVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::libmesh_assert(), libMesh::DistributedVector< T >::local_size(), and libMesh::DistributedVector< T >::size().

776 {
777  libmesh_assert (this->initialized());
778  libmesh_assert_equal_to (_values.size(), _local_size);
779  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
780  libmesh_assert_less (i, size());
781  libmesh_assert_less (i-first_local_index(), local_size());
782 
783  _values[i - _first_local_index] += value;
784 }
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type local_size() const libmesh_override
virtual numeric_index_type first_local_index() const libmesh_override
numeric_index_type _local_size
virtual numeric_index_type size() const libmesh_override
numeric_index_type _first_local_index
template<typename T >
void libMesh::DistributedVector< T >::add ( const T  s)
virtual

$U(0-LIBMESH_DIM)+=s$. Addition of s to all components. Note that s is a scalar and not a vector.

Implements libMesh::NumericVector< T >.

Definition at line 211 of file distributed_vector.C.

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

212 {
213  libmesh_assert (this->initialized());
214  libmesh_assert_equal_to (_values.size(), _local_size);
215  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
216 
217  for (numeric_index_type i=0; i<local_size(); i++)
218  _values[i] += v;
219 }
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type local_size() const libmesh_override
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T >
void libMesh::DistributedVector< T >::add ( const NumericVector< T > &  V)
virtual

$U+=V$. Simple vector addition, equal to the operator +=.

Implements libMesh::NumericVector< T >.

Definition at line 224 of file distributed_vector.C.

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

225 {
226  libmesh_assert (this->initialized());
227  libmesh_assert_equal_to (_values.size(), _local_size);
228  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
229 
230  add (1., v);
231 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual void add(const numeric_index_type i, const T value) libmesh_override
libmesh_assert(j)
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T >
void libMesh::DistributedVector< T >::add ( const T  a,
const NumericVector< T > &  v 
)
virtual

$U+=a*V$. Simple vector addition, equal to the operator +=.

Implements libMesh::NumericVector< T >.

Definition at line 236 of file distributed_vector.C.

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

237 {
238  libmesh_assert (this->initialized());
239  libmesh_assert_equal_to (_values.size(), _local_size);
240  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
241 
242  add(a, v);
243 }
numeric_index_type _last_local_index
virtual bool initialized() const
virtual void add(const numeric_index_type i, const T value) libmesh_override
libmesh_assert(j)
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T>
virtual void libMesh::DistributedVector< T >::add_vector ( const NumericVector< T > &  ,
const SparseMatrix< T > &   
)
inlinevirtual

$U+=A*V$. Add the product of a Sparse matrix A and a Numeric vector V to this Numeric vector. Not implemented.

Implements libMesh::NumericVector< T >.

Definition at line 326 of file distributed_vector.h.

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

$ U+=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 383 of file numeric_vector.C.

Referenced by libMesh::NumericVector< T >::add_vector(), libMesh::LinearImplicitSystem::assembly(), libMesh::NumericVector< Number >::operator/=(), libMesh::NewmarkSystem::update_rhs(), and libMesh::SparseMatrix< T >::vector_mult_add().

385 {
386  int n = dof_indices.size();
387  for (int i=0; i<n; i++)
388  this->add (dof_indices[i], v[i]);
389 }
virtual void add(const numeric_index_type i, const T value)=0
template<typename T>
void libMesh::NumericVector< T >::add_vector ( const std::vector< T > &  v,
const std::vector< numeric_index_type > &  dof_indices 
)
inlineinherited

$ U+=v $ where v is a std::vector and each dof_indices[i] specifies where to add value v[i]

Definition at line 840 of file numeric_vector.h.

References libMesh::NumericVector< T >::add_vector(), and libMesh::libmesh_assert().

842 {
843  libmesh_assert(v.size() == dof_indices.size());
844  if (!v.empty())
845  this->add_vector(&v[0], dof_indices);
846 }
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
libmesh_assert(j)
template<typename T>
void libMesh::NumericVector< T >::add_vector ( const NumericVector< T > &  V,
const std::vector< numeric_index_type > &  dof_indices 
)
virtualinherited

$ U+=v $ where v is a NumericVector and each dof_indices[i] specifies where to add value v(i)

Definition at line 394 of file numeric_vector.C.

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

396 {
397  int n = dof_indices.size();
398  libmesh_assert_equal_to(v.size(), static_cast<unsigned>(n));
399  for (int i=0; i<n; i++)
400  this->add (dof_indices[i], v(i));
401 }
virtual void add(const numeric_index_type i, const T value)=0
template<typename T>
void libMesh::NumericVector< T >::add_vector ( const DenseVector< T > &  V,
const std::vector< numeric_index_type > &  dof_indices 
)
inlineinherited

$ U+=v $ where v is a DenseVector and each dof_indices[i] specifies where to add value v(i)

Definition at line 852 of file numeric_vector.h.

References libMesh::NumericVector< T >::add_vector(), libMesh::DenseVector< T >::empty(), libMesh::libmesh_assert(), and libMesh::DenseVector< T >::size().

854 {
855  libmesh_assert(v.size() == dof_indices.size());
856  if (!v.empty())
857  this->add_vector(&v(0), dof_indices);
858 }
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
libmesh_assert(j)
template<typename T>
void libMesh::NumericVector< T >::add_vector ( const NumericVector< T > &  v,
const ShellMatrix< T > &  a 
)
inherited

$U+=A*V$, add the product of a ShellMatrix A and a NumericVector V to this, where this=U.

Definition at line 406 of file numeric_vector.C.

References libMesh::ShellMatrix< T >::vector_mult_add().

408 {
409  a.vector_mult_add(*this,v);
410 }
template<typename T>
virtual void libMesh::DistributedVector< T >::add_vector_transpose ( const NumericVector< T > &  ,
const SparseMatrix< T > &   
)
inlinevirtual

$U+=A^T*V$. Add the product of the transpose of a Sparse matrix A_trans and a Numeric vector V to this Numeric vector. Not implemented.

Implements libMesh::NumericVector< T >.

Definition at line 336 of file distributed_vector.h.

References libMesh::DistributedVector< T >::abs(), libMesh::DistributedVector< T >::dot(), libMesh::DistributedVector< T >::localize(), libMesh::DistributedVector< T >::localize_to_one(), libMesh::DistributedVector< T >::pointwise_mult(), libMesh::DistributedVector< T >::scale(), and libMesh::DistributedVector< T >::swap().

338  { libmesh_not_implemented(); }
template<typename T >
UniquePtr< 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 46 of file numeric_vector.C.

References libMesh::AUTOMATIC, libMesh::EIGEN_SOLVERS, libMesh::LASPACK_SOLVERS, libMesh::PETSC_SOLVERS, and libMesh::TRILINOS_SOLVERS.

Referenced by libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::System::add_vector(), libMesh::NumericVector< T >::build(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::System::calculate_norm(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::EquationSystems::get_solution(), and libMesh::System::project_vector().

47 {
48  // Build the appropriate vector
49  switch (solver_package)
50  {
51 
52 #ifdef LIBMESH_HAVE_LASPACK
53  case LASPACK_SOLVERS:
54  return UniquePtr<NumericVector<T> >(new LaspackVector<T>(comm, AUTOMATIC));
55 #endif
56 
57 #ifdef LIBMESH_HAVE_PETSC
58  case PETSC_SOLVERS:
59  return UniquePtr<NumericVector<T> >(new PetscVector<T>(comm, AUTOMATIC));
60 #endif
61 
62 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA
63  case TRILINOS_SOLVERS:
64  return UniquePtr<NumericVector<T> >(new EpetraVector<T>(comm, AUTOMATIC));
65 #endif
66 
67 #ifdef LIBMESH_HAVE_EIGEN
68  case EIGEN_SOLVERS:
69  return UniquePtr<NumericVector<T> >(new EigenSparseVector<T>(comm, AUTOMATIC));
70 #endif
71 
72  default:
73  return UniquePtr<NumericVector<T> >(new DistributedVector<T>(comm, AUTOMATIC));
74  }
75 
76  libmesh_error_msg("We'll never get here!");
77  return UniquePtr<NumericVector<T> >();
78 }
EIGEN_SOLVERS
Definition: libmesh.C:260
TRILINOS_SOLVERS
Definition: libmesh.C:258
LASPACK_SOLVERS
Definition: libmesh.C:262
const Parallel::Communicator & comm() const
template<typename T >
UniquePtr< NumericVector< T > > libMesh::NumericVector< T >::build ( const SolverPackage  solver_package = libMesh::default_solver_package())
staticinherited

Builds a NumericVector on the processors in communicator CommWorld using the linear solver package specified by solver_package. Deprecated.

Definition at line 84 of file numeric_vector.C.

References libMesh::NumericVector< T >::build(), and libMesh::CommWorld.

85 {
86  libmesh_deprecated();
87  return NumericVector<T>::build(CommWorld, solver_package);
88 }
static UniquePtr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
Parallel::FakeCommunicator CommWorld
Definition: libmesh.C:206
template<typename T >
void libMesh::DistributedVector< T >::clear ( )
inlinevirtual
Returns
the DistributedVector to a pristine state.

Reimplemented from libMesh::NumericVector< T >.

Definition at line 638 of file distributed_vector.h.

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, and libMesh::DistributedVector< T >::_values.

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

639 {
640  _values.clear();
641 
642  _global_size =
643  _local_size =
645  _last_local_index = 0;
646 
647 
648  this->_is_closed = this->_is_initialized = false;
649 }
numeric_index_type _last_local_index
numeric_index_type _global_size
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T >
UniquePtr< NumericVector< T > > libMesh::DistributedVector< T >::clone ( ) const
inlinevirtual

Creates a copy of this vector and returns it in an UniquePtr.

Implements libMesh::NumericVector< T >.

Definition at line 681 of file distributed_vector.h.

References libMesh::ParallelObject::comm(), and libMesh::NumericVector< T >::init().

682 {
683  NumericVector<T> * cloned_vector = new DistributedVector<T>(this->comm());
684  cloned_vector->init(*this, true);
685  *cloned_vector = *this;
686  return UniquePtr<NumericVector<T> >(cloned_vector);
687 }
const Parallel::Communicator & comm() const
template<typename T >
void libMesh::DistributedVector< T >::close ( )
inlinevirtual

Call the assemble functions

Implements libMesh::NumericVector< T >.

Definition at line 627 of file distributed_vector.h.

References libMesh::NumericVector< T >::_is_closed, libMesh::NumericVector< T >::initialized(), and libMesh::libmesh_assert().

628 {
629  libmesh_assert (this->initialized());
630 
631  this->_is_closed = true;
632 }
virtual bool initialized() const
libmesh_assert(j)
const Parallel::Communicator& libMesh::ParallelObject::comm ( ) const
inlineinherited
Returns
a reference to the Parallel::Communicator object used by this mesh.

Definition at line 87 of file parallel_object.h.

References libMesh::ParallelObject::_communicator.

Referenced by libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_petsc_snes_jacobian(), libMesh::__libmesh_petsc_snes_postcheck(), libMesh::__libmesh_petsc_snes_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::ParmetisPartitioner::_do_repartition(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult(), libMesh::SlepcEigenSolver< 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::ImplicitSystem::add_matrix(), libMesh::System::add_vector(), libMesh::EigenSparseLinearSolver< T >::adjoint_solve(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assemble_qoi(), libMesh::MeshCommunication::assign_global_indices(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::DofMap::attach_matrix(), libMesh::Parallel::BinSorter< KeyType, IdxType >::binsort(), libMesh::Parallel::Sort< KeyType, IdxType >::binsort(), libMesh::MeshTools::bounding_box(), libMesh::MeshCommunication::broadcast(), libMesh::SparseMatrix< T >::build(), libMesh::MeshTools::Generation::build_extrusion(), libMesh::Parallel::Histogram< KeyType, IdxType >::build_histogram(), libMesh::PetscNonlinearSolver< T >::build_mat_null_space(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::MeshBase::cache_elem_dims(), libMesh::System::calculate_norm(), libMesh::DofMap::check_dirichlet_bcid_consistency(), libMesh::DistributedVector< T >::clone(), libMesh::EigenSparseVector< T >::clone(), libMesh::LaspackVector< T >::clone(), libMesh::EpetraVector< T >::clone(), libMesh::PetscVector< T >::clone(), libMesh::EpetraVector< T >::close(), libMesh::Parallel::Sort< KeyType, IdxType >::communicate_bins(), 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::MeshTools::correct_node_proc_ids(), libMesh::MeshRefinement::create_parent_error_vector(), 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::MeshCommunication::gather(), libMesh::MeshCommunication::gather_neighboring_elements(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::DofMap::get_info(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::EquationSystems::get_solution(), libMesh::LocationMap< T >::init(), libMesh::PetscDiffSolver::init(), libMesh::TimeSolver::init(), libMesh::TopologyMap::init(), libMesh::TaoOptimizationSolver< T >::init(), libMesh::PetscNonlinearSolver< T >::init(), libMesh::DistributedVector< T >::init(), libMesh::EpetraVector< T >::init(), libMesh::PetscVector< T >::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::EigenSystem::init_data(), libMesh::EigenSystem::init_matrices(), libMesh::ParmetisPartitioner::initialize(), libMesh::OptimizationSystem::initialize_equality_constraints_storage(), libMesh::OptimizationSystem::initialize_inequality_constraints_storage(), 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::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::DistributedVector< T >::max(), libMesh::FEMSystem::mesh_position_set(), libMesh::MeshSerializer::MeshSerializer(), libMesh::DistributedVector< T >::min(), 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::Partitioner::partition(), libMesh::MetisPartitioner::partition_range(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::SparseMatrix< T >::print(), libMesh::MeshTools::processor_bounding_box(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshCommunication::redistribute(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::MeshCommunication::send_coarse_ghosts(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::Partitioner::set_parent_processor_ids(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::Parallel::Sort< KeyType, IdxType >::sort(), libMesh::MeshTools::subdomain_bounding_box(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), libMesh::Parallel::sync_element_data_by_parent_id(), libMesh::Parallel::sync_node_data_by_element_id(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::DistributedVector< T >::zero_clone(), libMesh::LaspackVector< T >::zero_clone(), libMesh::EigenSparseVector< T >::zero_clone(), libMesh::EpetraVector< T >::zero_clone(), and libMesh::PetscVector< T >::zero_clone().

88  { return _communicator; }
const Parallel::Communicator & _communicator
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. When differences occur, the return value contains the first index i where the difference (a[i]-b[i]) exceeded the threshold. When no threshold is given, the libMesh TOLERANCE is used.

Definition at line 116 of file numeric_vector.C.

References std::abs(), libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::initialized(), libMesh::NumericVector< T >::last_local_index(), libMesh::libmesh_assert(), and std::max().

Referenced by libMesh::NumericVector< Number >::operator/=().

118 {
119  libmesh_assert (this->initialized());
120  libmesh_assert (other_vector.initialized());
121  libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
122  libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
123 
124  int first_different_i = std::numeric_limits<int>::max();
126 
127  do
128  {
129  if ( std::abs( (*this)(i) - other_vector(i) ) > threshold )
130  first_different_i = i;
131  else
132  i++;
133  }
134  while (first_different_i==std::numeric_limits<int>::max()
135  && i<last_local_index());
136 
137  // Find the correct first differing index in parallel
138  this->comm().min(first_different_i);
139 
140  if (first_different_i == std::numeric_limits<int>::max())
141  return -1;
142 
143  return first_different_i;
144 }
double abs(double a)
virtual numeric_index_type last_local_index() const =0
virtual bool initialized() const
long double max(long double a, double b)
libmesh_assert(j)
dof_id_type numeric_index_type
Definition: id_types.h:92
const Parallel::Communicator & comm() const
virtual numeric_index_type first_local_index() const =0
template<typename T >
void libMesh::DistributedVector< T >::conjugate ( )
virtual

Replace each entry v_i = real(v_i) + imag(v_i) of this vector by its complex conjugate, real(v_i) - imag(v_i)

Implements libMesh::NumericVector< T >.

Definition at line 197 of file distributed_vector.C.

References libMesh::libmesh_conj().

198 {
199  for (numeric_index_type i=0; i<local_size(); i++)
200  {
201  // Replace values by complex conjugate
202  _values[i] = libmesh_conj(_values[i]);
203  }
204 }
T libmesh_conj(T a)
virtual numeric_index_type local_size() const libmesh_override
dof_id_type numeric_index_type
Definition: id_types.h:92
template<typename T>
virtual void libMesh::NumericVector< T >::create_subvector ( NumericVector< T > &  ,
const std::vector< numeric_index_type > &   
) const
inlinevirtualinherited

Creates the subvector "subvector" from the indices in the "rows" array. 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 688 of file numeric_vector.h.

690  {
691  libmesh_not_implemented();
692  }
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited
template<typename T >
T libMesh::DistributedVector< T >::dot ( const NumericVector< T > &  V) const
virtual

Computes the dot product, p = U.V

Implements libMesh::NumericVector< T >.

Definition at line 273 of file distributed_vector.C.

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

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

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

Definition at line 345 of file numeric_vector.h.

345 { return (*this)(i); }
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.

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

101 {
102  _enable_print_counter = true;
103  return;
104 }
template<typename T >
numeric_index_type libMesh::DistributedVector< T >::first_local_index ( ) const
inlinevirtual
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 811 of file numeric_vector.h.

Referenced by libMesh::DofMap::enforce_adjoint_constraints_exactly(), libMesh::DofMap::enforce_constraints_exactly(), libMesh::FEMContext::pre_fe_reinit(), and libMesh::System::project_vector().

813 {
814  const std::size_t num = index.size();
815  for(std::size_t i=0; i<num; i++)
816  {
817  values[i] = (*this)(index[i]);
818  }
819 }
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 825 of file numeric_vector.h.

827 {
828  const std::size_t num = index.size();
829  values.resize(num);
830  if (!num)
831  return;
832 
833  this->get(index, &values[0]);
834 }
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 (Counts::iterator it = _counts.begin();
59  it != _counts.end(); ++it)
60  {
61  const std::string name(it->first);
62  const unsigned int creations = it->second.first;
63  const unsigned int destructions = it->second.second;
64 
65  oss << "| " << name << " reference count information:\n"
66  << "| Creations: " << creations << '\n'
67  << "| Destructions: " << destructions << '\n';
68  }
69 
70  oss << " ---------------------------------------------------------------------------- \n";
71 
72  return oss.str();
73 
74 #else
75 
76  return "";
77 
78 #endif
79 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
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 local relative threshold. When differences occur, the return value contains the first index where the difference (a[i]-b[i])/max_j(a[j],b[j]) exceeded the threshold. When no threshold is given, the libMesh TOLERANCE is used.

Definition at line 181 of file numeric_vector.C.

References std::abs(), libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::initialized(), libMesh::NumericVector< T >::last_local_index(), libMesh::libmesh_assert(), libMesh::NumericVector< T >::linfty_norm(), std::max(), and libMesh::Real.

Referenced by libMesh::NumericVector< Number >::operator/=().

183 {
184  libmesh_assert (this->initialized());
185  libmesh_assert (other_vector.initialized());
186  libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
187  libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
188 
189  int first_different_i = std::numeric_limits<int>::max();
191 
192  const Real my_norm = this->linfty_norm();
193  const Real other_norm = other_vector.linfty_norm();
194  const Real abs_threshold = std::max(my_norm, other_norm) * threshold;
195 
196  do
197  {
198  if ( std::abs( (*this)(i) - other_vector(i) ) > abs_threshold )
199  first_different_i = i;
200  else
201  i++;
202  }
203  while (first_different_i==std::numeric_limits<int>::max()
204  && i<last_local_index());
205 
206  // Find the correct first differing index in parallel
207  this->comm().min(first_different_i);
208 
209  if (first_different_i == std::numeric_limits<int>::max())
210  return -1;
211 
212  return first_different_i;
213 }
double abs(double a)
virtual numeric_index_type last_local_index() const =0
virtual bool initialized() const
long double max(long double a, double b)
libmesh_assert(j)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual Real linfty_norm() const =0
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Parallel::Communicator & comm() const
virtual numeric_index_type first_local_index() const =0
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 160 of file reference_counter.h.

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

Referenced by libMesh::ReferenceCounter::n_objects(), and libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

161 {
162  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
163  std::pair<unsigned int, unsigned int> & p = _counts[name];
164 
165  p.first++;
166 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
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 173 of file reference_counter.h.

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

Referenced by libMesh::ReferenceCounter::n_objects(), and libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

174 {
175  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
176  std::pair<unsigned int, unsigned int> & p = _counts[name];
177 
178  p.second++;
179 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
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 
)
inlinevirtual

Change the dimension of the vector to N. The reserved memory for this vector remains unchanged if possible, to make things faster, but this may waste some memory, so take this in the back of your head. However, if N==0 all memory is freed, i.e. 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 cited 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 510 of file distributed_vector.h.

References libMesh::DistributedVector< T >::_first_local_index, libMesh::DistributedVector< T >::_global_size, libMesh::NumericVector< T >::_is_initialized, libMesh::DistributedVector< T >::_last_local_index, libMesh::DistributedVector< T >::_local_size, libMesh::NumericVector< T >::_type, libMesh::DistributedVector< T >::_values, libMesh::AUTOMATIC, libMesh::DistributedVector< T >::clear(), libMesh::ParallelObject::comm(), libMesh::NumericVector< T >::initialized(), libMesh::libmesh_assert(), libMesh::ParallelObject::n_processors(), libMesh::PARALLEL, libMesh::ParallelObject::processor_id(), libMesh::SERIAL, libMesh::Parallel::Communicator::sum(), and libMesh::DistributedVector< T >::zero().

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

514 {
515  // This function must be run on all processors at once
516  parallel_object_only();
517 
518  libmesh_assert_less_equal (n_local, n);
519 
520  if (ptype == AUTOMATIC)
521  {
522  if (n == n_local)
523  this->_type = SERIAL;
524  else
525  this->_type = PARALLEL;
526  }
527  else
528  this->_type = ptype;
529 
530  libmesh_assert ((this->_type==SERIAL && n==n_local) ||
531  this->_type==PARALLEL);
532 
533  // Clear the data structures if already initialized
534  if (this->initialized())
535  this->clear();
536 
537  // Initialize data structures
538  _values.resize(n_local);
539  _local_size = n_local;
540  _global_size = n;
541 
542  _first_local_index = 0;
543 
544 #ifdef LIBMESH_HAVE_MPI
545 
546  std::vector<numeric_index_type> local_sizes (this->n_processors(), 0);
547 
548  local_sizes[this->processor_id()] = n_local;
549 
550  this->comm().sum(local_sizes);
551 
552  // _first_local_index is the sum of _local_size
553  // for all processor ids less than ours
554  for (processor_id_type p=0; p!=this->processor_id(); p++)
555  _first_local_index += local_sizes[p];
556 
557 
558 # ifdef DEBUG
559  // Make sure all the local sizes sum up to the global
560  // size, otherwise there is big trouble!
561  numeric_index_type dbg_sum=0;
562 
563  for (processor_id_type p=0; p!=this->n_processors(); p++)
564  dbg_sum += local_sizes[p];
565 
566  libmesh_assert_equal_to (dbg_sum, n);
567 
568 # endif
569 
570 #else
571 
572  // No other options without MPI!
573  if (n != n_local)
574  libmesh_error_msg("ERROR: MPI is required for n != n_local!");
575 
576 #endif
577 
579 
580  // Set the initialized flag
581  this->_is_initialized = true;
582 
583  // Zero the components unless directed otherwise
584  if (!fast)
585  this->zero();
586 }
numeric_index_type _last_local_index
processor_id_type n_processors() const
uint8_t processor_id_type
Definition: id_types.h:99
virtual bool initialized() const
numeric_index_type _global_size
libmesh_assert(j)
virtual void zero() libmesh_override
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual void clear() libmesh_override
numeric_index_type _local_size
const Parallel::Communicator & comm() const
numeric_index_type _first_local_index
processor_id_type processor_id() const
template<typename T >
void libMesh::DistributedVector< T >::init ( const numeric_index_type  N,
const bool  fast = false,
const ParallelType  ptype = AUTOMATIC 
)
inlinevirtual

call init with n_local = N,

Implements libMesh::NumericVector< T >.

Definition at line 616 of file distributed_vector.h.

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

619 {
620  this->init(n,n,fast,ptype);
621 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
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 > &  ,
const bool  fast = false,
const ParallelType  ptype = AUTOMATIC 
)
inlinevirtual

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

Implements libMesh::NumericVector< T >.

Definition at line 591 of file distributed_vector.h.

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

596 {
597  // TODO: we shouldn't ignore the ghost sparsity pattern
598  this->init(n, n_local, fast, ptype);
599 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
template<class T >
void libMesh::DistributedVector< T >::init ( const NumericVector< T > &  other,
const bool  fast = false 
)
virtual

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

Implements libMesh::NumericVector< T >.

Definition at line 606 of file distributed_vector.h.

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

608 {
609  this->init(other.size(),other.local_size(),fast,other.type());
610 }
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) libmesh_override
template<typename T>
virtual bool libMesh::NumericVector< T >::initialized ( ) const
inlinevirtualinherited
Returns
true if the vector has been initialized, false otherwise.

Definition at line 132 of file numeric_vector.h.

Referenced by libMesh::DistributedVector< T >::add(), libMesh::EigenSparseVector< T >::add(), libMesh::LaspackVector< T >::add(), libMesh::ImplicitSystem::assemble(), libMesh::LaspackVector< T >::clear(), libMesh::EpetraVector< T >::clear(), libMesh::PetscVector< T >::clear(), libMesh::DistributedVector< T >::close(), libMesh::EigenSparseVector< T >::close(), libMesh::LaspackVector< T >::close(), libMesh::EpetraVector< T >::close(), libMesh::NumericVector< T >::compare(), libMesh::PetscVector< T >::create_subvector(), libMesh::DistributedVector< T >::first_local_index(), libMesh::EigenSparseVector< T >::first_local_index(), libMesh::LaspackVector< T >::first_local_index(), libMesh::EpetraVector< T >::first_local_index(), libMesh::PetscVector< T >::first_local_index(), libMesh::NumericVector< T >::global_relative_compare(), libMesh::DistributedVector< T >::init(), libMesh::EigenSparseVector< T >::init(), libMesh::PetscVector< T >::init(), libMesh::DistributedVector< T >::last_local_index(), libMesh::EigenSparseVector< T >::last_local_index(), libMesh::LaspackVector< T >::last_local_index(), libMesh::EpetraVector< T >::last_local_index(), libMesh::PetscVector< T >::last_local_index(), libMesh::NumericVector< T >::local_relative_compare(), libMesh::DistributedVector< T >::local_size(), libMesh::EigenSparseVector< T >::local_size(), libMesh::LaspackVector< T >::local_size(), libMesh::EpetraVector< T >::local_size(), libMesh::PetscVector< T >::local_size(), libMesh::PetscVector< T >::map_global_to_local_index(), libMesh::DistributedVector< T >::max(), libMesh::EpetraVector< T >::max(), libMesh::DistributedVector< T >::min(), libMesh::EpetraVector< T >::min(), libMesh::DistributedVector< T >::operator()(), libMesh::LaspackVector< T >::operator()(), libMesh::EigenSparseVector< T >::operator()(), libMesh::EpetraVector< T >::operator()(), libMesh::NumericVector< T >::print(), libMesh::NumericVector< T >::print_global(), libMesh::DistributedVector< T >::set(), libMesh::EigenSparseVector< T >::set(), libMesh::LaspackVector< T >::set(), libMesh::DistributedVector< T >::size(), libMesh::LaspackVector< T >::size(), libMesh::EigenSparseVector< T >::size(), libMesh::EpetraVector< T >::size(), libMesh::PetscVector< T >::size(), libMesh::DistributedVector< T >::zero(), libMesh::LaspackVector< T >::zero(), libMesh::EigenSparseVector< T >::zero(), libMesh::EpetraVector< T >::zero(), and libMesh::LaspackVector< T >::~LaspackVector().

132 { return _is_initialized; }
template<typename T>
void libMesh::NumericVector< T >::insert ( const T *  v,
const std::vector< numeric_index_type > &  dof_indices 
)
virtualinherited

$ U=v $ where v is a T[] or T * and you want to specify WHERE to insert it

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

Definition at line 94 of file numeric_vector.C.

Referenced by libMesh::NumericVector< T >::insert(), and libMesh::NumericVector< Number >::operator/=().

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
template<typename T>
void libMesh::NumericVector< T >::insert ( const std::vector< T > &  v,
const std::vector< numeric_index_type > &  dof_indices 
)
inlineinherited

$ U=v $ where v is a std::vector<T> and you want to specify WHERE to insert it

Definition at line 864 of file numeric_vector.h.

References libMesh::NumericVector< T >::insert(), and libMesh::libmesh_assert().

866 {
867  libmesh_assert(v.size() == dof_indices.size());
868  if (!v.empty())
869  this->insert(&v[0], dof_indices);
870 }
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices)
libmesh_assert(j)
template<typename T>
void libMesh::NumericVector< T >::insert ( const NumericVector< T > &  V,
const std::vector< numeric_index_type > &  dof_indices 
)
virtualinherited

$U=V$, where U and V are type NumericVector<T> and you want to specify WHERE to insert the NumericVector<T> V

Definition at line 104 of file numeric_vector.C.

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

106 {
107  libmesh_assert_equal_to (V.size(), dof_indices.size());
108 
109  for (numeric_index_type i=0; i<dof_indices.size(); i++)
110  this->set (dof_indices[i], V(i));
111 }
dof_id_type numeric_index_type
Definition: id_types.h:92
template<typename T>
void libMesh::NumericVector< T >::insert ( const DenseVector< T > &  V,
const std::vector< numeric_index_type > &  dof_indices 
)
inlineinherited

$ U=V $ where U and V are type DenseVector<T> and you want to specify WHERE to insert the DenseVector<T> V

Definition at line 876 of file numeric_vector.h.

References libMesh::DenseVector< T >::empty(), libMesh::NumericVector< T >::insert(), libMesh::libmesh_assert(), and libMesh::DenseVector< T >::size().

878 {
879  libmesh_assert(v.size() == dof_indices.size());
880  if (!v.empty())
881  this->insert(&v(0), dof_indices);
882 }
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices)
libmesh_assert(j)
template<typename T>
void libMesh::NumericVector< T >::insert ( const DenseSubVector< T > &  V,
const std::vector< numeric_index_type > &  dof_indices 
)
inlineinherited

$ U=V $ where V is a DenseSubVector<T> and you want to specify WHERE to insert it

Definition at line 888 of file numeric_vector.h.

References libMesh::DenseSubVector< T >::empty(), libMesh::NumericVector< T >::insert(), libMesh::libmesh_assert(), and libMesh::DenseSubVector< T >::size().

890 {
891  libmesh_assert(v.size() == dof_indices.size());
892  if (!v.empty())
893  this->insert(&v(0), dof_indices);
894 }
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices)
libmesh_assert(j)
template<typename T >
Real libMesh::DistributedVector< T >::l1_norm ( ) const
virtual
Returns
the $l_1$-norm of the vector, i.e. the sum of the absolute values.

Implements libMesh::NumericVector< T >.

Definition at line 64 of file distributed_vector.C.

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

65 {
66  // This function must be run on all processors at once
67  parallel_object_only();
68 
69  libmesh_assert (this->initialized());
70  libmesh_assert_equal_to (_values.size(), _local_size);
71  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
72 
73  double local_l1 = 0.;
74 
75  for (numeric_index_type i=0; i<local_size(); i++)
76  local_l1 += std::abs(_values[i]);
77 
78  this->comm().sum(local_l1);
79 
80  return local_l1;
81 }
double abs(double a)
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type local_size() const libmesh_override
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
const Parallel::Communicator & comm() const
numeric_index_type _first_local_index
template<typename T >
Real libMesh::DistributedVector< T >::l2_norm ( ) const
virtual
Returns
the $l_2$-norm of the vector, i.e. the square root of the sum of the squares of the elements.

Implements libMesh::NumericVector< T >.

Definition at line 86 of file distributed_vector.C.

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

87 {
88  // This function must be run on all processors at once
89  parallel_object_only();
90 
91  libmesh_assert (this->initialized());
92  libmesh_assert_equal_to (_values.size(), _local_size);
93  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
94 
95  double local_l2 = 0.;
96 
97  for (numeric_index_type i=0; i<local_size(); i++)
98  local_l2 += TensorTools::norm_sq(_values[i]);
99 
100  this->comm().sum(local_l2);
101 
102  return std::sqrt(local_l2);
103 }
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type local_size() const libmesh_override
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
const Parallel::Communicator & comm() const
numeric_index_type _first_local_index
template<typename T >
numeric_index_type libMesh::DistributedVector< T >::last_local_index ( ) const
inlinevirtual
template<typename T >
Real libMesh::DistributedVector< T >::linfty_norm ( ) const
virtual
Returns
the maximum absolute value of the elements of this vector, which is the $l_\infty$-norm of a vector.

Implements libMesh::NumericVector< T >.

Definition at line 108 of file distributed_vector.C.

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

109 {
110  // This function must be run on all processors at once
111  parallel_object_only();
112 
113  libmesh_assert (this->initialized());
114  libmesh_assert_equal_to (_values.size(), _local_size);
115  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
116 
117  Real local_linfty = 0.;
118 
119  for (numeric_index_type i=0; i<local_size(); i++)
120  local_linfty = std::max(local_linfty,
121  static_cast<Real>(std::abs(_values[i]))
122  ); // Note we static_cast so that both
123  // types are the same, as required
124  // by std::max
125 
126  this->comm().max(local_linfty);
127 
128  return local_linfty;
129 }
double abs(double a)
numeric_index_type _last_local_index
virtual bool initialized() const
long double max(long double a, double b)
libmesh_assert(j)
virtual numeric_index_type local_size() const libmesh_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
const Parallel::Communicator & comm() const
numeric_index_type _first_local_index
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. When differences occur, the return value contains the first index where the difference (a[i]-b[i])/max(a[i],b[i]) exceeded the threshold. When no threshold is given, the libMesh TOLERANCE is used.

Definition at line 148 of file numeric_vector.C.

References std::abs(), libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::initialized(), libMesh::NumericVector< T >::last_local_index(), libMesh::libmesh_assert(), and std::max().

Referenced by libMesh::NumericVector< Number >::operator/=().

150 {
151  libmesh_assert (this->initialized());
152  libmesh_assert (other_vector.initialized());
153  libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
154  libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
155 
156  int first_different_i = std::numeric_limits<int>::max();
158 
159  do
160  {
161  if ( std::abs( (*this)(i) - other_vector(i) ) > threshold *
162  std::max(std::abs((*this)(i)), std::abs(other_vector(i))))
163  first_different_i = i;
164  else
165  i++;
166  }
167  while (first_different_i==std::numeric_limits<int>::max()
168  && i<last_local_index());
169 
170  // Find the correct first differing index in parallel
171  this->comm().min(first_different_i);
172 
173  if (first_different_i == std::numeric_limits<int>::max())
174  return -1;
175 
176  return first_different_i;
177 }
double abs(double a)
virtual numeric_index_type last_local_index() const =0
virtual bool initialized() const
long double max(long double a, double b)
libmesh_assert(j)
dof_id_type numeric_index_type
Definition: id_types.h:92
const Parallel::Communicator & comm() const
virtual numeric_index_type first_local_index() const =0
template<typename T >
numeric_index_type libMesh::DistributedVector< T >::local_size ( ) const
inlinevirtual
template<typename T >
void libMesh::DistributedVector< T >::localize ( std::vector< T > &  v_local) const
virtual

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

Implements libMesh::NumericVector< T >.

Definition at line 553 of file distributed_vector.C.

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

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

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

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

Implements libMesh::NumericVector< T >.

Definition at line 376 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, libMesh::initialized(), and libMesh::libmesh_assert().

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

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 408 of file distributed_vector.C.

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

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

Fill in the local std::vector "v_local" with the global indices given in "indices". See numeric_vector.h for more details.

Implements libMesh::NumericVector< T >.

Definition at line 422 of file distributed_vector.C.

References libMesh::n_processors(), and libMesh::processor_id().

424 {
425  // Resize v_local so there is enough room to hold all the local values.
426  v_local.resize(indices.size());
427 
428  // We need to know who has the values we want, so get everyone's _local_size
429  std::vector<numeric_index_type> local_sizes;
430  this->comm().allgather (_local_size, local_sizes);
431 
432  // Make a vector of partial sums of local sizes
433  std::vector<numeric_index_type> local_size_sums(this->n_processors());
434  local_size_sums[0] = local_sizes[0];
435  for (numeric_index_type i=1; i<local_sizes.size(); i++)
436  local_size_sums[i] = local_size_sums[i-1] + local_sizes[i];
437 
438  // We now fill in 'requested_ids' based on the indices. Also keep
439  // track of the local index (in the indices vector) for each of
440  // these, since we need that when unpacking.
441  std::vector<std::vector<numeric_index_type> >
442  requested_ids(this->n_processors()),
443  local_requested_ids(this->n_processors());
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  std::iterator_traits<iter_t>::difference_type on_proc
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  // First trade indices, then trade values with all other processors.
466  for (processor_id_type pid=0; pid<this->n_processors(); pid++)
467  {
468  // Another data structure which holds the ids that other processors request us to fill.
469  std::vector<numeric_index_type> requests_for_me_to_fill;
470 
471  // Trade my requests with processor procup and procdown
472  const processor_id_type procup = (this->processor_id() + pid) % this->n_processors();
473  const processor_id_type procdown = (this->n_processors() + this->processor_id() - pid) % this->n_processors();
474 
475  this->comm().send_receive (procup, requested_ids[procup],
476  procdown, requests_for_me_to_fill);
477 
478  // The first send/receive we did was for indices, the second one will be
479  // for corresponding floating point values, so create storage for that now...
480  std::vector<T> values_to_send(requests_for_me_to_fill.size());
481 
482  for (std::size_t i=0; i<requests_for_me_to_fill.size(); i++)
483  {
484  // The index of the requested value
485  const numeric_index_type requested_index = requests_for_me_to_fill[i];
486 
487  // Transform into local numbering, and get requested value.
488  values_to_send[i] = _values[requested_index - _first_local_index];
489  }
490 
491  // Send the values we were supposed to send, receive the values
492  // we were supposed to receive.
493  std::vector<T> values_to_receive(requested_ids[procup].size());
494  this->comm().send_receive (procdown, values_to_send,
495  procup, values_to_receive);
496 
497  // Error checking: make sure that requested_ids[procup] is same
498  // length as values_to_receive.
499  if (requested_ids[procup].size() != values_to_receive.size())
500  libmesh_error_msg("Requested " << requested_ids[procup].size() << " values, but received " << values_to_receive.size() << " values.");
501 
502  // Now write the received values to the appropriate place(s) in v_local
503  for (std::size_t i=0; i<values_to_receive.size(); i++)
504  {
505  // Get the index in v_local where this value needs to be inserted.
506  numeric_index_type local_requested_index = local_requested_ids[procup][i];
507 
508  // Actually set the value in v_local
509  v_local[local_requested_index] = values_to_receive[i];
510  }
511  }
512 }
processor_id_type n_processors() const
uint8_t processor_id_type
Definition: id_types.h:99
void send_receive(const unsigned int dest_processor_id, const T1 &send, const unsigned int source_processor_id, T2 &recv, const MessageTag &send_tag=no_tag, const MessageTag &recv_tag=any_tag) const
dof_id_type numeric_index_type
Definition: id_types.h:92
numeric_index_type _local_size
const Parallel::Communicator & comm() const
virtual numeric_index_type size() const libmesh_override
numeric_index_type _first_local_index
processor_id_type processor_id() const
void allgather(const T &send, std::vector< T > &recv) const
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 
)
virtual

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

Implements libMesh::NumericVector< T >.

Definition at line 517 of file distributed_vector.C.

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

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

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 574 of file distributed_vector.C.

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

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

576 {
577  // This function must be run on all processors at once
578  parallel_object_only();
579 
580  libmesh_assert (this->initialized());
581  libmesh_assert_equal_to (_values.size(), _local_size);
582  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
583 
584  v_local = this->_values;
585 
586  this->comm().gather (pid, v_local);
587 
588 #ifndef LIBMESH_HAVE_MPI
589  libmesh_assert_equal_to (local_size(), size());
590 #endif
591 }
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type local_size() const libmesh_override
numeric_index_type _local_size
void gather(const unsigned int root_id, const T &send, std::vector< T > &recv) const
const Parallel::Communicator & comm() const
virtual numeric_index_type size() const libmesh_override
numeric_index_type _first_local_index
template<typename T >
Real libMesh::DistributedVector< T >::max ( ) const
inlinevirtual
Returns
the maximum element in the vector. In case of complex numbers, this returns the maximum Real part.

Implements libMesh::NumericVector< T >.

Definition at line 813 of file distributed_vector.h.

References libMesh::DistributedVector< T >::_first_local_index, libMesh::DistributedVector< T >::_last_local_index, libMesh::DistributedVector< T >::_local_size, libMesh::DistributedVector< T >::_values, libMesh::ParallelObject::comm(), libMesh::NumericVector< T >::initialized(), libMesh::libmesh_assert(), libMesh::libmesh_real(), std::max(), libMesh::Parallel::Communicator::max(), and libMesh::Real.

814 {
815  // This function must be run on all processors at once
816  parallel_object_only();
817 
818  libmesh_assert (this->initialized());
819  libmesh_assert_equal_to (_values.size(), _local_size);
820  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
821 
822  Real local_max = _values.size() ?
823  libmesh_real(_values[0]) : -std::numeric_limits<Real>::max();
824  for (numeric_index_type i = 1; i < _values.size(); ++i)
825  local_max = std::max(libmesh_real(_values[i]), local_max);
826 
827  this->comm().max(local_max);
828 
829  return local_max;
830 }
T libmesh_real(T a)
numeric_index_type _last_local_index
virtual bool initialized() const
long double max(long double a, double b)
libmesh_assert(j)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual Real max() const libmesh_override
numeric_index_type _local_size
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Parallel::Communicator & comm() const
numeric_index_type _first_local_index
template<typename T >
Real libMesh::DistributedVector< T >::min ( ) const
inlinevirtual
Returns
the minimum element in the vector. In case of complex numbers, this returns the minimum Real part.

Implements libMesh::NumericVector< T >.

Definition at line 790 of file distributed_vector.h.

References libMesh::DistributedVector< T >::_first_local_index, libMesh::DistributedVector< T >::_last_local_index, libMesh::DistributedVector< T >::_local_size, libMesh::DistributedVector< T >::_values, libMesh::ParallelObject::comm(), libMesh::NumericVector< T >::initialized(), libMesh::libmesh_assert(), libMesh::libmesh_real(), std::max(), std::min(), libMesh::Parallel::Communicator::min(), and libMesh::Real.

791 {
792  // This function must be run on all processors at once
793  parallel_object_only();
794 
795  libmesh_assert (this->initialized());
796  libmesh_assert_equal_to (_values.size(), _local_size);
797  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
798 
799  Real local_min = _values.size() ?
800  libmesh_real(_values[0]) : std::numeric_limits<Real>::max();
801  for (numeric_index_type i = 1; i < _values.size(); ++i)
802  local_min = std::min(libmesh_real(_values[i]), local_min);
803 
804  this->comm().min(local_min);
805 
806  return local_min;
807 }
T libmesh_real(T a)
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual Real max() const libmesh_override
numeric_index_type _local_size
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Parallel::Communicator & comm() const
long double min(long double a, double b)
numeric_index_type _first_local_index
static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited
processor_id_type libMesh::ParallelObject::n_processors ( ) const
inlineinherited
Returns
the number of processors in the group.

Definition at line 93 of file parallel_object.h.

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

Referenced by libMesh::ParmetisPartitioner::_do_repartition(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::DistributedMesh::add_elem(), libMesh::DistributedMesh::add_node(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assembly(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::MeshCommunication::broadcast(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::DistributedMesh::clear(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DistributedMesh::DistributedMesh(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshCommunication::gather(), libMesh::MeshCommunication::gather_neighboring_elements(), libMesh::MeshBase::get_info(), libMesh::EquationSystems::get_solution(), libMesh::DistributedVector< T >::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::DistributedMesh::insert_elem(), 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::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::Partitioner::partition(), libMesh::MeshBase::partition(), libMesh::Partitioner::partition_unpartitioned_elements(), libMesh::PetscLinearSolver< T >::PetscLinearSolver(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::SparseMatrix< T >::print(), libMesh::MeshTools::processor_bounding_box(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::MeshCommunication::redistribute(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::Partitioner::repartition(), libMesh::MeshCommunication::send_coarse_ghosts(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::Parallel::Sort< KeyType, IdxType >::sort(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

94  { return cast_int<processor_id_type>(_communicator.size()); }
unsigned int size() const
Definition: parallel.h:679
const Parallel::Communicator & _communicator
template<typename T >
T libMesh::DistributedVector< T >::operator() ( const numeric_index_type  i) const
inlinevirtual

Access components, returns U(i).

Implements libMesh::NumericVector< T >.

Definition at line 745 of file distributed_vector.h.

References libMesh::DistributedVector< T >::_first_local_index, libMesh::DistributedVector< T >::_last_local_index, libMesh::DistributedVector< T >::_local_size, libMesh::DistributedVector< T >::_values, libMesh::DistributedVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::DistributedVector< T >::last_local_index(), and libMesh::libmesh_assert().

746 {
747  libmesh_assert (this->initialized());
748  libmesh_assert_equal_to (_values.size(), _local_size);
749  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
750  libmesh_assert ( ((i >= first_local_index()) &&
751  (i < last_local_index())) );
752 
753  return _values[i - _first_local_index];
754 }
virtual numeric_index_type last_local_index() const libmesh_override
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type first_local_index() const libmesh_override
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T>
NumericVector<T>& libMesh::NumericVector< T >::operator*= ( const T  a)
inlineinherited

Multiplication operator. Equivalent to U.scale(a)

Definition at line 381 of file numeric_vector.h.

381 { this->scale(a); return *this; }
virtual void scale(const T factor)=0
template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator+= ( const NumericVector< T > &  v)
virtual

Addition operator. Fast equivalent to U.add(1, V).

Implements libMesh::NumericVector< T >.

Definition at line 134 of file distributed_vector.C.

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

135 {
136  libmesh_assert (this->closed());
137  libmesh_assert (this->initialized());
138  libmesh_assert_equal_to (_values.size(), _local_size);
139  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
140 
141  add(1., v);
142 
143  return *this;
144 }
virtual bool closed() const
numeric_index_type _last_local_index
virtual bool initialized() const
virtual void add(const numeric_index_type i, const T value) libmesh_override
libmesh_assert(j)
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator-= ( const NumericVector< T > &  v)
virtual

Subtraction operator. Fast equivalent to U.add(-1, V).

Implements libMesh::NumericVector< T >.

Definition at line 149 of file distributed_vector.C.

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

150 {
151  libmesh_assert (this->closed());
152  libmesh_assert (this->initialized());
153  libmesh_assert_equal_to (_values.size(), _local_size);
154  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
155 
156  add(-1., v);
157 
158  return *this;
159 }
virtual bool closed() const
numeric_index_type _last_local_index
virtual bool initialized() const
virtual void add(const numeric_index_type i, const T value) libmesh_override
libmesh_assert(j)
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator/= ( NumericVector< T > &  v)
virtual

Pointwise Division operator. ie divide every entry in this vector by the entry in v

Implements libMesh::NumericVector< T >.

Definition at line 164 of file distributed_vector.C.

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

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

Division operator. Equivalent to U.scale(1./a)

Definition at line 387 of file numeric_vector.h.

Referenced by libMesh::NumericVector< Number >::operator/=().

387 { this->scale(1./a); return *this; }
virtual void scale(const T factor)=0
template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator= ( const T  s)
virtual

$U(0-N) = s$: fill all components.

Implements libMesh::NumericVector< T >.

Definition at line 301 of file distributed_vector.C.

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

302 {
303  libmesh_assert (this->initialized());
304  libmesh_assert_equal_to (_values.size(), _local_size);
305  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
306 
307  for (std::size_t i=0; i<local_size(); i++)
308  _values[i] = s;
309 
310  return *this;
311 }
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type local_size() const libmesh_override
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator= ( const NumericVector< T > &  v)
virtual

$U = V$: copy all components.

Implements libMesh::NumericVector< T >.

Definition at line 317 of file distributed_vector.C.

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

$U = V$: copy all components.

Definition at line 331 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().

332 {
333  this->_is_initialized = v._is_initialized;
334  this->_is_closed = v._is_closed;
335 
336  _global_size = v._global_size;
337  _local_size = v._local_size;
338  _first_local_index = v._first_local_index;
339  _last_local_index = v._last_local_index;
340 
341  if (v.local_size() == this->local_size())
342  _values = v._values;
343 
344  else
345  libmesh_error_msg("v.local_size() = " << v.local_size() << " must be equal to this->local_size() = " << this->local_size());
346 
347  return *this;
348 }
numeric_index_type _last_local_index
numeric_index_type _global_size
virtual numeric_index_type local_size() const libmesh_override
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T >
NumericVector< T > & libMesh::DistributedVector< T >::operator= ( const std::vector< T > &  v)
virtual

$U = V$: copy all components.

Implements libMesh::NumericVector< T >.

Definition at line 354 of file distributed_vector.C.

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

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

Computes the pointwise (i.e. component-wise) product of vec1 and vec2 and stores the result in *this.

Implements libMesh::NumericVector< T >.

Definition at line 596 of file distributed_vector.C.

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

600 {
601  libmesh_not_implemented();
602 }
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 921 of file numeric_vector.h.

References libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::NumericVector< T >::last_local_index(), libMesh::libmesh_assert(), libMesh::NumericVector< T >::local_size(), and libMesh::NumericVector< T >::size().

Referenced by libMesh::NumericVector< Number >::operator/=().

922 {
923  libmesh_assert (this->initialized());
924  os << "Size\tglobal = " << this->size()
925  << "\t\tlocal = " << this->local_size() << std::endl;
926 
927  os << "#\tValue" << std::endl;
928  for (numeric_index_type i=this->first_local_index(); i<this->last_local_index(); i++)
929  os << i << "\t" << (*this)(i) << std::endl;
930 }
virtual numeric_index_type size() const =0
virtual numeric_index_type last_local_index() const =0
virtual bool initialized() const
libmesh_assert(j)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type local_size() const =0
virtual numeric_index_type first_local_index() const =0
template<>
void libMesh::NumericVector< Complex >::print ( std::ostream &  os) const
inlineinherited

Definition at line 903 of file numeric_vector.h.

References libMesh::NumericVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::NumericVector< T >::last_local_index(), libMesh::libmesh_assert(), libMesh::NumericVector< T >::local_size(), and libMesh::NumericVector< T >::size().

904 {
905  libmesh_assert (this->initialized());
906  os << "Size\tglobal = " << this->size()
907  << "\t\tlocal = " << this->local_size() << std::endl;
908 
909  // std::complex<>::operator<<() is defined, but use this form
910  os << "#\tReal part\t\tImaginary part" << std::endl;
911  for (numeric_index_type i=this->first_local_index(); i<this->last_local_index(); i++)
912  os << i << "\t"
913  << (*this)(i).real() << "\t\t"
914  << (*this)(i).imag() << std::endl;
915 }
virtual numeric_index_type size() const =0
virtual numeric_index_type last_local_index() const =0
virtual bool initialized() const
libmesh_assert(j)
dof_id_type numeric_index_type
Definition: id_types.h:92
virtual numeric_index_type local_size() const =0
virtual numeric_index_type first_local_index() const =0
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 958 of file numeric_vector.h.

References libMesh::NumericVector< T >::initialized(), libMesh::libmesh_assert(), libMesh::NumericVector< T >::localize(), libMesh::ParallelObject::processor_id(), and libMesh::NumericVector< T >::size().

Referenced by libMesh::NumericVector< Number >::operator/=().

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

Definition at line 936 of file numeric_vector.h.

References libMesh::NumericVector< T >::initialized(), libMesh::libmesh_assert(), libMesh::NumericVector< T >::localize(), libMesh::ParallelObject::processor_id(), and libMesh::NumericVector< T >::size().

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

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

Definition at line 88 of file reference_counter.C.

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

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

89 {
91 }
static std::string get_info()
template<typename T>
virtual void libMesh::NumericVector< T >::print_matlab ( const std::string &  = "") const
inlinevirtualinherited

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

Reimplemented in libMesh::PetscVector< T >.

Definition at line 677 of file numeric_vector.h.

678  {
679  libmesh_not_implemented();
680  }
processor_id_type libMesh::ParallelObject::processor_id ( ) const
inlineinherited
Returns
the rank of this processor in the group.

Definition at line 99 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::DistributedMesh::add_elem(), libMesh::BoundaryInfo::add_elements(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::MeshRefinement::add_node(), libMesh::MeshTools::Modification::all_tri(), libMesh::FEMSystem::assembly(), libMesh::ParmetisPartitioner::assign_partitioning(), libMesh::MeshCommunication::broadcast(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::ParmetisPartitioner::build_graph(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::DofMap::build_sparsity(), 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::EnsightIO::EnsightIO(), libMesh::MeshFunction::find_element(), libMesh::MeshFunction::find_elements(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::MeshCommunication::gather(), libMesh::MeshCommunication::gather_neighboring_elements(), 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::EquationSystems::get_solution(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SparsityPattern::Build::handle_vi_vj(), libMesh::DistributedVector< T >::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ParmetisPartitioner::initialize(), 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::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), 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::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::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::SparsityPattern::Build::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::MetisPartitioner::partition_range(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::SparseMatrix< T >::print(), libMesh::NumericVector< T >::print_global(), 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::ExodusII_IO_Helper::read_elem_num_map(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::MeshCommunication::redistribute(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::MeshCommunication::send_coarse_ghosts(), libMesh::Partitioner::set_node_processor_ids(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::MeshTools::total_weight(), libMesh::Parallel::Packing< Node * >::unpack(), libMesh::Parallel::Packing< Elem * >::unpack(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::MeshTools::weight(), libMesh::NameBasedIO::write(), libMesh::XdrIO::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::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::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), and libMesh::ExodusII_IO_Helper::write_timestep().

100  { return cast_int<processor_id_type>(_communicator.rank()); }
const Parallel::Communicator & _communicator
unsigned int rank() const
Definition: parallel.h:677
template<typename T >
void libMesh::DistributedVector< T >::reciprocal ( )
virtual

Replace each entry v_i of this vector by its reciprocal, 1/v_i.

Implements libMesh::NumericVector< T >.

Definition at line 182 of file distributed_vector.C.

183 {
184  for (numeric_index_type i=0; i<local_size(); i++)
185  {
186  // Don't divide by zero
187  libmesh_assert_not_equal_to (_values[i], T(0));
188 
189  _values[i] = 1. / _values[i];
190  }
191 }
virtual numeric_index_type local_size() const libmesh_override
dof_id_type numeric_index_type
Definition: id_types.h:92
template<typename T >
void libMesh::DistributedVector< T >::scale ( const T  factor)
virtual

Scale each element of the vector by the given factor.

Implements libMesh::NumericVector< T >.

Definition at line 248 of file distributed_vector.C.

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

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

249 {
250  libmesh_assert (this->initialized());
251  libmesh_assert_equal_to (_values.size(), _local_size);
252  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
253 
254  for (std::size_t i=0; i<local_size(); i++)
255  _values[i] *= factor;
256 }
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type local_size() const libmesh_override
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T >
void libMesh::DistributedVector< T >::set ( const numeric_index_type  i,
const T  value 
)
inlinevirtual

v(i) = value

Implements libMesh::NumericVector< T >.

Definition at line 760 of file distributed_vector.h.

References libMesh::DistributedVector< T >::_first_local_index, libMesh::DistributedVector< T >::_last_local_index, libMesh::DistributedVector< T >::_local_size, libMesh::DistributedVector< T >::_values, libMesh::DistributedVector< T >::first_local_index(), libMesh::NumericVector< T >::initialized(), libMesh::libmesh_assert(), libMesh::DistributedVector< T >::local_size(), and libMesh::DistributedVector< T >::size().

761 {
762  libmesh_assert (this->initialized());
763  libmesh_assert_equal_to (_values.size(), _local_size);
764  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
765  libmesh_assert_less (i, size());
766  libmesh_assert_less (i-first_local_index(), local_size());
767 
768  _values[i - _first_local_index] = value;
769 }
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
virtual numeric_index_type local_size() const libmesh_override
virtual numeric_index_type first_local_index() const libmesh_override
numeric_index_type _local_size
virtual numeric_index_type size() const libmesh_override
numeric_index_type _first_local_index
template<typename T >
numeric_index_type libMesh::DistributedVector< T >::size ( ) const
inlinevirtual
Returns
dimension of the vector. This function was formerly called n(), but was renamed to get the DistributedVector class closer to the C++ standard library's std::vector container.

Implements libMesh::NumericVector< T >.

Definition at line 693 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, libMesh::NumericVector< T >::initialized(), and libMesh::libmesh_assert().

Referenced by libMesh::DistributedVector< T >::add(), libMesh::DistributedVector< T >::operator/=(), and libMesh::DistributedVector< T >::set().

694 {
695  libmesh_assert (this->initialized());
696  libmesh_assert_equal_to (_values.size(), _local_size);
697  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
698 
699  return _global_size;
700 }
numeric_index_type _last_local_index
virtual bool initialized() const
numeric_index_type _global_size
libmesh_assert(j)
numeric_index_type _local_size
numeric_index_type _first_local_index
template<class T >
Real libMesh::NumericVector< T >::subset_l1_norm ( const std::set< numeric_index_type > &  indices) const
virtualinherited
Returns
the $l_1$-norm of the vector, i.e. the sum of the absolute values for the specified entries in the vector.

Note that the indices must necessary live on this processor.

Definition at line 323 of file numeric_vector.C.

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

Referenced by libMesh::NumericVector< Number >::closed(), and libMesh::System::discrete_var_norm().

324 {
325  const NumericVector<T> & v = *this;
326 
327  std::set<numeric_index_type>::const_iterator it = indices.begin();
328  const std::set<numeric_index_type>::const_iterator it_end = indices.end();
329 
330  Real norm = 0;
331 
332  for(; it!=it_end; ++it)
333  norm += std::abs(v(*it));
334 
335  this->comm().sum(norm);
336 
337  return norm;
338 }
double abs(double a)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Parallel::Communicator & comm() const
template<class T >
Real libMesh::NumericVector< T >::subset_l2_norm ( const std::set< numeric_index_type > &  indices) const
virtualinherited
Returns
the $l_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 that the indices must necessary live on this processor.

Definition at line 341 of file numeric_vector.C.

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

Referenced by libMesh::NumericVector< Number >::closed(), and libMesh::System::discrete_var_norm().

342 {
343  const NumericVector<T> & v = *this;
344 
345  std::set<numeric_index_type>::const_iterator it = indices.begin();
346  const std::set<numeric_index_type>::const_iterator it_end = indices.end();
347 
348  Real norm = 0;
349 
350  for(; it!=it_end; ++it)
351  norm += TensorTools::norm_sq(v(*it));
352 
353  this->comm().sum(norm);
354 
355  return std::sqrt(norm);
356 }
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Parallel::Communicator & comm() const
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 $l_\infty$-norm of a vector.

Note that the indices must necessary live on this processor.

Definition at line 359 of file numeric_vector.C.

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

Referenced by libMesh::NumericVector< Number >::closed(), and libMesh::System::discrete_var_norm().

360 {
361  const NumericVector<T> & v = *this;
362 
363  std::set<numeric_index_type>::const_iterator it = indices.begin();
364  const std::set<numeric_index_type>::const_iterator it_end = indices.end();
365 
366  Real norm = 0;
367 
368  for(; it!=it_end; ++it)
369  {
370  Real value = std::abs(v(*it));
371  if(value > norm)
372  norm = value;
373  }
374 
375  this->comm().max(norm);
376 
377  return norm;
378 }
double abs(double a)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Parallel::Communicator & comm() const
template<typename T >
T libMesh::DistributedVector< T >::sum ( ) const
virtual
Returns
the sum of all values in the vector

Implements libMesh::NumericVector< T >.

Definition at line 42 of file distributed_vector.C.

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

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

Swaps the vector data and metadata

Reimplemented from libMesh::NumericVector< T >.

Definition at line 835 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().

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

836 {
837  DistributedVector<T> & v = cast_ref<DistributedVector<T> &>(other);
838 
839  std::swap(_global_size, v._global_size);
840  std::swap(_local_size, v._local_size);
841  std::swap(_first_local_index, v._first_local_index);
842  std::swap(_last_local_index, v._last_local_index);
843 
844  // This should be O(1) with any reasonable STL implementation
845  std::swap(_values, v._values);
846 }
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
template<typename T>
ParallelType& libMesh::NumericVector< T >::type ( )
inlineinherited
Returns
the type (SERIAL, PARALLEL, GHOSTED) of the vector.

Definition at line 142 of file numeric_vector.h.

142 { return _type; }
template<typename T >
void libMesh::DistributedVector< T >::zero ( )
inlinevirtual

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

Implements libMesh::NumericVector< T >.

Definition at line 655 of file distributed_vector.h.

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

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

656 {
657  libmesh_assert (this->initialized());
658  libmesh_assert_equal_to (_values.size(), _local_size);
659  libmesh_assert_equal_to ((_last_local_index - _first_local_index), _local_size);
660 
661  std::fill (_values.begin(),
662  _values.end(),
663  0.);
664 }
numeric_index_type _last_local_index
virtual bool initialized() const
libmesh_assert(j)
numeric_index_type _local_size
numeric_index_type _first_local_index
template<typename T >
UniquePtr< NumericVector< T > > libMesh::DistributedVector< T >::zero_clone ( ) const
inlinevirtual

Creates a vector which has the same type, size and partitioning as this vector, but whose data is all zero. Returns it in an UniquePtr.

Implements libMesh::NumericVector< T >.

Definition at line 670 of file distributed_vector.h.

References libMesh::ParallelObject::comm(), and libMesh::NumericVector< T >::init().

671 {
672  NumericVector<T> * cloned_vector = new DistributedVector<T>(this->comm());
673  cloned_vector->init(*this);
674  return UniquePtr<NumericVector<T> >(cloned_vector);
675 }
const Parallel::Communicator & comm() const

Member Data Documentation

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited
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 134 of file reference_counter.h.

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

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 128 of file reference_counter.h.

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 123 of file reference_counter.h.

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


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